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
HHO_brinkman

Implementation of the HHO scheme for the Stokes and Brinkman equations. More...

Classes

struct  HArDCore3D::BrinkmanParameters
 Structure to store physical parameter (and compute derived parameters) More...
 
struct  HArDCore3D::BrinkmanNorms
 Structure to store norm components (energy for velocity, L2 pressure, and global) More...
 
struct  HArDCore3D::Brinkman
 Structure for Brinkman model. More...
 

Typedefs

typedef IntegralWeight HArDCore3D::BrinkmanParameters::ViscosityType
 
typedef IntegralWeight HArDCore3D::BrinkmanParameters::PermeabilityInvType
 
typedef Eigen::SparseMatrix< doubleHArDCore3D::Brinkman::SystemMatrixType
 
typedef std::function< VectorRd(const VectorRd &)> HArDCore3D::Brinkman::MomentumForcingTermType
 
typedef std::function< double(const VectorRd &)> HArDCore3D::Brinkman::CompressibilityForcingTermType
 
typedef std::function< VectorRd(const VectorRd &)> HArDCore3D::Brinkman::VelocityType
 
typedef std::function< MatrixRd(const VectorRd &)> HArDCore3D::Brinkman::VelocityGradientType
 
typedef std::function< double(const VectorRd &)> HArDCore3D::Brinkman::PressureType
 
typedef std::function< VectorRd(const VectorRd &)> HArDCore3D::Brinkman::PressureGradientType
 

Functions

 HArDCore3D::BrinkmanParameters::BrinkmanParameters (const ViscosityType &_mu, const PermeabilityInvType &_kappainv)
 Constructor.
 
double HArDCore3D::BrinkmanParameters::Cf (const Cell &T) const
 Returns the friction coefficient.
 
 HArDCore3D::BrinkmanNorms::BrinkmanNorms (double norm_u, double norm_p)
 Constructor.
 
 HArDCore3D::Brinkman::Brinkman (const VHHOSpace &vhho_space, const GlobalDOFSpace &p_space, const BoundaryConditions &BC, bool use_threads, std::ostream &output=std::cout)
 Constructor.
 
void HArDCore3D::Brinkman::assembleLinearSystem (const MomentumForcingTermType &f, const CompressibilityForcingTermType &g, const BrinkmanParameters &para, const VelocityType &u, Eigen::VectorXd &UDir)
 Assemble the global system

 
size_t HArDCore3D::Brinkman::nloc_sc_u () const
 Returns the local number of velocity statically condensed DOFs.
 
size_t HArDCore3D::Brinkman::nloc_sc_p () const
 Returns the local number of pressure statically condensed DOFs.
 
size_t HArDCore3D::Brinkman::numSCDOFs_u () const
 Returns the number of velocity statically condensed DOFs.
 
size_t HArDCore3D::Brinkman::numSCDOFs_p () const
 Returns the number of pressure statically condensed DOFs.
 
size_t HArDCore3D::Brinkman::numSCDOFs () const
 Returns the number of statically condensed DOFs.
 
size_t HArDCore3D::Brinkman::numDirDOFs () const
 Returns the number of Dirichlet DOFs.
 
size_t HArDCore3D::Brinkman::dimVelocity () const
 Returns the dimension of velocity space.
 
size_t HArDCore3D::Brinkman::dimPressure () const
 Returns the dimension of pressure space.
 
size_t HArDCore3D::Brinkman::numNonSCDOFs () const
 Returns the number of DOFs after SC and with Lagrange multiplier, but before eliminating Dirichlet DOFs.
 
size_t HArDCore3D::Brinkman::sizeSystem () const
 Returns the size of the final system with Lagrange multiplier, after application of SC and removal of Dirichlet BCs.
 
const VHHOSpaceHArDCore3D::Brinkman::vhhospace () const
 Returns the velocity space.
 
const GlobalDOFSpaceHArDCore3D::Brinkman::pspace () const
 Returns the pressure space.
 
const MeshHArDCore3D::Brinkman::mesh () const
 Returns the mesh.
 
const SystemMatrixTypeHArDCore3D::Brinkman::systemMatrix () const
 Returns the linear system matrix.
 
SystemMatrixTypeHArDCore3D::Brinkman::systemMatrix ()
 Returns the linear system matrix.
 
const Eigen::VectorXd & HArDCore3D::Brinkman::systemVector () const
 Returns the linear system right-hand side vector.
 
Eigen::VectorXd & HArDCore3D::Brinkman::systemVector ()
 Returns the linear system right-hand side vector.
 
const std::pair< double, double > & HArDCore3D::Brinkman::stabilizationParameter () const
 Returns the stabilization parameter (scaling)
 
std::pair< double, double > & HArDCore3D::Brinkman::stabilizationParameter ()
 Returns the stabilization parameter.
 
const SystemMatrixTypeHArDCore3D::Brinkman::scMatrix () const
 Returns the static condensation recovery operator.
 
Eigen::VectorXd & HArDCore3D::Brinkman::scVector ()
 Returns the static condensation rhs.
 
Eigen::VectorXd HArDCore3D::Brinkman::interpolate (const VelocityType &u, const PressureType &p, const int doe_cell=-1, const int doe_face=-1) const
 Interpolates velocity and pressure.
 
std::vector< BrinkmanNormsHArDCore3D::Brinkman::computeBrinkmanNorms (const std::vector< Eigen::VectorXd > &list_dofs, const BrinkmanParameters &para) const
 Compute the discrete energy norm of a family of vectors representing the dofs.
 
std::vector< doubleHArDCore3D::Brinkman::pressureVertexValues (const Eigen::VectorXd &p) const
 Create vertex values for the pressure (from the element values), for plotting.
 
std::pair< double, doubleHArDCore3D::Brinkman::computeConditionNum () const
 Computes the condition number of the matrix.
 
std::pair< double, doubleHArDCore3D::Brinkman::computeFlux (const Eigen::VectorXd &u, const std::pair< std::vector< VectorRd >, VectorRd > &surf) const
 Compute the velocity flux (integral of u.n) across a given surface, and the area of the surface, assumed to be star-shaped with respect the barycenter of its vertices.
 

Variables

ViscosityType HArDCore3D::BrinkmanParameters::mu
 
PermeabilityInvType HArDCore3D::BrinkmanParameters::kappainv
 
