HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
sddr-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 // Author: Jerome Droniou (jerome.droniou@monash.edu)
4 //
5 
6 
7 #ifndef SDDR_STOKES_HPP
8 #define SDDR_STOKES_HPP
9 
10 #include <iostream>
11 
12 #include <boost/math/constants/constants.hpp>
13 
14 #include <Eigen/Sparse>
15 #include <unsupported/Eigen/SparseExtra>
16 
17 #include <mesh.hpp>
18 #include <mesh_builder.hpp>
20 
21 #include <sxgrad.hpp>
22 #include <sxcurl.hpp>
23 #include <sxdiv.hpp>
24 
30 namespace HArDCore3D
31 {
32 
39  struct StokesNorms
40  {
42  StokesNorms(double norm_u, double norm_curl_u, double norm_p, double norm_grad_p):
43  u(norm_u),
44  curl_u(norm_curl_u),
45  p(norm_p),
46  grad_p(norm_grad_p),
47  hcurl_u( std::sqrt( std::pow(norm_u, 2) + std::pow(norm_curl_u, 2) ) ),
48  hgrad_p( std::sqrt( std::pow(norm_p, 2) + std::pow(norm_grad_p, 2) ) )
49  {
50  // do nothing
51  };
52 
53  double u;
54  double curl_u;
55  double p;
56  double grad_p;
57  double hcurl_u;
58  double hgrad_p;
59 
60  };
61 
63  struct Stokes
64  {
65  typedef Eigen::SparseMatrix<double> SystemMatrixType;
66 
67  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> ForcingTermType;
68  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> VelocityType;
69  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> VorticityType;
70  typedef std::function<double(const Eigen::Vector3d &)> PressureType;
71  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> PressureGradientType;
73 
76  const DDRCore & ddrcore,
77  bool use_threads,
78  std::ostream & output = std::cout
79  );
80 
83  const ForcingTermType & f,
84  const VelocityType & u,
85  const VorticityType & omega,
86  const ViscosityType & nu
87  );
88 
90  inline size_t dimensionSpace() const
91  {
92  return m_sxcurl.dimension() + m_sxgrad.dimension();
93  }
94 
96  inline size_t nbSCDOFs_u() const
97  {
98  return m_nloc_sc_u.sum();
99  }
100 
102  inline size_t nbSCDOFs_p() const
103  {
104  return m_nloc_sc_p.sum();
105  }
106 
108  inline size_t nbSCDOFs() const
109  {
110  return nbSCDOFs_u() + nbSCDOFs_p();
111  }
112 
114  inline size_t sizeSystem() const
115  {
116  return dimensionSpace() + 1 - nbSCDOFs();
117  }
118 
120  inline const SXGrad & sxGrad() const
121  {
122  return m_sxgrad;
123  }
124 
126  inline const SXCurl & sxCurl() const
127  {
128  return m_sxcurl;
129  }
130 
132  inline const SXDiv & sxDiv() const
133  {
134  return m_sxdiv;
135  }
136 
138  inline const SystemMatrixType & systemMatrix() const {
139  return m_A;
140  }
141 
144  return m_A;
145  }
146 
148  inline const Eigen::VectorXd & systemVector() const {
149  return m_b;
150  }
151 
153  inline Eigen::VectorXd & systemVector() {
154  return m_b;
155  }
156 
158  inline const SystemMatrixType & scMatrix() const {
159  return m_sc_A;
160  }
161 
163  inline Eigen::VectorXd & scVector() {
164  return m_sc_b;
165  }
166 
168  inline const double & stabilizationParameter() const {
169  return m_stab_par;
170  }
171 
173  inline double & stabilizationParameter() {
174  return m_stab_par;
175  }
176 
178  std::vector<StokesNorms> computeStokesNorms(
179  const std::vector<Eigen::VectorXd> & list_dofs
180  ) const;
181 
183  std::pair<StokesNorms, StokesNorms> computeContinuousErrorsNorms(
184  const Eigen::VectorXd & v,
185  const VelocityType & u,
186  const VorticityType & curl_u,
187  const PressureType & p,
188  const PressureGradientType & grad_p
189  ) const;
190 
191 
192  private:
194  std::pair<Eigen::MatrixXd, Eigen::VectorXd>
195  _compute_local_contribution(
196  size_t iT,
197  const Eigen::VectorXd & interp_f,
198  const ViscosityType & nu
199  );
200 
202  LocalStaticCondensation _compute_static_condensation(const size_t & iT) const;
203 
205  void _assemble_local_contribution(
206  size_t iT,
207  const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT,
208  std::list<Eigen::Triplet<double> > & A1,
209  Eigen::VectorXd & b1,
210  std::list<Eigen::Triplet<double> > & A2,
211  Eigen::VectorXd & b2
212  );
213 
214  const DDRCore & m_ddrcore;
215  bool m_use_threads;
216  std::ostream & m_output;
217  SerendipityProblem m_ser_pro;
218  const SXGrad m_sxgrad;
219  const SXCurl m_sxcurl;
220  const SXDiv m_sxdiv;
221  const Eigen::VectorXd m_nloc_sc_u; // Nb of statically condensed DOFs for velocity in each cell (cell unknowns)
222  const Eigen::VectorXd m_nloc_sc_p; // Nb of statically condensed DOFs for pressure in each cell (cell unknowns)
223  SystemMatrixType m_A; // Matrix and RHS of statically condensed system
224  Eigen::VectorXd m_b;
225  SystemMatrixType m_sc_A; // Static condensation operator and RHS (to recover statically condensed DOFs)
226  Eigen::VectorXd m_sc_b;
227  double m_stab_par;
228  };
229 
230  //------------------------------------------------------------------------------
231  // Exact solutions
232  //------------------------------------------------------------------------------
233 
234  static const double PI = boost::math::constants::pi<double>();
235  using std::sin;
236  using std::cos;
237  using std::pow;
238 
239  //------------------------------------------------------------------------------
240  // Trigonometric solution
241  //------------------------------------------------------------------------------
242 
243  double pressure_scaling = 1.;
244  static Stokes::VelocityType
245  trigonometric_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
246  return Eigen::Vector3d(
247  0.5 * sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2)),
248  0.5 * cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
249  -cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2))
250  );
251  };
252 
253  static Stokes::VorticityType
254  trigonometric_curl_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
255  return 3. * PI * Eigen::Vector3d(
256  cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
257  -sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
258  0.
259  );
260  };
261 
262  static Stokes::PressureType
263  trigonometric_p = [](const Eigen::Vector3d & x) -> double {
264  return pressure_scaling * sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2));
265  };
266 
268  trigonometric_grad_p = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
269  return pressure_scaling * 2. * PI * Eigen::Vector3d(
270  cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
271  sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
272  sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2))
273  );
274  };
275 
277  trigonometric_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
278  return 6. * std::pow(PI, 2) * Eigen::Vector3d(
279  sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2)),
280  cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
281  -2. * cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2))
282  )
284  };
285 
286  static Stokes::ViscosityType
288 
289  //------------------------------------------------------------------------------
290  // Linear solution
291  //------------------------------------------------------------------------------
292 
293  static Stokes::VelocityType
294  linear_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
295  return Eigen::Vector3d(x(0), -x(1), 0.);
296  };
297 
298  static Stokes::VorticityType
299  linear_curl_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
300  return Eigen::Vector3d::Zero();
301  };
302 
303  static Stokes::PressureType
304  linear_p = [](const Eigen::Vector3d & x) -> double {
305  return pressure_scaling * (x(0) - x(1));
306  };
307 
309  linear_grad_p = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
310  return pressure_scaling * Eigen::Vector3d(1., -1., 0.);
311  };
312 
314  linear_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
315  return linear_grad_p(x);
316  };
317 
318  static Stokes::ViscosityType
320 
321 
322  //------------------------------------------------------------------------------
323  // Solution derived from a field
324  //------------------------------------------------------------------------------
325  // Setting f(x,y,z)= (sin(pi x) sin(pi y) sin(pi z))^3, we take
326  // U=(-d_y f, d_x f, 0), p=0
327  // Calculations are made symbolically in matlab
328 
329  static Stokes::VelocityType
330  field_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
331  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;
332  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);
333  return Eigen::Vector3d(ux, uy, 0.);
334  };
335 
336  static Stokes::VorticityType
337  field_curl_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
338  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;
339  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;
340  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;
341 
342  return Eigen::Vector3d(cux, cuy, cuz);
343  };
344 
345  static Stokes::PressureType
346  field_p = [](const Eigen::Vector3d & x) -> double {
347  return 0.;
348  };
349 
351  field_grad_p = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
352  return Eigen::Vector3d(0., 0., 0.);
353  };
354 
356  field_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
357  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;
358 
359  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;
360 
361  double fz = 0.;
362  return Eigen::Vector3d(fx, fy, fz);
363  };
364 
365  static Stokes::ViscosityType
367 
368  //------------------------------------------------------------------------------
369 
370 
371 } // end of namespace HArDCore3D
372 
373 #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
Discrete Serendipity Hgrad space: local operators, L2 product and global interpolator.
Definition: sxgrad.hpp:20
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: variabledofspace.hpp:158
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
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
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
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
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
size_t nbSCDOFs() const
Returns the number of statically condensed DOFs (both velocity and pressure)
Definition: ddr-stokes.hpp:109
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> PressureGradientType
Definition: ddr-stokes.hpp:72
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
static Stokes::PressureType trigonometric_p
Definition: ddr-stokes.hpp:269
IntegralWeight ViscosityType
Definition: ddr-stokes.hpp:73
double p
Norm of pressure.
Definition: ddr-stokes.hpp:56
static Stokes::VorticityType field_curl_u
Definition: ddr-stokes.hpp:343
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
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
size_t dimensionSpace() const
Returns the global problem dimension (without Lagrange multiplier, just velocity & pressure)
Definition: sddr-stokes.hpp:90
const SXCurl & sxCurl() const
Returns the space SXCurl.
Definition: sddr-stokes.hpp:126
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...
size_t nbSCDOFs_u() const
Returns the number of statically condensed DOFs (velocity)
Definition: sddr-stokes.hpp:96
Eigen::SparseMatrix< double > SystemMatrixType
Definition: sddr-stokes.hpp:65
const SXDiv & sxDiv() const
Returns the space XDiv.
Definition: sddr-stokes.hpp:132
size_t nbSCDOFs_p() const
Returns the number of statically condensed DOFs (pressure)
Definition: sddr-stokes.hpp:102
const double & stabilizationParameter() const
Returns the stabilization parameter.
Definition: sddr-stokes.hpp:168
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> VelocityType
Definition: sddr-stokes.hpp:68
size_t nbSCDOFs() const
Returns the number of statically condensed DOFs (both velocity and pressure)
Definition: sddr-stokes.hpp:108
size_t sizeSystem() const
Returns the size of the statically condensed system (with Lagrange multiplier)
Definition: sddr-stokes.hpp:114
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> PressureGradientType
Definition: sddr-stokes.hpp:71
Eigen::VectorXd & systemVector()
Returns the linear system right-hand side vector.
Definition: sddr-stokes.hpp:153
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> VorticityType
Definition: sddr-stokes.hpp:69
std::function< double(const Eigen::Vector3d &)> PressureType
Definition: sddr-stokes.hpp:70
double & stabilizationParameter()
Returns the stabilization parameter.
Definition: sddr-stokes.hpp:173
SystemMatrixType & systemMatrix()
Returns the linear system matrix.
Definition: sddr-stokes.hpp:143
const Eigen::VectorXd & systemVector() const
Returns the linear system right-hand side vector.
Definition: sddr-stokes.hpp:148
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,...
const SystemMatrixType & scMatrix() const
Returns the static condensation recovery operator.
Definition: sddr-stokes.hpp:158
const SXGrad & sxGrad() const
Returns the space SXGrad.
Definition: sddr-stokes.hpp:120
void assembleLinearSystem(const ForcingTermType &f, const VelocityType &u, const VorticityType &omega, const ViscosityType &nu)
Assemble the global system
const SystemMatrixType & systemMatrix() const
Returns the linear system matrix.
Definition: sddr-stokes.hpp:138
IntegralWeight ViscosityType
Definition: sddr-stokes.hpp:72
StokesNorms(double norm_u, double norm_curl_u, double norm_p, double norm_grad_p)
Constructor.
Definition: sddr-stokes.hpp:42
Eigen::VectorXd & scVector()
Returns the static condensation rhs.
Definition: sddr-stokes.hpp:163
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> ForcingTermType
Definition: sddr-stokes.hpp:67
Stokes(const DDRCore &ddrcore, bool use_threads, std::ostream &output=std::cout)
Constructor.
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