HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
vem-stokes.hpp
Go to the documentation of this file.
1 // Implementation of the Virtual Element complex for the Stokes problem in curl-curl formulation.
2 //
3 // Author: Jerome Droniou (jerome.droniou@monash.edu)
4 //
5 
6 #ifndef STOKES_HPP
7 #define STOKES_HPP
8 
9 #include <iostream>
10 
11 #include <boost/math/constants/constants.hpp>
12 
13 #include <Eigen/Sparse>
14 #include <unsupported/Eigen/SparseExtra>
15 
16 #include <mesh.hpp>
17 #include <mesh_builder.hpp>
19 
20 #include <vdiv.hpp>
21 #include <vgrad.hpp>
22 #include <vcurl.hpp>
23 
29 namespace HArDCore3D
30 {
31 
38  struct StokesNorms
39  {
41  StokesNorms(double norm_u, double norm_curl_u, double norm_p, double norm_grad_p):
42  u(norm_u),
43  curl_u(norm_curl_u),
44  p(norm_p),
45  grad_p(norm_grad_p),
46  hcurl_u( std::sqrt( std::pow(norm_u, 2) + std::pow(norm_curl_u, 2) ) ),
47  hgrad_p( std::sqrt( std::pow(norm_p, 2) + std::pow(norm_grad_p, 2) ) )
48  {
49  // do nothing
50  };
51 
52  double u;
53  double curl_u;
54  double p;
55  double grad_p;
56  double hcurl_u;
57  double hgrad_p;
58 
59  };
60 
62  struct Stokes
63  {
64  typedef Eigen::SparseMatrix<double> SystemMatrixType;
65 
66  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> ForcingTermType;
67  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> VelocityType;
68  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> VorticityType;
69  typedef std::function<double(const Eigen::Vector3d &)> PressureType;
70  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> PressureGradientType;
72 
74  Stokes(
75  const VEMCore & vemcore,
76  bool use_threads,
77  std::ostream & output = std::cout
78  );
79 
82  const ForcingTermType & f,
83  const ForcingTermType & curl_f,
84  const VelocityType & u,
85  const VorticityType & omega,
86  const ViscosityType & nu
87  );
88 
90  inline size_t dimension() const
91  {
92  return m_vcurl.dimension() + m_vgrad.dimension();
93  }
94 
96  inline size_t nbSCDOFs_u() const
97  {
98  return m_vemcore.mesh().n_cells() * m_nloc_sc_u;
99  }
100 
102  inline size_t nbSCDOFs_p() const
103  {
104  return m_vemcore.mesh().n_cells() * m_nloc_sc_p;
105  }
106 
108  inline size_t nbSCDOFs() const
109  {
110  return m_vemcore.mesh().n_cells() * (m_nloc_sc_u + m_nloc_sc_p);
111  }
112 
114  inline size_t sizeSystem() const
115  {
116  return dimension() + 1 - nbSCDOFs();
117  }
119  inline const VGrad & vGrad() const
120  {
121  return m_vgrad;
122  }
123 
125  inline const VCurl & vCurl() const
126  {
127  return m_vcurl;
128  }
129 
131  inline const VDiv & vDiv() const
132  {
133  return m_vdiv;
134  }
135 
137  inline const SystemMatrixType & systemMatrix() const {
138  return m_A;
139  }
140 
143  return m_A;
144  }
145 
147  inline const Eigen::VectorXd & systemVector() const {
148  return m_b;
149  }
150 
152  inline Eigen::VectorXd & systemVector() {
153  return m_b;
154  }
155 
157  inline const SystemMatrixType & scMatrix() const {
158  return m_sc_A;
159  }
160 
162  inline Eigen::VectorXd & scVector() {
163  return m_sc_b;
164  }
165 
167  inline const double & stabilizationParameter() const {
168  return m_stab_par;
169  }
170 
172  inline double & stabilizationParameter() {
173  return m_stab_par;
174  }
175 
177  std::vector<StokesNorms> computeStokesNorms(
178  const std::vector<Eigen::VectorXd> & list_dofs
179  ) const;
180 
182  std::pair<StokesNorms, StokesNorms> computeContinuousErrorsNorms(
183  const Eigen::VectorXd & v,
184  const VelocityType & u,
185  const VorticityType & curl_u,
186  const PressureType & p,
187  const PressureGradientType & grad_p
188  ) const;
189 
190 
191  private:
193  std::pair<Eigen::MatrixXd, Eigen::VectorXd>
194  _compute_local_contribution(
195  size_t iT,
196  const Eigen::VectorXd & interp_f,
197  const ViscosityType & nu
198  );
199 
201  LocalStaticCondensation _compute_static_condensation(const size_t & iT) const;
202 
204  void _assemble_local_contribution(
205  size_t iT,
206  const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT,
207  std::list<Eigen::Triplet<double> > & A1,
208  Eigen::VectorXd & b1,
209  std::list<Eigen::Triplet<double> > & A2,
210  Eigen::VectorXd & b2
211  );
212 
213  const VEMCore & m_vemcore;
214  bool m_use_threads;
215  std::ostream & m_output;
216  VGrad m_vgrad;
217  VCurl m_vcurl;
218  VDiv m_vdiv;
219  const size_t m_nloc_sc_u; // Nb of statically condensed DOFs for velocity in each cell (cell unknowns)
220  const size_t m_nloc_sc_p; // Nb of statically condensed DOFs for pressure in each cell (cell unknowns)
221  SystemMatrixType m_A; // Matrix and RHS of statically condensed system
222  Eigen::VectorXd m_b;
223  SystemMatrixType m_sc_A; // Static condensation operator and RHS (to recover statically condensed DOFs)
224  Eigen::VectorXd m_sc_b;
225  double m_stab_par;
226  };
227 
228  //------------------------------------------------------------------------------
229  // Exact solutions
230  //------------------------------------------------------------------------------
231 
232  static const double PI = boost::math::constants::pi<double>();
233  using std::sin;
234  using std::cos;
235  using std::pow;
236 
237  //------------------------------------------------------------------------------
238  // Trigonometric solution
239  //------------------------------------------------------------------------------
240 
241  double pressure_scaling = 1.;
242  static Stokes::VelocityType
243  trigonometric_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
244  return Eigen::Vector3d(
245  0.5 * sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2)),
246  0.5 * cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
247  -cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2))
248  );
249  };
250 
251  static Stokes::VorticityType
252  trigonometric_curl_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
253  return 3. * PI * Eigen::Vector3d(
254  cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
255  -sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
256  0.
257  );
258  };
259 
260  static Stokes::PressureType
261  trigonometric_p = [](const Eigen::Vector3d & x) -> double {
262  return pressure_scaling * sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2));
263  };
264 
266  trigonometric_grad_p = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
267  return pressure_scaling * 2. * PI * Eigen::Vector3d(
268  cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
269  sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
270  sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2))
271  );
272  };
273 
275  trigonometric_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
276  return 6. * std::pow(PI, 2) * Eigen::Vector3d(
277  sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2)),
278  cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
279  -2. * cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2))
280  )
282  };
283 
285  trigonometric_curl_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
286  return 36. * std::pow(PI, 3) * Eigen::Vector3d(
287  cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
288  -sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
289  0.
290  );
291  };
292 
293  static Stokes::ViscosityType
295 
296  //------------------------------------------------------------------------------
297  // Linear solution
298  //------------------------------------------------------------------------------
299 
300  static Stokes::VelocityType
301  linear_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
302  return Eigen::Vector3d(x(0), -x(1), 0.);
303  };
304 
305  static Stokes::VorticityType
306  linear_curl_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
307  return Eigen::Vector3d::Zero();
308  };
309 
310  static Stokes::PressureType
311  linear_p = [](const Eigen::Vector3d & x) -> double {
312  return pressure_scaling * (x(0) - x(1));
313  };
314 
316  linear_grad_p = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
317  return pressure_scaling * Eigen::Vector3d(1., -1., 0.);
318  };
319 
321  linear_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
322  return linear_grad_p(x);
323  };
324 
326  linear_curl_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
327  return Eigen::Vector3d::Zero();
328  };
329 
330  static Stokes::ViscosityType
332 
333 
334  //------------------------------------------------------------------------------
335  // Solution derived from a field
336  //------------------------------------------------------------------------------
337  // Setting f(x,y,z)= (sin(pi x) sin(pi y) sin(pi z))^3, we take
338  // U=(-d_y f, d_x f, 0), p=0
339  // Calculations are made symbolically in matlab
340 
341  static Stokes::VelocityType
342  field_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
343  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;
344  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);
345  return Eigen::Vector3d(ux, uy, 0.);
346  };
347 
348  static Stokes::VorticityType
349  field_curl_u = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
350  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;
351  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;
352  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;
353 
354  return Eigen::Vector3d(cux, cuy, cuz);
355  };
356 
357  static Stokes::PressureType
358  field_p = [](const Eigen::Vector3d & x) -> double {
359  return 0.;
360  };
361 
363  field_grad_p = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
364  return Eigen::Vector3d(0., 0., 0.);
365  };
366 
368  field_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
369  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;
370 
371  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;
372 
373  double fz = 0.;
374  return Eigen::Vector3d(fx, fy, fz);
375  };
376 
378  field_curl_f = [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
379  double curl_fx = (PI*PI*PI*PI)*cos(PI*x(0))*pow(cos(PI*x(2)),3.0)*pow(sin(PI*x(0)),2.0)*pow(sin(PI*x(1)),3.0)*-1.8E+1-(PI*PI*PI*PI)*pow(cos(PI*x(0)),3.0)*cos(PI*x(2))*pow(sin(PI*x(1)),3.0)*pow(sin(PI*x(2)),2.0)*1.8E+1+(PI*PI*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)*1.53E+2-(PI*PI*PI*PI)*cos(PI*x(0))*pow(cos(PI*x(1)),2.0)*cos(PI*x(2))*pow(sin(PI*x(0)),2.0)*sin(PI*x(1))*pow(sin(PI*x(2)),2.0)*5.4E+1;
380 ;
381 
382  double curl_fy = (PI*PI*PI*PI)*cos(PI*x(1))*pow(cos(PI*x(2)),3.0)*pow(sin(PI*x(0)),3.0)*pow(sin(PI*x(1)),2.0)*-1.8E+1-(PI*PI*PI*PI)*pow(cos(PI*x(1)),3.0)*cos(PI*x(2))*pow(sin(PI*x(0)),3.0)*pow(sin(PI*x(2)),2.0)*1.8E+1+(PI*PI*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)*1.53E+2-(PI*PI*PI*PI)*pow(cos(PI*x(0)),2.0)*cos(PI*x(1))*cos(PI*x(2))*sin(PI*x(0))*pow(sin(PI*x(1)),2.0)*pow(sin(PI*x(2)),2.0)*5.4E+1;
383 ;
384 
385  double curl_fz = (PI*PI*PI*PI)*pow(sin(PI*x(0)),3.0)*pow(sin(PI*x(1)),3.0)*pow(sin(PI*x(2)),3.0)*7.8E+1-(PI*PI*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)*1.14E+2-(PI*PI*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)*1.14E+2-(PI*PI*PI*PI)*pow(cos(PI*x(2)),2.0)*pow(sin(PI*x(0)),3.0)*pow(sin(PI*x(1)),3.0)*sin(PI*x(2))*3.6E+1+(PI*PI*PI*PI)*pow(cos(PI*x(0)),2.0)*pow(cos(PI*x(1)),2.0)*sin(PI*x(0))*sin(PI*x(1))*pow(sin(PI*x(2)),3.0)*7.2E+1+(PI*PI*PI*PI)*pow(cos(PI*x(0)),2.0)*pow(cos(PI*x(2)),2.0)*sin(PI*x(0))*pow(sin(PI*x(1)),3.0)*sin(PI*x(2))*3.6E+1+(PI*PI*PI*PI)*pow(cos(PI*x(1)),2.0)*pow(cos(PI*x(2)),2.0)*pow(sin(PI*x(0)),3.0)*sin(PI*x(1))*sin(PI*x(2))*3.6E+1;
386 ;
387  return Eigen::Vector3d(curl_fx, curl_fy, curl_fz);
388  };
389 
390  static Stokes::ViscosityType
392 
393  //------------------------------------------------------------------------------
394 
395 
396 } // end of namespace HArDCore3D
397 
398 #endif
Virtual Hcurl space: local operators, L2 product and global interpolator.
Definition: vcurl.hpp:28
Virtual Hdiv space: local operators, L2 product and global interpolator.
Definition: vdiv.hpp:26
Construct all polynomial spaces for the VEM sequence.
Definition: vemcore.hpp:39
Virtual H1 space: local operators, L2 product and global interpolator.
Definition: vgrad.hpp:27
size_t dimension() const
Returns the dimension of the global space (all DOFs for all geometric entities)
Definition: localdofspace.hpp:80
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
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
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
Stokes(const DDRCore &ddrcore, bool use_threads, std::ostream &output=std::cout)
Constructor.
Definition: ddr-stokes.cpp:312
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
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
std::size_t n_cells() const
number of cells in the mesh.
Definition: MeshND.hpp:60
const Mesh & mesh() const
Return a const reference to the mesh.
Definition: vemcore.hpp:117
const VDiv & vDiv() const
Returns the space VDiv.
Definition: vem-stokes.hpp:131
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: vem-stokes.hpp:96
Eigen::SparseMatrix< double > SystemMatrixType
Definition: vem-stokes.hpp:64
size_t nbSCDOFs_p() const
Returns the number of statically condensed DOFs (pressure)
Definition: vem-stokes.hpp:102
const double & stabilizationParameter() const
Returns the stabilization parameter.
Definition: vem-stokes.hpp:167
static Stokes::ForcingTermType linear_curl_f
Definition: vem-stokes.hpp:326
static Stokes::ForcingTermType field_curl_f
Definition: vem-stokes.hpp:378
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> VelocityType
Definition: vem-stokes.hpp:67
const VCurl & vCurl() const
Returns the space VCurl.
Definition: vem-stokes.hpp:125
size_t nbSCDOFs() const
Returns the number of statically condensed DOFs (both velocity and pressure)
Definition: vem-stokes.hpp:108
size_t sizeSystem() const
Returns the size of the statically condensed system (with Lagrange multiplier)
Definition: vem-stokes.hpp:114
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> PressureGradientType
Definition: vem-stokes.hpp:70
Eigen::VectorXd & systemVector()
Returns the linear system right-hand side vector.
Definition: vem-stokes.hpp:152
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> VorticityType
Definition: vem-stokes.hpp:68
std::function< double(const Eigen::Vector3d &)> PressureType
Definition: vem-stokes.hpp:69
double & stabilizationParameter()
Returns the stabilization parameter.
Definition: vem-stokes.hpp:172
SystemMatrixType & systemMatrix()
Returns the linear system matrix.
Definition: vem-stokes.hpp:142
const VGrad & vGrad() const
Returns the space VGrad.
Definition: vem-stokes.hpp:119
const Eigen::VectorXd & systemVector() const
Returns the linear system right-hand side vector.
Definition: vem-stokes.hpp:147
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: vem-stokes.hpp:157
const SystemMatrixType & systemMatrix() const
Returns the linear system matrix.
Definition: vem-stokes.hpp:137
size_t dimension() const
Returns the global problem dimension (without Lagrange multiplier, just velocity & pressure)
Definition: vem-stokes.hpp:90
IntegralWeight ViscosityType
Definition: vem-stokes.hpp:71
StokesNorms(double norm_u, double norm_curl_u, double norm_p, double norm_grad_p)
Constructor.
Definition: vem-stokes.hpp:41
static Stokes::ForcingTermType trigonometric_curl_f
Definition: vem-stokes.hpp:285
Eigen::VectorXd & scVector()
Returns the static condensation rhs.
Definition: vem-stokes.hpp:162
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> ForcingTermType
Definition: vem-stokes.hpp:66
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