double HArDCore3D::BrinkmanNorms::u
 Energy norm of velocity.
 
double HArDCore3D::BrinkmanNorms::p
 L2 norm of p.
 
double HArDCore3D::BrinkmanNorms::energy
 Global energy norm

 
static const double HArDCore3D::PI = boost::math::constants::pi<double>()
 
double HArDCore3D::scaling_mu = 1.
 
double HArDCore3D::scaling_kappainv = 1.
 
double constexpr HArDCore3D::eps =1e-14
 
static BrinkmanParameters::ViscosityType HArDCore3D::linear_mu = BrinkmanParameters::ViscosityType(1.)
 
static BrinkmanParameters::PermeabilityInvType HArDCore3D::linear_kappainv = BrinkmanParameters::PermeabilityInvType(1.)
 
static const MatrixRd HArDCore3D::mat_u = (MatrixRd() << 1.,2.,3., 0.,-1.,1. , 0.,1.,1.).finished()
 
static const VectorRd HArDCore3D::vec_p = VectorRd(-1., 2., -5.)
 
static Brinkman::VelocityType HArDCore3D::linear_u
 
static Brinkman::PressureType HArDCore3D::linear_p
 
static Brinkman::VelocityGradientType HArDCore3D::linear_gradu
 
static Brinkman::PressureGradientType HArDCore3D::linear_gradp
 
static Brinkman::MomentumForcingTermType HArDCore3D::linear_f
 
static Brinkman::CompressibilityForcingTermType HArDCore3D::linear_g = [](const VectorRd & x)->double { return linear_gradu(x).trace();}
 
static BrinkmanParameters::ViscosityType HArDCore3D::quadratic_mu = BrinkmanParameters::ViscosityType(1.)
 
static BrinkmanParameters::PermeabilityInvType HArDCore3D::quadratic_kappainv = BrinkmanParameters::PermeabilityInvType(1.)
 
static Brinkman::VelocityType HArDCore3D::quadratic_u
 
static Brinkman::VelocityGradientType HArDCore3D::quadratic_gradu
 
static Brinkman::PressureType HArDCore3D::quadratic_p
 
static Brinkman::PressureGradientType HArDCore3D::quadratic_gradp
 
static Brinkman::MomentumForcingTermType HArDCore3D::quadratic_f
 
static Brinkman::CompressibilityForcingTermType HArDCore3D::quadratic_g = [](const VectorRd & x)->double { return quadratic_gradu(x).trace();}
 
static BrinkmanParameters::ViscosityType HArDCore3D::cubic_mu = BrinkmanParameters::ViscosityType(1.)
 
static BrinkmanParameters::PermeabilityInvType HArDCore3D::cubic_kappainv = BrinkmanParameters::PermeabilityInvType(1.)
 
static Brinkman::VelocityType HArDCore3D::cubic_u
 
static Brinkman::VelocityGradientType HArDCore3D::cubic_gradu
 
static Brinkman::PressureType HArDCore3D::cubic_p
 
static Brinkman::PressureGradientType HArDCore3D::cubic_gradp
 
static Brinkman::MomentumForcingTermType HArDCore3D::cubic_f
 
static Brinkman::CompressibilityForcingTermType HArDCore3D::cubic_g = [](const VectorRd & x)->double { return cubic_gradu(x).trace();}
 
static BrinkmanParameters::ViscosityType HArDCore3D::trigonometric_mu = BrinkmanParameters::ViscosityType(1.)
 
static BrinkmanParameters::PermeabilityInvType HArDCore3D::trigonometric_kappainv = BrinkmanParameters::PermeabilityInvType(1.)
 
static Brinkman::VelocityType HArDCore3D::trigonometric_u
 
static Brinkman::VelocityGradientType HArDCore3D::trigonometric_gradu
 
static Brinkman::PressureType HArDCore3D::trigonometric_p
 
static Brinkman::PressureGradientType HArDCore3D::trigonometric_gradp
 
static Brinkman::MomentumForcingTermType HArDCore3D::trigonometric_f
 
static Brinkman::CompressibilityForcingTermType HArDCore3D::trigonometric_g = [](const VectorRd & x)->double { return trigonometric_gradu(x).trace();}
 
static BrinkmanParameters::ViscosityType HArDCore3D::regimes_mu = BrinkmanParameters::ViscosityType(1.)
 
static BrinkmanParameters::PermeabilityInvType HArDCore3D::regimes_kappainv = BrinkmanParameters::PermeabilityInvType(1.)
 
static Brinkman::PressureType HArDCore3D::regimes_p = trigonometric_p
 
static Brinkman::PressureGradientType HArDCore3D::regimes_gradp = trigonometric_gradp
 
static Brinkman::VelocityType HArDCore3D::regimes_uD
 
static Brinkman::VelocityGradientType HArDCore3D::regimes_graduD
 
static Brinkman::VelocityType HArDCore3D::regimes_uS = trigonometric_u
 
static Brinkman::VelocityGradientType HArDCore3D::regimes_graduS = trigonometric_gradu
 
static std::function< double(const VectorRd &)> HArDCore3D::XiS
 
static Brinkman::VelocityType HArDCore3D::regimes_u = [](const VectorRd & x) -> VectorRd { return XiS(x) * regimes_uS(x) + (1.-XiS(x)) * regimes_uD(x);}
 
static Brinkman::VelocityGradientType HArDCore3D::regimes_gradu = [](const VectorRd & x) -> MatrixRd { return XiS(x) * regimes_graduS(x) + (1.-XiS(x)) * regimes_graduD(x);}
 
static Brinkman::MomentumForcingTermType HArDCore3D::regimes_f
 
static Brinkman::CompressibilityForcingTermType HArDCore3D::regimes_g = [](const VectorRd & x)->double { return regimes_gradu(x).trace() ;}
 
static std::function< double(const VectorRd &)> HArDCore3D::XiV = [](const VectorRd & x)->double { return (x.z() > 2*std::abs(x.x()-1.)+0.5-eps); }
 
static BrinkmanParameters::ViscosityType HArDCore3D::vcracked_mu = BrinkmanParameters::ViscosityType([](const VectorRd & x)->double { return 0.01 * XiV(x);})
 
