HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Classes | Typedefs | Functions | Variables
Quadratures

Classes providing quadratures on edges and in cells. More...

Classes

struct  HArDCore2D::VecHash
 Hash function for VectorZd type. More...
 
class  HArDCore2D::QuadRuleEdge
 
class  HArDCore2D::LegendreGauss
 
struct  HArDCore2D::QuadratureNode
 Description of one node and one weight from a quadrature rule. More...
 

Typedefs

typedef std::unordered_map< VectorZd, double, VecHashHArDCore2D::MonomialCellIntegralsType
 Type for list of integrals of monomials. More...
 
typedef std::vector< double > HArDCore2D::MonomialEdgeIntegralsType
 Type for list of edge integrals of monomials. More...
 
typedef std::vector< QuadratureNodeHArDCore2D::QuadratureRule
 

Functions

size_t HArDCore2D::VecHash::operator() (const VectorZd &p) const
 
MonomialCellIntegralsType HArDCore2D::IntegrateCellMonomials (const Cell &T, const size_t maxdeg)
 Compute all integrals on a cell of monomials up to a total degree, using vertex values. More...
 
MonomialCellIntegralsType HArDCore2D::CheckIntegralsDegree (const Cell &T, const size_t degree, const MonomialCellIntegralsType &mono_int_map={})
 Checks if the degree of an existing list of monomial integrals is sufficient, other re-compute and return a proper list. More...
 
template<typename BasisType >
Eigen::MatrixXd HArDCore2D::transformGM (const Family< BasisType > &family_basis, const char RC, const Eigen::MatrixXd &anc_GM)
 Transforms a Gram Matrix from an ancestor to a family basis. More...
 
template<typename BasisType >
Eigen::MatrixXd HArDCore2D::transformGM (const RestrictedBasis< BasisType > &restr_basis, const char RC, const Eigen::MatrixXd &anc_GM)
 Transforms a Gram Matrix from an ancestor to a restricted basis. More...
 
template<typename BasisType >
Eigen::MatrixXd HArDCore2D::transformGM (const ShiftedBasis< BasisType > &shifted_basis, const char RC, const Eigen::MatrixXd &anc_GM)
 Transforms a Gram Matrix from an ancestor to a shifted basis. More...
 
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const MonomialScalarBasisCell &basis1, const MonomialScalarBasisCell &basis2, MonomialCellIntegralsType mono_int_map={})
 Computes the Gram Matrix of a pair of local scalar monomial bases. More...
 
template<typename BasisType >
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const BasisType &basis, MonomialCellIntegralsType mono_int_map={})
 This overload to simplify the call to GramMatrix in case the two bases are the same. More...
 
template<typename BasisType1 , typename BasisType2 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const TensorizedVectorFamily< BasisType1, N > &basis1, const TensorizedVectorFamily< BasisType2, N > &basis2, MonomialCellIntegralsType mono_int_map={})
 Template to compute the Gram Matrix of any pair of tensorized scalar bases. More...
 
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const RolyComplBasisCell &basis1, const RolyComplBasisCell &basis2, MonomialCellIntegralsType mono_int_map={})
 Computes the Gram Matrix of a pair of RolyCompl bases. More...
 
template<typename BasisType1 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const RolyComplBasisCell &rolycompl_basis, const TensorizedVectorFamily< BasisType1, N > &tens_family, MonomialCellIntegralsType mono_int_map={})
 Template to compute the Gram Matrix of a RolyCompl basis and a tensorized scalar basis. More...
 
template<typename BasisType1 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const TensorizedVectorFamily< BasisType1, N > &tens_family, const RolyComplBasisCell &rolycompl_basis, MonomialCellIntegralsType mono_int_map={})
 Template to compute the Gram Matrix of a tensorized scalar basis and a RolyCompl basis. More...
 
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const GolyComplBasisCell &basis1, const GolyComplBasisCell &basis2, MonomialCellIntegralsType mono_int_map={})
 Gram Matrix of a pair of GolyCompl bases. More...
 
Eigen::MatrixXd HArDCore2D::GMRolyComplScalar (const Cell &T, const RolyComplBasisCell &rolycompl_basis, const MonomialScalarBasisCell &mono_basis, const size_t m, MonomialCellIntegralsType mono_int_map)
 Computes the Gram Matrix of the mth component of a RolyCompl Basis and a monomial basis. More...
 
template<typename BasisType >
Eigen::MatrixXd HArDCore2D::GMRolyComplScalar (const Cell &T, const RolyComplBasisCell &basis1, const BasisType &basis2, const size_t m, MonomialCellIntegralsType mono_int_map)
 Generic template to compute the Gram Matrix of the mth component of a RolyCompl Basis and any basis. More...
 
template<typename BasisType >
constexpr bool HArDCore2D::useAncestor ()
 Determines if the ancestor of a basis will be used to compute a Gram matrix for this basis. More...
 
template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const BasisType1 &basis1, const BasisType2 &basis2, MonomialCellIntegralsType mono_int_map={})
 Generic template to compute the Gram Matrix of any pair of bases. More...
 
