HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
Loading...
Searching...
No Matches
sddr-magnetostatics.hpp
Go to the documentation of this file.
1// Implementation of the serendipity discrete de Rham sequence for the magnetostatic problem.
2//
3// Author: Jerome Droniou (jerome.droniou@monash.edu)
4//
5
6/*
7 * The DDR scheme is described in
8 *
9 * An arbitrary-order method for magnetostatics on polyhedral meshes based on a discrete de Rham
10 sequence.
11 * D. A. Di Pietro and J. Droniou, 31p, 2020. url: https://arxiv.org/abs/2005.06890.
12 *
13 * and its serendipity version in
14 *
15 * TBC ////
16 *
17 * If you use this code in a scientific publication, please mention the above articles.
18 *
19 */
20
21#ifndef S_MAGNETOSTATICS_HPP
22#define S_MAGNETOSTATICS_HPP
23
24#include <iostream>
25
26#include <boost/math/constants/constants.hpp>
27
28#include <Eigen/Sparse>
29#include <unsupported/Eigen/SparseExtra>
30
31#include <mesh.hpp>
32#include <mesh_builder.hpp>
34#include <linearsolver.hpp>
35
36#include <sxdiv.hpp>
37#include <sxcurl.hpp>
38
44namespace HArDCore3D
45{
46
53 struct Magnetostatics
54 {
55 typedef Eigen::SparseMatrix<double> SystemMatrixType;
56
57 typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> ForcingTermType;
58 typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> SolutionPotentialType;
59 typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> SolutionCurlType;
61
64 const DDRCore & ddrcore,
65 bool use_threads,
66 std::ostream & output = std::cout
67 );
68
71 const ForcingTermType & f,
72 const PermeabilityType & mu,
73 const SolutionPotentialType & u
74 );
75
76
78 inline size_t dimensionSpace() const
79 {
80 return m_sxcurl.dimension() + m_sxdiv.dimension();
81 }
82
84 inline size_t nbSCDOFs() const
85 {
86 return m_nloc_sc_H.sum();
87 }
88
90 inline size_t sizeSystem() const
91 {
92 return dimensionSpace() - nbSCDOFs();
93 }
94
96 inline const SXCurl & sxCurl() const
97 {
98 return m_sxcurl;
99 }
100
102 inline const SXDiv & sxDiv() const
103 {
104 return m_sxdiv;
105 }
106
108 inline const SystemMatrixType & systemMatrix() const {
109 return m_A;
110 }
111
114 return m_A;
115 }
116
118 inline const Eigen::VectorXd & systemVector() const {
119 return m_b;
120 }
121
123 inline Eigen::VectorXd & systemVector() {
124 return m_b;
125 }
126
128 inline const double & stabilizationParameter() const {
129 return m_stab_par;
130 }
131
133 inline double & stabilizationParameter() {
134 return m_stab_par;
135 }
136
138 inline const SystemMatrixType & scMatrix() const {
139 return m_sc_A;
140 }
141
143 inline Eigen::VectorXd & scVector() {
144 return m_sc_b;
145 }
146
148 std::vector<double> computeNorms(
149 const std::vector<Eigen::VectorXd> & list_dofs
150 ) const;
151
152 private:
154 std::pair<Eigen::MatrixXd, Eigen::VectorXd>
155 _compute_local_contribution(
156 size_t iT,
157 const ForcingTermType & f,
158 const PermeabilityType & mu
159 );
160
162 LocalStaticCondensation _compute_static_condensation(const size_t & iT) const;
163
165 void _assemble_local_contribution(
166 size_t iT,
167 const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT,
168 std::list<Eigen::Triplet<double> > & A1,
169 Eigen::VectorXd & b1,
170 std::list<Eigen::Triplet<double> > & A2,
171 Eigen::VectorXd & b2
172 );
173
174 const DDRCore & m_ddrcore;
175 bool m_use_threads;
176 std::ostream & m_output;
177 SerendipityProblem m_ser_pro;
178 SXCurl m_sxcurl;
179 SXDiv m_sxdiv;
180 Eigen::VectorXi m_nloc_sc_H; // Number of statically condensed DOF in each cell
181 SystemMatrixType m_A; // Matrix and RHS of statically condensed system
182 Eigen::VectorXd m_b;
183 SystemMatrixType m_sc_A; // Static condensation operator and RHS (to recover statically condensed DOFs)
184 Eigen::VectorXd m_sc_b;
185 double m_stab_par;
186 };
187
188 //------------------------------------------------------------------------------
189 // Exact solutions
190 //------------------------------------------------------------------------------
191
192 static const double PI = boost::math::constants::pi<double>();
193 using std::sin;
194 using std::cos;
195
196 //------------------------------------------------------------------------------
197
199 constant_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
200 return Eigen::Vector3d(1., 2., 3.);
201 };
202
204 constant_sigma = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
205 return Eigen::Vector3d::Zero();
206 };
207
209 constant_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
210 return Eigen::Vector3d::Zero();
211 };
212
215
216 //------------------------------------------------------------------------------
217
219 linear_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
220 return Eigen::Vector3d(
221 x(0) + x(1) + x(2),
222 2.*(x(0)-x(1)+x(2)),
223 -x(0) - x(1) + x(2)
224 );
225 };
226
228 linear_sigma = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
229 return Eigen::Vector3d(-3., 2., 1.);
230 };
231
233 linear_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
234 return Eigen::Vector3d::Zero();
235 };
236
239
240 //------------------------------------------------------------------------------
241
243 trigonometric_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
244 return Eigen::Vector3d(
245 cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
246 -2. * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
247 sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
248 );
249 };
250
252 trigonometric_sigma = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
253 return 3. * PI * Eigen::Vector3d(
254 sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2)),
255 0.,
256 -cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
257 );
258 };
259
261 trigonometric_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
262 return 3. * std::pow(PI, 2) * Eigen::Vector3d(
263 cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
264 -2. * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
265 sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
266 );
267 };
268
271
272 //------------------------------------------------------------------------------
273
274
276 variable_permeability_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
277 return Eigen::Vector3d(
278 cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
279 -2. * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
280 sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
281 );
282 };
283
285 variable_permeability_sigma = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
286 return 3 * PI / (1. + x(0) + x(1) + x(2))
287 * Eigen::Vector3d(
288 sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2)),
289 0.,
290 -cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
291 );
292 };
293
295 variable_permeability_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
296 double a = 3 * PI / std::pow(1. + x(0) + x(1) + x(2), 2);
297 double b = 3 * std::pow(PI, 2) / (1. + x(0) + x(1) + x(2));
298 double dy_sigmax = -a * sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2))
299 - b * sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2));
300 double dz_sigmax = -a * sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2))
301 - b * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2));
302 ;
303 double dx_sigmaz = a * cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
304 + b * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2));
305
306 double dy_sigmaz = a * cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
307 + b * cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2));
308
309 return Eigen::Vector3d(
310 dy_sigmaz,
312 -dy_sigmax
313 );
314 };
315
318 [](const Cell & T, const Eigen::Vector3d & x) -> double { return 1. + x(0) + x(1) + x(2); },
319 [](const Cell & T) -> size_t { return 1; }
320 );
321
322} // end of namespace HArDCore3D
323
324#endif
Construct all polynomial spaces for the DDR sequence.
Definition ddrcore.hpp:62
Discrete Serendipity Hcurl space: local operators, L2 product and global interpolator.
Definition sxcurl.hpp:21
Discrete Serendipity Hdiv space: local operators, L2 product and global interpolator.
Definition sxdiv.hpp:18
Construct all polynomial spaces for the DDR sequence.
Definition serendipity_problem.hpp:40
@ Matrix
Definition basis.hpp:67
size_t dimension() const
Returns the dimension of the global space (all DOFs for all geometric entities)
Definition localdofspace.hpp:80
size_t dimension() const
Returns the dimension of the global space (all DOFs for all geometric entities)
Definition variabledofspace.hpp:158
static const double PI
Definition ddr-magnetostatics.hpp:187
static Magnetostatics::SolutionPotentialType linear_u
Definition ddr-magnetostatics.hpp:214
static Magnetostatics::SolutionCurlType linear_sigma
Definition ddr-magnetostatics.hpp:223
static Magnetostatics::SolutionCurlType trigonometric_sigma
Definition ddr-magnetostatics.hpp:247
Eigen::SparseMatrix< double > SystemMatrixType
Definition sddr-magnetostatics.hpp:55
static Magnetostatics::PermeabilityType constant_mu
Definition ddr-magnetostatics.hpp:209
static Magnetostatics::PermeabilityType linear_mu
Definition ddr-magnetostatics.hpp:233
size_t nbSCDOFs() const
Returns the number of statically condensed DOFs (here, the cell magnetic field DOFs)
Definition ddr-magnetostatics.hpp:80
static Magnetostatics::SolutionPotentialType constant_u
Definition ddr-magnetostatics.hpp:194
static Magnetostatics::SolutionPotentialType trigonometric_u
Definition ddr-magnetostatics.hpp:238
static Magnetostatics::PermeabilityType variable_permeability_mu
Definition ddr-magnetostatics.hpp:312
static Magnetostatics::ForcingTermType constant_f
Definition ddr-magnetostatics.hpp:204
const SXCurl & sxCurl() const
Returns the space SXCurl.
Definition sddr-magnetostatics.hpp:96
size_t dimensionSpace() const
Returns the dimension of the magnetic field + potential space.
Definition ddr-magnetostatics.hpp:74
static Magnetostatics::SolutionCurlType constant_sigma
Definition ddr-magnetostatics.hpp:199
static Magnetostatics::PermeabilityType trigonometric_mu
Definition ddr-magnetostatics.hpp:265
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> SolutionCurlType
Definition sddr-magnetostatics.hpp:59
IntegralWeight PermeabilityType
Definition sddr-magnetostatics.hpp:60
static Magnetostatics::SolutionCurlType variable_permeability_sigma
Definition ddr-magnetostatics.hpp:280
const SXDiv & sxDiv() const
Returns the space SXDiv.
Definition sddr-magnetostatics.hpp:102
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> SolutionPotentialType
Definition sddr-magnetostatics.hpp:58
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> ForcingTermType
Definition sddr-magnetostatics.hpp:57
static Magnetostatics::SolutionPotentialType variable_permeability_u
Definition ddr-magnetostatics.hpp:271
static Magnetostatics::ForcingTermType trigonometric_f
Definition ddr-magnetostatics.hpp:256
static Magnetostatics::ForcingTermType linear_f
Definition ddr-magnetostatics.hpp:228
static Magnetostatics::ForcingTermType variable_permeability_f
Definition ddr-magnetostatics.hpp:290
bool use_threads
Definition HHO_DiffAdvecReac.hpp:47
Definition ddr-magnetostatics.hpp:41
Structure for weights (scalar, at the moment) in integral.
Definition integralweight.hpp:36
Structure to store information for, and perform, local static condensation.
Definition local_static_condensation.hpp:25
Eigen::VectorXd & scVector()
Returns the static condensation rhs.
Definition sddr-magnetostatics.hpp:143
double & stabilizationParameter()
Returns the stabilization parameter.
Definition sddr-magnetostatics.hpp:133
const double & stabilizationParameter() const
Returns the stabilization parameter.
Definition sddr-magnetostatics.hpp:128
SystemMatrixType & systemMatrix()
Returns the linear system matrix.
Definition sddr-magnetostatics.hpp:113
size_t nbSCDOFs() const
Returns the number of statically condensed DOFs (here, the cell magnetic field DOFs)
Definition sddr-magnetostatics.hpp:84
const SystemMatrixType & systemMatrix() const
Returns the linear system matrix.
Definition sddr-magnetostatics.hpp:108
size_t dimensionSpace() const
Returns the dimension of the magnetic field + potential space.
Definition sddr-magnetostatics.hpp:78
Magnetostatics(const DDRCore &ddrcore, bool use_threads, std::ostream &output=std::cout)
Constructor.
std::vector< double > computeNorms(const std::vector< Eigen::VectorXd > &list_dofs) const
Compute the discrete Hcurl \times Hdiv norm of a family of vectors representing the dofs.
size_t sizeSystem() const
Returns the size of the statically condensed system.
Definition sddr-magnetostatics.hpp:90
const Eigen::VectorXd & systemVector() const
Returns the linear system right-hand side vector.
Definition sddr-magnetostatics.hpp:118
Eigen::VectorXd & systemVector()
Returns the linear system right-hand side vector.
Definition sddr-magnetostatics.hpp:123
void assembleLinearSystem(const ForcingTermType &f, const PermeabilityType &mu, const SolutionPotentialType &u)
Assemble the global system
const SystemMatrixType & scMatrix() const
Returns the static condensation recovery operator.
Definition sddr-magnetostatics.hpp:138