static BrinkmanParameters::PermeabilityInvType HArDCore3D::vcracked_kappainv = BrinkmanParameters::PermeabilityInvType([](const VectorRd & x)->double { return 1e-5*(1.-XiV(x)) + 0.001*XiV(x);})
 
static Brinkman::VelocityType HArDCore3D::vcracked_u
 
static Brinkman::VelocityGradientType HArDCore3D::vcracked_gradu
 
static Brinkman::PressureType HArDCore3D::vcracked_p
 
static Brinkman::PressureGradientType HArDCore3D::vcracked_gradp
 
static Brinkman::MomentumForcingTermType HArDCore3D::vcracked_f
 
static Brinkman::CompressibilityForcingTermType HArDCore3D::vcracked_g
 
static std::function< double(const VectorRd &)> HArDCore3D::XiCavity
 
static std::function< double(const VectorRd &)> HArDCore3D::XiWedge
 
static std::function< double(const VectorRd &)> HArDCore3D::XiBox = [](const VectorRd & x)->double { return 1. - XiCavity(x) - XiWedge(x); }
 
constexpr double HArDCore3D::viscosity_in_cavity = 0.01
 
static BrinkmanParameters::ViscosityType HArDCore3D::cavity_mu = viscosity_in_cavity * BrinkmanParameters::ViscosityType(XiCavity)
 
constexpr double HArDCore3D::permeability_in_wedge = 1e-2
 
constexpr double HArDCore3D::permeability_in_box = 1e-7
 
static BrinkmanParameters::PermeabilityInvType HArDCore3D::cavity_kappainv
 
static Brinkman::VelocityType HArDCore3D::cavity_u
 
static Brinkman::VelocityGradientType HArDCore3D::cavity_gradu
 
static Brinkman::PressureType HArDCore3D::cavity_p
 
static Brinkman::PressureGradientType HArDCore3D::cavity_gradp
 
static Brinkman::MomentumForcingTermType HArDCore3D::cavity_f
 
static Brinkman::CompressibilityForcingTermType HArDCore3D::cavity_g
 
static std::function< double(const VectorRd &)> HArDCore3D::Xi_S
 
static std::function< double(const VectorRd &)> HArDCore3D::Xi_D
 
constexpr double HArDCore3D::viscosity_S = 1.
 
constexpr double HArDCore3D::viscosity_D = 0.
 
static BrinkmanParameters::ViscosityType HArDCore3D::tworegions_mu
 
constexpr double HArDCore3D::permeability_S = 1e-7
 
constexpr double HArDCore3D::permeability_D = 1e-2
 
static BrinkmanParameters::PermeabilityInvType HArDCore3D::tworegions_kappainv
 
static std::function< VectorRd(const VectorRd &)> HArDCore3D::alpha
 
static std::function< VectorRd(const VectorRd &)> HArDCore3D::beta_S
 
static std::function< VectorRd(const VectorRd &)> HArDCore3D::beta_D
 
static Brinkman::VelocityType HArDCore3D::tworegions_u
 
static Brinkman::VelocityGradientType HArDCore3D::grad_alpha
 
static Brinkman::VelocityGradientType HArDCore3D::grad_beta_S
 
static Brinkman::VelocityGradientType HArDCore3D::grad_beta_D
 
static Brinkman::VelocityGradientType HArDCore3D::tworegions_gradu
 
static Brinkman::PressureType HArDCore3D::tworegions_p = trigonometric_p
 
static Brinkman::PressureGradientType HArDCore3D::tworegions_gradp = trigonometric_gradp
 
static Brinkman::MomentumForcingTermType HArDCore3D::DivGrad_alpha
 
static Brinkman::MomentumForcingTermType HArDCore3D::DivGrad_beta_S
 
static Brinkman::MomentumForcingTermType HArDCore3D::DivGrad_beta_D
 
static Brinkman::MomentumForcingTermType HArDCore3D::tworegions_f
 
static Brinkman::CompressibilityForcingTermType HArDCore3D::tworegions_g = [](const VectorRd & x)->double { return tworegions_gradu(x).trace();}
 

Detailed Description

Implementation of the HHO scheme for the Stokes and Brinkman equations.

Typedef Documentation

◆ CompressibilityForcingTermType

◆ MomentumForcingTermType

typedef std::function<VectorRd(const VectorRd &)> HArDCore3D::Brinkman::MomentumForcingTermType

◆ PermeabilityInvType

◆ PressureGradientType

typedef std::function<VectorRd(const VectorRd &)> HArDCore3D::Brinkman::PressureGradientType

◆ PressureType

◆ SystemMatrixType

◆ VelocityGradientType

◆ VelocityType

typedef std::function<VectorRd(const VectorRd &)> HArDCore3D::Brinkman::VelocityType

◆ ViscosityType

Function Documentation

◆ assembleLinearSystem()

void Brinkman::assembleLinearSystem ( const MomentumForcingTermType f,
const CompressibilityForcingTermType g,
const BrinkmanParameters para,
const VelocityType u,
Eigen::VectorXd &  UDir 
)

Assemble the global system

Parameters
fForcing term
gForcing term for the divergence equation
paraViscosity and permeability
uExact solution for boundary condition
UDirVector filled in by Dirichlets BCs

◆ Brinkman()

Brinkman::Brinkman ( const VHHOSpace vhho_space,
const GlobalDOFSpace p_space,
const BoundaryConditions BC,
bool  use_threads,
std::ostream &  output = std::cout 
)

Constructor.

Parameters
vhho_spaceVelocity space and operators
p_spaceSpace for the pressure
BCBoundary conditions
use_threadsTrue for parallel execution, false for sequential execution
outputOutput stream to print status messages

◆ BrinkmanNorms()

HArDCore3D::BrinkmanNorms::BrinkmanNorms ( double  norm_u,
double  norm_p 
)
inline

Constructor.

◆ BrinkmanParameters()

HArDCore3D::BrinkmanParameters::BrinkmanParameters ( const ViscosityType _mu,
const PermeabilityInvType _kappainv 
)
inline

Constructor.

Parameters
_muViscosity
_kappainvInverse of permeability

◆ Cf()

double HArDCore3D::BrinkmanParameters::Cf ( const Cell &  T) const
inline

Returns the friction coefficient.

◆ computeBrinkmanNorms()