Eigen::MatrixXd HArDCore2D::GMScalarDerivative (const Cell &T, const MonomialScalarBasisCell &basis1, const MonomialScalarBasisCell &basis2, const size_t m, MonomialCellIntegralsType mono_int_map={})
 Computes the Gram Matrix of a pair of local scalar monomial bases, taking a partial derivative of the first (w.r.t. the homogeneous coordinates, without scaling) More...
 
Eigen::MatrixXd HArDCore2D::GMScalarDerivative (const Cell &T, const MonomialScalarBasisCell &basis1, const MonomialScalarBasisCell &basis2, const size_t m, const size_t l, MonomialCellIntegralsType mono_int_map={})
 Computes the Gram Matrix of a pair of local scalar monomial bases, taking partial derivatives of each of them (w.r.t. the homogeneous coordinates, without scaling) More...
 
template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GMScalarDerivative (const Cell &T, const BasisType1 &basis1, const BasisType2 &basis2, const size_t m, MonomialCellIntegralsType mono_int_map={})
 Generic template to compute the Gram Matrix of any pair of scalar bases, taking a partial derivative of the first (w.r.t. the homogeneous coordinates, without scaling) More...
 
template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GMScalarDerivative (const Cell &T, const BasisType1 &basis1, const BasisType2 &basis2, const size_t m, const size_t l, MonomialCellIntegralsType mono_int_map={})
 Generic template to compute the Gram Matrix of any pair of scalar bases, taking partial derivatives of each of them (w.r.t. the homogeneous coordinates, without scaling) More...
 
template<typename BasisType1 , typename BasisType2 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const GradientBasis< BasisType1 > &grad_basis, const TensorizedVectorFamily< BasisType2, N > &tens_family, MonomialCellIntegralsType mono_int_map={})
 Template to compute the Gram Matrix of a gradient basis and a tensorized scalar basis. More...
 
template<typename BasisType1 , typename BasisType2 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const TensorizedVectorFamily< BasisType1, N > &tens_family, const GradientBasis< BasisType2 > &grad_basis, MonomialCellIntegralsType mono_int_map={})
 Template to compute the Gram Matrix of a tensorized scalar basis and a gradient basis. More...
 
template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const GradientBasis< BasisType1 > &grad_basis1, const GradientBasis< BasisType2 > &grad_basis2, MonomialCellIntegralsType mono_int_map={})
 Template to compute the Gram Matrix of a gradient basis and another gradient basis. More...
 
template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const CurlBasis< BasisType1 > &basis1, const CurlBasis< BasisType2 > &basis2, MonomialCellIntegralsType mono_int_map={})
 Generic template to compute the Gram Matrix of a pair of Curl bases. More...
 
template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const CurlBasis< BasisType1 > &curl_basis, const TensorizedVectorFamily< BasisType2, 2 > &tens_family, MonomialCellIntegralsType mono_int_map={})
 Template to compute the Gram Matrix of a curl basis and a tensorized scalar basis. More...
 
template<typename BasisType1 , typename BasisType2 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const TensorizedVectorFamily< BasisType1, N > &tens_family, const CurlBasis< BasisType2 > &curl_basis, MonomialCellIntegralsType mono_int_map={})
 Template to compute the Gram Matrix of a tensorized scalar basis and a gradient basis. More...
 
template<typename BasisType1 , typename BasisType2 >
boost::disable_if< boost::is_same< BasisType2, MonomialCellIntegralsType >, Eigen::MatrixXd >::type HArDCore2D::GramMatrix (const Cell &T, const DivergenceBasis< BasisType1 > &basis1, const BasisType2 &basis2, MonomialCellIntegralsType mono_int_map={})
 Generic template to compute the Gram Matrix of a Divergence basis and any other basis. More...
 
template<typename BasisType1 , typename Basis2 >
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const BasisType1 &basis1, const DivergenceBasis< Basis2 > &basis2, MonomialCellIntegralsType mono_int_map={})
 Template to compute the Gram Matrix of any basis and a Divergence basis. More...
 
template<typename BasisType1 >
Eigen::MatrixXd HArDCore2D::GramMatrixDiv (const Cell &T, const TensorizedVectorFamily< BasisType1, 2 > &basis1, const MonomialScalarBasisCell &basis2, MonomialCellIntegralsType mono_int_map={})
 Template to compute the Gram Matrix of a Divergence<Tensorized> basis and a monomial scalar basis. More...
 
Eigen::MatrixXd HArDCore2D::GramMatrixDiv (const Cell &T, const RolyComplBasisCell &basis1, const MonomialScalarBasisCell &basis2, MonomialCellIntegralsType mono_int_map={})
 Computes the Gram Matrix of a Divergence<RolyCompl> basis and a monomial scalar basis. More...
 
template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrixDiv (const Cell &T, const BasisType1 &basis1, const BasisType2 &basis2, MonomialCellIntegralsType mono_int_map={})
 Template to compute the Gram Matrix of the divergence of any basis and any other basis. More...
 
template<typename BasisType1 , typename BasisType2 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const MatrixFamily< BasisType1, N > &basis1, const MatrixFamily< BasisType2, N > &basis2, MonomialCellIntegralsType mono_int_map={})
 Gram Matrix of any pair of MatrixFamily bases. More...
 
