HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
ddr-stokes.hpp
Go to the documentation of this file.
1 // Implementation of the discrete de Rham sequence for the Stokes problem in curl-curl formulation.
2 //
3 // Authors: Daniele Di Pietro (daniele.di-pietro@umontpellier.fr)
4 // Jerome Droniou (jerome.droniou@monash.edu)
5 //
6 
7 
8 #ifndef STOKES_HPP
9 #define STOKES_HPP
10 
11 #include <iostream>
12 
13 #include <boost/math/constants/constants.hpp>
14 
15 #include <Eigen/Sparse>
16 #include <unsupported/Eigen/SparseExtra>
17 
18 #include <mesh.hpp>
19 #include <mesh_builder.hpp>
21 
22 #include <xgrad.hpp>
23 #include <xcurl.hpp>
24 #include <xdiv.hpp>
25 
31 namespace HArDCore3D
32 {
33 
40  struct StokesNorms
41  {
43  StokesNorms(double norm_u, double norm_curl_u, double norm_p, double norm_grad_p):
44  u(norm_u),
45  curl_u(norm_curl_u),
46  p(norm_p),
47  grad_p(norm_grad_p),
48  hcurl_u( std::sqrt( std::pow(norm_u, 2) + std::pow(norm_curl_u, 2) ) ),
49  hgrad_p( std::sqrt( std::pow(norm_p, 2) + std::pow(norm_grad_p, 2) ) )
50  {
51  // do nothing
52  };
53 
54  double u;
55  double curl_u;
56  double p;
57  double grad_p;
58  double hcurl_u;
59  double hgrad_p;
60 
61  };
62 
64  struct Stokes
65  {
66  typedef Eigen::SparseMatrix<double> SystemMatrixType;
67 
68  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> ForcingTermType;
69  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> VelocityType;
70  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> VorticityType;
71  typedef std::function<double(const Eigen::Vector3d &)> PressureType;
72  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> PressureGradientType;
74 
76  Stokes(
77  const DDRCore & ddrcore,
78  bool use_threads,
79  std::ostream & output = std::cout
80  );
81 
84  const ForcingTermType & f,
85  const VelocityType & u,
86  const VorticityType & omega,
87  const ViscosityType & nu
88  );
89 
91  inline size_t dimension() const
92  {
93  return m_xcurl.dimension() + m_xgrad.dimension();
94  }
95 
97  inline size_t nbSCDOFs_u() const
98  {
99  return m_ddrcore.mesh().n_cells() * m_nloc_sc_u;
100  }
101 
103  inline size_t nbSCDOFs_p() const
104  {
105  return m_ddrcore.mesh().n_cells() * m_nloc_sc_p;
106  }
107 
109  inline size_t nbSCDOFs() const
110  {
111  return m_ddrcore.mesh().n_cells() * (m_nloc_sc_u + m_nloc_sc_p);
112  }
113 
115  inline size_t sizeSystem() const
116  {
117  return dimension() + 1 - nbSCDOFs();
118  }
119 
121  inline const XGrad & xGrad() const
122  {
123  return m_xgrad;
124  }
125 
127  inline const XCurl & xCurl() const
128  {
129  return m_xcurl;
130  }
131 
133  inline const XDiv & xDiv() const
134  {
135  return m_xdiv;
136  }
137 
139  inline const SystemMatrixType & systemMatrix() const {
140  return m_A;
141  }
142 
145  return m_A;
146  }
147 
149  inline const Eigen::VectorXd & systemVector() const {
150  return m_b;
151  }
152 
154  inline Eigen::VectorXd & systemVector() {
155  return m_b;
156  }
157 
159  inline const SystemMatrixType & scMatrix() const {
160  return m_sc_A;
161  }
162 
164  inline Eigen::VectorXd & scVector() {
165  return m_sc_b;
166  }
167 
169  inline const double & stabilizationParameter() const {
170  return m_stab_par;
171  }
172 
174  inline double & stabilizationParameter() {
175  return m_stab_par;
176  }
177 
179  std::vector<StokesNorms> computeStokesNorms(
180  const std::vector<Eigen::VectorXd> & list_dofs
181  ) const;
182 
184  std::pair<StokesNorms, StokesNorms> computeContinuousErrorsNorms(
185  const Eigen::VectorXd & v,
186  const VelocityType & u,
187  const VorticityType & curl_u,
188  const PressureType & p,
189  const PressureGradientType & grad_p
190  ) const;
191 
194  const PressureType & p,
195  const PressureGradientType & grad_p,
196  const size_t deg_quad_interpolate
197  ) const;
198 
199  private:
201  std::pair<Eigen::MatrixXd, Eigen::VectorXd>
202  _compute_local_contribution(
203  size_t iT,
204  const Eigen::VectorXd & interp_f,
205  const ViscosityType & nu
206  );
207 
209  LocalStaticCondensation _compute_static_condensation(const size_t & iT) const;
210 
212  void _assemble_local_contribution(
213  size_t iT,
214  const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT,
215  std::list<Eigen::Triplet<double> > & A1,
216  Eigen::VectorXd & b1,
217  std::list<Eigen::Triplet<double> > & A2,
218  Eigen::VectorXd & b2
219  );
220 
221  const DDRCore & m_ddrcore;
222  bool m_use_threads;
223  std::ostream & m_output;
224  const XGrad m_xgrad;
225  const XCurl m_xcurl;
226  const XDiv m_xdiv;
227  const size_t m_nloc_sc_u; // Nb of statically condensed DOFs for velocity in each cell (cell unknowns)
228  const size_t m_nloc_sc_p; // Nb of statically condensed DOFs for pressure in each cell (cell unknowns)
229  SystemMatrixType m_A; // Matrix and RHS of statically condensed system
230  Eigen::VectorXd m_b;
231  SystemMatrixType m_sc_A; // Static condensation operator and RHS (to recover statically condensed DOFs)
232  Eigen::VectorXd m_sc_b;
233  double m_stab_par;
234  };
235 
236  //------------------------------------------------------------------------------
237  // Exact solutions
238  //------------------------------------------------------------------------------
239 
240  static const double PI = boost::math::constants::pi<double>();
241  using std::sin;
242  using std::cos;
243  using std::pow;
244 
245  //------------------------------------------------------------------------------
246  // Trigonometric solution
247  //------------------------------------------------------------------------------
248 
249  double pressure_scaling = 1.;
250  static Stokes::VelocityType
251  trigonometric_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
252  return Eigen::Vector3d(
253  0.5 * sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2)),
254  0.5 * cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
255  -cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2))
256  );
257  };
258 
259  static Stokes::VorticityType
260  trigonometric_curl_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
261  return 3. * PI * Eigen::Vector3d(
262  cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
263  -sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
264  0.
265  );
266  };
267 
268  static Stokes::PressureType
269  trigonometric_p = [](const Eigen::Vector3d & x) -> double {
270  return pressure_scaling * sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2));
271  };
272 
274  trigonometric_grad_p = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
275  return pressure_scaling * 2. * PI * Eigen::Vector3d(
276  cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
277  sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
278  sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2))
279  );
280  };
281 
283  trigonometric_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
284  return 6. * std::pow(PI, 2) * Eigen::Vector3d(
285  sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2)),
286  cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
287  -2. * cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2))
288  )
290  };
291 
292  static Stokes::ViscosityType
294 
295  //------------------------------------------------------------------------------
296  // Linear solution
297  //------------------------------------------------------------------------------
298 
299  static Stokes::VelocityType
300  linear_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
301  return Eigen::Vector3d(x(0), -x(1), 0.);
302  };
303 
304  static Stokes::VorticityType
305  linear_curl_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
306  return Eigen::Vector3d::Zero();
307  };
308 
309  static Stokes::PressureType
310  linear_p = [](const Eigen::Vector3d & x) -> double {
311  return pressure_scaling * (x(0) - x(1));
312  };
313 
315  linear_grad_p = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
316  return pressure_scaling * Eigen::Vector3d(1., -1., 0.);
317  };
318 
320  linear_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
321  return linear_grad_p(x);
322  };
323 
324  static Stokes::ViscosityType
326 
327 
328  //------------------------------------------------------------------------------
329  // Solution derived from a field
330  //------------------------------------------------------------------------------
331  // Setting f(x,y,z)= (sin(pi x) sin(pi y) sin(pi z))^3, we take
332  // U=(-d_y f, d_x f, 0), p=0
333  // Calculations are made symbolically in matlab
334 
335  static Stokes::VelocityType
336  field_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
337  double ux = PI*cos(PI*x(1))*pow(sin(PI*x(0)),3.0)*pow(sin(PI*x(1)),2.0)*pow(sin(PI*x(2)),3.0)*3.0;
338  double uy = PI*cos(PI*x(0))*pow(sin(PI*x(0)),2.0)*pow(sin(PI*x(1)),3.0)*pow(sin(PI*x(2)),3.0)*(-3.0);
339  return Eigen::Vector3d(ux, uy, 0.);
340  };
341 
342  static Stokes::VorticityType
343  field_curl_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
344  double cux = (PI*PI)*cos(PI*x(0))*cos(PI*x(2))*pow(sin(PI*x(0)),2.0)*pow(sin(PI*x(1)),3.0)*pow(sin(PI*x(2)),2.0)*9.0;
345  double cuy = (PI*PI)*cos(PI*x(1))*cos(PI*x(2))*pow(sin(PI*x(0)),3.0)*pow(sin(PI*x(1)),2.0)*pow(sin(PI*x(2)),2.0)*9.0;
346  double cuz = (PI*PI)*pow(sin(PI*x(0)),3.0)*pow(sin(PI*x(1)),3.0)*pow(sin(PI*x(2)),3.0)*6.0-(PI*PI)*pow(cos(PI*x(0)),2.0)*sin(PI*x(0))*pow(sin(PI*x(1)),3.0)*pow(sin(PI*x(2)),3.0)*6.0-(PI*PI)*pow(cos(PI*x(1)),2.0)*pow(sin(PI*x(0)),3.0)*sin(PI*x(1))*pow(sin(PI*x(2)),3.0)*6.0;
347 
348  return Eigen::Vector3d(cux, cuy, cuz);
349  };
350 
351  static Stokes::PressureType
352  field_p = [](const Eigen::Vector3d & x) -> double {
353  return 0.;
354  };
355 
357  field_grad_p = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
358  return Eigen::Vector3d(0., 0., 0.);
359  };
360 
362  field_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
363  double fx = (PI*PI*PI)*pow(cos(PI*x(1)),3.0)*pow(sin(PI*x(0)),3.0)*pow(sin(PI*x(2)),3.0)*-6.0+(PI*PI*PI)*cos(PI*x(1))*pow(sin(PI*x(0)),3.0)*pow(sin(PI*x(1)),2.0)*pow(sin(PI*x(2)),3.0)*3.9E+1-(PI*PI*PI)*pow(cos(PI*x(0)),2.0)*cos(PI*x(1))*sin(PI*x(0))*pow(sin(PI*x(1)),2.0)*pow(sin(PI*x(2)),3.0)*1.8E+1-(PI*PI*PI)*cos(PI*x(1))*pow(cos(PI*x(2)),2.0)*pow(sin(PI*x(0)),3.0)*pow(sin(PI*x(1)),2.0)*sin(PI*x(2))*1.8E+1;
364 
365  double fy = (PI*PI*PI)*pow(cos(PI*x(0)),3.0)*pow(sin(PI*x(1)),3.0)*pow(sin(PI*x(2)),3.0)*6.0-(PI*PI*PI)*cos(PI*x(0))*pow(sin(PI*x(0)),2.0)*pow(sin(PI*x(1)),3.0)*pow(sin(PI*x(2)),3.0)*3.9E+1+(PI*PI*PI)*cos(PI*x(0))*pow(cos(PI*x(1)),2.0)*pow(sin(PI*x(0)),2.0)*sin(PI*x(1))*pow(sin(PI*x(2)),3.0)*1.8E+1+(PI*PI*PI)*cos(PI*x(0))*pow(cos(PI*x(2)),2.0)*pow(sin(PI*x(0)),2.0)*pow(sin(PI*x(1)),3.0)*sin(PI*x(2))*1.8E+1;
366 
367  double fz = 0.;
368  return Eigen::Vector3d(fx, fy, fz);
369  };
370 
371  static Stokes::ViscosityType
373 
374  //------------------------------------------------------------------------------
375 
376 
377 } // end of namespace HArDCore3D
378 
379 #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
Discrete H1 space: local operators, L2 product and global interpolator.
Definition: xgrad.hpp:18
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
static Magnetostatics::SolutionPotentialType linear_u
Definition: ddr-magnetostatics.hpp:213
static Magnetostatics::SolutionPotentialType trigonometric_u
Definition: ddr-magnetostatics.hpp:237
static Magnetostatics::ForcingTermType trigonometric_f
Definition: ddr-magnetostatics.hpp:255
static Magnetostatics::ForcingTermType linear_f
Definition: ddr-magnetostatics.hpp:227
double grad_p
Norm of grad of pressure.
Definition: ddr-stokes.hpp:57
static Stokes::ForcingTermType field_f
Definition: ddr-stokes.hpp:362
void assembleLinearSystem(const ForcingTermType &f, const VelocityType &u, const VorticityType &omega, const ViscosityType &nu)
Assemble the global system
Definition: ddr-stokes.cpp:338
static Stokes::ViscosityType field_nu
Definition: ddr-stokes.hpp:372
size_t nbSCDOFs_u() const
Returns the number of statically condensed DOFs (velocity)
Definition: ddr-stokes.hpp:97
std::pair< StokesNorms, StokesNorms > computeContinuousErrorsNorms(const Eigen::VectorXd &v, const VelocityType &u, const VorticityType &curl_u, const PressureType &p, const PressureGradientType &grad_p) const
Compute the continuous Hcurl errors in u and Hgrad error in p (1st component), and same with continuo...
Definition: ddr-stokes.cpp:612
static Stokes::ViscosityType trigonometric_nu
Definition: ddr-stokes.hpp:293
Eigen::SparseMatrix< double > SystemMatrixType
Definition: ddr-stokes.hpp:66
double hgrad_p
Hgrad norm of p.
Definition: ddr-stokes.hpp:59
std::vector< StokesNorms > computeStokesNorms(const std::vector< Eigen::VectorXd > &list_dofs) const
Compute the discrete L2 norms, for a family of Eigen::VectorXd representing velocities & pressures,...
Definition: ddr-stokes.cpp:556
static Stokes::VelocityType field_u
Definition: ddr-stokes.hpp:336
size_t nbSCDOFs_p() const
Returns the number of statically condensed DOFs (pressure)
Definition: ddr-stokes.hpp:103
static Stokes::PressureGradientType field_grad_p
Definition: ddr-stokes.hpp:357
double hcurl_u
Hcurl norm of u.
Definition: ddr-stokes.hpp:58
const double & stabilizationParameter() const
Returns the stabilization parameter.
Definition: ddr-stokes.hpp:169
double evaluateDefectCommutationGradInterp(const PressureType &p, const PressureGradientType &grad_p, const size_t deg_quad_interpolate) const
Evaluate the defect of commutation G_h (Igrad p) = Icurl (grad p) by computing the norm of the differ...
Definition: ddr-stokes.cpp:706
static Stokes::PressureGradientType linear_grad_p
Definition: ddr-stokes.hpp:315
static Stokes::ViscosityType linear_nu
Definition: ddr-stokes.hpp:325
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> VelocityType
Definition: ddr-stokes.hpp:69
const XGrad & xGrad() const
Returns the space XGrad.
Definition: ddr-stokes.hpp:121
size_t nbSCDOFs() const
Returns the number of statically condensed DOFs (both velocity and pressure)
Definition: ddr-stokes.hpp:109
size_t sizeSystem() const
Returns the size of the statically condensed system (with Lagrange multiplier)
Definition: ddr-stokes.hpp:115
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> PressureGradientType
Definition: ddr-stokes.hpp:72
Eigen::VectorXd & systemVector()
Returns the linear system right-hand side vector.
Definition: ddr-stokes.hpp:154
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> VorticityType
Definition: ddr-stokes.hpp:70
double u
Norm of velocity.
Definition: ddr-stokes.hpp:52
std::function< double(const Eigen::Vector3d &)> PressureType
Definition: ddr-stokes.hpp:71
double & stabilizationParameter()
Returns the stabilization parameter.
Definition: ddr-stokes.hpp:174
SystemMatrixType & systemMatrix()
Returns the linear system matrix.
Definition: ddr-stokes.hpp:144
const Eigen::VectorXd & systemVector() const
Returns the linear system right-hand side vector.
Definition: ddr-stokes.hpp:149
const SystemMatrixType & scMatrix() const
Returns the static condensation recovery operator.
Definition: ddr-stokes.hpp:159
const XDiv & xDiv() const
Returns the space XDiv.
Definition: ddr-stokes.hpp:133
Stokes(const DDRCore &ddrcore, bool use_threads, std::ostream &output=std::cout)
Constructor.
Definition: ddr-stokes.cpp:312
const SystemMatrixType & systemMatrix() const
Returns the linear system matrix.
Definition: ddr-stokes.hpp:139
size_t dimension() const
Returns the global problem dimension (without Lagrange multiplier, just velocity & pressure)
Definition: ddr-stokes.hpp:91
static Stokes::PressureType trigonometric_p
Definition: ddr-stokes.hpp:269
IntegralWeight ViscosityType
Definition: ddr-stokes.hpp:73
StokesNorms(double norm_u, double norm_curl_u, double norm_p, double norm_grad_p)
Constructor.
Definition: ddr-stokes.hpp:43
double p
Norm of pressure.
Definition: ddr-stokes.hpp:56
static Stokes::VorticityType field_curl_u
Definition: ddr-stokes.hpp:343
Eigen::VectorXd & scVector()
Returns the static condensation rhs.
Definition: ddr-stokes.hpp:164
static Stokes::VorticityType linear_curl_u
Definition: ddr-stokes.hpp:305
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> ForcingTermType
Definition: ddr-stokes.hpp:68
static Stokes::VorticityType trigonometric_curl_u
Definition: ddr-stokes.hpp:260
double curl_u
Norm of curl of velocity (vorticity)
Definition: ddr-stokes.hpp:55
static Stokes::PressureType linear_p
Definition: ddr-stokes.hpp:310
static Stokes::PressureType field_p
Definition: ddr-stokes.hpp:352
const XCurl & xCurl() const
Returns the space XCurl.
Definition: ddr-stokes.hpp:127
static Stokes::PressureGradientType trigonometric_grad_p
Definition: ddr-stokes.hpp:274
double pressure_scaling
Definition: ddr-stokes.hpp:249
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
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
Structure to store norm components (for velocity, its curl, pressure, and its gradient),...
Definition: ddr-stokes.hpp:41
Assemble a Stokes problem.
Definition: ddr-stokes.hpp:65