std::vector< BrinkmanNorms > Brinkman::computeBrinkmanNorms ( const std::vector< Eigen::VectorXd > &  list_dofs,
const BrinkmanParameters para 
) const

Compute the discrete energy norm of a family of vectors representing the dofs.

Parameters
list_dofsThe list of vectors representing the dofs
paraViscosity and permeability

◆ computeConditionNum()

std::pair< double, double > Brinkman::computeConditionNum ( ) const

Computes the condition number of the matrix.

◆ computeFlux()

std::pair< double, double > Brinkman::computeFlux ( const Eigen::VectorXd &  u,
const std::pair< std::vector< VectorRd >, VectorRd > &  surf 
) const

Compute the velocity flux (integral of u.n) across a given surface, and the area of the surface, assumed to be star-shaped with respect the barycenter of its vertices.

Parameters
uVector of velocity DOFs
surfSurface described by coplanar vertices (in order around the boundary) and one unit normal vector to indicate the flux to compute

◆ dimPressure()

size_t HArDCore3D::Brinkman::dimPressure ( ) const
inline

Returns the dimension of pressure space.

◆ dimVelocity()

size_t HArDCore3D::Brinkman::dimVelocity ( ) const
inline

Returns the dimension of velocity space.

◆ interpolate()

Eigen::VectorXd Brinkman::interpolate ( const VelocityType u,
const PressureType p,
const int  doe_cell = -1,
const int  doe_face = -1 
) const

Interpolates velocity and pressure.

Parameters
uContinuous velocity to interpolate
pContinuous pressure to interpolate
doe_cellThe optional degre of cell quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used.
doe_faceThe optional degre of face quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used.

◆ mesh()

const Mesh & HArDCore3D::Brinkman::mesh ( ) const
inline

Returns the mesh.

◆ nloc_sc_p()

size_t HArDCore3D::Brinkman::nloc_sc_p ( ) const
inline

Returns the local number of pressure statically condensed DOFs.

◆ nloc_sc_u()

size_t HArDCore3D::Brinkman::nloc_sc_u ( ) const
inline

Returns the local number of velocity statically condensed DOFs.

◆ numDirDOFs()

size_t HArDCore3D::Brinkman::numDirDOFs ( ) const
inline

Returns the number of Dirichlet DOFs.

◆ numNonSCDOFs()

size_t HArDCore3D::Brinkman::numNonSCDOFs ( ) const
inline

Returns the number of DOFs after SC and with Lagrange multiplier, but before eliminating Dirichlet DOFs.

◆ numSCDOFs()

size_t HArDCore3D::Brinkman::numSCDOFs ( ) const
inline

Returns the number of statically condensed DOFs.

◆ numSCDOFs_p()

size_t HArDCore3D::Brinkman::numSCDOFs_p ( ) const
inline

Returns the number of pressure statically condensed DOFs.

◆ numSCDOFs_u()

size_t HArDCore3D::Brinkman::numSCDOFs_u ( ) const
inline

Returns the number of velocity statically condensed DOFs.

◆ pressureVertexValues()

std::vector< double > Brinkman::pressureVertexValues ( const Eigen::VectorXd &  p) const

Create vertex values for the pressure (from the element values), for plotting.

Parameters
pVector of pressure DOFs

◆ pspace()

const GlobalDOFSpace & HArDCore3D::Brinkman::pspace ( ) const
inline

Returns the pressure space.

◆ scMatrix()

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

Returns the static condensation recovery operator.

◆ scVector()

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

Returns the static condensation rhs.

◆ sizeSystem()

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

Returns the size of the final system with Lagrange multiplier, after application of SC and removal of Dirichlet BCs.

◆ stabilizationParameter() [1/2]

std::pair< double, double > & HArDCore3D::Brinkman::stabilizationParameter ( )
inline

Returns the stabilization parameter.

◆ stabilizationParameter() [2/2]

const std::pair< double, double > & HArDCore3D::Brinkman::stabilizationParameter ( ) const
inline

Returns the stabilization parameter (scaling)

◆ systemMatrix() [1/2]

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

Returns the linear system matrix.

◆ systemMatrix() [2/2]

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

Returns the linear system matrix.

◆ systemVector() [1/2]

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

Returns the linear system right-hand side vector.

◆ systemVector() [2/2]

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

Returns the linear system right-hand side vector.

◆ vhhospace()

const VHHOSpace & HArDCore3D::Brinkman::vhhospace ( ) const
inline

Returns the velocity space.

Variable Documentation

◆ alpha

std::function<VectorRd(const VectorRd &)> HArDCore3D::alpha
static
Initial value:
= [](const VectorRd & x)->VectorRd {
return VectorRd(
exp(-x.y()-x.z()),
sin(PI * x.y()) * sin(PI * x.z()),
x.y() * x.z()
);
}
@ Matrix
Definition basis.hpp:67
static const double PI
Definition ddr-magnetostatics.hpp:187

◆ beta_D

std::function<VectorRd(const VectorRd &)> HArDCore3D::beta_D
static
Initial value:
= [](const VectorRd & x)->VectorRd {
return cos(PI * x.x()) * (x.x() - 0.5)
* VectorRd(
sin(PI * x.y()) + cos(PI * x.z()),
std::pow(x.z(), 3),
std::pow(x.y(), 2) * std::pow(x.z(), 2)
);
}

◆ beta_S

std::function<VectorRd(const VectorRd &)> HArDCore3D::beta_S
static
Initial value:
= [](const VectorRd & x)->VectorRd {
return cos(PI * x.x()) * (x.x() - 0.5)
* VectorRd(
x.y() + x.z(),
x.y() + cos(PI * x.z()),
sin(PI * x.y())
);
}

◆ cavity_f

Brinkman::MomentumForcingTermType HArDCore3D::cavity_f
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return VectorRd(0., 0., -.98);
}

◆ cavity_g

Brinkman::CompressibilityForcingTermType HArDCore3D::cavity_g
static
Initial value:
= [](const VectorRd & x) -> double {
return 0.;
}

◆ cavity_gradp

Brinkman::PressureGradientType HArDCore3D::cavity_gradp
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return VectorRd::Zero();
}

◆ cavity_gradu

Brinkman::VelocityGradientType HArDCore3D::cavity_gradu
static
Initial value:
= [](const VectorRd & x) -> MatrixRd {
return MatrixRd::Zero();
}
Eigen::Matrix3d MatrixRd
Definition basis.hpp:51

