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
DDRCore

Classes defining the Discrete De Rham sequence. More...

Classes

class  HArDCore3D::DDRCore
 Construct all polynomial spaces for the DDR sequence. More...
 
struct  HArDCore3D::DDRCore::CellBases
 Structure to store element bases. More...
 
struct  HArDCore3D::DDRCore::FaceBases
 Structure to store face bases. More...
 
struct  HArDCore3D::DDRCore::EdgeBases
 Structure to store edge bases. More...
 
class  HArDCore3D::SerendipityProblem
 Construct all polynomial spaces for the DDR sequence. More...
 
class  HArDCore3D::SXCurl
 Discrete Serendipity Hcurl space: local operators, L2 product and global interpolator. More...
 
struct  HArDCore3D::SXCurl::TransferOperators
 A structure to store the serendipity, extension and reduction operators. More...
 
class  HArDCore3D::SXDiv
 Discrete Serendipity Hdiv space: local operators, L2 product and global interpolator. More...
 
class  HArDCore3D::SXGrad
 Discrete Serendipity Hgrad space: local operators, L2 product and global interpolator. More...
 
struct  HArDCore3D::SXGrad::TransferOperators
 A structure to store the serendipity, extension and reduction operators. More...
 
class  HArDCore3D::XCurl
 Discrete Hcurl space: local operators, L2 product and global interpolator. More...
 
struct  HArDCore3D::XCurl::LocalOperators
 A structure to store the local operators (curl and potential) More...
 
class  HArDCore3D::XDiv
 Discrete Hdiv space: local operators, L2 product and global interpolator. More...
 
struct  HArDCore3D::XDiv::LocalOperators
 A structure to store local operators (divergence and potential) More...
 
class  HArDCore3D::XGrad
 Discrete H1 space: local operators, L2 product and global interpolator. More...
 
struct  HArDCore3D::XGrad::LocalOperators
 A structure to store local operators (gradient and potential) More...
 

Typedefs

typedef Family< MonomialScalarBasisCellHArDCore3D::DDRCore::PolyBasisCellType
 
typedef TensorizedVectorFamily< PolyBasisCellType, 3 > HArDCore3D::DDRCore::Poly3BasisCellType
 
typedef Family< GradientBasis< ShiftedBasis< MonomialScalarBasisCell > > > HArDCore3D::DDRCore::GolyBasisCellType
 
typedef Family< GolyComplBasisCellHArDCore3D::DDRCore::GolyComplBasisCellType
 
typedef Family< CurlBasis< GolyComplBasisCell > > HArDCore3D::DDRCore::RolyBasisCellType
 
typedef Family< RolyComplBasisCellHArDCore3D::DDRCore::RolyComplBasisCellType
 
typedef Family< MonomialScalarBasisFaceHArDCore3D::DDRCore::PolyBasisFaceType
 
typedef TangentFamily< PolyBasisFaceTypeHArDCore3D::DDRCore::Poly2BasisFaceType
 
typedef Family< CurlBasis< ShiftedBasis< MonomialScalarBasisFace > > > HArDCore3D::DDRCore::RolyBasisFaceType
 
typedef Family< RolyComplBasisFaceHArDCore3D::DDRCore::RolyComplBasisFaceType
 
typedef Family< MonomialScalarBasisEdgeHArDCore3D::DDRCore::PolyBasisEdgeType
 
typedef Cell HArDCore3D::DDRCore::CellBases::GeometricSupport
 Geometric support.
 
typedef Face HArDCore3D::DDRCore::FaceBases::GeometricSupport
 Geometric support.
 
typedef Edge HArDCore3D::DDRCore::EdgeBases::GeometricSupport
 Geometric support.
 
typedef RestrictedBasis< DDRCore::PolyBasisFaceTypeHArDCore3D::SerendipityProblem::PolylBasisFaceType
 
typedef RestrictedBasis< DDRCore::PolyBasisCellTypeHArDCore3D::SerendipityProblem::PolylBasisCellType
 
typedef RestrictedBasis< DDRCore::RolyComplBasisFaceTypeHArDCore3D::SerendipityProblem::RolyCompllpoBasisFaceType
 
typedef RestrictedBasis< DDRCore::RolyComplBasisCellTypeHArDCore3D::SerendipityProblem::RolyCompllpoBasisCellType
 
typedef Eigen::FullPivLU< Eigen::MatrixXd > HArDCore3D::SerendipityProblem::InverseProblem
 Type for inverses of matrix for serendipity problem.
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::SXCurl::FunctionType
 
typedef std::function< double(const Eigen::Vector3d &)> HArDCore3D::SXGrad::FunctionType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::XCurl::FunctionType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::XDiv::FunctionType
 
typedef std::function< double(const Eigen::Vector3d &)> HArDCore3D::XGrad::FunctionType
 

Functions

 HArDCore3D::DDRCore::DDRCore (const Mesh &mesh, size_t K, bool use_threads=true, std::ostream &output=std::cout)
 Constructor.
 
const MeshHArDCore3D::DDRCore::mesh () const
 Return a const reference to the mesh.
 
const size_tHArDCore3D::DDRCore::degree () const
 Return the polynomial degree.
 
const CellBasesHArDCore3D::DDRCore::cellBases (size_t iT) const
 Return cell bases for element iT.
 
const FaceBasesHArDCore3D::DDRCore::faceBases (size_t iF) const
 Return face bases for face iF.
 
const EdgeBasesHArDCore3D::DDRCore::edgeBases (size_t iE) const
 Return edge bases for edge iE.
 
 HArDCore3D::SerendipityProblem::SerendipityProblem (const DDRCore &ddrcore, bool use_threads=true, std::ostream &output=std::cout, bool use_serendipity=true)
 Constructor.
 
const std::vector< size_t > & HArDCore3D::SerendipityProblem::serendipityEdges (size_t iF) const
 Return the list of serendipity edges in a face.
 
const int HArDCore3D::SerendipityProblem::n_serendipityEdges (size_t iF) const
 Return the number of serendipity edges in a face.
 
const int HArDCore3D::SerendipityProblem::serDegreeFace (size_t iF) const
 Return the serendipity degree ell_F in a face.
 
const std::vector< size_t > & HArDCore3D::SerendipityProblem::serendipityFaces (size_t iT) const
 Return the list of serendipity face in a cell.
 
const int HArDCore3D::SerendipityProblem::n_serendipityFaces (size_t iT) const
 Return the number of serendipity faces in a cell.
 
const int HArDCore3D::SerendipityProblem::serDegreeCell (size_t iT) const
 Return the serendipity degree ell_T in a cell.
 
size_t HArDCore3D::SerendipityProblem::dimFacePolyl (size_t iF) const
 Return the dimension of P^l on face of index iF.
 
size_t HArDCore3D::SerendipityProblem::dimFacePolyl (const Face &F) const
 Return the dimension of P^{l+1} on face F.
 
const PolylBasisFaceTypeHArDCore3D::SerendipityProblem::faceBasisPolyl (size_t iF) const
 Return the basis of P^l on face of index iF.
 
const PolylBasisFaceTypeHArDCore3D::SerendipityProblem::faceBasisPolyl (const Face &F) const
 Return the basis of P^l on face F.
 
size_t HArDCore3D::SerendipityProblem::dimCellPolyl (size_t iT) const
 Return the dimension of P^l on cell of index iT.
 
size_t HArDCore3D::SerendipityProblem::dimCellPolyl (const Cell &T) const
 Return the dimension of P^l on cell T.
 
const PolylBasisCellTypeHArDCore3D::SerendipityProblem::cellBasisPolyl (size_t iT) const
 Return the basis of P^l on cell of index iT.
 
const PolylBasisCellTypeHArDCore3D::SerendipityProblem::cellBasisPolyl (const Cell &T) const
 Return the basis of P^l on cell T.
 
Eigen::VectorXi HArDCore3D::SerendipityProblem::nDOFs_faces_SXGrad () const
 Number of DOFs on faces for serendipity XGrad space.
 
Eigen::VectorXi HArDCore3D::SerendipityProblem::nDOFs_cells_SXGrad () const
 Number of DOFs on cells for serendipity XGrad space.
 
size_t HArDCore3D::SerendipityProblem::dimFaceRolyCompllpo (size_t iF) const
 Return the dimension of R^{c,l+1} on face of index iF.
 
size_t HArDCore3D::SerendipityProblem::dimFaceRolyCompllpo (const Face &F) const
 Return the dimension of R^{c,l+1} on face F.
 
const RolyCompllpoBasisFaceTypeHArDCore3D::SerendipityProblem::faceBasisRolyCompllpo (size_t iF) const
 Return the basis of R^{c,l+1} on face of index iF.
 
const RolyCompllpoBasisFaceTypeHArDCore3D::SerendipityProblem::faceBasisRolyCompllpo (const Face &F) const
 Return the basis of R^{c,l+1} on face F.
 