template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const DivergenceBasis< MatrixFamily< BasisType1, dimspace >> &basis1, const TensorizedVectorFamily< BasisType2, dimspace > &basis2, MonomialCellIntegralsType mono_int_map={})
 
template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const GradientBasis< TensorizedVectorFamily< BasisType1, dimspace >> &basis1, const MatrixFamily< BasisType2, dimspace > &basis2, MonomialCellIntegralsType mono_int_map={})
 Gram Matrix of the gradient basis of a tensorized family and a matrix family (only valid if N=dimspace in Tensorized and MatrixFamily). More...
 
template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const MatrixFamily< BasisType1, dimspace > &basis1, const GradientBasis< TensorizedVectorFamily< BasisType2, dimspace >> &basis2, MonomialCellIntegralsType mono_int_map={})
 Gram Matrix of a Matrix family and the gradient of a tensorized family. More...
 
template<typename BasisType1 , typename BasisType2 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix (const Cell &T, const GradientBasis< TensorizedVectorFamily< BasisType1, N >> &basis1, const GradientBasis< TensorizedVectorFamily< BasisType2, N >> &basis2, MonomialCellIntegralsType mono_int_map={})
 Gram Matrix of two gradient bases of tensorized families. More...
 
MonomialEdgeIntegralsType HArDCore2D::IntegrateEdgeMonomials (const Edge &E, const size_t maxdeg)
 Compute all integrals of edge monomials up to a total degree. More...
 
MonomialEdgeIntegralsType HArDCore2D::CheckIntegralsDegree (const Edge &E, const size_t degree, const MonomialEdgeIntegralsType &mono_int_map={})
 Checks if the degree of an existing list of monomial integrals is sufficient, other re-compute and return a proper list. More...
 
Eigen::MatrixXd HArDCore2D::GramMatrix (const Edge &E, const MonomialScalarBasisEdge &basis1, const MonomialScalarBasisEdge &basis2, MonomialEdgeIntegralsType mono_int_map={})
 Computes the Gram Matrix of a pair of local scalar monomial bases. More...
 
template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix (const Edge &E, const BasisType1 &basis1, const BasisType2 &basis2, MonomialEdgeIntegralsType mono_int_map={})
 Generic template to compute the Gram Matrix of any pair of bases. More...
 
template<typename BasisType1 , typename BasisType2 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix (const Edge &E, const TensorizedVectorFamily< BasisType1, N > &basis1, const TensorizedVectorFamily< BasisType2, N > &basis2, MonomialEdgeIntegralsType mono_int_map={})
 
template<typename BasisType >
Eigen::MatrixXd HArDCore2D::GramMatrix (const Edge &E, const BasisType &basis, MonomialEdgeIntegralsType mono_int_map={})
 This overload to simplify the call to GramMatrix in case the two bases are the same. More...
 
Eigen::MatrixXd HArDCore2D::GMDer (const Edge &E, const MonomialScalarBasisEdge &basis1, const MonomialScalarBasisEdge &basis2, MonomialEdgeIntegralsType mono_int_map={})
 Computes the Gram Matrix of the derivative of a monomial basis with another monomial basis. More...
 
template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GMDer (const Edge &E, const BasisType1 &basis1, const BasisType2 &basis2, MonomialEdgeIntegralsType mono_int_map={})
 Generic template for GMDer with derived bases. More...
 
template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix (const Edge &E, const GradientBasis< BasisType1 > &basis1, const BasisType2 &basis2, MonomialEdgeIntegralsType mono_int_map={})
 Computes the Gram Matrix of a gradient basis (considering the tangential gradient as a scalar) and a scalar basis. More...
 
template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix (const Edge &E, const BasisType1 &basis1, const GradientBasis< BasisType2 > &basis2, MonomialEdgeIntegralsType mono_int_map={})
 Computes the Gram Matrix of a scalar basis and a gradient basis (considering the tangential gradient as a scalar) More...
 
 HArDCore2D::QuadRuleEdge::QuadRuleEdge (size_t doe, bool warn)
 
 HArDCore2D::QuadRuleEdge::~QuadRuleEdge ()
 
size_t HArDCore2D::QuadRuleEdge::nq ()
 
double HArDCore2D::QuadRuleEdge::xq (size_t i)
 
double HArDCore2D::QuadRuleEdge::yq (size_t i)
 
double HArDCore2D::QuadRuleEdge::wq (size_t i)
 
void HArDCore2D::QuadRuleEdge::setup (double xV[], double yV[])
 
 HArDCore2D::LegendreGauss::LegendreGauss (size_t doe)
 
 HArDCore2D::LegendreGauss::~LegendreGauss ()
 
void HArDCore2D::LegendreGauss::sub_rule_01 ()
 
void HArDCore2D::LegendreGauss::sub_rule_02 ()
 
void HArDCore2D::LegendreGauss::sub_rule_03 ()
 
void HArDCore2D::LegendreGauss::sub_rule_04 ()
 
void HArDCore2D::LegendreGauss::sub_rule_05 ()
 