◆ cavity_kappainv

BrinkmanParameters::PermeabilityInvType HArDCore3D::cavity_kappainv
static
Initial value:
static std::function< double(const VectorRd &)> XiBox
Definition hho-brinkman.hpp:708
constexpr double permeability_in_wedge
Definition hho-brinkman.hpp:715
static std::function< double(const VectorRd &)> XiWedge
Definition hho-brinkman.hpp:704
IntegralWeight ViscosityType
Definition hho-brinkman.hpp:61
constexpr double permeability_in_box
Definition hho-brinkman.hpp:716

◆ cavity_mu

◆ cavity_p

Brinkman::PressureType HArDCore3D::cavity_p
static
Initial value:
= [](const VectorRd & x) -> double {
return 0.;
}

◆ cavity_u

Brinkman::VelocityType HArDCore3D::cavity_u
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return XiCavity(x) * VectorRd(x.x()*(1.-x.x()), 0., 0.);
}
static std::function< double(const VectorRd &)> XiCavity
Definition hho-brinkman.hpp:700

◆ cubic_f

Brinkman::MomentumForcingTermType HArDCore3D::cubic_f
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return - scaling_mu * VectorRd( 6 * x(0) + 2 * x(1), 2., 6 * x(2)) + scaling_kappainv * cubic_u(x) + cubic_gradp(x);
}
static Brinkman::PressureGradientType cubic_gradp
Definition hho-brinkman.hpp:480
static Brinkman::VelocityType cubic_u
Definition hho-brinkman.hpp:450
double scaling_mu
Definition hho-brinkman.hpp:341
double scaling_kappainv
Definition hho-brinkman.hpp:342

◆ cubic_g

Brinkman::CompressibilityForcingTermType HArDCore3D::cubic_g = [](const VectorRd & x)->double { return cubic_gradu(x).trace();}
static

◆ cubic_gradp

Brinkman::PressureGradientType HArDCore3D::cubic_gradp
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return VectorRd( 3 * std::pow(x(0), 2),
x(2),
x(1)
);
}

◆ cubic_gradu

Brinkman::VelocityGradientType HArDCore3D::cubic_gradu
static
Initial value:
= [](const VectorRd & x) -> MatrixRd {
MatrixRd M = MatrixRd::Zero();
M.row(0) << 3. * std::pow(x(0), 2) + 2 * x(0) * x(1),
std::pow(x(0), 2),
0.;
M.row(1) << x(1) * x(2),
2 * x(1) + x(0) * x(2),
x(0) * x(1);
M.row(2) << 0.,
0.,
3 * std::pow(x(2), 2);
return M;
}

◆ cubic_kappainv

◆ cubic_mu

◆ cubic_p

Brinkman::PressureType HArDCore3D::cubic_p
static
Initial value:
= [](const VectorRd & x) -> double {
return std::pow(x(0), 3) + x(1)*x(2) - 1./4. - 1./4.;
}

◆ cubic_u

Brinkman::VelocityType HArDCore3D::cubic_u
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return VectorRd(
std::pow(x(0), 3) + std::pow(x(0), 2)*x(1),
std::pow(x(1), 2) + x(0) * x(1) * x(2),
std::pow(x(2), 3)
);
}

◆ DivGrad_alpha

Brinkman::MomentumForcingTermType HArDCore3D::DivGrad_alpha
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return VectorRd(
2*exp(-x.y() - x.z()),
-2*pow(PI, 2)*sin(PI*x.y())*sin(PI*x.z()),
0
);
}

◆ DivGrad_beta_D

Brinkman::MomentumForcingTermType HArDCore3D::DivGrad_beta_D
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return VectorRd(
-pow(PI, 2)*(x.x() - 0.5)*(sin(PI*x.y()) + cos(PI*x.z()))*cos(PI*x.x()) - pow(PI, 2)*(x.x() - 0.5)*sin(PI*x.y())*cos(PI*x.x()) - pow(PI, 2)*(x.x() - 0.5)*cos(PI*x.x())*cos(PI*x.z()) - 2*PI*(sin(PI*x.y()) + cos(PI*x.z()))*sin(PI*x.x()),
-pow(PI, 2)*pow(x.z(), 3)*(x.x() - 0.5)*cos(PI*x.x()) - 2*PI*pow(x.z(), 3)*sin(PI*x.x()) + 6*x.z()*(x.x() - 0.5)*cos(PI*x.x()),
-pow(PI, 2)*pow(x.y(), 2)*pow(x.z(), 2)*(x.x() - 0.5)*cos(PI*x.x()) - 2*PI*pow(x.y(), 2)*pow(x.z(), 2)*sin(PI*x.x()) + 2*pow(x.y(), 2)*(x.x() - 0.5)*cos(PI*x.x()) + 2*pow(x.z(), 2)*(x.x() - 0.5)*cos(PI*x.x())
);
}

◆ DivGrad_beta_S

Brinkman::MomentumForcingTermType HArDCore3D::DivGrad_beta_S
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return VectorRd(
-pow(PI, 2)*(x.x() - 0.5)*(x.y() + x.z())*cos(PI*x.x()) - 2*PI*(x.y() + x.z())*sin(PI*x.x()),
-pow(PI, 2)*(x.x() - 0.5)*(x.y() + cos(PI*x.z()))*cos(PI*x.x()) - pow(PI, 2)*(x.x() - 0.5)*cos(PI*x.x())*cos(PI*x.z()) - 2*PI*(x.y() + cos(PI*x.z()))*sin(PI*x.x()),
-2*pow(PI, 2)*(x.x() - 0.5)*sin(PI*x.y())*cos(PI*x.x()) - 2*PI*sin(PI*x.x())*sin(PI*x.y())
);
}

◆ energy

double HArDCore3D::BrinkmanNorms::energy

Global energy norm

◆ eps

double constexpr HArDCore3D::eps =1e-14
constexpr

◆ grad_alpha

Brinkman::VelocityGradientType HArDCore3D::grad_alpha
static
Initial value:
= [](const VectorRd & x) -> MatrixRd {
MatrixRd M = MatrixRd::Zero();
M.row(0) << 0, -exp(-x.y() - x.z()), -exp(-x.y() - x.z());
M.row(1) << 0, PI*sin(PI*x.z())*cos(PI*x.y()), PI*sin(PI*x.y())*cos(PI*x.z());
M.row(2) << 0, x.z(), x.y();
return M;
}