size_t HArDCore3D::SerendipityProblem::dimCellRolyCompllpo (size_t iT) const
 Return the dimension of R^{c,l+1} on cell of index iT.
 
size_t HArDCore3D::SerendipityProblem::dimCellRolyCompllpo (const Cell &T) const
 Return the dimension of R^{c,l+1} on cell T.
 
const RolyCompllpoBasisCellTypeHArDCore3D::SerendipityProblem::cellBasisRolyCompllpo (size_t iT) const
 Return the basis of R^{c,l+1} on cell of index iT.
 
const RolyCompllpoBasisCellTypeHArDCore3D::SerendipityProblem::cellBasisRolyCompllpo (const Cell &T) const
 Return the basis of R^{c,l+1} on cell T.
 
Eigen::VectorXi HArDCore3D::SerendipityProblem::nDOFs_faces_SXCurl () const
 Number of DOFs on faces for serendipity XCurl space.
 
Eigen::VectorXi HArDCore3D::SerendipityProblem::nDOFs_cells_SXCurl () const
 Number of DOFs on cells for serendipity XCurl space.
 
const Eigen::MatrixXd HArDCore3D::SerendipityProblem::SerendipityOperatorFace (const size_t iF, const Eigen::MatrixXd &LF) const
 Compute the serendipity operator on the face of index iF.
 
const Eigen::MatrixXd HArDCore3D::SerendipityProblem::SerendipityOperatorFace (const Face &F, const Eigen::MatrixXd &LF) const
 Compute the serendipity operator on the face F.
 
const Eigen::MatrixXd HArDCore3D::SerendipityProblem::SerendipityOperatorCell (const size_t iT, const Eigen::MatrixXd &LT) const
 Compute the serendipity operator on the cell of index iT.
 
const Eigen::MatrixXd HArDCore3D::SerendipityProblem::SerendipityOperatorCell (const Cell &T, const Eigen::MatrixXd &LT) const
 Compute the serendipity operator on the Cell T.
 
const MeshHArDCore3D::SerendipityProblem::mesh () const
 Return a const reference to the mesh.
 
const DDRCoreHArDCore3D::SerendipityProblem::ddrCore () const
 Return a const reference to the underlying DDR core.
 
 HArDCore3D::SXCurl::TransferOperators::TransferOperators (const Eigen::MatrixXd &_serendipity, const Eigen::MatrixXd &_extension, const Eigen::MatrixXd &_reduction)
 
 HArDCore3D::SXCurl::SXCurl (const DDRCore &ddr_core, const SerendipityProblem &ser_pro, bool use_threads=true, std::ostream &output=std::cout)
 Constructor.
 
const MeshHArDCore3D::SXCurl::mesh () const
 Return the mesh.
 
const size_tHArDCore3D::SXCurl::degree () const
 Return the polynomial degree.
 
const SerendipityProblemHArDCore3D::SXCurl::serPro () const
 
Eigen::VectorXd HArDCore3D::SXCurl::interpolate (const FunctionType &v, const int doe_cell=-1, const int doe_face=-1, const int doe_edge=-1) const
 Interpolator of a continuous function.
 
const Eigen::MatrixXd & HArDCore3D::SXCurl::ScurlFace (size_t iF) const
 Return the serendipity reconstruction for the face of index iF.
 
const Eigen::MatrixXd & HArDCore3D::SXCurl::EcurlFace (size_t iF) const
 Return the extension for the face of index iF.
 
const Eigen::MatrixXd & HArDCore3D::SXCurl::RcurlFace (size_t iF) const
 Return the reduction for the face of index iF.
 
const Eigen::MatrixXd & HArDCore3D::SXCurl::ScurlFace (const Face &F) const
 Return the serendipity reconstruction for face F.
 
const Eigen::MatrixXd & HArDCore3D::SXCurl::EcurlFace (const Face &F) const
 Return the extension for face F.
 
const Eigen::MatrixXd & HArDCore3D::SXCurl::RcurlFace (const Face &F) const
 Return cell reduction for cell T.
 
const Eigen::MatrixXd & HArDCore3D::SXCurl::ScurlCell (size_t iT) const
 Return the serendipity reconstruction for the cell of index iT.
 
const Eigen::MatrixXd & HArDCore3D::SXCurl::EcurlCell (size_t iT) const
 Return the extension for the cell of index iT.
 
const Eigen::MatrixXd & HArDCore3D::SXCurl::RcurlCell (size_t iT) const
 Return the reduction for the cell of index iT.
 
const Eigen::MatrixXd & HArDCore3D::SXCurl::ScurlCell (const Cell &T) const
 Return the serendipity reconstruction for cell T.
 
const Eigen::MatrixXd & HArDCore3D::SXCurl::EcurlCell (const Cell &T) const
 Return the extension for cell T.
 
const Eigen::MatrixXd & HArDCore3D::SXCurl::RcurlCell (const Cell &T) const
 Return the reduction for cell T.
 
const Eigen::MatrixXd HArDCore3D::SXCurl::faceCurl (size_t iF) const
 Return the full curl operator on the face of index iF.
 
const Eigen::MatrixXd HArDCore3D::SXCurl::faceCurl (const Face &F) const
 Return the full curl operator on face F.
 
const Eigen::MatrixXd HArDCore3D::SXCurl::facePotential (size_t iF) const
 Return the potential operator on the face of index iF.
 
const Eigen::MatrixXd HArDCore3D::SXCurl::facePotential (const Face &F) const
 Return the potential operator on face F.
 
const Eigen::MatrixXd HArDCore3D::SXCurl::cellCurl (size_t iT) const
 Return the full curl operator on the cell of index iT.
 
const Eigen::MatrixXd HArDCore3D::SXCurl::cellCurl (const Cell &T) const
 Return the full curl operator on cell T.
 
const Eigen::MatrixXd HArDCore3D::SXCurl::cellPotential (size_t iT) const
 Return the potential operator on the cell of index iT.
 
const Eigen::MatrixXd HArDCore3D::SXCurl::cellPotential (const Cell &T) const
 Return the potential operator on cell T.
 
