HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Loading...
Searching...
No Matches
Classes | Typedefs | Functions | Variables
DDRCore

Classes providing tools to the Discrete De Rham sequence. More...

Classes

class  HArDCore2D::DDRCore
 Construct all polynomial spaces for the DDR sequence. More...
 
struct  HArDCore2D::DDRCore::CellBases
 Structure to store element bases. More...
 
struct  HArDCore2D::DDRCore::EdgeBases
 Structure to store edge bases. More...
 
class  HArDCore2D::EXCurl
 Extended XCurl space, with vector-valued polynomials on the edges. More...
 
struct  HArDCore2D::EXCurl::hhoLocalOperators
 A structure to store the local HHO operators. More...
 
class  HArDCore2D::XCurl
 Discrete Hcurl space: local operators, L2 product and global interpolator. More...
 
struct  HArDCore2D::XCurl::LocalOperators
 A structure to store the local operators (curl and potential) More...
 
class  HArDCore2D::XGrad
 Discrete H1 space: local operators, L2 product and global interpolator. More...
 
struct  HArDCore2D::XGrad::LocalOperators
 A structure to store local operators (gradient and potential) More...
 
class  HArDCore2D::XRot
 Discrete H1 space: local operators, L2 product and global interpolator. More...
 
struct  HArDCore2D::XRot::LocalOperators
 A structure to store local operators (vector rotor and potential) More...
 
class  HArDCore2D::XRotRot
 Discrete HRotRot space: local operators, L2 product and global interpolator. More...
 
struct  HArDCore2D::XRotRot::LocalOperators
 A structure to store the local operators (scalar rotor and potential) More...
 

Typedefs

typedef Family< MonomialScalarBasisCellHArDCore2D::DDRCore::PolyBasisCellType
 
typedef TensorizedVectorFamily< PolyBasisCellType, 2 > HArDCore2D::DDRCore::Poly2BasisCellType
 
typedef Family< GradientBasis< ShiftedBasis< MonomialScalarBasisCell > > > HArDCore2D::DDRCore::GolyBasisCellType
 
typedef Family< GolyComplBasisCellHArDCore2D::DDRCore::GolyComplBasisCellType
 
typedef Family< CurlBasis< ShiftedBasis< MonomialScalarBasisCell > > > HArDCore2D::DDRCore::RolyBasisCellType
 
typedef Family< RolyComplBasisCellHArDCore2D::DDRCore::RolyComplBasisCellType
 
typedef Family< MonomialScalarBasisEdgeHArDCore2D::DDRCore::PolyBasisEdgeType
 
typedef Cell HArDCore2D::DDRCore::CellBases::GeometricSupport
 Geometric support.
 
typedef Edge HArDCore2D::DDRCore::EdgeBases::GeometricSupport
 Geometric support.
 
typedef std::function< Eigen::Vector2d(const Eigen::Vector2d &)> HArDCore2D::EXCurl::FunctionType
 
typedef MatrixFamily< DDRCore::PolyBasisCellType, dimspaceHArDCore2D::EXCurl::Polyk2x2Type
 
typedef TensorizedVectorFamily< DDRCore::PolyBasisCellType, dimspaceHArDCore2D::EXCurl::Polykpo2Type
 
typedef std::function< Eigen::Vector2d(const Eigen::Vector2d &)> HArDCore2D::XCurl::FunctionType
 
typedef std::function< double(const Eigen::Vector2d &)> HArDCore2D::XGrad::FunctionType
 
typedef std::function< double(const Eigen::Vector2d &)> HArDCore2D::XRot::FunctionType
 
typedef std::function< Eigen::Vector2d(const Eigen::Vector2d &)> HArDCore2D::XRotRot::FunctionType
 
typedef std::function< double(const Eigen::Vector2d &)> HArDCore2D::XRotRot::RotType
 

Functions

 HArDCore2D::DDRCore::DDRCore (const Mesh &mesh, size_t K, bool use_threads=true, std::ostream &output=std::cout)
 Constructor.
 
const MeshHArDCore2D::DDRCore::mesh () const
 Return a const reference to the mesh.
 
const size_t & HArDCore2D::DDRCore::degree () const
 Return the polynomial degree.
 
const CellBasesHArDCore2D::DDRCore::cellBases (size_t iT) const
 Return cell bases for element iT.
 
const EdgeBasesHArDCore2D::DDRCore::edgeBases (size_t iE) const
 Return edge bases for edge iE.
 
 HArDCore2D::EXCurl::hhoLocalOperators::hhoLocalOperators (const Eigen::MatrixXd &_gradient, const Eigen::MatrixXd &_symmetric_gradient, const Eigen::MatrixXd &_potential, const Eigen::MatrixXd &_stabilisation)
 
 HArDCore2D::EXCurl::EXCurl (const DDRCore &ddr_core, bool use_threads=true, std::ostream &output=std::cout)
 Constructor.
 
const MeshHArDCore2D::EXCurl::mesh () const
 Return the mesh.
 
const size_t & HArDCore2D::EXCurl::degree () const
 Return the polynomial degree.
 
Eigen::VectorXd HArDCore2D::EXCurl::interpolate (const FunctionType &v, const int deg_quad=-1) const
 Interpolator of a continuous function.
 
const Eigen::MatrixXd HArDCore2D::EXCurl::ddrPotential (size_t iT) const
 Return cell potential for the cell of index iT.
 
const Eigen::MatrixXd HArDCore2D::EXCurl::ddrCurl (size_t iT) const
 Return cell (scalar) curl for the cell of index iT.
 
const Eigen::MatrixXd HArDCore2D::EXCurl::ddrPotential (const Cell &T) const
 Return cell potential for cell T.
 
const Eigen::MatrixXd HArDCore2D::EXCurl::ddrCurl (const Cell &T) const
 Return cell (scalar) curl for cell T.
 
const hhoLocalOperatorsHArDCore2D::EXCurl::hhoOperators (size_t iT) const
 Return hho operators for the cell of index iT.
 
const hhoLocalOperatorsHArDCore2D::EXCurl::hhoOperators (const Cell &T) const
 Return cell operators for cell T.
 
const DDRCore::CellBasesHArDCore2D::EXCurl::cellBases (size_t iT) const
 Return ddrcore cell bases for the cell of index iT.
 
