7 #ifndef SDDR_STOKES_HPP
8 #define SDDR_STOKES_HPP
12 #include <boost/math/constants/constants.hpp>
14 #include <Eigen/Sparse>
15 #include <unsupported/Eigen/SparseExtra>
42 StokesNorms(
double norm_u,
double norm_curl_u,
double norm_p,
double 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) ) )
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;
78 std::ostream & output = std::cout
98 return m_nloc_sc_u.sum();
104 return m_nloc_sc_p.sum();
179 const std::vector<Eigen::VectorXd> & list_dofs
184 const Eigen::VectorXd & v,
194 std::pair<Eigen::MatrixXd, Eigen::VectorXd>
195 _compute_local_contribution(
197 const Eigen::VectorXd & interp_f,
205 void _assemble_local_contribution(
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,
216 std::ostream & m_output;
221 const Eigen::VectorXd m_nloc_sc_u;
222 const Eigen::VectorXd m_nloc_sc_p;
226 Eigen::VectorXd m_sc_b;
234 static const double PI = boost::math::constants::pi<double>();
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))
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)),
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))
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))
294 linear_u = [](
const Eigen::Vector3d & x) -> Eigen::Vector3d {
295 return Eigen::Vector3d(x(0), -x(1), 0.);
300 return Eigen::Vector3d::Zero();
304 linear_p = [](
const Eigen::Vector3d & x) ->
double {
314 linear_f = [](
const Eigen::Vector3d & x) -> Eigen::Vector3d {
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.);
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;
342 return Eigen::Vector3d(cux, cuy, cuz);
346 field_p = [](
const Eigen::Vector3d & x) ->
double {
352 return Eigen::Vector3d(0., 0., 0.);
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;
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;
362 return Eigen::Vector3d(fx, fy, fz);
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