Eigen::MatrixXd HArDCore3D::SXCurl::computeL2Product (const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk3_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const
 Compute the matrix of the (weighted) L2-product.
 
Eigen::MatrixXd HArDCore3D::SXCurl::computeL2ProductGradient (const size_t iT, const SXGrad &sx_grad, const std::string &side, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk3_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").
 
std::vector< VectorRd > HArDCore3D::SXCurl::computeVertexValues (const Eigen::VectorXd &u) const
 Computes the values of the potential reconstruction at the mesh vertices.
 
const DDRCore::CellBasesHArDCore3D::SXCurl::cellBases (size_t iT) const
 Return cell bases for the face of index iT.
 
const DDRCore::CellBasesHArDCore3D::SXCurl::cellBases (const Cell &T) const
 Return cell bases for cell T.
 
const DDRCore::FaceBasesHArDCore3D::SXCurl::faceBases (size_t iF) const
 Return face bases for the face of index iF.
 
const DDRCore::FaceBasesHArDCore3D::SXCurl::faceBases (const Face &F) const
 Return cell bases for face F.
 
const DDRCore::EdgeBasesHArDCore3D::SXCurl::edgeBases (size_t iE) const
 Return edge bases for the edge of index iE.
 
const DDRCore::EdgeBasesHArDCore3D::SXCurl::edgeBases (const Edge &E) const
 Return edge bases for edge E.
 
 HArDCore3D::SXDiv::SXDiv (const DDRCore &ddr_core, bool use_threads=true, std::ostream &output=std::cout)
 Constructor.
 
Eigen::MatrixXd HArDCore3D::SXDiv::computeL2ProductCurl (const size_t iT, const SXCurl &sx_curl, const std::string &side, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk3_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 serendipity curl on the left/right/both sides (depending on argument "side").
 
 HArDCore3D::SXGrad::TransferOperators::TransferOperators (const Eigen::MatrixXd &_serendipity, const Eigen::MatrixXd &_extension, const Eigen::MatrixXd &_reduction)
 
 HArDCore3D::SXGrad::SXGrad (const DDRCore &ddr_core, const SerendipityProblem &ser_pro, bool use_threads=true, std::ostream &output=std::cout)
 Constructor.
 
const MeshHArDCore3D::SXGrad::mesh () const
 Return the mesh.
 
const size_tHArDCore3D::SXGrad::degree () const
 Return the polynomial degree.
 
const SerendipityProblemHArDCore3D::SXGrad::serPro () const
 
Eigen::VectorXd HArDCore3D::SXGrad::interpolate (const FunctionType &q, const int doe_cell=-1, const int doe_face=-1, const int doe_edge=-1) const
 Interpolator of a continuous function.
 
const Eigen::MatrixXd & HArDCore3D::SXGrad::SgradFace (size_t iF) const
 Return the serendipity reconstruction for the face of index iF.
 
const Eigen::MatrixXd & HArDCore3D::SXGrad::EgradFace (size_t iF) const
 Return the extension for the face of index iF.
 
const Eigen::MatrixXd & HArDCore3D::SXGrad::RgradFace (size_t iF) const
 Return the reduction for the face of index iF.
 
const Eigen::MatrixXd & HArDCore3D::SXGrad::SgradFace (const Face &F) const
 Return the serendipity reconstruction for face F.
 
const Eigen::MatrixXd & HArDCore3D::SXGrad::EgradFace (const Face &F) const
 Return the extension for face F.
 
const Eigen::MatrixXd & HArDCore3D::SXGrad::RgradFace (const Face &F) const
 Return cell reduction for cell T.
 
const Eigen::MatrixXd & HArDCore3D::SXGrad::SgradCell (size_t iT) const
 Return the serendipity reconstruction for the cell of index iT.
 
const Eigen::MatrixXd & HArDCore3D::SXGrad::EgradCell (size_t iT) const
 Return the extension for the cell of index iT.
 
const Eigen::MatrixXd & HArDCore3D::SXGrad::RgradCell (size_t iT) const
 Return the reduction for the cell of index iT.
 
const Eigen::MatrixXd & HArDCore3D::SXGrad::SgradCell (const Cell &T) const
 Return the serendipity reconstruction for cell T.
 
const Eigen::MatrixXd & HArDCore3D::SXGrad::EgradCell (const Cell &T) const
 Return the extension for cell T.
 
const Eigen::MatrixXd & HArDCore3D::SXGrad::RgradCell (const Cell &T) const
 Return the reduction for cell T.
 
const Eigen::MatrixXd HArDCore3D::SXGrad::edgeGradient (size_t iE) const
 Return the full gradient operator on the edge of index iE.
 
const Eigen::MatrixXd HArDCore3D::SXGrad::edgeGradient (const Edge &E) const
 Return the full gradient operator on edge E.
 
const Eigen::MatrixXd HArDCore3D::SXGrad::edgePotential (size_t iE) const
 Return the potential operator on the edge of index iE.
 
const Eigen::MatrixXd HArDCore3D::SXGrad::edgePotential (const Edge &E) const
 Return the potential operator on edge E.
 
const Eigen::MatrixXd HArDCore3D::SXGrad::faceGradient (size_t iF) const
 Return the full gradient operator on the face of index iF.
 
const Eigen::MatrixXd HArDCore3D::SXGrad::faceGradient (const Face &F) const
 Return the full gradient operator on face F.
 
const Eigen::MatrixXd HArDCore3D::SXGrad::facePotential (size_t iF) const
 Return the potential operator on the face of index iF.
 
const Eigen::MatrixXd HArDCore3D::SXGrad::facePotential (const Face &F) const
 Return the potential operator on face F.
 
const Eigen::MatrixXd HArDCore3D::SXGrad::cellGradient (size_t iT) const
 Return the full gradient operator on the cell of index iT.
 
const Eigen::MatrixXd HArDCore3D::SXGrad::cellGradient (const Cell &T) const
 Return the full gradient operator on cell T.
 
const Eigen::MatrixXd HArDCore3D::SXGrad::cellPotential (size_t iT) const
 Return the potential operator on the cell of index iT.
 
const Eigen::MatrixXd HArDCore3D::SXGrad::cellPotential (const Cell &T) const
 Return the potential operator on cell T.
 
Eigen::MatrixXd HArDCore3D::SXGrad::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.
 
Eigen::MatrixXd HArDCore3D::SXGrad::computeStabilisation (const size_t iT, const IntegralWeight &weight=IntegralWeight(1.)) const
 Computes only the stabilisation matrix of the (weighted) L2-product for the cell of index iT.
 
std::vector< doubleHArDCore3D::SXGrad::computeVertexValues (const Eigen::VectorXd &u) const
 Computes the values of the potential reconstruction at the mesh vertices.
 
const DDRCore::CellBasesHArDCore3D::SXGrad::cellBases (size_t iT) const
 Return cell bases for the face of index iT.
 
const DDRCore::CellBasesHArDCore3D::SXGrad::cellBases (const Cell &T) const
 Return cell bases for cell T.
 
const DDRCore::FaceBasesHArDCore3D::SXGrad::faceBases (size_t iF) const
 Return face bases for the face of index iF.
 
const DDRCore::FaceBasesHArDCore3D::SXGrad::faceBases (const Face &F) const
 Return face bases for face F.
 
const DDRCore::EdgeBasesHArDCore3D::SXGrad::edgeBases (size_t iE) const
 Return edge bases for the edge of index iE.
 
const DDRCore::EdgeBasesHArDCore3D::SXGrad::edgeBases (const Edge &E) const
 Return edge bases for edge E.
 
 HArDCore3D::XCurl::LocalOperators::LocalOperators (const Eigen::MatrixXd &_curl, const Eigen::MatrixXd &_potential)
 
 HArDCore3D::XCurl::XCurl (const DDRCore &ddr_core, bool use_threads=true, std::ostream &output=std::cout)
 Constructor.
 
const MeshHArDCore3D::XCurl::mesh () const
 Return the mesh.
 
const size_tHArDCore3D::XCurl::degree () const
 Return the polynomial degree.
 
Eigen::VectorXd HArDCore3D::XCurl::interpolate (const FunctionType &v, const int doe_cell=-1, const int doe_face=-1, const int doe_edge=-1) const
 Interpolator of a continuous function.
 
const LocalOperatorsHArDCore3D::XCurl::cellOperators (size_t iT) const
 Return cell operators for the cell of index iT.
 
const LocalOperatorsHArDCore3D::XCurl::cellOperators (const Cell &T) const
 Return cell operators for cell T.
 
const LocalOperatorsHArDCore3D::XCurl::faceOperators (size_t iF) const
 Return face operators for the face of index iF.
 
const LocalOperatorsHArDCore3D::XCurl::faceOperators (const Face &F) const
 Return face operators for face F.
 
const DDRCore::CellBasesHArDCore3D::XCurl::cellBases (size_t iT) const
 Return cell bases for the face of index iT.
 
const DDRCore::CellBasesHArDCore3D::XCurl::cellBases (const Cell &T) const
 Return cell bases for cell T.
 
const DDRCore::FaceBasesHArDCore3D::XCurl::faceBases (size_t iF) const
 Return face bases for the face of index iF.
 
const DDRCore::FaceBasesHArDCore3D::XCurl::faceBases (const Face &F) const
 Return cell bases for face F.
 
const DDRCore::EdgeBasesHArDCore3D::XCurl::edgeBases (size_t iE) const
 Return edge bases for the edge of index iE.
 
const DDRCore::EdgeBasesHArDCore3D::XCurl::edgeBases (const Edge &E) const
 Return edge bases for edge E.
 
Eigen::MatrixXd HArDCore3D::XCurl::computeL2Product (const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk3_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. The stabilisation here is based on cell and face potentials.
 
Eigen::MatrixXd HArDCore3D::XCurl::computeL2ProductGradient (const size_t iT, const XGrad &x_grad, const std::string &side, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk3_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 HArDCore3D::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_Pk3_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.
 
Eigen::MatrixXd HArDCore3D::XCurl::computeL2ProductDOFs (size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk3_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.))
 Alternate version of L2 product, using only cell potentials and 'DOFs' (as in Di-Pietro & Droniou, JCP 2020)
 
 HArDCore3D::XDiv::LocalOperators::LocalOperators (const Eigen::MatrixXd &_divergence, const Eigen::MatrixXd &_divergence_rhs, const Eigen::MatrixXd &_potential)
 
 HArDCore3D::XDiv::XDiv (const DDRCore &ddr_core, bool use_threads=true, std::ostream &output=std::cout)
 Constructor.
 
const MeshHArDCore3D::XDiv::mesh () const
 Return the mesh.
 
const size_tHArDCore3D::XDiv::degree () const
 Return the polynomial degree.
 
Eigen::VectorXd HArDCore3D::XDiv::interpolate (const FunctionType &v, const int doe_cell=-1, const int doe_face=-1) const
 Interpolator of a continuous function.
 
const LocalOperatorsHArDCore3D::XDiv::cellOperators (size_t iT) const
 Return cell operators for the cell of index iT.
 
const LocalOperatorsHArDCore3D::XDiv::cellOperators (const Cell &T) const
 Return cell operators for cell T.
 
const DDRCore::CellBasesHArDCore3D::XDiv::cellBases (size_t iT) const
 Return cell bases for the face of index iT.
 
const DDRCore::CellBasesHArDCore3D::XDiv::cellBases (const Cell &T) const
 Return cell bases for cell T.
 
const DDRCore::FaceBasesHArDCore3D::XDiv::faceBases (size_t iF) const
 Return face bases for the face of index iF.
 
const DDRCore::FaceBasesHArDCore3D::XDiv::faceBases (const Face &F) const
 Return cell bases for face F.
 
Eigen::MatrixXd HArDCore3D::XDiv::computeL2Product (const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk3_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 HArDCore3D::XDiv::computeL2ProductCurl (const size_t iT, const XCurl &x_curl, const std::string &side, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk3_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 curl on the left/right/both sides (depending on argument "side").
 
Eigen::MatrixXd HArDCore3D::XDiv::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_Pk3_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 computeL2ProductCurl.
 
LocalOperators HArDCore3D::XDiv::_compute_cell_divergence_potential (size_t iT)
 
 HArDCore3D::XGrad::LocalOperators::LocalOperators (const Eigen::MatrixXd &_gradient, const Eigen::MatrixXd &_potential)
 
 HArDCore3D::XGrad::XGrad (const DDRCore &ddr_core, bool use_threads=true, std::ostream &output=std::cout)
 Constructor.
 
const MeshHArDCore3D::XGrad::mesh () const
 Return the mesh.
 
const size_tHArDCore3D::XGrad::degree () const
 Return the polynomial degree.
 
Eigen::VectorXd HArDCore3D::XGrad::interpolate (const FunctionType &q, const int doe_cell=-1, const int doe_face=-1, const int doe_edge=-1) const
 Interpolator of a continuous function.
 
const LocalOperatorsHArDCore3D::XGrad::edgeOperators (size_t iE) const
 Return edge operators for the edge of index iE.
 
const LocalOperatorsHArDCore3D::XGrad::edgeOperators (const Edge &E) const
 Return edge operators for edge E.
 
const LocalOperatorsHArDCore3D::XGrad::faceOperators (size_t iF) const
 Return face operators for the face of index iF.
 
const LocalOperatorsHArDCore3D::XGrad::faceOperators (const Face &F) const
 Return face operators for face F.
 
const LocalOperatorsHArDCore3D::XGrad::cellOperators (size_t iT) const
 Return cell operators for the cell of index iT.
 
const LocalOperatorsHArDCore3D::XGrad::cellOperators (const Cell &T) const
 Return cell operators for cell T.
 
const DDRCore::CellBasesHArDCore3D::XGrad::cellBases (size_t iT) const
 Return cell bases for the cell of index iT.
 
const DDRCore::CellBasesHArDCore3D::XGrad::cellBases (const Cell &T) const
 Return cell bases for cell T.
 
const DDRCore::FaceBasesHArDCore3D::XGrad::faceBases (size_t iF) const
 Return face bases for the face of index iF.
 
const DDRCore::FaceBasesHArDCore3D::XGrad::faceBases (const Face &F) const
 Return cell bases for face F.
 
const DDRCore::EdgeBasesHArDCore3D::XGrad::edgeBases (size_t iE) const
 Return edge bases for the edge of index iE.
 
const DDRCore::EdgeBasesHArDCore3D::XGrad::edgeBases (const Edge &E) const
 Return edge bases for edge E.
 
Eigen::MatrixXd HArDCore3D::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. The stabilisation here is based on cell and face potentials.
 
Eigen::MatrixXd HArDCore3D::XGrad::computeStabilisation (const size_t iT, const IntegralWeight &weight=IntegralWeight(1.)) const
 Computes only the stabilisation matrix of the (weighted) L2-product for the cell of index iT. This stabilisation is based on cell and face potentials.
 
Eigen::MatrixXd HArDCore3D::XGrad::buildGradientComponentsCell (size_t iT) const
 Build the components of the gradient operator (probably not useful in practice to implement schemes)
 

Variables

std::unique_ptr< PolyBasisCellTypeHArDCore3D::DDRCore::CellBases::Polykpo
 
std::unique_ptr< PolyBasisCellTypeHArDCore3D::DDRCore::CellBases::Polyk
 
std::unique_ptr< PolyBasisCellTypeHArDCore3D::DDRCore::CellBases::Polykmo
 
std::unique_ptr< Poly3BasisCellTypeHArDCore3D::DDRCore::CellBases::Polyk3
 
std::unique_ptr< GolyBasisCellTypeHArDCore3D::DDRCore::CellBases::Golykmo
 
std::unique_ptr< GolyComplBasisCellTypeHArDCore3D::DDRCore::CellBases::GolyComplk
 
std::unique_ptr< GolyComplBasisCellTypeHArDCore3D::DDRCore::CellBases::GolyComplkpo
 
std::unique_ptr< RolyBasisCellTypeHArDCore3D::DDRCore::CellBases::Rolykmo
 
std::unique_ptr< RolyComplBasisCellTypeHArDCore3D::DDRCore::CellBases::RolyComplk
 
std::unique_ptr< RolyComplBasisCellTypeHArDCore3D::DDRCore::CellBases::RolyComplkp2
 
std::unique_ptr< PolyBasisFaceTypeHArDCore3D::DDRCore::FaceBases::Polykpo
 
std::unique_ptr< PolyBasisFaceTypeHArDCore3D::DDRCore::FaceBases::Polyk
 
std::unique_ptr< PolyBasisFaceTypeHArDCore3D::DDRCore::FaceBases::Polykmo
 
std::unique_ptr< Poly2BasisFaceTypeHArDCore3D::DDRCore::FaceBases::Polyk2
 
std::unique_ptr< RolyBasisFaceTypeHArDCore3D::DDRCore::FaceBases::Rolykmo
 
std::unique_ptr< RolyComplBasisFaceTypeHArDCore3D::DDRCore::FaceBases::RolyComplk
 
std::unique_ptr< RolyComplBasisFaceTypeHArDCore3D::DDRCore::FaceBases::RolyComplkp2
 
std::unique_ptr< PolyBasisEdgeTypeHArDCore3D::DDRCore::EdgeBases::Polykpo
 
std::unique_ptr< PolyBasisEdgeTypeHArDCore3D::DDRCore::EdgeBases::Polyk
 
std::unique_ptr< PolyBasisEdgeTypeHArDCore3D::DDRCore::EdgeBases::Polykmo
 
Eigen::MatrixXd HArDCore3D::SXCurl::TransferOperators::serendipity
 
Eigen::MatrixXd HArDCore3D::SXCurl::TransferOperators::extension
 
Eigen::MatrixXd HArDCore3D::SXCurl::TransferOperators::reduction
 
Eigen::MatrixXd HArDCore3D::SXGrad::TransferOperators::serendipity
 
Eigen::MatrixXd HArDCore3D::SXGrad::TransferOperators::extension
 
Eigen::MatrixXd HArDCore3D::SXGrad::TransferOperators::reduction
 
Eigen::MatrixXd HArDCore3D::XCurl::LocalOperators::curl
 
Eigen::MatrixXd HArDCore3D::XCurl::LocalOperators::potential
 
Eigen::MatrixXd HArDCore3D::XDiv::LocalOperators::divergence
 
Eigen::MatrixXd HArDCore3D::XDiv::LocalOperators::divergence_rhs
 
Eigen::MatrixXd HArDCore3D::XDiv::LocalOperators::potential
 
const DDRCoreHArDCore3D::XDiv::m_ddr_core
 
bool HArDCore3D::XDiv::m_use_threads
 
std::ostream & HArDCore3D::XDiv::m_output
 
std::vector< std::unique_ptr< LocalOperators > > HArDCore3D::XDiv::m_cell_operators
 
Eigen::MatrixXd HArDCore3D::XGrad::LocalOperators::gradient
 
Eigen::MatrixXd HArDCore3D::XGrad::LocalOperators::potential
 

Detailed Description

Classes defining the Discrete De Rham sequence.

Typedef Documentation

◆ FunctionType [1/5]

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::SXCurl::FunctionType

◆ FunctionType [2/5]

typedef std::function<double(const Eigen::Vector3d &)> HArDCore3D::SXGrad::FunctionType

◆ FunctionType [3/5]

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::XCurl::FunctionType

◆ FunctionType [4/5]

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::XDiv::FunctionType

◆ FunctionType [5/5]

typedef std::function<double(const Eigen::Vector3d &)> HArDCore3D::XGrad::FunctionType

◆ GeometricSupport [1/3]

Geometric support.

◆ GeometricSupport [2/3]

Geometric support.

◆ GeometricSupport [3/3]

Geometric support.

◆ GolyBasisCellType

◆ GolyComplBasisCellType

◆ InverseProblem

typedef Eigen::FullPivLU<Eigen::MatrixXd> HArDCore3D::SerendipityProblem::InverseProblem

Type for inverses of matrix for serendipity problem.

◆ Poly2BasisFaceType

◆ Poly3BasisCellType

◆ PolyBasisCellType

◆ PolyBasisEdgeType

◆ PolyBasisFaceType

◆ PolylBasisCellType

◆ PolylBasisFaceType

◆ RolyBasisCellType

◆ RolyBasisFaceType

◆ RolyComplBasisCellType

◆ RolyComplBasisFaceType

◆ RolyCompllpoBasisCellType

◆ RolyCompllpoBasisFaceType

Function Documentation

◆ _compute_cell_divergence_potential()

XDiv::LocalOperators XDiv::_compute_cell_divergence_potential ( size_t  iT)
protected

◆ buildGradientComponentsCell()

Eigen::MatrixXd XGrad::buildGradientComponentsCell ( size_t  iT) const

Build the components of the gradient operator (probably not useful in practice to implement schemes)

◆ cellBases() [1/11]

const DDRCore::CellBases & HArDCore3D::SXCurl::cellBases ( const Cell &  T) const
inline

Return cell bases for cell T.

◆ cellBases() [2/11]

const DDRCore::CellBases & HArDCore3D::SXGrad::cellBases ( const Cell &  T) const
inline

Return cell bases for cell T.

◆ cellBases() [3/11]

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

Return cell bases for cell T.

◆ cellBases() [4/11]

const DDRCore::CellBases & HArDCore3D::XDiv::cellBases ( const Cell &  T) const
inline

Return cell bases for cell T.

◆ cellBases() [5/11]

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

Return cell bases for cell T.

◆ cellBases() [6/11]

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

Return cell bases for element iT.

◆ cellBases() [7/11]

const DDRCore::CellBases & HArDCore3D::SXCurl::cellBases ( size_t  iT) const
inline

Return cell bases for the face of index iT.

◆ cellBases() [8/11]

const DDRCore::CellBases & HArDCore3D::SXGrad::cellBases ( size_t  iT) const
inline

Return cell bases for the face of index iT.

◆ cellBases() [9/11]

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

Return cell bases for the face of index iT.

◆ cellBases() [10/11]

const DDRCore::CellBases & HArDCore3D::XDiv::cellBases ( size_t  iT) const
inline

Return cell bases for the face of index iT.

◆ cellBases() [11/11]

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

Return cell bases for the cell of index iT.

◆ cellBasisPolyl() [1/2]

const PolylBasisCellType & HArDCore3D::SerendipityProblem::cellBasisPolyl ( const Cell &  T) const
inline

Return the basis of P^l on cell T.

◆ cellBasisPolyl() [2/2]

const PolylBasisCellType & HArDCore3D::SerendipityProblem::cellBasisPolyl ( size_t  iT) const
inline

Return the basis of P^l on cell of index iT.

◆ cellBasisRolyCompllpo() [1/2]

const RolyCompllpoBasisCellType & HArDCore3D::SerendipityProblem::cellBasisRolyCompllpo ( const Cell &  T) const
inline

Return the basis of R^{c,l+1} on cell T.

◆ cellBasisRolyCompllpo() [2/2]

const RolyCompllpoBasisCellType & HArDCore3D::SerendipityProblem::cellBasisRolyCompllpo ( size_t  iT) const
inline

Return the basis of R^{c,l+1} on cell of index iT.

◆ cellCurl() [1/2]

const Eigen::MatrixXd HArDCore3D::SXCurl::cellCurl ( const Cell &  T) const
inline

Return the full curl operator on cell T.

◆ cellCurl() [2/2]

const Eigen::MatrixXd HArDCore3D::SXCurl::cellCurl ( size_t  iT) const
inline

Return the full curl operator on the cell of index iT.

◆ cellGradient() [1/2]

const Eigen::MatrixXd HArDCore3D::SXGrad::cellGradient ( const Cell &  T) const
inline

Return the full gradient operator on cell T.

◆ cellGradient() [2/2]

const Eigen::MatrixXd HArDCore3D::SXGrad::cellGradient ( size_t  iT) const
inline

Return the full gradient operator on the cell of index iT.

◆ cellOperators() [1/6]

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

Return cell operators for cell T.

◆ cellOperators() [2/6]

const LocalOperators & HArDCore3D::XDiv::cellOperators ( const Cell &  T) const
inline

Return cell operators for cell T.

◆ cellOperators() [3/6]

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

Return cell operators for cell T.

◆ cellOperators() [4/6]

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

Return cell operators for the cell of index iT.

◆ cellOperators() [5/6]

const LocalOperators & HArDCore3D::XDiv::cellOperators ( size_t  iT) const
inline

Return cell operators for the cell of index iT.

◆ cellOperators() [6/6]

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

Return cell operators for the cell of index iT.

◆ cellPotential() [1/4]

const Eigen::MatrixXd HArDCore3D::SXCurl::cellPotential ( const Cell &  T) const
inline

Return the potential operator on cell T.

◆ cellPotential() [2/4]

const Eigen::MatrixXd HArDCore3D::SXGrad::cellPotential ( const Cell &  T) const
inline

Return the potential operator on cell T.

◆ cellPotential() [3/4]

const Eigen::MatrixXd HArDCore3D::SXCurl::cellPotential ( size_t  iT) const
inline

Return the potential operator on the cell of index iT.

◆ cellPotential() [4/4]

const Eigen::MatrixXd HArDCore3D::SXGrad::cellPotential ( size_t  iT) const
inline

Return the potential operator on the cell of index iT.

◆ computeL2Product() [1/5]

Eigen::MatrixXd HArDCore3D::SXCurl::computeL2Product ( const size_t  iT,
const double penalty_factor = 1.,
const Eigen::MatrixXd &  mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
const IntegralWeight weight = IntegralWeight(1.) 
) const
inline

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

Parameters
iTindex of the cell
penalty_factorpre-factor for stabilisation term
mass_Pk3_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/5]

Eigen::MatrixXd XCurl::computeL2Product ( const size_t  iT,
const double penalty_factor = 1.,
const Eigen::MatrixXd &  mass_Pk3_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. The stabilisation here is based on cell and face potentials.

Parameters
iTindex of the cell
penalty_factorpre-factor for stabilisation term
mass_Pk3_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/5]

Eigen::MatrixXd XDiv::computeL2Product ( const size_t  iT,
const double penalty_factor = 1.,
const Eigen::MatrixXd &  mass_Pk3_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_Pk3_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() [4/5]

Eigen::MatrixXd HArDCore3D::SXGrad::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
inline

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

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() [5/5]

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. The stabilisation here is based on cell and face potentials.

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_with_Ops() [1/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_Pk3_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, face and element operators to apply on the left
rightOpedge, face and element operators to apply on the right
penalty_factorpre-factor for stabilisation term
w_mass_Pk3_Tmass matrix of (P^k(T))^3 weighted by weight
weightweight function in the L2 product

◆ computeL2Product_with_Ops() [2/2]

Eigen::MatrixXd XDiv::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_Pk3_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 computeL2ProductCurl.

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

◆ computeL2ProductCurl() [1/2]

Eigen::MatrixXd SXDiv::computeL2ProductCurl ( const size_t  iT,
const SXCurl sx_curl,
const std::string &  side,
const double penalty_factor = 1.,
const Eigen::MatrixXd &  mass_Pk3_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 serendipity curl on the left/right/both sides (depending on argument "side").

Parameters
iTindex of the cell
sx_curlinstance of SXCurl to access the serendipity curls
sidewhich side (left, right, both) we apply the gradient to
penalty_factorpre-factor for stabilisation term
mass_Pk3_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.

◆ computeL2ProductCurl() [2/2]

Eigen::MatrixXd XDiv::computeL2ProductCurl ( const size_t  iT,
const XCurl x_curl,
const std::string &  side,
const double penalty_factor = 1.,
const Eigen::MatrixXd &  mass_Pk3_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 curl on the left/right/both sides (depending on argument "side").

Parameters
iTindex of the cell
x_curlinstance of XCurl to access the full curls
sidewhich side (left, right, both) we apply the gradient to
penalty_factorpre-factor for stabilisation term
mass_Pk3_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.

◆ computeL2ProductDOFs()

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

Alternate version of L2 product, using only cell potentials and 'DOFs' (as in Di-Pietro & Droniou, JCP 2020)

Parameters
iTindex of the cell
penalty_factorpre-factor for stabilisation term
mass_Pk3_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, default to constant equal to 1.

◆ computeL2ProductGradient() [1/2]

Eigen::MatrixXd SXCurl::computeL2ProductGradient ( const size_t  iT,
const SXGrad sx_grad,
const std::string &  side,
const double penalty_factor = 1.,
const Eigen::MatrixXd &  mass_Pk3_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
sx_gradinstance of SXGrad to access the full gradients
sidewhich side (left, right, both) we apply the gradient to
penalty_factorpre-factor for stabilisation term
mass_Pk3_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_Pk3_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_Pk3_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.

◆ computeStabilisation() [1/2]

Eigen::MatrixXd HArDCore3D::SXGrad::computeStabilisation ( const size_t  iT,
const IntegralWeight weight = IntegralWeight(1.) 
) const
inline

Computes only the stabilisation matrix of the (weighted) L2-product for the cell of index iT.

Parameters
iTindex of the cell
weightweight function in the stabilisation, defaults to 1

◆ computeStabilisation() [2/2]

Eigen::MatrixXd XGrad::computeStabilisation ( const size_t  iT,
const IntegralWeight weight = IntegralWeight(1.) 
) const

Computes only the stabilisation matrix of the (weighted) L2-product for the cell of index iT. This stabilisation is based on cell and face potentials.

Parameters
iTindex of the cell
weightweight function in the stabilisation, defaults to 1

◆ computeVertexValues() [1/2]

std::vector< VectorRd > SXCurl::computeVertexValues ( const Eigen::VectorXd &  u) const

Computes the values of the potential reconstruction at the mesh vertices.

Parameters
uDOFs in the discrete space

◆ computeVertexValues() [2/2]

std::vector< double > SXGrad::computeVertexValues ( const Eigen::VectorXd &  u) const

Computes the values of the potential reconstruction at the mesh vertices.

Parameters
uDOFs in the discrete space

◆ ddrCore()

const DDRCore & HArDCore3D::SerendipityProblem::ddrCore ( ) const
inline

Return a const reference to the underlying DDR core.

◆ DDRCore()

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

Constructor.

◆ degree() [1/6]

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

Return the polynomial degree.

◆ degree() [2/6]

const size_t & HArDCore3D::SXCurl::degree ( ) const
inline

Return the polynomial degree.

◆ degree() [3/6]

const size_t & HArDCore3D::SXGrad::degree ( ) const
inline

Return the polynomial degree.

◆ degree() [4/6]

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

Return the polynomial degree.

◆ degree() [5/6]

const size_t & HArDCore3D::XDiv::degree ( ) const
inline

Return the polynomial degree.

◆ degree() [6/6]

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

Return the polynomial degree.

◆ dimCellPolyl() [1/2]

size_t HArDCore3D::SerendipityProblem::dimCellPolyl ( const Cell &  T) const
inline

Return the dimension of P^l on cell T.

◆ dimCellPolyl() [2/2]

size_t HArDCore3D::SerendipityProblem::dimCellPolyl ( size_t  iT) const
inline

Return the dimension of P^l on cell of index iT.

◆ dimCellRolyCompllpo() [1/2]

size_t HArDCore3D::SerendipityProblem::dimCellRolyCompllpo ( const Cell &  T) const
inline

Return the dimension of R^{c,l+1} on cell T.

◆ dimCellRolyCompllpo() [2/2]

size_t HArDCore3D::SerendipityProblem::dimCellRolyCompllpo ( size_t  iT) const
inline

Return the dimension of R^{c,l+1} on cell of index iT.

◆ dimFacePolyl() [1/2]

size_t HArDCore3D::SerendipityProblem::dimFacePolyl ( const Face &  F) const
inline

Return the dimension of P^{l+1} on face F.

◆ dimFacePolyl() [2/2]

size_t HArDCore3D::SerendipityProblem::dimFacePolyl ( size_t  iF) const
inline

Return the dimension of P^l on face of index iF.

◆ dimFaceRolyCompllpo() [1/2]

size_t HArDCore3D::SerendipityProblem::dimFaceRolyCompllpo ( const Face &  F) const
inline

Return the dimension of R^{c,l+1} on face F.

◆ dimFaceRolyCompllpo() [2/2]

size_t HArDCore3D::SerendipityProblem::dimFaceRolyCompllpo ( size_t  iF) const
inline

Return the dimension of R^{c,l+1} on face of index iF.

◆ EcurlCell() [1/2]

const Eigen::MatrixXd & HArDCore3D::SXCurl::EcurlCell ( const Cell &  T) const
inline

Return the extension for cell T.

◆ EcurlCell() [2/2]

const Eigen::MatrixXd & HArDCore3D::SXCurl::EcurlCell ( size_t  iT) const
inline

Return the extension for the cell of index iT.

◆ EcurlFace() [1/2]

const Eigen::MatrixXd & HArDCore3D::SXCurl::EcurlFace ( const Face &  F) const
inline

Return the extension for face F.

◆ EcurlFace() [2/2]

const Eigen::MatrixXd & HArDCore3D::SXCurl::EcurlFace ( size_t  iF) const
inline

Return the extension for the face of index iF.

◆ edgeBases() [1/9]

const DDRCore::EdgeBases & HArDCore3D::SXCurl::edgeBases ( const Edge &  E) const
inline

Return edge bases for edge E.

◆ edgeBases() [2/9]

const DDRCore::EdgeBases & HArDCore3D::SXGrad::edgeBases ( const Edge &  E) const
inline

Return edge bases for edge E.

◆ edgeBases() [3/9]

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

Return edge bases for edge E.

◆ edgeBases() [4/9]

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

Return edge bases for edge E.

◆ edgeBases() [5/9]

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

Return edge bases for edge iE.

◆ edgeBases() [6/9]

const DDRCore::EdgeBases & HArDCore3D::SXCurl::edgeBases ( size_t  iE) const
inline

Return edge bases for the edge of index iE.

◆ edgeBases() [7/9]

const DDRCore::EdgeBases & HArDCore3D::SXGrad::edgeBases ( size_t  iE) const
inline

Return edge bases for the edge of index iE.

◆ edgeBases() [8/9]

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

Return edge bases for the edge of index iE.

◆ edgeBases() [9/9]

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

Return edge bases for the edge of index iE.

◆ edgeGradient() [1/2]

const Eigen::MatrixXd HArDCore3D::SXGrad::edgeGradient ( const Edge &  E) const
inline

Return the full gradient operator on edge E.

◆ edgeGradient() [2/2]

const Eigen::MatrixXd HArDCore3D::SXGrad::edgeGradient ( size_t  iE) const
inline

Return the full gradient operator on the edge of index iE.

◆ edgeOperators() [1/2]

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

Return edge operators for edge E.

◆ edgeOperators() [2/2]

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

Return edge operators for the edge of index iE.

◆ edgePotential() [1/2]

const Eigen::MatrixXd HArDCore3D::SXGrad::edgePotential ( const Edge &  E) const
inline

Return the potential operator on edge E.

◆ edgePotential() [2/2]

const Eigen::MatrixXd HArDCore3D::SXGrad::edgePotential ( size_t  iE) const
inline

Return the potential operator on the edge of index iE.

◆ EgradCell() [1/2]

const Eigen::MatrixXd & HArDCore3D::SXGrad::EgradCell ( const Cell &  T) const
inline

Return the extension for cell T.

◆ EgradCell() [2/2]

const Eigen::MatrixXd & HArDCore3D::SXGrad::EgradCell ( size_t  iT) const
inline

Return the extension for the cell of index iT.

◆ EgradFace() [1/2]

const Eigen::MatrixXd & HArDCore3D::SXGrad::EgradFace ( const Face &  F) const
inline

Return the extension for face F.

◆ EgradFace() [2/2]

const Eigen::MatrixXd & HArDCore3D::SXGrad::EgradFace ( size_t  iF) const
inline

Return the extension for the face of index iF.

◆ faceBases() [1/11]

const DDRCore::FaceBases & HArDCore3D::SXCurl::faceBases ( const Face &  F) const
inline

Return cell bases for face F.

◆ faceBases() [2/11]

const DDRCore::FaceBases & HArDCore3D::SXGrad::faceBases ( const Face &  F) const
inline

Return face bases for face F.

◆ faceBases() [3/11]

const DDRCore::FaceBases & HArDCore3D::XCurl::faceBases ( const Face &  F) const
inline

Return cell bases for face F.

◆ faceBases() [4/11]

const DDRCore::FaceBases & HArDCore3D::XDiv::faceBases ( const Face &  F) const
inline

Return cell bases for face F.

◆ faceBases() [5/11]

const DDRCore::FaceBases & HArDCore3D::XGrad::faceBases ( const Face &  F) const
inline

Return cell bases for face F.

◆ faceBases() [6/11]

const FaceBases & HArDCore3D::DDRCore::faceBases ( size_t  iF) const
inline

Return face bases for face iF.

◆ faceBases() [7/11]

const DDRCore::FaceBases & HArDCore3D::SXCurl::faceBases ( size_t  iF) const
inline

Return face bases for the face of index iF.

◆ faceBases() [8/11]

const DDRCore::FaceBases & HArDCore3D::SXGrad::faceBases ( size_t  iF) const
inline

Return face bases for the face of index iF.

◆ faceBases() [9/11]

const DDRCore::FaceBases & HArDCore3D::XCurl::faceBases ( size_t  iF) const
inline

Return face bases for the face of index iF.

◆ faceBases() [10/11]

const DDRCore::FaceBases & HArDCore3D::XDiv::faceBases ( size_t  iF) const
inline

Return face bases for the face of index iF.

◆ faceBases() [11/11]

const DDRCore::FaceBases & HArDCore3D::XGrad::faceBases ( size_t  iF) const
inline

Return face bases for the face of index iF.

◆ faceBasisPolyl() [1/2]

const PolylBasisFaceType & HArDCore3D::SerendipityProblem::faceBasisPolyl ( const Face &  F) const
inline

Return the basis of P^l on face F.

◆ faceBasisPolyl() [2/2]

const PolylBasisFaceType & HArDCore3D::SerendipityProblem::faceBasisPolyl ( size_t  iF) const
inline

Return the basis of P^l on face of index iF.

◆ faceBasisRolyCompllpo() [1/2]

const RolyCompllpoBasisFaceType & HArDCore3D::SerendipityProblem::faceBasisRolyCompllpo ( const Face &  F) const
inline

Return the basis of R^{c,l+1} on face F.

◆ faceBasisRolyCompllpo() [2/2]

const RolyCompllpoBasisFaceType & HArDCore3D::SerendipityProblem::faceBasisRolyCompllpo ( size_t  iF) const
inline

Return the basis of R^{c,l+1} on face of index iF.

◆ faceCurl() [1/2]

const Eigen::MatrixXd HArDCore3D::SXCurl::faceCurl ( const Face &  F) const
inline

Return the full curl operator on face F.

◆ faceCurl() [2/2]

const Eigen::MatrixXd HArDCore3D::SXCurl::faceCurl ( size_t  iF) const
inline

Return the full curl operator on the face of index iF.

◆ faceGradient() [1/2]

const Eigen::MatrixXd HArDCore3D::SXGrad::faceGradient ( const Face &  F) const
inline

Return the full gradient operator on face F.

◆ faceGradient() [2/2]

const Eigen::MatrixXd HArDCore3D::SXGrad::faceGradient ( size_t  iF) const
inline

Return the full gradient operator on the face of index iF.

◆ faceOperators() [1/4]

const LocalOperators & HArDCore3D::XCurl::faceOperators ( const Face &  F) const
inline

Return face operators for face F.

◆ faceOperators() [2/4]

const LocalOperators & HArDCore3D::XGrad::faceOperators ( const Face &  F) const
inline

Return face operators for face F.

◆ faceOperators() [3/4]

const LocalOperators & HArDCore3D::XCurl::faceOperators ( size_t  iF) const
inline

Return face operators for the face of index iF.

◆ faceOperators() [4/4]

const LocalOperators & HArDCore3D::XGrad::faceOperators ( size_t  iF) const
inline

Return face operators for the face of index iF.

◆ facePotential() [1/4]

const Eigen::MatrixXd HArDCore3D::SXCurl::facePotential ( const Face &  F) const
inline

Return the potential operator on face F.

◆ facePotential() [2/4]

const Eigen::MatrixXd HArDCore3D::SXGrad::facePotential ( const Face &  F) const
inline

Return the potential operator on face F.

◆ facePotential() [3/4]

const Eigen::MatrixXd HArDCore3D::SXCurl::facePotential ( size_t  iF) const
inline

Return the potential operator on the face of index iF.

◆ facePotential() [4/4]

const Eigen::MatrixXd HArDCore3D::SXGrad::facePotential ( size_t  iF) const
inline

Return the potential operator on the face of index iF.

◆ interpolate() [1/5]

Eigen::VectorXd SXGrad::interpolate ( const FunctionType q,
const int  doe_cell = -1,
const int  doe_face = -1,
const int  doe_edge = -1 
) const

Interpolator of a continuous function.

Parameters
qThe function 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.
doe_edgeThe optional degre of edge quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used.

◆ interpolate() [2/5]

Eigen::VectorXd XGrad::interpolate ( const FunctionType q,
const int  doe_cell = -1,
const int  doe_face = -1,
const int  doe_edge = -1 
) const

Interpolator of a continuous function.

Parameters
qThe function 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.
doe_edgeThe optional degre of edge quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used.

◆ interpolate() [3/5]

Eigen::VectorXd XDiv::interpolate ( const FunctionType v,
const int  doe_cell = -1,
const int  doe_face = -1 
) const

Interpolator of a continuous function.

Parameters
vThe function 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.

◆ interpolate() [4/5]

Eigen::VectorXd SXCurl::interpolate ( const FunctionType v,
const int  doe_cell = -1,
const int  doe_face = -1,
const int  doe_edge = -1 
) const

Interpolator of a continuous function.

Parameters
vThe function 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.
doe_edgeThe optional degre of edge quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used.

◆ interpolate() [5/5]

Eigen::VectorXd XCurl::interpolate ( const FunctionType v,
const int  doe_cell = -1,
const int  doe_face = -1,
const int  doe_edge = -1 
) const

Interpolator of a continuous function.

Parameters
vThe function 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.
doe_edgeThe optional degre of edge quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used.

◆ LocalOperators() [1/3]

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

◆ LocalOperators() [2/3]

HArDCore3D::XDiv::LocalOperators::LocalOperators ( const Eigen::MatrixXd &  _divergence,
const Eigen::MatrixXd &  _divergence_rhs,
const Eigen::MatrixXd &  _potential 
)
inline
Parameters
_divergenceDivergence operator
_divergence_rhsRHS for the divergence operator problem
_potentialPotential operator

◆ LocalOperators() [3/3]

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

◆ mesh() [1/7]

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

Return a const reference to the mesh.

◆ mesh() [2/7]

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

Return a const reference to the mesh.

◆ mesh() [3/7]

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

Return the mesh.

◆ mesh() [4/7]

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

Return the mesh.

◆ mesh() [5/7]

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

Return the mesh.

◆ mesh() [6/7]

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

Return the mesh.

◆ mesh() [7/7]

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

Return the mesh.

◆ n_serendipityEdges()

const int HArDCore3D::SerendipityProblem::n_serendipityEdges ( size_t  iF) const
inline

Return the number of serendipity edges in a face.

◆ n_serendipityFaces()

const int HArDCore3D::SerendipityProblem::n_serendipityFaces ( size_t  iT) const
inline

Return the number of serendipity faces in a cell.

◆ nDOFs_cells_SXCurl()

Eigen::VectorXi SerendipityProblem::nDOFs_cells_SXCurl ( ) const

Number of DOFs on cells for serendipity XCurl space.

◆ nDOFs_cells_SXGrad()

Eigen::VectorXi SerendipityProblem::nDOFs_cells_SXGrad ( ) const

Number of DOFs on cells for serendipity XGrad space.

◆ nDOFs_faces_SXCurl()

Eigen::VectorXi SerendipityProblem::nDOFs_faces_SXCurl ( ) const

Number of DOFs on faces for serendipity XCurl space.

◆ nDOFs_faces_SXGrad()

Eigen::VectorXi SerendipityProblem::nDOFs_faces_SXGrad ( ) const

Number of DOFs on faces for serendipity XGrad space.

◆ RcurlCell() [1/2]

const Eigen::MatrixXd & HArDCore3D::SXCurl::RcurlCell ( const Cell &  T) const
inline

Return the reduction for cell T.

◆ RcurlCell() [2/2]

const Eigen::MatrixXd & HArDCore3D::SXCurl::RcurlCell ( size_t  iT) const
inline

Return the reduction for the cell of index iT.

◆ RcurlFace() [1/2]

const Eigen::MatrixXd & HArDCore3D::SXCurl::RcurlFace ( const Face &  F) const
inline

Return cell reduction for cell T.

◆ RcurlFace() [2/2]

const Eigen::MatrixXd & HArDCore3D::SXCurl::RcurlFace ( size_t  iF) const
inline

Return the reduction for the face of index iF.

◆ RgradCell() [1/2]

const Eigen::MatrixXd & HArDCore3D::SXGrad::RgradCell ( const Cell &  T) const
inline

Return the reduction for cell T.

◆ RgradCell() [2/2]

const Eigen::MatrixXd & HArDCore3D::SXGrad::RgradCell ( size_t  iT) const
inline

Return the reduction for the cell of index iT.

◆ RgradFace() [1/2]

const Eigen::MatrixXd & HArDCore3D::SXGrad::RgradFace ( const Face &  F) const
inline

Return cell reduction for cell T.

◆ RgradFace() [2/2]

const Eigen::MatrixXd & HArDCore3D::SXGrad::RgradFace ( size_t  iF) const
inline

Return the reduction for the face of index iF.

◆ ScurlCell() [1/2]

const Eigen::MatrixXd & HArDCore3D::SXCurl::ScurlCell ( const Cell &  T) const
inline

Return the serendipity reconstruction for cell T.

◆ ScurlCell() [2/2]

const Eigen::MatrixXd & HArDCore3D::SXCurl::ScurlCell ( size_t  iT) const
inline

Return the serendipity reconstruction for the cell of index iT.

◆ ScurlFace() [1/2]

const Eigen::MatrixXd & HArDCore3D::SXCurl::ScurlFace ( const Face &  F) const
inline

Return the serendipity reconstruction for face F.

◆ ScurlFace() [2/2]

const Eigen::MatrixXd & HArDCore3D::SXCurl::ScurlFace ( size_t  iF) const
inline

Return the serendipity reconstruction for the face of index iF.

◆ serDegreeCell()

const int HArDCore3D::SerendipityProblem::serDegreeCell ( size_t  iT) const
inline

Return the serendipity degree ell_T in a cell.

◆ serDegreeFace()

const int HArDCore3D::SerendipityProblem::serDegreeFace ( size_t  iF) const
inline

Return the serendipity degree ell_F in a face.

◆ serendipityEdges()

const std::vector< size_t > & HArDCore3D::SerendipityProblem::serendipityEdges ( size_t  iF) const
inline

Return the list of serendipity edges in a face.

◆ serendipityFaces()

const std::vector< size_t > & HArDCore3D::SerendipityProblem::serendipityFaces ( size_t  iT) const
inline

Return the list of serendipity face in a cell.

◆ SerendipityOperatorCell() [1/2]

const Eigen::MatrixXd HArDCore3D::SerendipityProblem::SerendipityOperatorCell ( const Cell &  T,
const Eigen::MatrixXd &  LT 
) const
inline

Compute the serendipity operator on the Cell T.

◆ SerendipityOperatorCell() [2/2]

const Eigen::MatrixXd SerendipityProblem::SerendipityOperatorCell ( const size_t  iT,
const Eigen::MatrixXd &  LT 
) const

Compute the serendipity operator on the cell of index iT.

◆ SerendipityOperatorFace() [1/2]

const Eigen::MatrixXd HArDCore3D::SerendipityProblem::SerendipityOperatorFace ( const Face &  F,
const Eigen::MatrixXd &  LF 
) const
inline

Compute the serendipity operator on the face F.

◆ SerendipityOperatorFace() [2/2]

const Eigen::MatrixXd SerendipityProblem::SerendipityOperatorFace ( const size_t  iF,
const Eigen::MatrixXd &  LF 
) const

Compute the serendipity operator on the face of index iF.

◆ SerendipityProblem()

SerendipityProblem::SerendipityProblem ( const DDRCore ddrcore,
bool  use_threads = true,
std::ostream &  output = std::cout,
bool  use_serendipity = true 
)

Constructor.

Parameters
ddrcoreDDR class from which this serendipity is constructed
use_threadsUsing multithreading for construction
outputOutput for messages
use_serendipityFalse means that we actually use the DDR class (no serendipity is applied)

◆ serPro() [1/2]

const SerendipityProblem & HArDCore3D::SXCurl::serPro ( ) const
inline

◆ serPro() [2/2]

const SerendipityProblem & HArDCore3D::SXGrad::serPro ( ) const
inline

◆ SgradCell() [1/2]

const Eigen::MatrixXd & HArDCore3D::SXGrad::SgradCell ( const Cell &  T) const
inline

Return the serendipity reconstruction for cell T.

◆ SgradCell() [2/2]

const Eigen::MatrixXd & HArDCore3D::SXGrad::SgradCell ( size_t  iT) const
inline

Return the serendipity reconstruction for the cell of index iT.

◆ SgradFace() [1/2]

const Eigen::MatrixXd & HArDCore3D::SXGrad::SgradFace ( const Face &  F) const
inline

Return the serendipity reconstruction for face F.

◆ SgradFace() [2/2]

const Eigen::MatrixXd & HArDCore3D::SXGrad::SgradFace ( size_t  iF) const
inline

Return the serendipity reconstruction for the face of index iF.

◆ SXCurl()

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

Constructor.

◆ SXDiv()

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

Constructor.

◆ SXGrad()

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

Constructor.

◆ TransferOperators() [1/2]

HArDCore3D::SXCurl::TransferOperators::TransferOperators ( const Eigen::MatrixXd &  _serendipity,
const Eigen::MatrixXd &  _extension,
const Eigen::MatrixXd &  _reduction 
)
inline
Parameters
_serendipitySerendipity reconstruction
_extensionExtension SXCurl->XCurl
_reductionReduction XCurl->SXCurl

◆ TransferOperators() [2/2]

HArDCore3D::SXGrad::TransferOperators::TransferOperators ( const Eigen::MatrixXd &  _serendipity,
const Eigen::MatrixXd &  _extension,
const Eigen::MatrixXd &  _reduction 
)
inline
Parameters
_serendipitySerendipity reconstruction
_extensionExtension SXCurl->XCurl
_reductionReduction XCurl->SXCurl

◆ XCurl()

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

Constructor.

◆ XDiv()

XDiv::XDiv ( 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.

Variable Documentation

◆ curl

Eigen::MatrixXd HArDCore3D::XCurl::LocalOperators::curl

◆ divergence

Eigen::MatrixXd HArDCore3D::XDiv::LocalOperators::divergence

◆ divergence_rhs

Eigen::MatrixXd HArDCore3D::XDiv::LocalOperators::divergence_rhs

◆ extension [1/2]

Eigen::MatrixXd HArDCore3D::SXCurl::TransferOperators::extension

◆ extension [2/2]

Eigen::MatrixXd HArDCore3D::SXGrad::TransferOperators::extension

◆ GolyComplk

std::unique_ptr<GolyComplBasisCellType> HArDCore3D::DDRCore::CellBases::GolyComplk

◆ GolyComplkpo

std::unique_ptr<GolyComplBasisCellType> HArDCore3D::DDRCore::CellBases::GolyComplkpo

◆ Golykmo

std::unique_ptr<GolyBasisCellType> HArDCore3D::DDRCore::CellBases::Golykmo

◆ gradient

Eigen::MatrixXd HArDCore3D::XGrad::LocalOperators::gradient

◆ m_cell_operators

std::vector<std::unique_ptr<LocalOperators> > HArDCore3D::XDiv::m_cell_operators
protected

◆ m_ddr_core

const DDRCore& HArDCore3D::XDiv::m_ddr_core
protected

◆ m_output

std::ostream& HArDCore3D::XDiv::m_output
protected

◆ m_use_threads

bool HArDCore3D::XDiv::m_use_threads
protected

◆ Polyk [1/3]

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

◆ Polyk [2/3]

std::unique_ptr<PolyBasisFaceType> HArDCore3D::DDRCore::FaceBases::Polyk

◆ Polyk [3/3]

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

◆ Polyk2

std::unique_ptr<Poly2BasisFaceType> HArDCore3D::DDRCore::FaceBases::Polyk2

◆ Polyk3

std::unique_ptr<Poly3BasisCellType> HArDCore3D::DDRCore::CellBases::Polyk3

◆ Polykmo [1/3]

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

◆ Polykmo [2/3]

std::unique_ptr<PolyBasisFaceType> HArDCore3D::DDRCore::FaceBases::Polykmo

◆ Polykmo [3/3]

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

◆ Polykpo [1/3]

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

◆ Polykpo [2/3]

std::unique_ptr<PolyBasisFaceType> HArDCore3D::DDRCore::FaceBases::Polykpo

◆ Polykpo [3/3]

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

◆ potential [1/3]

Eigen::MatrixXd HArDCore3D::XCurl::LocalOperators::potential

◆ potential [2/3]

Eigen::MatrixXd HArDCore3D::XDiv::LocalOperators::potential

◆ potential [3/3]

Eigen::MatrixXd HArDCore3D::XGrad::LocalOperators::potential

◆ reduction [1/2]

Eigen::MatrixXd HArDCore3D::SXCurl::TransferOperators::reduction

◆ reduction [2/2]

Eigen::MatrixXd HArDCore3D::SXGrad::TransferOperators::reduction

◆ RolyComplk [1/2]

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

◆ RolyComplk [2/2]

std::unique_ptr<RolyComplBasisFaceType> HArDCore3D::DDRCore::FaceBases::RolyComplk

◆ RolyComplkp2 [1/2]

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

◆ RolyComplkp2 [2/2]

std::unique_ptr<RolyComplBasisFaceType> HArDCore3D::DDRCore::FaceBases::RolyComplkp2

◆ Rolykmo [1/2]

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

◆ Rolykmo [2/2]

std::unique_ptr<RolyBasisFaceType> HArDCore3D::DDRCore::FaceBases::Rolykmo

◆ serendipity [1/2]

Eigen::MatrixXd HArDCore3D::SXCurl::TransferOperators::serendipity

◆ serendipity [2/2]

Eigen::MatrixXd HArDCore3D::SXGrad::TransferOperators::serendipity