◆ grad_beta_D

Brinkman::VelocityGradientType HArDCore3D::grad_beta_D
static
Initial value:
= [](const VectorRd & x) -> MatrixRd {
MatrixRd M = MatrixRd::Zero();
M.row(0) << -PI*(x.x() - 0.5)*(sin(PI*x.y()) + cos(PI*x.z()))*sin(PI*x.x()) + (sin(PI*x.y()) + cos(PI*x.z()))*cos(PI*x.x()),
PI*(x.x() - 0.5)*cos(PI*x.x())*cos(PI*x.y()),
-PI*(x.x() - 0.5)*sin(PI*x.z())*cos(PI*x.x());
M.row(1) << -PI*pow(x.z(), 3)*(x.x() - 0.5)*sin(PI*x.x()) + pow(x.z(), 3)*cos(PI*x.x()),
0,
3*pow(x.z(), 2)*(x.x() - 0.5)*cos(PI*x.x());
M.row(2) << -PI*pow(x.y(), 2)*pow(x.z(), 2)*(x.x() - 0.5)*sin(PI*x.x()) + pow(x.y(), 2)*pow(x.z(), 2)*cos(PI*x.x()),
2*x.y()*pow(x.z(), 2)*(x.x() - 0.5)*cos(PI*x.x()),
2*pow(x.y(), 2)*x.z()*(x.x() - 0.5)*cos(PI*x.x());
return M;
}

◆ grad_beta_S

Brinkman::VelocityGradientType HArDCore3D::grad_beta_S
static
Initial value:
= [](const VectorRd & x) -> MatrixRd {
MatrixRd M = MatrixRd::Zero();
M.row(0) << -PI*(x.x() - 0.5)*(x.y() + x.z())*sin(PI*x.x()) + (x.y() + x.z())*cos(PI*x.x()),
(x.x() - 0.5)*cos(PI*x.x()),
(x.x() - 0.5)*cos(PI*x.x());
M.row(1) << -PI*(x.x() - 0.5)*(x.y() + cos(PI*x.z()))*sin(PI*x.x()) + (x.y() + cos(PI*x.z()))*cos(PI*x.x()),
(x.x() - 0.5)*cos(PI*x.x()),
-PI*(x.x() - 0.5)*sin(PI*x.z())*cos(PI*x.x());
M.row(2) << -PI*(x.x() - 0.5)*sin(PI*x.x())*sin(PI*x.y()) + sin(PI*x.y())*cos(PI*x.x()),
PI*(x.x() - 0.5)*cos(PI*x.x())*cos(PI*x.y()),
0;
return M;
}

◆ kappainv

PermeabilityInvType HArDCore3D::BrinkmanParameters::kappainv

◆ linear_f

Brinkman::MomentumForcingTermType HArDCore3D::linear_f
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return VectorRd::Zero() + scaling_kappainv*linear_u(x) + linear_gradp(x);
}
static Magnetostatics::SolutionPotentialType linear_u
Definition ddr-magnetostatics.hpp:214
static Brinkman::PressureGradientType linear_gradp
Definition hho-brinkman.hpp:375

◆ linear_g

Brinkman::CompressibilityForcingTermType HArDCore3D::linear_g = [](const VectorRd & x)->double { return linear_gradu(x).trace();}
static

◆ linear_gradp

Brinkman::PressureGradientType HArDCore3D::linear_gradp
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return vec_p;
}
static const VectorRd vec_p
Definition hho-brinkman.hpp:357

◆ linear_gradu

Brinkman::VelocityGradientType HArDCore3D::linear_gradu
static
Initial value:
= [](const VectorRd & x) -> MatrixRd {
return mat_u;
}
static const MatrixRd mat_u
Definition hho-brinkman.hpp:356

◆ linear_kappainv

◆ linear_mu

◆ linear_p

Brinkman::PressureType HArDCore3D::linear_p
static
Initial value:
= [](const VectorRd & x) -> double {
return vec_p.dot(x - VectorRd(.5,.5,.5));
}

◆ linear_u

Brinkman::VelocityType HArDCore3D::linear_u
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return mat_u * x;
}

◆ mat_u

const MatrixRd HArDCore3D::mat_u = (MatrixRd() << 1.,2.,3., 0.,-1.,1. , 0.,1.,1.).finished()
static

◆ mu

ViscosityType HArDCore3D::BrinkmanParameters::mu

◆ p

double HArDCore3D::BrinkmanNorms::p

L2 norm of p.

◆ permeability_D

constexpr double HArDCore3D::permeability_D = 1e-2
constexpr

◆ permeability_in_box

constexpr double HArDCore3D::permeability_in_box = 1e-7
constexpr

◆ permeability_in_wedge

constexpr double HArDCore3D::permeability_in_wedge = 1e-2
constexpr

◆ permeability_S

constexpr double HArDCore3D::permeability_S = 1e-7
constexpr

◆ PI

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

◆ quadratic_f

Brinkman::MomentumForcingTermType HArDCore3D::quadratic_f
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return - scaling_mu * VectorRd( 2., 2., 2.) + scaling_kappainv * quadratic_u(x) + quadratic_gradp(x);
}
static Brinkman::PressureGradientType quadratic_gradp
Definition hho-brinkman.hpp:426
static Brinkman::VelocityType quadratic_u
Definition hho-brinkman.hpp:396

◆ quadratic_g

Brinkman::CompressibilityForcingTermType HArDCore3D::quadratic_g = [](const VectorRd & x)->double { return quadratic_gradu(x).trace();}
static

◆ quadratic_gradp

Brinkman::PressureGradientType HArDCore3D::quadratic_gradp
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return VectorRd( 2 * x(0),
x(2),
x(1)
);
}

◆ quadratic_gradu

Brinkman::VelocityGradientType HArDCore3D::quadratic_gradu
static
Initial value:
= [](const VectorRd & x) -> MatrixRd {
MatrixRd M = MatrixRd::Zero();
M.row(0) << 2. * x(0),
0.,
0.;
M.row(1) << 0.,
2 * x(1),
0.;
M.row(2) << 0.,
0.,
2 * x(2);
return M;
}

