HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
ddr-magnetostatics.hpp
Go to the documentation of this file.
1 // Implementation of the discrete de Rham sequence for the magnetostatic problem.
2 //
3 // Author: Daniele Di Pietro (daniele.di-pietro@umontpellier.fr)
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  * If you use this code in a scientific publication, please mention the above article.
14  *
15  */
16 
17 #ifndef MAGNETOSTATICS_HPP
18 #define MAGNETOSTATICS_HPP
19 
20 #include <iostream>
21 
22 #include <boost/math/constants/constants.hpp>
23 
24 #include <Eigen/Sparse>
25 #include <unsupported/Eigen/SparseExtra>
26 
27 #include <mesh.hpp>
28 #include <mesh_builder.hpp>
30 
31 #include <xdiv.hpp>
32 #include <xcurl.hpp>
33 
39 namespace HArDCore3D
40 {
41 
49  {
50  typedef Eigen::SparseMatrix<double> SystemMatrixType;
51 
52  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> ForcingTermType;
53  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> SolutionPotentialType;
54  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> SolutionCurlType;
56 
59  const DDRCore & ddrcore,
60  bool use_threads,
61  std::ostream & output = std::cout
62  );
63 
66  const ForcingTermType & f,
67  const PermeabilityType & mu,
68  const SolutionPotentialType & u
69  );
70 
71 
73  inline size_t dimensionSpace() const
74  {
75  return m_xcurl.dimension() + m_xdiv.dimension();
76  }
77 
79  inline size_t nbSCDOFs() const
80  {
81  return m_ddrcore.mesh().n_cells() * m_nloc_sc_H;
82  }
83 
85  inline size_t sizeSystem() const
86  {
87  return dimensionSpace() - nbSCDOFs();
88  }
89 
91  inline const XCurl & xCurl() const
92  {
93  return m_xcurl;
94  }
95 
97  inline const XDiv & xDiv() const
98  {
99  return m_xdiv;
100  }
101 
103  inline const SystemMatrixType & systemMatrix() const {
104  return m_A;
105  }
106 
109  return m_A;
110  }
111 
113  inline const Eigen::VectorXd & systemVector() const {
114  return m_b;
115  }
116 
118  inline Eigen::VectorXd & systemVector() {
119  return m_b;
120  }
121 
123  inline const double & stabilizationParameter() const {
124  return m_stab_par;
125  }
126 
128  inline double & stabilizationParameter() {
129  return m_stab_par;
130  }
131 
133  inline const SystemMatrixType & scMatrix() const {
134  return m_sc_A;
135  }
136 
138  inline Eigen::VectorXd & scVector() {
139  return m_sc_b;
140  }
141 
143  std::vector<double> computeNorms(
144  const std::vector<Eigen::VectorXd> & list_dofs
145  ) const;
146 
147  private:
149  std::pair<Eigen::MatrixXd, Eigen::VectorXd>
150  _compute_local_contribution(
151  size_t iT,
152  const ForcingTermType & f,
153  const PermeabilityType & mu
154  );
155 
157  LocalStaticCondensation _compute_static_condensation(const size_t & iT) const;
158 
160  void _assemble_local_contribution(
161  size_t iT,
162  const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT,
163  std::list<Eigen::Triplet<double> > & A1,
164  Eigen::VectorXd & b1,
165  std::list<Eigen::Triplet<double> > & A2,
166  Eigen::VectorXd & b2
167  );
168 
169  const DDRCore & m_ddrcore;
170  bool m_use_threads;
171  std::ostream & m_output;
172  XCurl m_xcurl;
173  XDiv m_xdiv;
174  const size_t m_nloc_sc_H; // Number of statically condensed DOFs of H in each cell
175  SystemMatrixType m_A; // Matrix and RHS of statically condensed system
176  Eigen::VectorXd m_b;
177  SystemMatrixType m_sc_A; // Static condensation operator and RHS (to recover statically condensed DOFs)
178  Eigen::VectorXd m_sc_b;
179  double m_stab_par;
180  };
181 
182  //------------------------------------------------------------------------------
183  // Exact solutions
184  //------------------------------------------------------------------------------
185 
186  static const double PI = boost::math::constants::pi<double>();
187  using std::sin;
188  using std::cos;
189 
190  //------------------------------------------------------------------------------
191 
193  constant_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
194  return Eigen::Vector3d(1., 2., 3.);
195  };
196 
198  constant_sigma = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
199  return Eigen::Vector3d::Zero();
200  };
201 
203  constant_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
204  return Eigen::Vector3d::Zero();
205  };
206 
209 
210  //------------------------------------------------------------------------------
211 
213  linear_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
214  return Eigen::Vector3d(
215  x(0) + x(1) + x(2),
216  2.*(x(0)-x(1)+x(2)),
217  -x(0) - x(1) + x(2)
218  );
219  };
220 
222  linear_sigma = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
223  return Eigen::Vector3d(-3., 2., 1.);
224  };
225 
227  linear_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
228  return Eigen::Vector3d::Zero();
229  };
230 
233 
234  //------------------------------------------------------------------------------
235 
237  trigonometric_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
238  return Eigen::Vector3d(
239  cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
240  -2. * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
241  sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
242  );
243  };
244 
246  trigonometric_sigma = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
247  return 3. * PI * Eigen::Vector3d(
248  sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2)),
249  0.,
250  -cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
251  );
252  };
253 
255  trigonometric_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
256  return 3. * std::pow(PI, 2) * Eigen::Vector3d(
257  cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
258  -2. * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
259  sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
260  );
261  };
262 
265 
266  //------------------------------------------------------------------------------
267 
268 
270  variable_permeability_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
271  return Eigen::Vector3d(
272  cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
273  -2. * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
274  sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
275  );
276  };
277 
279  variable_permeability_sigma = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
280  return 3 * PI / (1. + x(0) + x(1) + x(2))
281  * Eigen::Vector3d(
282  sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2)),
283  0.,
284  -cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
285  );
286  };
287 
289  variable_permeability_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
290  double a = 3 * PI / std::pow(1. + x(0) + x(1) + x(2), 2);
291  double b = 3 * std::pow(PI, 2) / (1. + x(0) + x(1) + x(2));
292  double dy_sigmax = -a * sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2))
293  - b * sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2));
294  double dz_sigmax = -a * sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2))
295  - b * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2));
296  ;
297  double dx_sigmaz = a * cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
298  + b * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2));
299 
300  double dy_sigmaz = a * cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
301  + b * cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2));
302 
303  return Eigen::Vector3d(
304  dy_sigmaz,
305  dz_sigmax - dx_sigmaz,
306  -dy_sigmax
307  );
308  };
309 
312  [](const Cell & T, const Eigen::Vector3d & x) -> double { return 1. + x(0) + x(1) + x(2); },
313  [](const Cell & T) -> size_t { return 1; }
314  );
315 
316 } // end of namespace HArDCore3D
317 
318 #endif
Construct all polynomial spaces for the DDR sequence.
Definition: ddrcore.hpp:62
Discrete Hcurl space: local operators, L2 product and global interpolator.
Definition: xcurl.hpp:19
Discrete Hdiv space: local operators, L2 product and global interpolator.
Definition: xdiv.hpp:20
size_t dimension() const
Returns the dimension of the global space (all DOFs for all geometric entities)
Definition: localdofspace.hpp:80
const Mesh & mesh() const
Return a const reference to the mesh.
Definition: ddrcore.hpp:134
static const double PI
Definition: ddr-magnetostatics.hpp:186
const XCurl & xCurl() const
Returns the space XCurl.
Definition: ddr-magnetostatics.hpp:91
static Magnetostatics::SolutionPotentialType linear_u
Definition: ddr-magnetostatics.hpp:213
static Magnetostatics::SolutionCurlType linear_sigma
Definition: ddr-magnetostatics.hpp:222
const Eigen::VectorXd & systemVector() const
Returns the linear system right-hand side vector.
Definition: ddr-magnetostatics.hpp:113
static Magnetostatics::SolutionCurlType trigonometric_sigma
Definition: ddr-magnetostatics.hpp:246
Eigen::SparseMatrix< double > SystemMatrixType
Definition: ddr-magnetostatics.hpp:50
const XDiv & xDiv() const
Returns the space XDiv.
Definition: ddr-magnetostatics.hpp:97
double & stabilizationParameter()
Returns the stabilization parameter.
Definition: ddr-magnetostatics.hpp:128
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
Eigen::VectorXd & systemVector()
Returns the linear system right-hand side vector.
Definition: ddr-magnetostatics.hpp:118
static Magnetostatics::SolutionPotentialType constant_u
Definition: ddr-magnetostatics.hpp:193
static Magnetostatics::SolutionPotentialType trigonometric_u
Definition: ddr-magnetostatics.hpp:237
const SystemMatrixType & scMatrix() const
Returns the static condensation recovery operator.
Definition: ddr-magnetostatics.hpp:133
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
Magnetostatics(const DDRCore &ddrcore, bool use_threads, std::ostream &output=std::cout)
Constructor.
Definition: ddr-magnetostatics.cpp:246
static Magnetostatics::SolutionCurlType constant_sigma
Definition: ddr-magnetostatics.hpp:198
const double & stabilizationParameter() const
Returns the stabilization parameter.
Definition: ddr-magnetostatics.hpp:123
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.
Definition: ddr-magnetostatics.cpp:451
static Magnetostatics::PermeabilityType trigonometric_mu
Definition: ddr-magnetostatics.hpp:264
SystemMatrixType & systemMatrix()
Returns the linear system matrix.
Definition: ddr-magnetostatics.hpp:108
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> SolutionCurlType
Definition: ddr-magnetostatics.hpp:54
IntegralWeight PermeabilityType
Definition: ddr-magnetostatics.hpp:55
size_t sizeSystem() const
Returns the size of the statically condensed system.
Definition: ddr-magnetostatics.hpp:85
static Magnetostatics::SolutionCurlType variable_permeability_sigma
Definition: ddr-magnetostatics.hpp:279
Eigen::VectorXd & scVector()
Returns the static condensation rhs.
Definition: ddr-magnetostatics.hpp:138
void assembleLinearSystem(const ForcingTermType &f, const PermeabilityType &mu, const SolutionPotentialType &u)
Assemble the global system
Definition: ddr-magnetostatics.cpp:271
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> SolutionPotentialType
Definition: ddr-magnetostatics.hpp:53
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> ForcingTermType
Definition: ddr-magnetostatics.hpp:52
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
const SystemMatrixType & systemMatrix() const
Returns the linear system matrix.
Definition: ddr-magnetostatics.hpp:103
static Magnetostatics::ForcingTermType variable_permeability_f
Definition: ddr-magnetostatics.hpp:289
bool use_threads
Definition: HHO_DiffAdvecReac.hpp:47
std::size_t n_cells() const
number of cells in the mesh.
Definition: MeshND.hpp:60
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
Assemble a magnetostatic problem.
Definition: ddr-magnetostatics.hpp:49