void HArDCore2D::LegendreGauss::sub_rule_06 ()
 
void HArDCore2D::LegendreGauss::sub_rule_07 ()
 
void HArDCore2D::LegendreGauss::sub_rule_08 ()
 
void HArDCore2D::LegendreGauss::sub_rule_09 ()
 
void HArDCore2D::LegendreGauss::sub_rule_10 ()
 
void HArDCore2D::LegendreGauss::sub_rule_11 ()
 
size_t HArDCore2D::LegendreGauss::npts ()
 
double HArDCore2D::LegendreGauss::wq (size_t i)
 
double HArDCore2D::LegendreGauss::tq (size_t i)
 
 HArDCore2D::QuadratureNode::QuadratureNode (double x, double y, double w)
 
Eigen::Vector2d HArDCore2D::QuadratureNode::vector () const
 Returns the quadrature point as an Eigen vector. More...
 
QuadratureRule HArDCore2D::generate_quadrature_rule (const Cell &T, const int doe, const bool force_split=false)
 Generate quadrature rule on mesh element. More...
 
QuadratureRule HArDCore2D::generate_quadrature_rule (const Edge &E, const int doe)
 Generate quadrature rule on mesh face. More...
 

Variables

double HArDCore2D::QuadratureNode::x
 
double HArDCore2D::QuadratureNode::y
 
double HArDCore2D::QuadratureNode::w
 

Detailed Description

Classes providing quadratures on edges and in cells.

Typedef Documentation

◆ MonomialCellIntegralsType

typedef std::unordered_map<VectorZd, double, VecHash> HArDCore2D::MonomialCellIntegralsType

Type for list of integrals of monomials.

◆ MonomialEdgeIntegralsType

typedef std::vector<double> HArDCore2D::MonomialEdgeIntegralsType

Type for list of edge integrals of monomials.

◆ QuadratureRule

Function Documentation

◆ CheckIntegralsDegree() [1/2]

MonomialCellIntegralsType HArDCore2D::CheckIntegralsDegree ( const Cell &  T,
const size_t  degree,
const MonomialCellIntegralsType mono_int_map = {} 
)

Checks if the degree of an existing list of monomial integrals is sufficient, other re-compute and return a proper list.

Parameters
TCell
degreeExpected degree
mono_int_mapExisting list, optional

◆ CheckIntegralsDegree() [2/2]

MonomialEdgeIntegralsType HArDCore2D::CheckIntegralsDegree ( const Edge &  E,
const size_t  degree,
const MonomialEdgeIntegralsType mono_int_map = {} 
)

Checks if the degree of an existing list of monomial integrals is sufficient, other re-compute and return a proper list.

Parameters
EEdge
degreeExpected degree
mono_int_mapExisting list, optional

◆ generate_quadrature_rule() [1/2]

QuadratureRule HArDCore2D::generate_quadrature_rule ( const Cell &  T,
const int  doe,
const bool  force_split = false 
)

Generate quadrature rule on mesh element.

Returns
list of quadrature nodes and weights
Parameters
TReference to the mesh cell
doeDegree of exactness
force_splitTRUE if we want the quadrature nodes to be computed by forcing the splitting of the cell into triangles based on its center of mass and edges (otherwise, for simple cells, quadrature nodes are computed by splitting in fewer triangles)

◆ generate_quadrature_rule() [2/2]

QuadratureRule HArDCore2D::generate_quadrature_rule ( const Edge &  E,
const int  doe 
)

Generate quadrature rule on mesh face.

Returns
list of quadrature nodes and weights
Parameters
EReference to the mesh face
doeDegree of exactness

◆ GMDer() [1/2]

template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GMDer ( const Edge &  E,
const BasisType1 &  basis1,
const BasisType2 &  basis2,
MonomialEdgeIntegralsType  mono_int_map = {} 
)

Generic template for GMDer with derived bases.

Parameters
EEdge to which the basis corresponds
basis1First basis (rows of the Gram matrix), to be differentiated
basis2Second basis (columns of the Gram matrix)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GMDer() [2/2]

Eigen::MatrixXd HArDCore2D::GMDer ( const Edge &  E,
const MonomialScalarBasisEdge basis1,
const MonomialScalarBasisEdge basis2,
MonomialEdgeIntegralsType  mono_int_map = {} 
)

Computes the Gram Matrix of the derivative of a monomial basis with another monomial basis.

Parameters
EEdge to which the basis corresponds
basis1First basis, to be differentiated
basis2Second basis
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GMRolyComplScalar() [1/2]

template<typename BasisType >
Eigen::MatrixXd HArDCore2D::GMRolyComplScalar ( const Cell &  T,
const RolyComplBasisCell basis1,
const BasisType &  basis2,
const size_t  m,
MonomialCellIntegralsType  mono_int_map 
)

Generic template to compute the Gram Matrix of the mth component of a RolyCompl Basis and any basis.

Parameters
TCell to which the basis corresponds
basis1First basis
basis2Second basis (columns of the Gram matrix)
mDifferentiate basis1 with respect to the mth variable
mono_int_maplist of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GMRolyComplScalar() [2/2]

