HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
Loading...
Searching...
No Matches
Classes | Typedefs | Functions | Variables
SDDR_stokes

Implementation of the serendipity DDR scheme for the Stokes problem in curl-curl formulation. More...

Collaboration diagram for SDDR_stokes:

Classes

struct  HArDCore3D::StokesNorms
 Structure to store norm components (for velocity, its curl, pressure, and its gradient), as well as Hcurl and Hgrad norms of velocity and pressure. More...
 
struct  HArDCore3D::Stokes
 Assemble a Stokes problem. More...
 

Typedefs

typedef Eigen::SparseMatrix< doubleHArDCore3D::Stokes::SystemMatrixType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Stokes::ForcingTermType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Stokes::VelocityType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Stokes::VorticityType
 
typedef std::function< double(const Eigen::Vector3d &)> HArDCore3D::Stokes::PressureType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Stokes::PressureGradientType
 
typedef IntegralWeight HArDCore3D::Stokes::ViscosityType
 

Functions

 HArDCore3D::StokesNorms::StokesNorms (double norm_u, double norm_curl_u, double norm_p, double norm_grad_p)
 Constructor.
 
 HArDCore3D::StokesNorms::StokesNorms ()
 
 HArDCore3D::Stokes::Stokes (const DDRCore &ddrcore, bool use_threads, std::ostream &output=std::cout)
 Constructor.
 
void HArDCore3D::Stokes::assembleLinearSystem (const ForcingTermType &f, const VelocityType &u, const VorticityType &omega, const ViscosityType &nu)
 Assemble the global system

 
size_t HArDCore3D::Stokes::dimensionSpace () const
 Returns the global problem dimension (without Lagrange multiplier, just velocity & pressure)
 
size_t HArDCore3D::Stokes::nbSCDOFs_u () const
 Returns the number of statically condensed DOFs (velocity)
 
size_t HArDCore3D::Stokes::nbSCDOFs_p () const
 Returns the number of statically condensed DOFs (pressure)
 
size_t HArDCore3D::Stokes::nbSCDOFs () const
 Returns the number of statically condensed DOFs (both velocity and pressure)
 
size_t HArDCore3D::Stokes::sizeSystem () const
 Returns the size of the statically condensed system (with Lagrange multiplier)
 
const SXGradHArDCore3D::Stokes::sxGrad () const
 Returns the space SXGrad.
 
const SXCurlHArDCore3D::Stokes::sxCurl () const
 Returns the space SXCurl.
 
const SXDivHArDCore3D::Stokes::sxDiv () const
 Returns the space XDiv.
 
const SystemMatrixTypeHArDCore3D::Stokes::systemMatrix () const
 Returns the linear system matrix.
 
SystemMatrixTypeHArDCore3D::Stokes::systemMatrix ()
 Returns the linear system matrix.
 
const Eigen::VectorXd & HArDCore3D::Stokes::systemVector () const
 Returns the linear system right-hand side vector.
 
Eigen::VectorXd & HArDCore3D::Stokes::systemVector ()
 Returns the linear system right-hand side vector.
 
const SystemMatrixTypeHArDCore3D::Stokes::scMatrix () const
 Returns the static condensation recovery operator.
 
Eigen::VectorXd & HArDCore3D::Stokes::scVector ()
 Returns the static condensation rhs.
 
const doubleHArDCore3D::Stokes::stabilizationParameter () const
 Returns the stabilization parameter.
 
doubleHArDCore3D::Stokes::stabilizationParameter ()
 Returns the stabilization parameter.
 
std::vector< StokesNormsHArDCore3D::Stokes::computeStokesNorms (const std::vector< Eigen::VectorXd > &list_dofs) const
 Compute the discrete L2 norms, for a family of Eigen::VectorXd representing velocities & pressures, of u, curl u, p and grad p.
 
std::pair< StokesNorms, StokesNormsHArDCore3D::Stokes::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 continuous norms (second component)
 

Variables

static const double HArDCore3D::PI = boost::math::constants::pi<double>()
 
static Stokes::VelocityType HArDCore3D::trigonometric_u
 
static Stokes::VorticityType HArDCore3D::trigonometric_curl_u
 
static Stokes::PressureType HArDCore3D::trigonometric_p
 
static Stokes::PressureGradientType HArDCore3D::trigonometric_grad_p
 
static Stokes::ForcingTermType HArDCore3D::trigonometric_f
 
static Stokes::ViscosityType HArDCore3D::trigonometric_nu = Stokes::ViscosityType(1.)
 
static Stokes::VelocityType HArDCore3D::linear_u
 
static Stokes::VorticityType HArDCore3D::linear_curl_u
 
static Stokes::PressureType HArDCore3D::linear_p
 
static Stokes::PressureGradientType HArDCore3D::linear_grad_p
 
static Stokes::ForcingTermType HArDCore3D::linear_f
 
static Stokes::ViscosityType HArDCore3D::linear_nu = Stokes::ViscosityType(1.)
 
static Stokes::VelocityType HArDCore3D::field_u
 
static Stokes::VorticityType HArDCore3D::field_curl_u
 