const DDRCore::CellBasesHArDCore2D::EXCurl::cellBases (const Cell &T) const
 Return ddrcore cell bases for cell T.
 
const Polyk2x2TypeHArDCore2D::EXCurl::Polyk2x2 (size_t iT) const
 Return basis for Matricial P^k space.
 
const Polykpo2TypeHArDCore2D::EXCurl::Polykpo2 (size_t iT) const
 Return basis for the (P^{k+1})^2 space.
 
const DDRCore::EdgeBasesHArDCore2D::EXCurl::edgeBases (size_t iE) const
 Return ddrcore edge bases for the edge of index iE.
 
const DDRCore::EdgeBasesHArDCore2D::EXCurl::edgeBases (const Edge &E) const
 Return ddrcore edge bases for edge E.
 
Eigen::MatrixXd HArDCore2D::EXCurl::computeL2Product (const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk2_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const
 Compute the matrix of the (weighted) L2-product for the cell of index iT.
 
Eigen::MatrixXd HArDCore2D::EXCurl::computeL2ProductGradient (const size_t iT, const XGrad &x_grad, const std::string &side, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk2_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const
 Compute the matrix of the (weighted) L2-product as 'computeL2Product', with application of the discrete gradient on the left/right/both sides (depending on argument "side").
 
Eigen::MatrixXd HArDCore2D::EXCurl::computeL2Product_with_Ops (const size_t iT, const std::vector< Eigen::MatrixXd > &leftOp, const std::vector< Eigen::MatrixXd > &rightOp, const double &penalty_factor, const Eigen::MatrixXd &w_mass_Pk2_T, const IntegralWeight &weight) const
 Compute the matrix of the L2 product, applying leftOp and rightOp to the variables. Probably not directly called, mostly invoked through the wrappers computeL2Product and computeL2ProductGradient.
 
 HArDCore2D::XCurl::LocalOperators::LocalOperators (const Eigen::MatrixXd &_curl, const Eigen::MatrixXd &_potential)
 
 HArDCore2D::XCurl::XCurl (const DDRCore &ddr_core, bool use_threads=true, std::ostream &output=std::cout)
 Constructor.
 
const MeshHArDCore2D::XCurl::mesh () const
 Return the mesh.
 
const size_t & HArDCore2D::XCurl::degree () const
 Return the polynomial degree.
 
Eigen::VectorXd HArDCore2D::XCurl::interpolate (const FunctionType &v, const int deg_quad=-1) const
 Interpolator of a continuous function.
 
const LocalOperatorsHArDCore2D::XCurl::cellOperators (size_t iT) const
 Return cell operators for the cell of index iT.
 
const LocalOperatorsHArDCore2D::XCurl::cellOperators (const Cell &T) const
 Return cell operators for cell T.
 
const DDRCore::CellBasesHArDCore2D::XCurl::cellBases (size_t iT) const
 Return cell bases for the cell of index iT.
 
const DDRCore::CellBasesHArDCore2D::XCurl::cellBases (const Cell &T) const
 Return cell bases for cell T.
 
const DDRCore::EdgeBasesHArDCore2D::XCurl::edgeBases (size_t iE) const
 Return edge bases for the edge of index iE.
 
const DDRCore::EdgeBasesHArDCore2D::XCurl::edgeBases (const Edge &E) const
 Return edge bases for edge E.
 
Eigen::MatrixXd HArDCore2D::XCurl::computeL2Product (const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk2_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const
 Compute the matrix of the (weighted) L2-product for the cell of index iT.
 
Eigen::MatrixXd HArDCore2D::XCurl::computeL2ProductGradient (const size_t iT, const XGrad &x_grad, const std::string &side, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk2_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const
 Compute the matrix of the (weighted) L2-product as 'computeL2Product', with application of the discrete gradient on the left/right/both sides (depending on argument "side").
 
Eigen::MatrixXd HArDCore2D::XCurl::computeL2Product_with_Ops (const size_t iT, const std::vector< Eigen::MatrixXd > &leftOp, const std::vector< Eigen::MatrixXd > &rightOp, const double &penalty_factor, const Eigen::MatrixXd &w_mass_Pk2_T, const IntegralWeight &weight) const
 Compute the matrix of the L2 product, applying leftOp and rightOp to the variables. Probably not directly called, mostly invoked through the wrappers computeL2Product and computeL2ProductGradient.
 
 HArDCore2D::XGrad::LocalOperators::LocalOperators (const Eigen::MatrixXd &_gradient, const Eigen::MatrixXd &_potential)
 
 HArDCore2D::XGrad::XGrad (const DDRCore &ddr_core, bool use_threads=true, std::ostream &output=std::cout)
 Constructor.
 
const MeshHArDCore2D::XGrad::mesh () const
 Return the mesh.
 
const size_t & HArDCore2D::XGrad::degree () const
 Return the polynomial degree.
 
Eigen::VectorXd HArDCore2D::XGrad::interpolate (const FunctionType &q, const int deg_quad=-1) const
 Interpolator of a continuous function.
 
const LocalOperatorsHArDCore2D::XGrad::edgeOperators (size_t iE) const
 Return edge operators for the edge of index iE.
 
const LocalOperatorsHArDCore2D::XGrad::edgeOperators (const Edge &E) const
 Return edge operators for edge E.
 
const LocalOperatorsHArDCore2D::XGrad::cellOperators (size_t iT) const
 Return cell operators for the cell of index iT.
 
const LocalOperatorsHArDCore2D::XGrad::cellOperators (const Cell &T) const
 Return cell operators for cell T.
 
const DDRCore::CellBasesHArDCore2D::XGrad::cellBases (size_t iT) const
 Return cell bases for the cell of index iT.
 
const DDRCore::CellBasesHArDCore2D::XGrad::cellBases (const Cell &T) const
 Return cell bases for cell T.
 
const DDRCore::EdgeBasesHArDCore2D::XGrad::edgeBases (size_t iE) const
 Return edge bases for the edge of index iE.
 
const DDRCore::EdgeBasesHArDCore2D::XGrad::edgeBases (const Edge &E) const
 Return edge bases for edge E.
 
Eigen::MatrixXd HArDCore2D::XGrad::computeL2Product (const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pkpo_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const
 Compute the matrix of the (weighted) L2-product for the cell of index iT.
 
double HArDCore2D::XGrad::evaluatePotential (const size_t iT, const Eigen::VectorXd &vT, const VectorRd &x) const
 Evaluate the value of the potential at a point x.
 
double HArDCore2D::XGrad::computeL2Norm (const Eigen::VectorXd &v) const
 Compute the L2-norm of a vector of the space.
 
 HArDCore2D::XRot::LocalOperators::LocalOperators (const Eigen::MatrixXd &_rotor, const Eigen::MatrixXd &_rotor_rhs, const Eigen::MatrixXd &_potential)
 
 HArDCore2D::XRot::XRot (const DDRCore &ddr_core, bool use_threads=true, std::ostream &output=std::cout)
 Constructor.
 
const MeshHArDCore2D::XRot::mesh () const
 Return the mesh.
 
const bool HArDCore2D::XRot::useThreads () const
 Return true if we use thread-based parallelism.
 
const size_t & HArDCore2D::XRot::degree () const
 Return the polynomial degree.
 
Eigen::VectorXd HArDCore2D::XRot::interpolate (const FunctionType &q, const int deg_quad=-1) const
 Interpolator of a continuous function.
 
const Eigen::MatrixXd & HArDCore2D::XRot::edgePotential (size_t iE) const
 Return edge potential for the edge of index iE.
 
const Eigen::MatrixXd & HArDCore2D::XRot::edgePotential (const Edge &E) const
 Return edge potential for the edge E.
 
const LocalOperatorsHArDCore2D::XRot::cellOperators (size_t iT) const
 Return cell operators for the cell of index iT.
 
const LocalOperatorsHArDCore2D::XRot::cellOperators (const Cell &T) const
 Return cell operators for cell T.
 
const DDRCore::CellBasesHArDCore2D::XRot::cellBases (size_t iT) const
 Return cell bases for the cell of index iT.
 
const DDRCore::CellBasesHArDCore2D::XRot::cellBases (const Cell &T) const
 Return cell bases for cell T.
 
const DDRCore::EdgeBasesHArDCore2D::XRot::edgeBases (size_t iE) const
 Return edge bases for the edge of index iE.
 
const DDRCore::EdgeBasesHArDCore2D::XRot::edgeBases (const Edge &E) const
 Return edge bases for edge E.
 
Eigen::MatrixXd HArDCore2D::XRot::computeRotorL2Product (const size_t iT, const double &penalty_factor_cell=1., const double &penalty_factor_edge=1.) const
 Compute rotor L2-product.
 
double HArDCore2D::XRot::computeRotorL2Norm (const Eigen::VectorXd &v) const
 Compute rotor L2-norm.
 
 HArDCore2D::XRotRot::LocalOperators::LocalOperators (const Eigen::MatrixXd &_rotor, const Eigen::MatrixXd &_potential)
 
 HArDCore2D::XRotRot::XRotRot (const DDRCore &ddr_core, bool use_threads=true, std::ostream &output=std::cout)
 Constructor.
 
const MeshHArDCore2D::XRotRot::mesh () const
 Return the mesh.
 
const size_t & HArDCore2D::XRotRot::degree () const
 Return the polynomial degree.
 
Eigen::VectorXd HArDCore2D::XRotRot::interpolate (const FunctionType &v, const RotType &rot_v, const int deg_quad=-1) const
 Interpolator of a continuous function.
 
const LocalOperatorsHArDCore2D::XRotRot::cellOperators (size_t iT) const
 Return cell operators for the cell of index iT.
 
const LocalOperatorsHArDCore2D::XRotRot::cellOperators (const Cell &T) const
 Return cell operators for cell T.
 
Eigen::MatrixXd HArDCore2D::XRotRot::cellRotor (size_t iT) const
 Return cell rotor for cell of index iT.
 
Eigen::MatrixXd HArDCore2D::XRotRot::cellRotor (const Cell &T) const
 Return cell rotor for cell T.
 
const DDRCore::CellBasesHArDCore2D::XRotRot::cellBases (size_t iT) const
 Return cell bases for the cell of index iT.
 
const DDRCore::CellBasesHArDCore2D::XRotRot::cellBases (const Cell &T) const
 Return cell bases for cell T.
 
const DDRCore::EdgeBasesHArDCore2D::XRotRot::edgeBases (size_t iE) const
 Return edge bases for the edge of index iE.
 
const DDRCore::EdgeBasesHArDCore2D::XRotRot::edgeBases (const Edge &E) const
 Return edge bases for edge E.
 
Eigen::MatrixXd HArDCore2D::XRotRot::computeL2Product (size_t iT, const double &penalty_factor=1., const IntegralWeight &weight=IntegralWeight(1.)) const
 Compute the matrix of the (weighted) L2-product for the cell of index iT.
 
Eigen::MatrixXd HArDCore2D::XRotRot::computeGradientPotentialL2Product (size_t iT, const XGrad *x_grad, const double &penalty_factor=1., const IntegralWeight &weight=IntegralWeight(1.)) const
 Compute the matrix of the (weighted) L2-product as 'computeL2Product', with application of the discrete gradient on the left side.
 
Eigen::MatrixXd HArDCore2D::XRotRot::computeGradientL2Product (size_t iT, const XGrad *x_grad, const double &penalty_factor=1., const IntegralWeight &weight=IntegralWeight(1.)) const
 Compute the matrix of the (weighted) gradient-gradient L2-product.
 
template<typename LeftOperatorFillerType , typename RightOperatorFillerType >
Eigen::MatrixXd HArDCore2D::XRotRot::computeL2ProductWithOperatorFillers (size_t iT, const double &penalty_factor, const IntegralWeight &weight, LeftOperatorFillerType fillLeftOp, RightOperatorFillerType fillRightOp) const
 Compute the matrix of the L2 product given operators filling functions.
 
double HArDCore2D::XRotRot::computeL2Norm (const Eigen::VectorXd &v) const
 Compute the L2-norm of a vector of the space.
 
double HArDCore2D::XRotRot::computeGradientL2Norm (const Eigen::VectorXd &v, const XGrad *x_grad) const
 Compute the L2-norm of the discrete gradient of a vector in XGrad.
 

Variables

std::unique_ptr< PolyBasisCellTypeHArDCore2D::DDRCore::CellBases::Polykpo
 
std::unique_ptr< PolyBasisCellTypeHArDCore2D::DDRCore::CellBases::Polyk
 
std::unique_ptr< PolyBasisCellTypeHArDCore2D::DDRCore::CellBases::Polykmo
 
std::unique_ptr< PolyBasisCellTypeHArDCore2D::DDRCore::CellBases::Polykmtwo
 
std::unique_ptr< Poly2BasisCellTypeHArDCore2D::DDRCore::CellBases::Polyk2
 
std::unique_ptr< Poly2BasisCellTypeHArDCore2D::DDRCore::CellBases::Polykmo2
 
std::unique_ptr< RolyBasisCellTypeHArDCore2D::DDRCore::CellBases::Rolykmo
 
std::unique_ptr< RolyComplBasisCellTypeHArDCore2D::DDRCore::CellBases::RolyComplk
 
std::unique_ptr< RolyComplBasisCellTypeHArDCore2D::DDRCore::CellBases::RolyComplkp2
 
std::unique_ptr< GolyComplBasisCellTypeHArDCore2D::DDRCore::CellBases::GolyComplkp2
 
std::unique_ptr< PolyBasisEdgeTypeHArDCore2D::DDRCore::EdgeBases::Polykpo
 
std::unique_ptr< PolyBasisEdgeTypeHArDCore2D::DDRCore::EdgeBases::Polyk
 
std::unique_ptr< PolyBasisEdgeTypeHArDCore2D::DDRCore::EdgeBases::Polykmo
 
Eigen::MatrixXd HArDCore2D::EXCurl::hhoLocalOperators::gradient
 
Eigen::MatrixXd HArDCore2D::EXCurl::hhoLocalOperators::symmetric_gradient
 
Eigen::MatrixXd HArDCore2D::EXCurl::hhoLocalOperators::potential
 
Eigen::MatrixXd HArDCore2D::EXCurl::hhoLocalOperators::stabilisation
 
Eigen::MatrixXd HArDCore2D::XCurl::LocalOperators::curl
 
Eigen::MatrixXd HArDCore2D::XCurl::LocalOperators::potential
 
Eigen::MatrixXd HArDCore2D::XGrad::LocalOperators::gradient
 
Eigen::MatrixXd HArDCore2D::XGrad::LocalOperators::potential
 
Eigen::MatrixXd HArDCore2D::XRot::LocalOperators::rotor
 
Eigen::MatrixXd HArDCore2D::XRot::LocalOperators::rotor_rhs
 
Eigen::MatrixXd HArDCore2D::XRot::LocalOperators::potential
 
Eigen::MatrixXd HArDCore2D::XRotRot::LocalOperators::rotor
 
Eigen::MatrixXd HArDCore2D::XRotRot::LocalOperators::potential
 

Detailed Description

Classes providing tools to the Discrete De Rham sequence.

Typedef Documentation

◆ FunctionType [1/5]

typedef std::function<Eigen::Vector2d(const Eigen::Vector2d &)> HArDCore2D::EXCurl::FunctionType

◆ FunctionType [2/5]

typedef std::function<Eigen::Vector2d(const Eigen::Vector2d &)> HArDCore2D::XCurl::FunctionType

◆ FunctionType [3/5]

typedef std::function<double(const Eigen::Vector2d &)> HArDCore2D::XGrad::FunctionType

◆ FunctionType [4/5]

typedef std::function<double(const Eigen::Vector2d &)> HArDCore2D::XRot::FunctionType

◆ FunctionType [5/5]

typedef std::function<Eigen::Vector2d(const Eigen::Vector2d &)> HArDCore2D::XRotRot::FunctionType

◆ GeometricSupport [1/2]

Geometric support.

◆ GeometricSupport [2/2]

Geometric support.

◆ GolyBasisCellType

◆ GolyComplBasisCellType

◆ Poly2BasisCellType

◆ PolyBasisCellType

◆ PolyBasisEdgeType

◆ Polyk2x2Type

◆ Polykpo2Type

◆ RolyBasisCellType

◆ RolyComplBasisCellType

◆ RotType

typedef std::function<double(const Eigen::Vector2d &)> HArDCore2D::XRotRot::RotType

Function Documentation

◆ cellBases() [1/11]

const DDRCore::CellBases & HArDCore2D::EXCurl::cellBases ( const Cell &  T) const
inline

Return ddrcore cell bases for cell T.

◆ cellBases() [2/11]

const DDRCore::CellBases & HArDCore2D::XCurl::cellBases ( const Cell &  T) const
inline

Return cell bases for cell T.

◆ cellBases() [3/11]

const DDRCore::CellBases & HArDCore2D::XGrad::cellBases ( const Cell &  T) const
inline

Return cell bases for cell T.

◆ cellBases() [4/11]

const DDRCore::CellBases & HArDCore2D::XRot::cellBases ( const Cell &  T) const
inline

Return cell bases for cell T.

◆ cellBases() [5/11]

const DDRCore::CellBases & HArDCore2D::XRotRot::cellBases ( const Cell &  T) const
inline

Return cell bases for cell T.

◆ cellBases() [6/11]

const CellBases & HArDCore2D::DDRCore::cellBases ( size_t  iT) const
inline

Return cell bases for element iT.

◆ cellBases() [7/11]

const DDRCore::CellBases & HArDCore2D::EXCurl::cellBases ( size_t  iT) const
inline

Return ddrcore cell bases for the cell of index iT.

◆ cellBases() [8/11]

const DDRCore::CellBases & HArDCore2D::XCurl::cellBases ( size_t  iT) const
inline

Return cell bases for the cell of index iT.

◆ cellBases() [9/11]

const DDRCore::CellBases & HArDCore2D::XGrad::cellBases ( size_t  iT) const
inline

Return cell bases for the cell of index iT.

◆ cellBases() [10/11]

const DDRCore::CellBases & HArDCore2D::XRot::cellBases ( size_t  iT) const
inline

Return cell bases for the cell of index iT.

◆ cellBases() [11/11]

const DDRCore::CellBases & HArDCore2D::XRotRot::cellBases ( size_t  iT) const
inline

Return cell bases for the cell of index iT.

◆ cellOperators() [1/8]

const LocalOperators & HArDCore2D::XCurl::cellOperators ( const Cell &  T) const
inline

Return cell operators for cell T.

◆ cellOperators() [2/8]

const LocalOperators & HArDCore2D::XGrad::cellOperators ( const Cell &  T) const
inline

Return cell operators for cell T.

◆ cellOperators() [3/8]

const LocalOperators & HArDCore2D::XRot::cellOperators ( const Cell &  T) const
inline

Return cell operators for cell T.

◆ cellOperators() [4/8]

const LocalOperators & HArDCore2D::XRotRot::cellOperators ( const Cell &  T) const
inline

Return cell operators for cell T.

◆ cellOperators() [5/8]

const LocalOperators & HArDCore2D::XCurl::cellOperators ( size_t  iT) const
inline

Return cell operators for the cell of index iT.

◆ cellOperators() [6/8]

const LocalOperators & HArDCore2D::XGrad::cellOperators ( size_t  iT) const
inline

Return cell operators for the cell of index iT.

◆ cellOperators() [7/8]

const LocalOperators & HArDCore2D::XRot::cellOperators ( size_t  iT) const
inline

Return cell operators for the cell of index iT.

◆ cellOperators() [8/8]

const LocalOperators & HArDCore2D::XRotRot::cellOperators ( size_t  iT) const
inline

Return cell operators for the cell of index iT.

◆ cellRotor() [1/2]

Eigen::MatrixXd HArDCore2D::XRotRot::cellRotor ( const Cell &  T) const
inline

Return cell rotor for cell T.

◆ cellRotor() [2/2]

Eigen::MatrixXd HArDCore2D::XRotRot::cellRotor ( size_t  iT) const
inline

Return cell rotor for cell of index iT.

◆ computeGradientL2Norm()

double XRotRot::computeGradientL2Norm ( const Eigen::VectorXd &  v,
const XGrad x_grad 
) const

Compute the L2-norm of the discrete gradient of a vector in XGrad.

◆ computeGradientL2Product()

Eigen::MatrixXd XRotRot::computeGradientL2Product ( size_t  iT,
const XGrad x_grad,
const double &  penalty_factor = 1.,
const IntegralWeight weight = IntegralWeight(1.) 
) const

Compute the matrix of the (weighted) gradient-gradient L2-product.

Parameters
iTIndex of the cell
x_gradInstance of XGrad to access the full gradients
penalty_factorPre-factor for stabilisation term
weightWeight function in the L2 product, defaults to constant 1.

◆ computeGradientPotentialL2Product()

Eigen::MatrixXd XRotRot::computeGradientPotentialL2Product ( size_t  iT,
const XGrad x_grad,
const double &  penalty_factor = 1.,
const IntegralWeight weight = IntegralWeight(1.) 
) const

Compute the matrix of the (weighted) L2-product as 'computeL2Product', with application of the discrete gradient on the left side.

Parameters
iTIndex of the cell
x_gradInstance of XGrad to access the full gradients
penalty_factorPre-factor for stabilisation term
weightWeight function in the L2 product, defaults to constant 1.

◆ computeL2Norm() [1/2]

double XGrad::computeL2Norm ( const Eigen::VectorXd &  v) const

Compute the L2-norm of a vector of the space.

◆ computeL2Norm() [2/2]

double XRotRot::computeL2Norm ( const Eigen::VectorXd &  v) const

Compute the L2-norm of a vector of the space.

◆ computeL2Product() [1/4]

Eigen::MatrixXd EXCurl::computeL2Product ( const size_t  iT,
const double &  penalty_factor = 1.,
const Eigen::MatrixXd &  mass_Pk2_T = Eigen::MatrixXd::Zero(1,1),
const IntegralWeight weight = IntegralWeight(1.) 
) const

Compute the matrix of the (weighted) L2-product for the cell of index iT.

Parameters
iTindex of the cell
penalty_factorpre-factor for stabilisation term
mass_Pk2_Tif pre-computed, the mass matrix of (P^k(T))^3; if none is pre-computed, passing Eigen::MatrixXd::Zero(1,1) will force the calculation
weightweight function in the L2 product, defaults to 1

◆ computeL2Product() [2/4]

Eigen::MatrixXd XCurl::computeL2Product ( const size_t  iT,
const double &  penalty_factor = 1.,
const Eigen::MatrixXd &  mass_Pk2_T = Eigen::MatrixXd::Zero(1,1),
const IntegralWeight weight = IntegralWeight(1.) 
) const

Compute the matrix of the (weighted) L2-product for the cell of index iT.

Parameters
iTindex of the cell
penalty_factorpre-factor for stabilisation term
mass_Pk2_Tif pre-computed, the mass matrix of (P^k(T))^3; if none is pre-computed, passing Eigen::MatrixXd::Zero(1,1) will force the calculation
weightweight function in the L2 product, defaults to 1

◆ computeL2Product() [3/4]

Eigen::MatrixXd XGrad::computeL2Product ( const size_t  iT,
const double &  penalty_factor = 1.,
const Eigen::MatrixXd &  mass_Pkpo_T = Eigen::MatrixXd::Zero(1,1),
const IntegralWeight weight = IntegralWeight(1.) 
) const

Compute the matrix of the (weighted) L2-product for the cell of index iT.

Parameters
iTindex of the cell
penalty_factorpre-factor for stabilisation term
mass_Pkpo_Tif pre-computed, the mass matrix of P^{k+1}(T); if none is pre-computed, passing Eigen::MatrixXd::Zero(1,1) will force the calculation
weightweight function in the L2 product, defaults to 1

◆ computeL2Product() [4/4]

Eigen::MatrixXd XRotRot::computeL2Product ( size_t  iT,
const double &  penalty_factor = 1.,
const IntegralWeight weight = IntegralWeight(1.) 
) const

Compute the matrix of the (weighted) L2-product for the cell of index iT.

Parameters
iTIndex of the cell
penalty_factorPre-factor for stabilisation term
weightWeight function in the L2 product, defaults to 1

◆ computeL2Product_with_Ops() [1/2]

Eigen::MatrixXd HArDCore2D::EXCurl::computeL2Product_with_Ops ( const size_t  iT,
const std::vector< Eigen::MatrixXd > &  leftOp,
const std::vector< Eigen::MatrixXd > &  rightOp,
const double &  penalty_factor,
const Eigen::MatrixXd &  w_mass_Pk2_T,
const IntegralWeight weight 
) const

Compute the matrix of the L2 product, applying leftOp and rightOp to the variables. Probably not directly called, mostly invoked through the wrappers computeL2Product and computeL2ProductGradient.

Parameters
iTindex of the cell
leftOpedge and element operators to apply on the left
rightOpedge and element operators to apply on the right
penalty_factorpre-factor for stabilisation term
w_mass_Pk2_Tmass matrix of (P^k(T))^3 weighted by weight
weightweight function in the L2 product

◆ computeL2Product_with_Ops() [2/2]

Eigen::MatrixXd XCurl::computeL2Product_with_Ops ( const size_t  iT,
const std::vector< Eigen::MatrixXd > &  leftOp,
const std::vector< Eigen::MatrixXd > &  rightOp,
const double &  penalty_factor,
const Eigen::MatrixXd &  w_mass_Pk2_T,
const IntegralWeight weight 
) const

Compute the matrix of the L2 product, applying leftOp and rightOp to the variables. Probably not directly called, mostly invoked through the wrappers computeL2Product and computeL2ProductGradient.

Parameters
iTindex of the cell
leftOpedge and element operators to apply on the left
rightOpedge and element operators to apply on the right
penalty_factorpre-factor for stabilisation term
w_mass_Pk2_Tmass matrix of (P^k(T))^3 weighted by weight
weightweight function in the L2 product

◆ computeL2ProductGradient() [1/2]

Eigen::MatrixXd EXCurl::computeL2ProductGradient ( const size_t  iT,
const XGrad x_grad,
const std::string &  side,
const double &  penalty_factor = 1.,
const Eigen::MatrixXd &  mass_Pk2_T = Eigen::MatrixXd::Zero(1,1),
const IntegralWeight weight = IntegralWeight(1.) 
) const

Compute the matrix of the (weighted) L2-product as 'computeL2Product', with application of the discrete gradient on the left/right/both sides (depending on argument "side").

Parameters
iTindex of the cell
x_gradinstance of XGrad to access the full gradients
sidewhich side (left, right, both) we apply the gradient to
penalty_factorpre-factor for stabilisation term
mass_Pk2_Tif pre-computed, the mass matrix of (P^k(T))^3; if none is pre-computed, passing Eigen::MatrixXd::Zero(1,1) will force the calculation
weightweight function in the L2 product, defaults to constant 1.

◆ computeL2ProductGradient() [2/2]

Eigen::MatrixXd XCurl::computeL2ProductGradient ( const size_t  iT,
const XGrad x_grad,
const std::string &  side,
const double &  penalty_factor = 1.,
const Eigen::MatrixXd &  mass_Pk2_T = Eigen::MatrixXd::Zero(1,1),
const IntegralWeight weight = IntegralWeight(1.) 
) const

Compute the matrix of the (weighted) L2-product as 'computeL2Product', with application of the discrete gradient on the left/right/both sides (depending on argument "side").

Parameters
iTindex of the cell
x_gradinstance of XGrad to access the full gradients
sidewhich side (left, right, both) we apply the gradient to
penalty_factorpre-factor for stabilisation term
mass_Pk2_Tif pre-computed, the mass matrix of (P^k(T))^3; if none is pre-computed, passing Eigen::MatrixXd::Zero(1,1) will force the calculation
weightweight function in the L2 product, defaults to constant 1.

◆ computeL2ProductWithOperatorFillers()

template<typename LeftOperatorFillerType , typename RightOperatorFillerType >
Eigen::MatrixXd HArDCore2D::XRotRot::computeL2ProductWithOperatorFillers ( size_t  iT,
const double &  penalty_factor,
const IntegralWeight weight,
LeftOperatorFillerType  fillLeftOp,
RightOperatorFillerType  fillRightOp 
) const

Compute the matrix of the L2 product given operators filling functions.

Parameters
iTIndex of the cell
penalty_factorPenalty factor for stabilisation term
weightWeight function in the L2 product
fillLeftOpFunction to fill the left-hand side operator
fillRightOpFunction to fill the right-hand side operator

◆ computeRotorL2Norm()

double XRot::computeRotorL2Norm ( const Eigen::VectorXd &  v) const

Compute rotor L2-norm.

◆ computeRotorL2Product()

Eigen::MatrixXd XRot::computeRotorL2Product ( const size_t  iT,
const double &  penalty_factor_cell = 1.,
const double &  penalty_factor_edge = 1. 
) const

Compute rotor L2-product.

Parameters
iTIndex of the cell
penalty_factor_cellCoefficient for the cell-based term
penalty_factor_edgeCoefficient for the edge-based term

◆ DDRCore()

DDRCore::DDRCore ( const Mesh mesh,
size_t  K,
bool  use_threads = true,
std::ostream &  output = std::cout 
)

Constructor.

◆ ddrCurl() [1/2]

const Eigen::MatrixXd HArDCore2D::EXCurl::ddrCurl ( const Cell &  T) const
inline

Return cell (scalar) curl for cell T.

◆ ddrCurl() [2/2]

const Eigen::MatrixXd HArDCore2D::EXCurl::ddrCurl ( size_t  iT) const
inline

Return cell (scalar) curl for the cell of index iT.

◆ ddrPotential() [1/2]

const Eigen::MatrixXd HArDCore2D::EXCurl::ddrPotential ( const Cell &  T) const
inline

Return cell potential for cell T.

◆ ddrPotential() [2/2]

const Eigen::MatrixXd HArDCore2D::EXCurl::ddrPotential ( size_t  iT) const
inline

Return cell potential for the cell of index iT.

◆ degree() [1/6]

const size_t & HArDCore2D::DDRCore::degree ( ) const
inline

Return the polynomial degree.

◆ degree() [2/6]

const size_t & HArDCore2D::EXCurl::degree ( ) const
inline

Return the polynomial degree.

◆ degree() [3/6]

const size_t & HArDCore2D::XCurl::degree ( ) const
inline

Return the polynomial degree.

◆ degree() [4/6]

const size_t & HArDCore2D::XGrad::degree ( ) const
inline

Return the polynomial degree.

◆ degree() [5/6]

const size_t & HArDCore2D::XRot::degree ( ) const
inline

Return the polynomial degree.

◆ degree() [6/6]

const size_t & HArDCore2D::XRotRot::degree ( ) const
inline

Return the polynomial degree.

◆ edgeBases() [1/11]

const DDRCore::EdgeBases & HArDCore2D::EXCurl::edgeBases ( const Edge &  E) const
inline

Return ddrcore edge bases for edge E.

◆ edgeBases() [2/11]

const DDRCore::EdgeBases & HArDCore2D::XCurl::edgeBases ( const Edge &  E) const
inline

Return edge bases for edge E.

◆ edgeBases() [3/11]

const DDRCore::EdgeBases & HArDCore2D::XGrad::edgeBases ( const Edge &  E) const
inline

Return edge bases for edge E.

◆ edgeBases() [4/11]

const DDRCore::EdgeBases & HArDCore2D::XRot::edgeBases ( const Edge &  E) const
inline

Return edge bases for edge E.

◆ edgeBases() [5/11]

const DDRCore::EdgeBases & HArDCore2D::XRotRot::edgeBases ( const Edge &  E) const
inline

Return edge bases for edge E.

◆ edgeBases() [6/11]

const EdgeBases & HArDCore2D::DDRCore::edgeBases ( size_t  iE) const
inline

Return edge bases for edge iE.

◆ edgeBases() [7/11]

const DDRCore::EdgeBases & HArDCore2D::EXCurl::edgeBases ( size_t  iE) const
inline

Return ddrcore edge bases for the edge of index iE.

◆ edgeBases() [8/11]

const DDRCore::EdgeBases & HArDCore2D::XCurl::edgeBases ( size_t  iE) const
inline

Return edge bases for the edge of index iE.

◆ edgeBases() [9/11]

const DDRCore::EdgeBases & HArDCore2D::XGrad::edgeBases ( size_t  iE) const
inline

Return edge bases for the edge of index iE.

◆ edgeBases() [10/11]

const DDRCore::EdgeBases & HArDCore2D::XRot::edgeBases ( size_t  iE) const
inline

Return edge bases for the edge of index iE.

◆ edgeBases() [11/11]

const DDRCore::EdgeBases & HArDCore2D::XRotRot::edgeBases ( size_t  iE) const
inline

Return edge bases for the edge of index iE.

◆ edgeOperators() [1/2]

const LocalOperators & HArDCore2D::XGrad::edgeOperators ( const Edge &  E) const
inline

Return edge operators for edge E.

◆ edgeOperators() [2/2]

const LocalOperators & HArDCore2D::XGrad::edgeOperators ( size_t  iE) const
inline

Return edge operators for the edge of index iE.

◆ edgePotential() [1/2]

const Eigen::MatrixXd & HArDCore2D::XRot::edgePotential ( const Edge &  E) const
inline

Return edge potential for the edge E.

◆ edgePotential() [2/2]

const Eigen::MatrixXd & HArDCore2D::XRot::edgePotential ( size_t  iE) const
inline

Return edge potential for the edge of index iE.

◆ evaluatePotential()

double XGrad::evaluatePotential ( const size_t  iT,
const Eigen::VectorXd &  vT,
const VectorRd x 
) const

Evaluate the value of the potential at a point x.

Parameters
iTindex of the cell in which to take the potential
vTvector of local DOFs
xpoint at which to evaluate the potential

◆ EXCurl()

EXCurl::EXCurl ( const DDRCore ddr_core,
bool  use_threads = true,
std::ostream &  output = std::cout 
)

Constructor.

◆ hhoLocalOperators()

HArDCore2D::EXCurl::hhoLocalOperators::hhoLocalOperators ( const Eigen::MatrixXd &  _gradient,
const Eigen::MatrixXd &  _symmetric_gradient,
const Eigen::MatrixXd &  _potential,
const Eigen::MatrixXd &  _stabilisation 
)
inline
Parameters
_gradientgradient operator
_symmetric_gradientsymmetric gradient operator
_potentialPotential operator
_stabilisationStabilisation

◆ hhoOperators() [1/2]

const hhoLocalOperators & HArDCore2D::EXCurl::hhoOperators ( const Cell &  T) const
inline

Return cell operators for cell T.

◆ hhoOperators() [2/2]

const hhoLocalOperators & HArDCore2D::EXCurl::hhoOperators ( size_t  iT) const
inline

Return hho operators for the cell of index iT.

◆ interpolate() [1/5]

Eigen::VectorXd XGrad::interpolate ( const FunctionType q,
const int  deg_quad = -1 
) const

Interpolator of a continuous function.

Parameters
qThe function to interpolate
deg_quadThe optional degre of quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used.

◆ interpolate() [2/5]

Eigen::VectorXd XRot::interpolate ( const FunctionType q,
const int  deg_quad = -1 
) const

Interpolator of a continuous function.

Parameters
qThe function to interpolate
deg_quadThe optional degre of quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used.

◆ interpolate() [3/5]

Eigen::VectorXd EXCurl::interpolate ( const FunctionType v,
const int  deg_quad = -1 
) const

Interpolator of a continuous function.

Parameters
vThe function to interpolate
deg_quadThe optional degre of quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used.

◆ interpolate() [4/5]

Eigen::VectorXd XCurl::interpolate ( const FunctionType v,
const int  deg_quad = -1 
) const

Interpolator of a continuous function.

Parameters
vThe function to interpolate
deg_quadThe optional degre of quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used.

◆ interpolate() [5/5]

Eigen::VectorXd XRotRot::interpolate ( const FunctionType v,
const RotType rot_v,
const int  deg_quad = -1 
) const

Interpolator of a continuous function.

Parameters
vThe function to interpolate
rot_vRotor of the function to interpolate
deg_quadThe optional degre of quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used.

◆ LocalOperators() [1/4]

HArDCore2D::XCurl::LocalOperators::LocalOperators ( const Eigen::MatrixXd &  _curl,
const Eigen::MatrixXd &  _potential 
)
inline
Parameters
_curlCurl operator
_potentialPotential operator

◆ LocalOperators() [2/4]

HArDCore2D::XGrad::LocalOperators::LocalOperators ( const Eigen::MatrixXd &  _gradient,
const Eigen::MatrixXd &  _potential 
)
inline
Parameters
_gradientGradient operator
_potentialPotential operator

◆ LocalOperators() [3/4]

HArDCore2D::XRotRot::LocalOperators::LocalOperators ( const Eigen::MatrixXd &  _rotor,
const Eigen::MatrixXd &  _potential 
)
inline
Parameters
_rotorRot operator
_potentialPotential operator

◆ LocalOperators() [4/4]

HArDCore2D::XRot::LocalOperators::LocalOperators ( const Eigen::MatrixXd &  _rotor,
const Eigen::MatrixXd &  _rotor_rhs,
const Eigen::MatrixXd &  _potential 
)
inline
Parameters
_rotorVector rotor operator
_rotor_rhsVector rotor right-hand side
_potentialPotential operator

◆ mesh() [1/6]

const Mesh & HArDCore2D::DDRCore::mesh ( ) const
inline

Return a const reference to the mesh.

◆ mesh() [2/6]

const Mesh & HArDCore2D::EXCurl::mesh ( ) const
inline

Return the mesh.

◆ mesh() [3/6]

const Mesh & HArDCore2D::XCurl::mesh ( ) const
inline

Return the mesh.

◆ mesh() [4/6]

const Mesh & HArDCore2D::XGrad::mesh ( ) const
inline

Return the mesh.

◆ mesh() [5/6]

const Mesh & HArDCore2D::XRot::mesh ( ) const
inline

Return the mesh.

◆ mesh() [6/6]

const Mesh & HArDCore2D::XRotRot::mesh ( ) const
inline

Return the mesh.

◆ Polyk2x2()

const Polyk2x2Type & HArDCore2D::EXCurl::Polyk2x2 ( size_t  iT) const
inline

Return basis for Matricial P^k space.

◆ Polykpo2()

const Polykpo2Type & HArDCore2D::EXCurl::Polykpo2 ( size_t  iT) const
inline

Return basis for the (P^{k+1})^2 space.

◆ useThreads()

const bool HArDCore2D::XRot::useThreads ( ) const
inline

Return true if we use thread-based parallelism.

◆ XCurl()

XCurl::XCurl ( const DDRCore ddr_core,
bool  use_threads = true,
std::ostream &  output = std::cout 
)

Constructor.

◆ XGrad()

XGrad::XGrad ( const DDRCore ddr_core,
bool  use_threads = true,
std::ostream &  output = std::cout 
)

Constructor.

◆ XRot()

XRot::XRot ( const DDRCore ddr_core,
bool  use_threads = true,
std::ostream &  output = std::cout 
)

Constructor.

◆ XRotRot()

XRotRot::XRotRot ( const DDRCore ddr_core,
bool  use_threads = true,
std::ostream &  output = std::cout 
)

Constructor.

Variable Documentation

◆ curl

Eigen::MatrixXd HArDCore2D::XCurl::LocalOperators::curl

◆ GolyComplkp2

std::unique_ptr<GolyComplBasisCellType> HArDCore2D::DDRCore::CellBases::GolyComplkp2

◆ gradient [1/2]

Eigen::MatrixXd HArDCore2D::EXCurl::hhoLocalOperators::gradient

◆ gradient [2/2]

Eigen::MatrixXd HArDCore2D::XGrad::LocalOperators::gradient

◆ Polyk [1/2]

std::unique_ptr<PolyBasisCellType> HArDCore2D::DDRCore::CellBases::Polyk

◆ Polyk [2/2]

std::unique_ptr<PolyBasisEdgeType> HArDCore2D::DDRCore::EdgeBases::Polyk

◆ Polyk2

std::unique_ptr<Poly2BasisCellType> HArDCore2D::DDRCore::CellBases::Polyk2

◆ Polykmo [1/2]

std::unique_ptr<PolyBasisCellType> HArDCore2D::DDRCore::CellBases::Polykmo

◆ Polykmo [2/2]

std::unique_ptr<PolyBasisEdgeType> HArDCore2D::DDRCore::EdgeBases::Polykmo

◆ Polykmo2

std::unique_ptr<Poly2BasisCellType> HArDCore2D::DDRCore::CellBases::Polykmo2

◆ Polykmtwo

std::unique_ptr<PolyBasisCellType> HArDCore2D::DDRCore::CellBases::Polykmtwo

◆ Polykpo [1/2]

std::unique_ptr<PolyBasisCellType> HArDCore2D::DDRCore::CellBases::Polykpo

◆ Polykpo [2/2]

std::unique_ptr<PolyBasisEdgeType> HArDCore2D::DDRCore::EdgeBases::Polykpo

◆ potential [1/5]

Eigen::MatrixXd HArDCore2D::EXCurl::hhoLocalOperators::potential

◆ potential [2/5]

Eigen::MatrixXd HArDCore2D::XCurl::LocalOperators::potential

◆ potential [3/5]

Eigen::MatrixXd HArDCore2D::XGrad::LocalOperators::potential

◆ potential [4/5]

Eigen::MatrixXd HArDCore2D::XRot::LocalOperators::potential

◆ potential [5/5]

Eigen::MatrixXd HArDCore2D::XRotRot::LocalOperators::potential

◆ RolyComplk

std::unique_ptr<RolyComplBasisCellType> HArDCore2D::DDRCore::CellBases::RolyComplk

◆ RolyComplkp2

std::unique_ptr<RolyComplBasisCellType> HArDCore2D::DDRCore::CellBases::RolyComplkp2

◆ Rolykmo

std::unique_ptr<RolyBasisCellType> HArDCore2D::DDRCore::CellBases::Rolykmo

◆ rotor [1/2]

Eigen::MatrixXd HArDCore2D::XRot::LocalOperators::rotor

◆ rotor [2/2]

Eigen::MatrixXd HArDCore2D::XRotRot::LocalOperators::rotor

◆ rotor_rhs

Eigen::MatrixXd HArDCore2D::XRot::LocalOperators::rotor_rhs

◆ stabilisation

Eigen::MatrixXd HArDCore2D::EXCurl::hhoLocalOperators::stabilisation

◆ symmetric_gradient

Eigen::MatrixXd HArDCore2D::EXCurl::hhoLocalOperators::symmetric_gradient