Eigen::MatrixXd HArDCore2D::GMRolyComplScalar ( const Cell &  T,
const RolyComplBasisCell rolycompl_basis,
const MonomialScalarBasisCell mono_basis,
const size_t  m,
MonomialCellIntegralsType  mono_int_map 
)

Computes the Gram Matrix of the mth component of a RolyCompl Basis and a monomial basis.

Parameters
TCell to which the basis corresponds
rolycompl_basisFirst basis
mono_basisSecond basis
mAdd one to the power of the mth variable
mono_int_maplist of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GMScalarDerivative() [1/4]

template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GMScalarDerivative ( const Cell &  T,
const BasisType1 &  basis1,
const BasisType2 &  basis2,
const size_t  m,
const size_t  l,
MonomialCellIntegralsType  mono_int_map = {} 
)

Generic template to compute the Gram Matrix of any pair of scalar bases, taking partial derivatives of each of them (w.r.t. the homogeneous coordinates, without scaling)

Parameters
TCell to which the basis corresponds
basis1First basis (rows of the Gram matrix)
basis2Second basis (columns of the Gram matrix)
mDifferentiate basis1 with respect to the mth variable
lDifferentiate basis2 with respect to the lth variable
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GMScalarDerivative() [2/4]

template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GMScalarDerivative ( const Cell &  T,
const BasisType1 &  basis1,
const BasisType2 &  basis2,
const size_t  m,
MonomialCellIntegralsType  mono_int_map = {} 
)

Generic template to compute the Gram Matrix of any pair of scalar bases, taking a partial derivative of the first (w.r.t. the homogeneous coordinates, without scaling)

Parameters
TCell to which the basis corresponds
basis1First basis (rows of the Gram matrix)
basis2Second basis (columns of the Gram matrix)
mDifferentiate basis1 with respect to the mth variable
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GMScalarDerivative() [3/4]

Eigen::MatrixXd HArDCore2D::GMScalarDerivative ( const Cell &  T,
const MonomialScalarBasisCell basis1,
const MonomialScalarBasisCell basis2,
const size_t  m,
const size_t  l,
MonomialCellIntegralsType  mono_int_map = {} 
)

Computes the Gram Matrix of a pair of local scalar monomial bases, taking partial derivatives of each of them (w.r.t. the homogeneous coordinates, without scaling)

Parameters
TCell to which the basis corresponds
basis1First basis
basis2Second basis
mDifferentiate basis1 with respect to the mth variable
lDifferentiate basis2 with respect to the lth variable
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GMScalarDerivative() [4/4]

Eigen::MatrixXd HArDCore2D::GMScalarDerivative ( const Cell &  T,
const MonomialScalarBasisCell basis1,
const MonomialScalarBasisCell basis2,
const size_t  m,
MonomialCellIntegralsType  mono_int_map = {} 
)

Computes the Gram Matrix of a pair of local scalar monomial bases, taking a partial derivative of the first (w.r.t. the homogeneous coordinates, without scaling)

Parameters
TCell to which the basis corresponds
basis1First basis
basis2Second basis
mDifferentiate basis1 with respect to the mth variable
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [1/27]

template<typename BasisType >
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const BasisType &  basis,
MonomialCellIntegralsType  mono_int_map = {} 
)

This overload to simplify the call to GramMatrix in case the two bases are the same.

◆ GramMatrix() [2/27]

template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const BasisType1 &  basis1,
const BasisType2 &  basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Generic template to compute the Gram Matrix of any pair of bases.

Parameters
TCell to which the basis corresponds
basis1First basis (rows of the Gram matrix)
basis2Second basis (columns of the Gram matrix)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [3/27]

template<typename BasisType1 , typename Basis2 >
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const BasisType1 &  basis1,
const DivergenceBasis< Basis2 > &  basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Template to compute the Gram Matrix of any basis and a Divergence basis.

Parameters
TCell to which the basis corresponds
basis1First basis (tensorized basis)
basis2Second basis (divergence basis)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [4/27]

template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const CurlBasis< BasisType1 > &  basis1,
const CurlBasis< BasisType2 > &  basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Generic template to compute the Gram Matrix of a pair of Curl bases.

Parameters
TCell to which the basis corresponds
basis1First basis
basis2Second basis
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [5/27]

template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const CurlBasis< BasisType1 > &  curl_basis,
const TensorizedVectorFamily< BasisType2, 2 > &  tens_family,
MonomialCellIntegralsType  mono_int_map = {} 
)

Template to compute the Gram Matrix of a curl basis and a tensorized scalar basis.

Parameters
TCell to which the basis corresponds
curl_basisFirst basis (rows of the Gram matrix) - curl basis
tens_familySecond basis (columns of the Gram matrix) - tensorized basis
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of the bases

◆ GramMatrix() [6/27]