◆ quadratic_kappainv

◆ quadratic_mu

◆ quadratic_p

Brinkman::PressureType HArDCore3D::quadratic_p
static
Initial value:
= [](const VectorRd & x) -> double {
return std::pow(x(0), 2) + x(1)*x(2) - 1./3. - 1./4.;
}

◆ quadratic_u

Brinkman::VelocityType HArDCore3D::quadratic_u
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return VectorRd(
std::pow(x(0), 2),
std::pow(x(1), 2),
std::pow(x(2), 2)
);
}

◆ regimes_f

Brinkman::MomentumForcingTermType HArDCore3D::regimes_f
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
VectorRd fvalue = XiS(x) * scaling_mu * 6. * std::pow(PI, 2) * VectorRd(
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))
)
- scaling_mu * (1-XiS(x))/scaling_kappainv * 12. * std::pow(PI, 2) * regimes_gradp(x);
}
return fvalue;
}
static Brinkman::VelocityType regimes_u
Definition hho-brinkman.hpp:613
double constexpr eps
Definition hho-brinkman.hpp:345
static Brinkman::PressureGradientType regimes_gradp
Definition hho-brinkman.hpp:571
static std::function< double(const VectorRd &)> XiS
Definition hho-brinkman.hpp:608

◆ regimes_g

Brinkman::CompressibilityForcingTermType HArDCore3D::regimes_g = [](const VectorRd & x)->double { return regimes_gradu(x).trace() ;}
static

◆ regimes_gradp

Brinkman::PressureGradientType HArDCore3D::regimes_gradp = trigonometric_gradp
static

◆ regimes_gradu

Brinkman::VelocityGradientType HArDCore3D::regimes_gradu = [](const VectorRd & x) -> MatrixRd { return XiS(x) * regimes_graduS(x) + (1.-XiS(x)) * regimes_graduD(x);}
static

◆ regimes_graduD

Brinkman::VelocityGradientType HArDCore3D::regimes_graduD
static
Initial value:
= [](const VectorRd & x) -> MatrixRd {
MatrixRd M = MatrixRd::Zero();
M.row(0) << sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
-cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
-cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2));
M.row(1) << -cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
-sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2));
M.row(2) << -cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
-sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2)),
sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2));
M *= 4. * std::pow(PI, 2) / scaling_kappainv;
}
return M;
}

◆ regimes_graduS

Brinkman::VelocityGradientType HArDCore3D::regimes_graduS = trigonometric_gradu
static

◆ regimes_kappainv

◆ regimes_mu

◆ regimes_p

Brinkman::PressureType HArDCore3D::regimes_p = trigonometric_p
static

◆ regimes_u

Brinkman::VelocityType HArDCore3D::regimes_u = [](const VectorRd & x) -> VectorRd { return XiS(x) * regimes_uS(x) + (1.-XiS(x)) * regimes_uD(x);}
static

◆ regimes_uD

Brinkman::VelocityType HArDCore3D::regimes_uD
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
VectorRd uDvalue = VectorRd::Zero();
uDvalue = - std::pow(scaling_kappainv, -1) * regimes_gradp(x) ;
}
return uDvalue;
}

◆ regimes_uS

Brinkman::VelocityType HArDCore3D::regimes_uS = trigonometric_u
static

◆ scaling_kappainv

double HArDCore3D::scaling_kappainv = 1.

◆ scaling_mu

double HArDCore3D::scaling_mu = 1.

◆ trigonometric_f

Brinkman::MomentumForcingTermType HArDCore3D::trigonometric_f
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return scaling_mu * 6. * std::pow(PI, 2) * VectorRd(
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 Magnetostatics::SolutionPotentialType trigonometric_u
Definition ddr-magnetostatics.hpp:238
static Brinkman::PressureGradientType trigonometric_gradp
Definition hho-brinkman.hpp:536

◆ trigonometric_g

Brinkman::CompressibilityForcingTermType HArDCore3D::trigonometric_g = [](const VectorRd & x)->double { return trigonometric_gradu(x).trace();}
static

◆ trigonometric_gradp

Brinkman::PressureGradientType HArDCore3D::trigonometric_gradp
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return pressure_scaling * 2. * PI * VectorRd(
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))
);
}
double pressure_scaling
Definition sddr-navier-stokes.hpp:464

◆ trigonometric_gradu

Brinkman::VelocityGradientType HArDCore3D::trigonometric_gradu
static
Initial value:
= [](const VectorRd & x) -> MatrixRd {
MatrixRd M = MatrixRd::Zero();
M.row(0) << PI * cos(2 * PI * x(0)) * cos(2 * PI * x(1)) * cos(2 * PI * x(2)),
-PI * sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
-PI * sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2));
M.row(1) << - PI * sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
PI * cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2)),
-PI * cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2));
M.row(2) << 2 * PI * sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
2 * PI * cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
- 2 * PI * cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2));
return M;
}

◆ trigonometric_kappainv

◆ trigonometric_mu

BrinkmanParameters::ViscosityType HArDCore3D::trigonometric_mu = BrinkmanParameters::ViscosityType(1.)
static

◆ trigonometric_p

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

◆ trigonometric_u

Brinkman::VelocityType HArDCore3D::trigonometric_u
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return VectorRd(
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))
);
}

◆ tworegions_f

Brinkman::MomentumForcingTermType HArDCore3D::tworegions_f
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return - scaling_mu * ( viscosity_S * Xi_S(x) * (DivGrad_alpha(x) + DivGrad_beta_S(x)) +
( std::pow(permeability_S, -1) * Xi_S(x) + std::pow(permeability_D, -1) * Xi_D(x) ) * tworegions_u(x)
}
static Brinkman::MomentumForcingTermType DivGrad_alpha
Definition hho-brinkman.hpp:881
static std::function< double(const VectorRd &)> Xi_D
Definition hho-brinkman.hpp:765
constexpr double permeability_D
Definition hho-brinkman.hpp:777
static Brinkman::MomentumForcingTermType DivGrad_beta_S
Definition hho-brinkman.hpp:890
static std::function< double(const VectorRd &)> Xi_S
Definition hho-brinkman.hpp:761
static Brinkman::PressureGradientType tworegions_gradp
Definition hho-brinkman.hpp:868
static Brinkman::MomentumForcingTermType DivGrad_beta_D
Definition hho-brinkman.hpp:900
constexpr double viscosity_S
Definition hho-brinkman.hpp:770
static Brinkman::VelocityType tworegions_u
Definition hho-brinkman.hpp:813
constexpr double viscosity_D
Definition hho-brinkman.hpp:771
constexpr double permeability_S
Definition hho-brinkman.hpp:776

