HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
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 
35 #include <sxdiv.hpp>
36 #include <sxcurl.hpp>
37 
43 namespace HArDCore3D
44 {
45 
52  struct Magnetostatics
53  {
54  typedef Eigen::SparseMatrix<double> SystemMatrixType;
55 
56  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> ForcingTermType;
57  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> SolutionPotentialType;
58  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> SolutionCurlType;
60 
63  const DDRCore & ddrcore,
64  bool use_threads,
65  std::ostream & output = std::cout
66  );
67 
70  const ForcingTermType & f,
71  const PermeabilityType & mu,
72  const SolutionPotentialType & u
73  );
74 
75 
77  inline size_t dimensionSpace() const
78  {
79  return m_sxcurl.dimension() + m_sxdiv.dimension();
80  }
81 
83  inline size_t nbSCDOFs() const
84  {
85  return m_nloc_sc_H.sum();
86  }
87 
89  inline size_t sizeSystem() const
90  {
91  return dimensionSpace() - nbSCDOFs();
92  }
93 
95  inline const SXCurl & sxCurl() const
96  {
97  return m_sxcurl;
98  }
99 
101  inline const SXDiv & sxDiv() const
102  {
103  return m_sxdiv;
104  }
105 
107  inline const SystemMatrixType & systemMatrix() const {
108  return m_A;
109  }
110 
113  return m_A;
114  }
115 
117  inline const Eigen::VectorXd & systemVector() const {
118  return m_b;
119  }
120 
122  inline Eigen::VectorXd & systemVector() {
123  return m_b;
124  }
125 
127  inline const double & stabilizationParameter() const {
128  return m_stab_par;
129  }
130 
132  inline double & stabilizationParameter() {
133  return m_stab_par;
134  }
135 
137  inline const SystemMatrixType & scMatrix() const {
138  return m_sc_A;
139  }
140 
142  inline Eigen::VectorXd & scVector() {
143  return m_sc_b;
144  }
145 
147  std::vector<double> computeNorms(
148  const std::vector<Eigen::VectorXd> & list_dofs
149  ) const;
150 
151  private:
153  std::pair<Eigen::MatrixXd, Eigen::VectorXd>
154  _compute_local_contribution(
155  size_t iT,
156  const ForcingTermType & f,
157  const PermeabilityType & mu
158  );
159 
161  LocalStaticCondensation _compute_static_condensation(const size_t & iT) const;
162 
164  void _assemble_local_contribution(
165  size_t iT,
166  const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT,
167  std::list<Eigen::Triplet<double> > & A1,
168  Eigen::VectorXd & b1,
169  std::list<Eigen::Triplet<double> > & A2,
170  Eigen::VectorXd & b2
171  );
172 
173  const DDRCore & m_ddrcore;
174  bool m_use_threads;
175  std::ostream & m_output;
176  SerendipityProblem m_ser_pro;
177  SXCurl m_sxcurl;
178  SXDiv m_sxdiv;
179  Eigen::VectorXd m_nloc_sc_H; // Number of statically condensed DOF in each cell
180  SystemMatrixType m_A; // Matrix and RHS of statically condensed system
181  Eigen::VectorXd m_b;
182  SystemMatrixType m_sc_A; // Static condensation operator and RHS (to recover statically condensed DOFs)
183  Eigen::VectorXd m_sc_b;
184  double m_stab_par;
185  };
186 
187  //------------------------------------------------------------------------------
188  // Exact solutions
189  //------------------------------------------------------------------------------
190 
191  static const double PI = boost::math::constants::pi<double>();
192  using std::sin;
193  using std::cos;
194 
195  //------------------------------------------------------------------------------
196 
198  constant_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
199  return Eigen::Vector3d(1., 2., 3.);
200  };
201 
203  constant_sigma = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
204  return Eigen::Vector3d::Zero();
205  };
206 
208  constant_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
209  return Eigen::Vector3d::Zero();
210  };
211 
214 
215  //------------------------------------------------------------------------------
216 
218  linear_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
219  return Eigen::Vector3d(
220  x(0) + x(1) + x(2),
221  2.*(x(0)-x(1)+x(2)),
222  -x(0) - x(1) + x(2)
223  );
224  };
225 
227  linear_sigma = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
228  return Eigen::Vector3d(-3., 2., 1.);
229  };
230 
232  linear_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
233  return Eigen::Vector3d::Zero();
234  };
235 
238 
239  //------------------------------------------------------------------------------
240 
242  trigonometric_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
243  return Eigen::Vector3d(
244  cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
245  -2. * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
246  sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
247  );
248  };
249 
251  trigonometric_sigma = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
252  return 3. * PI * Eigen::Vector3d(
253  sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2)),
254  0.,
255  -cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
256  );
257  };
258 
260  trigonometric_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
261  return 3. * std::pow(PI, 2) * Eigen::Vector3d(
262  cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
263  -2. * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
264  sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
265  );
266  };
267 
270 
271  //------------------------------------------------------------------------------
272 
273 
275  variable_permeability_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
276  return Eigen::Vector3d(
277  cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
278  -2. * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
279  sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
280  );
281  };
282 
284  variable_permeability_sigma = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
285  return 3 * PI / (1. + x(0) + x(1) + x(2))
286  * Eigen::Vector3d(
287  sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2)),
288  0.,
289  -cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
290  );
291  };
292 
294  variable_permeability_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
295  double a = 3 * PI / std::pow(1. + x(0) + x(1) + x(2), 2);
296  double b = 3 * std::pow(PI, 2) / (1. + x(0) + x(1) + x(2));
297  double dy_sigmax = -a * sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2))
298  - b * sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2));
299  double dz_sigmax = -a * sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2))
300  - b * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2));
301  ;
302  double dx_sigmaz = a * cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
303  + b * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2));
304 
305  double dy_sigmaz = a * cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
306  + b * cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2));
307 
308  return Eigen::Vector3d(
309  dy_sigmaz,
310  dz_sigmax - dx_sigmaz,
311  -dy_sigmax
312  );
313  };
314 
317  [](const Cell & T, const Eigen::Vector3d & x) -> double { return 1. + x(0) + x(1) + x(2); },
318  [](const Cell & T) -> size_t { return 1; }
319  );
320 
321 } // end of namespace HArDCore3D
322 
323 #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
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:186
static Magnetostatics::SolutionPotentialType linear_u
Definition: ddr-magnetostatics.hpp:213
static Magnetostatics::SolutionCurlType linear_sigma
Definition: ddr-magnetostatics.hpp:222
static Magnetostatics::SolutionCurlType trigonometric_sigma
Definition: ddr-magnetostatics.hpp:246
Eigen::SparseMatrix< double > SystemMatrixType
Definition: sddr-magnetostatics.hpp:54
const SXCurl & sxCurl() const
Returns the space SXCurl.
Definition: sddr-magnetostatics.hpp:95
static Magnetostatics::PermeabilityType constant_mu
Definition: ddr-magnetostatics.hpp:208
static Magnetostatics::PermeabilityType linear_mu
Definition: ddr-magnetostatics.hpp:232
size_t nbSCDOFs() const
Returns the number of statically condensed DOFs (here, the cell magnetic field DOFs)
Definition: ddr-magnetostatics.hpp:79
static Magnetostatics::SolutionPotentialType constant_u
Definition: ddr-magnetostatics.hpp:193
static Magnetostatics::SolutionPotentialType trigonometric_u
Definition: ddr-magnetostatics.hpp:237
static Magnetostatics::PermeabilityType variable_permeability_mu
Definition: ddr-magnetostatics.hpp:311
static Magnetostatics::ForcingTermType constant_f
Definition: ddr-magnetostatics.hpp:203
size_t dimensionSpace() const
Returns the dimension of the magnetic field + potential space.
Definition: ddr-magnetostatics.hpp:73
static Magnetostatics::SolutionCurlType constant_sigma
Definition: ddr-magnetostatics.hpp:198
static Magnetostatics::PermeabilityType trigonometric_mu
Definition: ddr-magnetostatics.hpp:264
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> SolutionCurlType
Definition: sddr-magnetostatics.hpp:58
IntegralWeight PermeabilityType
Definition: sddr-magnetostatics.hpp:59
const SXDiv & sxDiv() const
Returns the space SXDiv.
Definition: sddr-magnetostatics.hpp:101
static Magnetostatics::SolutionCurlType variable_permeability_sigma
Definition: ddr-magnetostatics.hpp:279
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> SolutionPotentialType
Definition: sddr-magnetostatics.hpp:57
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> ForcingTermType
Definition: sddr-magnetostatics.hpp:56
static Magnetostatics::SolutionPotentialType variable_permeability_u
Definition: ddr-magnetostatics.hpp:270
static Magnetostatics::ForcingTermType trigonometric_f
Definition: ddr-magnetostatics.hpp:255
static Magnetostatics::ForcingTermType linear_f
Definition: ddr-magnetostatics.hpp:227
static Magnetostatics::ForcingTermType variable_permeability_f
Definition: ddr-magnetostatics.hpp:289
bool use_threads
Definition: HHO_DiffAdvecReac.hpp:47
Definition: ddr-magnetostatics.hpp:40
MeshND::Cell< 2 > Cell
Definition: Mesh2D.hpp:13
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
const Eigen::VectorXd & systemVector() const
Returns the linear system right-hand side vector.
Definition: sddr-magnetostatics.hpp:117
double & stabilizationParameter()
Returns the stabilization parameter.
Definition: sddr-magnetostatics.hpp:132
size_t nbSCDOFs() const
Returns the number of statically condensed DOFs (here, the cell magnetic field DOFs)
Definition: sddr-magnetostatics.hpp:83
Eigen::VectorXd & systemVector()
Returns the linear system right-hand side vector.
Definition: sddr-magnetostatics.hpp:122
const SystemMatrixType & scMatrix() const
Returns the static condensation recovery operator.
Definition: sddr-magnetostatics.hpp:137
size_t dimensionSpace() const
Returns the dimension of the magnetic field + potential space.
Definition: sddr-magnetostatics.hpp:77
Magnetostatics(const DDRCore &ddrcore, bool use_threads, std::ostream &output=std::cout)
Constructor.
const double & stabilizationParameter() const
Returns the stabilization parameter.
Definition: sddr-magnetostatics.hpp:127
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.
SystemMatrixType & systemMatrix()
Returns the linear system matrix.
Definition: sddr-magnetostatics.hpp:112
size_t sizeSystem() const
Returns the size of the statically condensed system.
Definition: sddr-magnetostatics.hpp:89
Eigen::VectorXd & scVector()
Returns the static condensation rhs.
Definition: sddr-magnetostatics.hpp:142
void assembleLinearSystem(const ForcingTermType &f, const PermeabilityType &mu, const SolutionPotentialType &u)
Assemble the global system
const SystemMatrixType & systemMatrix() const
Returns the linear system matrix.
Definition: sddr-magnetostatics.hpp:107