template<typename BasisType1 , typename BasisType2 >
boost::disable_if<boost::is_same<BasisType2, MonomialCellIntegralsType>, Eigen::MatrixXd>::type HArDCore2D::GramMatrix ( const Cell &  T,
const DivergenceBasis< BasisType1 > &  basis1,
const BasisType2 &  basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Generic template to compute the Gram Matrix of a Divergence basis and any other basis.

Parameters
TCell to which the basis corresponds
basis1First basis (rows of the Gram matrix)
basis2Second basis (columns of the Gram matrix)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [7/27]

template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const DivergenceBasis< MatrixFamily< BasisType1, dimspace >> &  basis1,
const TensorizedVectorFamily< BasisType2, dimspace > &  basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)
Parameters
TCell to which the basis corresponds
basis1Divergence of MatrixFamily (rows of the Gram matrix)
basis2Tensorized family (columns of the Gram matrix)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [8/27]

Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const GolyComplBasisCell basis1,
const GolyComplBasisCell basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Gram Matrix of a pair of GolyCompl bases.

Parameters
TCell to which the basis corresponds
basis1First basis
basis2Second basis
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [9/27]

template<typename BasisType1 , typename BasisType2 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const GradientBasis< BasisType1 > &  grad_basis,
const TensorizedVectorFamily< BasisType2, N > &  tens_family,
MonomialCellIntegralsType  mono_int_map = {} 
)

Template to compute the Gram Matrix of a gradient basis and a tensorized scalar basis.

Parameters
TCell to which the basis corresponds
grad_basisFirst basis (rows of the Gram matrix) - gradient basis
tens_familySecond basis (columns of the Gram matrix) - tensorized basis
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of the bases

◆ GramMatrix() [10/27]

template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const GradientBasis< BasisType1 > &  grad_basis1,
const GradientBasis< BasisType2 > &  grad_basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Template to compute the Gram Matrix of a gradient basis and another gradient basis.

Parameters
TCell to which the basis corresponds
grad_basis1First basis (rows of the Gram matrix) - gradient basis
grad_basis2Second basis (columns of the Gram matrix) - gradient basis
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of the bases

◆ GramMatrix() [11/27]

template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const GradientBasis< TensorizedVectorFamily< BasisType1, dimspace >> &  basis1,
const MatrixFamily< BasisType2, dimspace > &  basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Gram Matrix of the gradient basis of a tensorized family and a matrix family (only valid if N=dimspace in Tensorized and MatrixFamily).

Parameters
TCell to which the basis corresponds
basis1Gradient of tensorized (rows of the Gram matrix)
basis2Matrix family (columns of the Gram matrix)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [12/27]

template<typename BasisType1 , typename BasisType2 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const GradientBasis< TensorizedVectorFamily< BasisType1, N >> &  basis1,
const GradientBasis< TensorizedVectorFamily< BasisType2, N >> &  basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Gram Matrix of two gradient bases of tensorized families.

Parameters
TCell to which the basis corresponds
basis1First family
basis2Second family (columns of the Gram matrix)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [13/27]

template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const MatrixFamily< BasisType1, dimspace > &  basis1,
const GradientBasis< TensorizedVectorFamily< BasisType2, dimspace >> &  basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Gram Matrix of a Matrix family and the gradient of a tensorized family.

Parameters
TCell to which the basis corresponds
basis1Matrix family (rows of the Gram matrix)
basis2Gradient of tensorized (columns of the Gram matrix)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [14/27]

template<typename BasisType1 , typename BasisType2 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const MatrixFamily< BasisType1, N > &  basis1,
const MatrixFamily< BasisType2, N > &  basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Gram Matrix of any pair of MatrixFamily bases.

Parameters
TCell to which the basis corresponds
basis1First basis (rows of the Gram matrix)
basis2Second basis (columns of the Gram matrix)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [15/27]

Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const MonomialScalarBasisCell basis1,
const MonomialScalarBasisCell basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Computes the Gram Matrix of a pair of local scalar monomial bases.

Parameters
TCell to which the basis corresponds
basis1First basis
basis2Second basis
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [16/27]

Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const RolyComplBasisCell basis1,
const RolyComplBasisCell basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Computes the Gram Matrix of a pair of RolyCompl bases.

Parameters
TCell to which the basis corresponds
basis1First basis
basis2Second basis
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [17/27]

template<typename BasisType1 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const RolyComplBasisCell rolycompl_basis,
const TensorizedVectorFamily< BasisType1, N > &  tens_family,
MonomialCellIntegralsType  mono_int_map = {} 
)

Template to compute the Gram Matrix of a RolyCompl basis and a tensorized scalar basis.

Parameters
TCell to which the basis corresponds
rolycompl_basisFirst basis (RolyCompl basis)
tens_familySecond basis (tensorized basis)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [18/27]

template<typename BasisType1 , typename BasisType2 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const TensorizedVectorFamily< BasisType1, N > &  basis1,
const TensorizedVectorFamily< BasisType2, N > &  basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Template to compute the Gram Matrix of any pair of tensorized scalar bases.

Parameters
TCell to which the basis corresponds
basis1First basis (rows of the Gram matrix)
basis2Second basis (columns of the Gram matrix)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [19/27]

template<typename BasisType1 , typename BasisType2 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const TensorizedVectorFamily< BasisType1, N > &  tens_family,
const CurlBasis< BasisType2 > &  curl_basis,
MonomialCellIntegralsType  mono_int_map = {} 
)