◆ tworegions_g

Brinkman::CompressibilityForcingTermType HArDCore3D::tworegions_g = [](const VectorRd & x)->double { return tworegions_gradu(x).trace();}
static

◆ tworegions_gradp

Brinkman::PressureGradientType HArDCore3D::tworegions_gradp = trigonometric_gradp
static

◆ tworegions_gradu

Brinkman::VelocityGradientType HArDCore3D::tworegions_gradu
static
Initial value:
= [](const VectorRd & x) -> MatrixRd {
return grad_alpha(x) + Xi_S(x) * grad_beta_S(x) + Xi_D(x) * grad_beta_D(x);
}
static Brinkman::VelocityGradientType grad_beta_S
Definition hho-brinkman.hpp:828
static Brinkman::VelocityGradientType grad_alpha
Definition hho-brinkman.hpp:818
static Brinkman::VelocityGradientType grad_beta_D
Definition hho-brinkman.hpp:844

◆ tworegions_kappainv

BrinkmanParameters::PermeabilityInvType HArDCore3D::tworegions_kappainv
static

◆ tworegions_mu

BrinkmanParameters::ViscosityType HArDCore3D::tworegions_mu
static

◆ tworegions_p

Brinkman::PressureType HArDCore3D::tworegions_p = trigonometric_p
static

◆ tworegions_u

Brinkman::VelocityType HArDCore3D::tworegions_u
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return alpha(x) + Xi_S(x) * beta_S(x) + Xi_D(x) * beta_D(x);
}
static std::function< VectorRd(const VectorRd &)> beta_D
Definition hho-brinkman.hpp:803
static std::function< VectorRd(const VectorRd &)> beta_S
Definition hho-brinkman.hpp:793
static std::function< VectorRd(const VectorRd &)> alpha
Definition hho-brinkman.hpp:784

◆ u

double HArDCore3D::BrinkmanNorms::u

Energy norm of velocity.

◆ vcracked_f

Brinkman::MomentumForcingTermType HArDCore3D::vcracked_f
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return VectorRd(0., 0., -.98);
}

◆ vcracked_g

Brinkman::CompressibilityForcingTermType HArDCore3D::vcracked_g
static
Initial value:
= [](const VectorRd & x) -> double {
return 0.;
}

◆ vcracked_gradp

Brinkman::PressureGradientType HArDCore3D::vcracked_gradp
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return VectorRd::Zero();
}

◆ vcracked_gradu

Brinkman::VelocityGradientType HArDCore3D::vcracked_gradu
static
Initial value:
= [](const VectorRd & x) -> MatrixRd {
return MatrixRd::Zero();
}

◆ vcracked_kappainv

BrinkmanParameters::PermeabilityInvType HArDCore3D::vcracked_kappainv = BrinkmanParameters::PermeabilityInvType([](const VectorRd & x)->double { return 1e-5*(1.-XiV(x)) + 0.001*XiV(x);})
static

◆ vcracked_mu

BrinkmanParameters::ViscosityType HArDCore3D::vcracked_mu = BrinkmanParameters::ViscosityType([](const VectorRd & x)->double { return 0.01 * XiV(x);})
static

◆ vcracked_p

Brinkman::PressureType HArDCore3D::vcracked_p
static
Initial value:
= [](const VectorRd & x) -> double {
return 0.;
}

◆ vcracked_u

Brinkman::VelocityType HArDCore3D::vcracked_u
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
VectorRd value = VectorRd::Zero();
if ( x.x() < .2) {
value = VectorRd(.3, 0., 0.);
}else if ( (x.z() > 0.8) && (x.x() > 0.75-eps) && (x.x() < 1.25+eps) ){
value = VectorRd(1, 0., 0.);
}
return value;
}

◆ vec_p

const VectorRd HArDCore3D::vec_p = VectorRd(-1., 2., -5.)
static

◆ viscosity_D

constexpr double HArDCore3D::viscosity_D = 0.
constexpr

◆ viscosity_in_cavity

constexpr double HArDCore3D::viscosity_in_cavity = 0.01
constexpr

◆ viscosity_S

constexpr double HArDCore3D::viscosity_S = 1.
constexpr

◆ Xi_D

std::function<double(const VectorRd &)> HArDCore3D::Xi_D
static
Initial value:
= [](const VectorRd & x)->double {
return 1.-Xi_S(x);
}

◆ Xi_S

std::function<double(const VectorRd &)> HArDCore3D::Xi_S
static
Initial value:
= [](const VectorRd & x)->double {
return ( (x.x()<0.5) ? 1. : 0.);
}

◆ XiBox

std::function<double(const VectorRd &)> HArDCore3D::XiBox = [](const VectorRd & x)->double { return 1. - XiCavity(x) - XiWedge(x); }
static

◆ XiCavity

std::function<double(const VectorRd &)> HArDCore3D::XiCavity
static
Initial value:
= [](const VectorRd & x)->double {
return ( (x.x()>0) && (x.x()<1) && (x.y()>0) && (x.y()<1) && (x.z()>-1) ? 1. : 0.);
}

◆ XiS

std::function<double(const VectorRd &)> HArDCore3D::XiS
static
Initial value:
= [](const VectorRd & x) -> double {
return (scaling_mu > eps ? std::exp(-scaling_kappainv/scaling_mu) : 0.);
}

◆ XiV

std::function<double(const VectorRd &)> HArDCore3D::XiV = [](const VectorRd & x)->double { return (x.z() > 2*std::abs(x.x()-1.)+0.5-eps); }
static

◆ XiWedge

std::function<double(const VectorRd &)> HArDCore3D::XiWedge
static
Initial value:
= [](const VectorRd & x)->double {
return ( (x.x()>=1) && (x.x()<2) && (x.y()>0) && (x.y()<1) && (x.z()> -.75 + .25*(x.x()-1)) ? 1. : 0.);
}