static Stokes::PressureType HArDCore3D::field_p
 
static Stokes::PressureGradientType HArDCore3D::field_grad_p
 
static Stokes::ForcingTermType HArDCore3D::field_f
 
static Stokes::ViscosityType HArDCore3D::field_nu = Stokes::ViscosityType(1.)
 
static Stokes::VelocityType HArDCore3D::vertical_u
 
static Stokes::VorticityType HArDCore3D::vertical_curl_u
 
static Stokes::PressureType HArDCore3D::vertical_p = trigonometric_p
 
static Stokes::PressureGradientType HArDCore3D::vertical_grad_p = trigonometric_grad_p
 
static Stokes::ForcingTermType HArDCore3D::vertical_curl_u_cross_u
 
static Stokes::ForcingTermType HArDCore3D::vertical_f
 

Detailed Description

Implementation of the serendipity DDR scheme for the Stokes problem in curl-curl formulation.

Typedef Documentation

◆ ForcingTermType

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Stokes::ForcingTermType

◆ PressureGradientType

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Stokes::PressureGradientType

◆ PressureType

typedef std::function<double(const Eigen::Vector3d &)> HArDCore3D::Stokes::PressureType

◆ SystemMatrixType

◆ VelocityType

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Stokes::VelocityType

◆ ViscosityType

◆ VorticityType

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Stokes::VorticityType

Function Documentation

◆ assembleLinearSystem()

void HArDCore3D::Stokes::assembleLinearSystem ( const ForcingTermType f,
const VelocityType u,
const VorticityType omega,
const ViscosityType nu 
)

Assemble the global system

Parameters
fForcing term
uBoundary condition on velocity
omegaBoundary condition on vorticity
nuViscosity

◆ computeContinuousErrorsNorms()

std::pair< StokesNorms, StokesNorms > HArDCore3D::Stokes::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 continuous norms (second component)

Parameters
vThe vector, contains both velocity and pressure
uContinuous velocity
curl_uCurl of continuous velocity
pContinuous pressure
grad_pGrad of continuous pressure

◆ computeStokesNorms()

std::vector< StokesNorms > HArDCore3D::Stokes::computeStokesNorms ( const std::vector< Eigen::VectorXd > &  list_dofs) const

Compute the discrete L2 norms, for a family of Eigen::VectorXd representing velocities & pressures, of u, curl u, p and grad p.

Parameters
list_dofsList of vectors representing velocities and pressures

◆ dimensionSpace()

size_t HArDCore3D::Stokes::dimensionSpace ( ) const
inline

Returns the global problem dimension (without Lagrange multiplier, just velocity & pressure)

◆ nbSCDOFs()

size_t HArDCore3D::Stokes::nbSCDOFs ( ) const
inline

Returns the number of statically condensed DOFs (both velocity and pressure)

◆ nbSCDOFs_p()

size_t HArDCore3D::Stokes::nbSCDOFs_p ( ) const
inline

Returns the number of statically condensed DOFs (pressure)

◆ nbSCDOFs_u()

size_t HArDCore3D::Stokes::nbSCDOFs_u ( ) const
inline

Returns the number of statically condensed DOFs (velocity)

◆ scMatrix()

const SystemMatrixType & HArDCore3D::Stokes::scMatrix ( ) const
inline

Returns the static condensation recovery operator.

◆ scVector()

Eigen::VectorXd & HArDCore3D::Stokes::scVector ( )
inline

Returns the static condensation rhs.

◆ sizeSystem()

size_t HArDCore3D::Stokes::sizeSystem ( ) const
inline

Returns the size of the statically condensed system (with Lagrange multiplier)

◆ stabilizationParameter() [1/2]

double & HArDCore3D::Stokes::stabilizationParameter ( )
inline

Returns the stabilization parameter.

◆ stabilizationParameter() [2/2]

const double & HArDCore3D::Stokes::stabilizationParameter ( ) const
inline

Returns the stabilization parameter.

◆ Stokes()

HArDCore3D::Stokes::Stokes ( const DDRCore ddrcore,
bool  use_threads,
std::ostream &  output = std::cout 
)

Constructor.

Parameters
ddrcoreCore for the DDR space sequence
use_threadsTrue for parallel execution, false for sequential execution
outputOutput stream to print status messages

◆ StokesNorms() [1/2]

HArDCore3D::StokesNorms::StokesNorms ( )
inline

◆ StokesNorms() [2/2]

HArDCore3D::StokesNorms::StokesNorms ( double  norm_u,
double  norm_curl_u,
double  norm_p,
double  norm_grad_p 
)
inline

Constructor.

◆ sxCurl()

const SXCurl & HArDCore3D::Stokes::sxCurl ( ) const
inline

Returns the space SXCurl.

◆ sxDiv()

const SXDiv & HArDCore3D::Stokes::sxDiv ( ) const
inline

Returns the space XDiv.

◆ sxGrad()

const SXGrad & HArDCore3D::Stokes::sxGrad ( ) const
inline

Returns the space SXGrad.

◆ systemMatrix() [1/2]

SystemMatrixType & HArDCore3D::Stokes::systemMatrix ( )
inline