Template to compute the Gram Matrix of a tensorized scalar basis and a gradient basis.

Parameters
TCell to which the basis corresponds
tens_familyFirst basis (rows of the Gram matrix) - gradient basis
curl_basisSecond basis (columns of the Gram matrix) - tensorized basis
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of the bases

◆ GramMatrix() [20/27]

template<typename BasisType1 , typename BasisType2 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const TensorizedVectorFamily< BasisType1, N > &  tens_family,
const GradientBasis< BasisType2 > &  grad_basis,
MonomialCellIntegralsType  mono_int_map = {} 
)

Template to compute the Gram Matrix of a tensorized scalar basis and a gradient basis.

Parameters
TCell to which the basis corresponds
tens_familyFirst basis (rows of the Gram matrix) - gradient basis
grad_basisSecond basis (columns of the Gram matrix) - tensorized basis
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of the bases

◆ GramMatrix() [21/27]

template<typename BasisType1 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Cell &  T,
const TensorizedVectorFamily< BasisType1, N > &  tens_family,
const RolyComplBasisCell rolycompl_basis,
MonomialCellIntegralsType  mono_int_map = {} 
)

Template to compute the Gram Matrix of a tensorized scalar basis and a RolyCompl basis.

Parameters
TCell to which the basis corresponds
tens_familyFirst basis (tensorized basis)
rolycompl_basisSecond basis (RolyCompl basis)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [22/27]

template<typename BasisType >
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Edge &  E,
const BasisType &  basis,
MonomialEdgeIntegralsType  mono_int_map = {} 
)

This overload to simplify the call to GramMatrix in case the two bases are the same.

◆ GramMatrix() [23/27]

template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Edge &  E,
const BasisType1 &  basis1,
const BasisType2 &  basis2,
MonomialEdgeIntegralsType  mono_int_map = {} 
)

Generic template to compute the Gram Matrix of any pair of bases.

Parameters
EEdge to which the basis corresponds
basis1First basis (rows of the Gram matrix)
basis2Second basis (columns of the Gram matrix)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [24/27]

template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Edge &  E,
const BasisType1 &  basis1,
const GradientBasis< BasisType2 > &  basis2,
MonomialEdgeIntegralsType  mono_int_map = {} 
)

Computes the Gram Matrix of a scalar basis and a gradient basis (considering the tangential gradient as a scalar)

Parameters
EEdge to which the basis corresponds
basis1Scalar basis
basis2Gradient basis
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [25/27]

template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Edge &  E,
const GradientBasis< BasisType1 > &  basis1,
const BasisType2 &  basis2,
MonomialEdgeIntegralsType  mono_int_map = {} 
)

Computes the Gram Matrix of a gradient basis (considering the tangential gradient as a scalar) and a scalar basis.

Parameters
EEdge to which the basis corresponds
basis1Gradient basis
basis2Scalar basis
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [26/27]

Eigen::MatrixXd HArDCore2D::GramMatrix ( const Edge &  E,
const MonomialScalarBasisEdge basis1,
const MonomialScalarBasisEdge basis2,
MonomialEdgeIntegralsType  mono_int_map = {} 
)

Computes the Gram Matrix of a pair of local scalar monomial bases.

Parameters
EEdge to which the basis corresponds
basis1First basis
basis2Second basis
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrix() [27/27]

template<typename BasisType1 , typename BasisType2 , size_t N>
Eigen::MatrixXd HArDCore2D::GramMatrix ( const Edge &  E,
const TensorizedVectorFamily< BasisType1, N > &  basis1,
const TensorizedVectorFamily< BasisType2, N > &  basis2,
MonomialEdgeIntegralsType  mono_int_map = {} 
)
Parameters
EEdge to which the basis corresponds
basis1First basis (rows of the Gram matrix)
basis2Second basis (columns of the Gram matrix)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrixDiv() [1/3]