Returns the linear system matrix.

◆ systemMatrix() [2/2]

const SystemMatrixType & HArDCore3D::Stokes::systemMatrix ( ) const
inline

Returns the linear system matrix.

◆ systemVector() [1/2]

Eigen::VectorXd & HArDCore3D::Stokes::systemVector ( )
inline

Returns the linear system right-hand side vector.

◆ systemVector() [2/2]

const Eigen::VectorXd & HArDCore3D::Stokes::systemVector ( ) const
inline

Returns the linear system right-hand side vector.

Variable Documentation

◆ field_curl_u

Stokes::VorticityType HArDCore3D::field_curl_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
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;
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;
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;
return Eigen::Vector3d(cux, cuy, cuz);
}

◆ field_f

Stokes::ForcingTermType HArDCore3D::field_f
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
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;
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;
double fz = 0.;
return Eigen::Vector3d(fx, fy, fz);
}

◆ field_grad_p

Stokes::PressureGradientType HArDCore3D::field_grad_p
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0., 0., 0.);
}

◆ field_nu

Stokes::ViscosityType HArDCore3D::field_nu = Stokes::ViscosityType(1.)
static

◆ field_p

Stokes::PressureType HArDCore3D::field_p
static
Initial value:
= [](const Eigen::Vector3d & x) -> double {
return 0.;
}

◆ field_u

Stokes::VelocityType HArDCore3D::field_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
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;
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);
return Eigen::Vector3d(ux, uy, 0.);
}

◆ linear_curl_u

Stokes::VorticityType HArDCore3D::linear_curl_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d::Zero();
}

◆ linear_f

Stokes::ForcingTermType HArDCore3D::linear_f
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return linear_grad_p(x);
}
static Stokes::PressureGradientType linear_grad_p
Definition sddr-stokes.hpp:312

◆ linear_grad_p

Stokes::PressureGradientType HArDCore3D::linear_grad_p
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return pressure_scaling * Eigen::Vector3d(1., -1., 0.);
}
double pressure_scaling
Definition sddr-navier-stokes.hpp:464

◆ linear_nu

Stokes::ViscosityType HArDCore3D::linear_nu = Stokes::ViscosityType(1.)
static

◆ linear_p

Stokes::PressureType HArDCore3D::linear_p
static
Initial value:
= [](const Eigen::Vector3d & x) -> double {
return pressure_scaling * (x(0) - x(1));
}

◆ linear_u

Stokes::VelocityType HArDCore3D::linear_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(x(0), -x(1), 0.);
}

◆ PI

const double HArDCore3D::PI = boost::math::constants::pi<double>()
static

◆ trigonometric_curl_u

Stokes::VorticityType HArDCore3D::trigonometric_curl_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return 3. * PI * Eigen::Vector3d(
cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
-sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
0.
);
}

◆ trigonometric_f

Stokes::ForcingTermType HArDCore3D::trigonometric_f
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return 6. * std::pow(PI, 2) * Eigen::Vector3d(
sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2)),
cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
-2. * cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2))
)
}
static Stokes::PressureGradientType trigonometric_grad_p
Definition sddr-stokes.hpp:271

◆ trigonometric_grad_p

Stokes::PressureGradientType HArDCore3D::trigonometric_grad_p
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return pressure_scaling * 2. * PI * Eigen::Vector3d(
cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2))
);
}

◆ trigonometric_nu

Stokes::ViscosityType HArDCore3D::trigonometric_nu = Stokes::ViscosityType(1.)
static

◆ trigonometric_p

Stokes::PressureType HArDCore3D::trigonometric_p
static
Initial value:
= [](const Eigen::Vector3d & x) -> double {
return pressure_scaling * sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2));
}

◆ trigonometric_u

Stokes::VelocityType HArDCore3D::trigonometric_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
0.5 * sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2)),
0.5 * cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
-cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2))
);
}

◆ vertical_curl_u

Stokes::VorticityType HArDCore3D::vertical_curl_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, -PI*sin(PI*x.x()) + PI*sin(PI*x.y()) );
}

◆ vertical_curl_u_cross_u

Stokes::ForcingTermType HArDCore3D::vertical_curl_u_cross_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d (
-(-PI*sin(PI*x.x()) + PI*sin(PI*x.y()))*cos(PI*x.x()),
(-PI*sin(PI*x.x()) + PI*sin(PI*x.y()))*cos(PI*x.y()),
0.
);
}

◆ vertical_f

Stokes::ForcingTermType HArDCore3D::vertical_f
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(pow(PI, 2)*cos(PI*x.y()), pow(PI, 2)*cos(PI*x.x()), 0.)
}
static Stokes::PressureGradientType vertical_grad_p
Definition sddr-stokes.hpp:391

◆ vertical_grad_p

Stokes::PressureGradientType HArDCore3D::vertical_grad_p = trigonometric_grad_p
static

◆ vertical_p

Stokes::PressureType HArDCore3D::vertical_p = trigonometric_p
static

◆ vertical_u

Stokes::VelocityType HArDCore3D::vertical_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d( cos(PI*x.y()), cos(PI*x.x()), 1.);
}