template<typename BasisType1 , typename BasisType2 >
Eigen::MatrixXd HArDCore2D::GramMatrixDiv ( const Cell &  T,
const BasisType1 &  basis1,
const BasisType2 &  basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Template to compute the Gram Matrix of the divergence of any basis and any other basis.

Parameters
TCell to which the basis corresponds
basis1First basis (vector basis)
basis2Second basis (scalar basis)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of the bases

◆ GramMatrixDiv() [2/3]

Eigen::MatrixXd HArDCore2D::GramMatrixDiv ( const Cell &  T,
const RolyComplBasisCell basis1,
const MonomialScalarBasisCell basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Computes the Gram Matrix of a Divergence<RolyCompl> basis and a monomial scalar basis.

Parameters
TCell to which the basis corresponds
basis1First basis (RolyCompl basis)
basis2Second basis (tensorized basis)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ GramMatrixDiv() [3/3]

template<typename BasisType1 >
Eigen::MatrixXd HArDCore2D::GramMatrixDiv ( const Cell &  T,
const TensorizedVectorFamily< BasisType1, 2 > &  basis1,
const MonomialScalarBasisCell basis2,
MonomialCellIntegralsType  mono_int_map = {} 
)

Template to compute the Gram Matrix of a Divergence<Tensorized> basis and a monomial scalar basis.

Parameters
TCell to which the basis corresponds
basis1First basis (divergence basis)
basis2Second basis (monomial scalar basis)
mono_int_mapOptional list of integrals of monomials up to the sum of max degree of basis1 and basis2

◆ IntegrateCellMonomials()

MonomialCellIntegralsType HArDCore2D::IntegrateCellMonomials ( const Cell &  T,
const size_t  maxdeg 
)

Compute all integrals on a cell of monomials up to a total degree, using vertex values.

Parameters
TCell
maxdegMaximal total degree

◆ IntegrateEdgeMonomials()

MonomialEdgeIntegralsType HArDCore2D::IntegrateEdgeMonomials ( const Edge &  E,
const size_t  maxdeg 
)

Compute all integrals of edge monomials up to a total degree.

Parameters
EEdge
maxdegMaximal total degree

◆ LegendreGauss()

LegendreGauss::LegendreGauss ( size_t  doe)

◆ npts()

size_t LegendreGauss::npts ( )

◆ nq()

size_t QuadRuleEdge::nq ( )

◆ operator()()

size_t HArDCore2D::VecHash::operator() ( const VectorZd p) const
inline

◆ QuadratureNode()

HArDCore2D::QuadratureNode::QuadratureNode ( double  x,
double  y,
double  w 
)
inline

◆ QuadRuleEdge()

QuadRuleEdge::QuadRuleEdge ( size_t  doe,
bool  warn 
)

◆ setup()

void QuadRuleEdge::setup ( double  xV[],
double  yV[] 
)

◆ sub_rule_01()

void LegendreGauss::sub_rule_01 ( )

◆ sub_rule_02()

void LegendreGauss::sub_rule_02 ( )

◆ sub_rule_03()

void LegendreGauss::sub_rule_03 ( )

◆ sub_rule_04()

void LegendreGauss::sub_rule_04 ( )

◆ sub_rule_05()

void LegendreGauss::sub_rule_05 ( )

◆ sub_rule_06()

void LegendreGauss::sub_rule_06 ( )

◆ sub_rule_07()

void LegendreGauss::sub_rule_07 ( )

◆ sub_rule_08()

void LegendreGauss::sub_rule_08 ( )

◆ sub_rule_09()

void LegendreGauss::sub_rule_09 ( )

◆ sub_rule_10()

void LegendreGauss::sub_rule_10 ( )

◆ sub_rule_11()

void LegendreGauss::sub_rule_11 ( )

◆ tq()

double LegendreGauss::tq ( size_t  i)

◆ transformGM() [1/3]

template<typename BasisType >
Eigen::MatrixXd HArDCore2D::transformGM ( const Family< BasisType > &  family_basis,
const char  RC,
const Eigen::MatrixXd &  anc_GM 
)
inline

Transforms a Gram Matrix from an ancestor to a family basis.

Parameters
family_basisFamily
RCR if transformation applied on rows (left), C if applied on columns (right)
anc_GMGram matrix of the ancestor basis

◆ transformGM() [2/3]

template<typename BasisType >
Eigen::MatrixXd HArDCore2D::transformGM ( const RestrictedBasis< BasisType > &  restr_basis,
const char  RC,
const Eigen::MatrixXd &  anc_GM 
)
inline

Transforms a Gram Matrix from an ancestor to a restricted basis.

Parameters
restr_basisRestricted basis
RCR if transformation applied on rows (left), C if applied on columns (right)
anc_GMGram matrix of the ancestor basis

◆ transformGM() [3/3]

template<typename BasisType >
Eigen::MatrixXd HArDCore2D::transformGM ( const ShiftedBasis< BasisType > &  shifted_basis,
const char  RC,
const Eigen::MatrixXd &  anc_GM 
)
inline

Transforms a Gram Matrix from an ancestor to a shifted basis.

Parameters
shifted_basisShifted basis
RCR if transformation applied on rows (left), C if applied on columns (right)
anc_GMGram matrix of the ancestor basis

◆ useAncestor()

template<typename BasisType >
constexpr bool HArDCore2D::useAncestor ( )
constexpr

Determines if the ancestor of a basis will be used to compute a Gram matrix for this basis.

◆ vector()

Eigen::Vector2d HArDCore2D::QuadratureNode::vector ( ) const
inline

Returns the quadrature point as an Eigen vector.

◆ wq() [1/2]

double QuadRuleEdge::wq ( size_t  i)

◆ wq() [2/2]

double LegendreGauss::wq ( size_t  i)

◆ xq()

double QuadRuleEdge::xq ( size_t  i)

◆ yq()

double QuadRuleEdge::yq ( size_t  i)

◆ ~LegendreGauss()

LegendreGauss::~LegendreGauss ( )

◆ ~QuadRuleEdge()

QuadRuleEdge::~QuadRuleEdge ( )

Variable Documentation

◆ w

double HArDCore2D::QuadratureNode::w

◆ x

double HArDCore2D::QuadratureNode::x

◆ y

double HArDCore2D::QuadratureNode::y