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

Classes providing cell and edge quadrature rules, and values of basis functions at the nodes. More...

Classes

class  HArDCore2D::ElementQuad
 
class  HArDCore2D::UVector
 
class  HArDCore2D::HybridCore
 

Typedefs

typedef Family< MonomialScalarBasisCellHArDCore2D::HybridCore::PolyCellBasisType
 type for cell basis More...
 
typedef Family< MonomialScalarBasisEdgeHArDCore2D::HybridCore::PolyEdgeBasisType
 type for edge basis More...
 

Functions

 HArDCore2D::ElementQuad::ElementQuad (const HybridCore &hho, const size_t iT, const size_t doeT, const size_t doeF)
 Class constructor: loads the quadrature rules and values of basis functions/gradients at these points. More...
 
const QuadratureRuleHArDCore2D::ElementQuad::get_quadT () const
 Returns quadrature rules in cell. More...
 
const QuadratureRuleHArDCore2D::ElementQuad::get_quadF (size_t ilF) const
 Returns quadrature rules on edge with local number ilF. More...
 
const boost::multi_array< double, 2 > & HArDCore2D::ElementQuad::get_phiT_quadT () const
 Returns values of cell basis functions at cell quadrature nodes. More...
 
const boost::multi_array< VectorRd, 2 > & HArDCore2D::ElementQuad::get_dphiT_quadT () const
 Returns values of gradients of cell basis functions at cell quadrature nodes. More...
 
const boost::multi_array< double, 2 > & HArDCore2D::ElementQuad::get_phiT_quadF (size_t ilF) const
 Returns values of cell basis functions at edge quadrature nodes, for edge with local number ilF. More...
 
const boost::multi_array< double, 2 > HArDCore2D::ElementQuad::get_phiF_quadF (size_t ilF) const
 Returns values of edge basis functions at edge quadrature nodes, for edge with local number ilF. More...
 
const boost::multi_array< VectorRd, 2 > & HArDCore2D::ElementQuad::get_dphiT_quadF (size_t ilF) const
 Returns values of gradients of cell basis functions at edge quadrature nodes, for edge with local number ilF. More...
 
boost::multi_array< VectorRd, 2 > HArDCore2D::ElementQuad::get_vec_phiT_quadT (size_t degree) const
 Builds on the fly the values of vector cell basis functions at cell quadrature nodes. The vector basis is obtained by tensorization of the scalar one: \((\phi_1,0), (\phi_2,0), ..., (\phi_N,0), (0,\phi_1), (0,\phi_2) ... (0,\phi_N)\). More...
 
boost::multi_array< VectorRd, 2 > HArDCore2D::ElementQuad::get_vec_phiT_quadF (size_t ilF, size_t degree) const
 Builds on the fly the values of vector cell basis functions at edge quadrature nodes. The vector basis is obtained by tensorization of the scalar one as in get_vec_phiT_quadT. More...
 
boost::multi_array< VectorRd, 2 > HArDCore2D::ElementQuad::get_vec_phiF_quadF (size_t ilF, size_t degree) const
 Builds on the fly the values of vector edge basis functions at edge quadrature nodes. The vector basis is obtained by tensorization of the scalar one as in get_vec_phiT_quadT. More...
 
template<typename GeometricSupport >
size_t const HArDCore2D::DimPoly (int m)
 
template<>
const size_t HArDCore2D::DimPoly< Cell > (const int m)
 Compute the size of the basis of 2-variate polynomials up to degree m. More...
 
template<>
const size_t HArDCore2D::DimPoly< Edge > (const int m)
 Compute the size of the basis of 1-variate polynomials up to degree m. More...
 
 HArDCore2D::UVector::UVector (const Eigen::VectorXd &values, const Mesh &mesh, const int cell_deg, const size_t edge_deg)
 
Eigen::VectorXd & HArDCore2D::UVector::asVectorXd () const
 Return the values as an Eigen vector. More...
 
const int HArDCore2D::UVector::get_cell_deg () const
 Return the cell degree. More...
 
const size_t HArDCore2D::UVector::get_edge_deg () const
 Return the edge degree. More...
 
const size_t HArDCore2D::UVector::n_cell_dofs () const
 Number of dofs in each cell. More...
 
const size_t HArDCore2D::UVector::n_edge_dofs () const
 Number of dofs on each edge. More...
 
const size_t HArDCore2D::UVector::n_total_cell_dofs () const
 Total number of cell dofs (in the vector, this is where the edge dofs start) More...
 
Eigen::VectorXd HArDCore2D::UVector::restr (size_t iT) const
 Extract the restriction of the unknowns corresponding to cell iT and its edges. More...
 
UVector HArDCore2D::UVector::operator+ (const UVector &b)
 Overloads the addition: adds the coefficients. More...
 
UVector HArDCore2D::UVector::operator- (const UVector &b)
 Overloads the subtraction: subtracts the coefficients. More...
 
void HArDCore2D::UVector::operator+= (const UVector &b)
 Overloads the increment: adds the coefficients. More...
 
void HArDCore2D::UVector::operator/= (const double &r)
 Overloads the increment: adds the coefficients. More...
 
double HArDCore2D::UVector::operator() (size_t index) const
 Overloads the (): returns the corresponding coefficient. More...
 
 HArDCore2D::HybridCore::HybridCore (const Mesh *mesh_ptr, const int cell_deg, const size_t edge_deg, const bool use_threads=true, std::ostream &output=std::cout, const bool ortho=true)
 Class constructor: initialises the data structure with the given mesh, and desired polynomial degrees of the basis functions. More...
 
const MeshHArDCore2D::HybridCore::get_mesh () const
 Returns a pointer to the mesh. More...
 
const int HArDCore2D::HybridCore::CellDegree () const
 Return the degree of cell polynomials. More...
 
const int HArDCore2D::HybridCore::CellDegreePos () const
 
const size_t HArDCore2D::HybridCore::EdgeDegree () const
 Return the degree of edge polynomials. More...
 
const PolyCellBasisTypeHArDCore2D::HybridCore::CellBasis (size_t iT) const
 Return cell basis for element with global index iT. More...
 
const PolyEdgeBasisTypeHArDCore2D::HybridCore::EdgeBasis (size_t iE) const
 Return edge basis for edge with global index iE. More...
 
double HArDCore2D::HybridCore::L2norm (const UVector &Xh) const
 Compute L2 norm of a discrete function (using cell values) More...
 
double HArDCore2D::HybridCore::H1norm (const UVector &Xh) const
 Compute discrete H1 norm of a discrete function. More...
 
template<typename ContinuousFunction >
UVector HArDCore2D::HybridCore::interpolate (const ContinuousFunction &f, const int deg_cell, const size_t deg_edge, size_t doe) const
 Compute the interpolant in the discrete space of a continuous function. More...
 
Eigen::VectorXd HArDCore2D::HybridCore::compute_weights (size_t iT) const
 Computes the weights to get cell values from edge values when l=-1. More...
 
double HArDCore2D::HybridCore::evaluate_in_cell (const UVector Xh, size_t iT, VectorRd x) const
 Evaluates a discrete function in the cell iT at point x. More...
 
double HArDCore2D::HybridCore::evaluate_in_edge (const UVector Xh, size_t iE, VectorRd x) const
 Evaluates a discrete function on the edge iE at point x. More...
 
Eigen::VectorXd HArDCore2D::HybridCore::VertexValues (const UVector Xh, const std::string from_dofs)
 From a hybrid function, computes a vector of values at the vertices of the mesh. More...
 

Detailed Description

Classes providing cell and edge quadrature rules, and values of basis functions at the nodes.

Classes providing tools to implement schemes having polynomial unknowns in the cells and on the edges.

Typedef Documentation

◆ PolyCellBasisType

type for cell basis

◆ PolyEdgeBasisType

type for edge basis

Function Documentation

◆ asVectorXd()

Eigen::VectorXd& HArDCore2D::UVector::asVectorXd ( ) const
inline

Return the values as an Eigen vector.

◆ CellBasis()

const PolyCellBasisType& HArDCore2D::HybridCore::CellBasis ( size_t  iT) const
inline

Return cell basis for element with global index iT.

◆ CellDegree()

const int HArDCore2D::HybridCore::CellDegree ( ) const
inline

Return the degree of cell polynomials.

◆ CellDegreePos()

const int HArDCore2D::HybridCore::CellDegreePos ( ) const
inline

◆ compute_weights()

Eigen::VectorXd HybridCore::compute_weights ( size_t  iT) const

Computes the weights to get cell values from edge values when l=-1.

◆ DimPoly()

template<typename GeometricSupport >
size_t const HArDCore2D::DimPoly ( int  m)
inline

◆ DimPoly< Cell >()

template<>
const size_t HArDCore2D::DimPoly< Cell > ( const int  m)
inline

Compute the size of the basis of 2-variate polynomials up to degree m.

Parameters
mPolynomial degree

◆ DimPoly< Edge >()

template<>
const size_t HArDCore2D::DimPoly< Edge > ( const int  m)
inline

Compute the size of the basis of 1-variate polynomials up to degree m.

Parameters
mPolynomial degree

◆ EdgeBasis()

const PolyEdgeBasisType& HArDCore2D::HybridCore::EdgeBasis ( size_t  iE) const
inline

Return edge basis for edge with global index iE.

◆ EdgeDegree()

const size_t HArDCore2D::HybridCore::EdgeDegree ( ) const
inline

Return the degree of edge polynomials.

◆ ElementQuad()

ElementQuad::ElementQuad ( const HybridCore hho,
const size_t  iT,
const size_t  doeT,
const size_t  doeF 
)

Class constructor: loads the quadrature rules and values of basis functions/gradients at these points.

Parameters
hhoA reference to the hybridcore instance
iTGlobal index of cell
doeTThe degree of exactness for cell quadratures
doeFThe degree of exactness of edge quadratures

◆ evaluate_in_cell()

double HybridCore::evaluate_in_cell ( const UVector  Xh,
size_t  iT,
VectorRd  x 
) const

Evaluates a discrete function in the cell iT at point x.

◆ evaluate_in_edge()

double HybridCore::evaluate_in_edge ( const UVector  Xh,
size_t  iE,
VectorRd  x 
) const

Evaluates a discrete function on the edge iE at point x.

◆ get_cell_deg()

const int HArDCore2D::UVector::get_cell_deg ( ) const
inline

Return the cell degree.

◆ get_dphiT_quadF()

const boost::multi_array<VectorRd, 2>& HArDCore2D::ElementQuad::get_dphiT_quadF ( size_t  ilF) const
inline

Returns values of gradients of cell basis functions at edge quadrature nodes, for edge with local number ilF.

◆ get_dphiT_quadT()

const boost::multi_array<VectorRd, 2>& HArDCore2D::ElementQuad::get_dphiT_quadT ( ) const
inline

Returns values of gradients of cell basis functions at cell quadrature nodes.

◆ get_edge_deg()

const size_t HArDCore2D::UVector::get_edge_deg ( ) const
inline

Return the edge degree.

◆ get_mesh()

const Mesh* HArDCore2D::HybridCore::get_mesh ( ) const
inline

Returns a pointer to the mesh.

◆ get_phiF_quadF()

const boost::multi_array<double, 2> HArDCore2D::ElementQuad::get_phiF_quadF ( size_t  ilF) const
inline

Returns values of edge basis functions at edge quadrature nodes, for edge with local number ilF.

◆ get_phiT_quadF()

const boost::multi_array<double, 2>& HArDCore2D::ElementQuad::get_phiT_quadF ( size_t  ilF) const
inline

Returns values of cell basis functions at edge quadrature nodes, for edge with local number ilF.

◆ get_phiT_quadT()

const boost::multi_array<double, 2>& HArDCore2D::ElementQuad::get_phiT_quadT ( ) const
inline

Returns values of cell basis functions at cell quadrature nodes.

◆ get_quadF()

const QuadratureRule& HArDCore2D::ElementQuad::get_quadF ( size_t  ilF) const
inline

Returns quadrature rules on edge with local number ilF.

◆ get_quadT()

const QuadratureRule& HArDCore2D::ElementQuad::get_quadT ( ) const
inline

Returns quadrature rules in cell.

◆ get_vec_phiF_quadF()

boost::multi_array< VectorRd, 2 > ElementQuad::get_vec_phiF_quadF ( size_t  ilF,
size_t  degree 
) const

Builds on the fly the values of vector edge basis functions at edge quadrature nodes. The vector basis is obtained by tensorization of the scalar one as in get_vec_phiT_quadT.

Returns
vec_phiT_quadT such that, for r=0,..,dim-1 and i=0,..,DimPoly<Cell>(degree)-1, vec_phiT_quadT[r*DimPoly<Cell>(degree)-1 + i] has, on its row r, the values of the i-th scalar basis function at the quadrature nodes, and 0 on its other rows.
Parameters
degreelocal number of edge required degree of basis function

◆ get_vec_phiT_quadF()

boost::multi_array< VectorRd, 2 > ElementQuad::get_vec_phiT_quadF ( size_t  ilF,
size_t  degree 
) const

Builds on the fly the values of vector cell basis functions at edge quadrature nodes. The vector basis is obtained by tensorization of the scalar one as in get_vec_phiT_quadT.

Returns
vec_phiT_quadT such that, for r=0,..,dim-1 and i=0,..,DimPoly<Cell>(degree)-1, vec_phiT_quadT[r*DimPoly<Cell>(degree)-1 + i] has, on its row r, the values of the i-th scalar basis function at the quadrature nodes, and 0 on its other rows.
Parameters
degreelocal number of edge maximum degree of basis functions required

◆ get_vec_phiT_quadT()

boost::multi_array< VectorRd, 2 > ElementQuad::get_vec_phiT_quadT ( size_t  degree) const

Builds on the fly the values of vector cell basis functions at cell quadrature nodes. The vector basis is obtained by tensorization of the scalar one: \((\phi_1,0), (\phi_2,0), ..., (\phi_N,0), (0,\phi_1), (0,\phi_2) ... (0,\phi_N)\).

Returns
vec_phiT_quadT such that, for r=0,..,dim-1 and i=0,..,DimPoly<Cell>(degree)-1, vec_phiT_quadT[r*DimPoly<Cell>(degree)-1 + i] has, on its row r, the values of the i-th scalar basis function at the quadrature nodes, and 0 on its other rows.
Parameters
degreemaximal degree of basis functions that is required

◆ H1norm()

double HybridCore::H1norm ( const UVector Xh) const

Compute discrete H1 norm of a discrete function.

Parameters
XhVector of unknowns

◆ HybridCore()

HybridCore::HybridCore ( const Mesh mesh_ptr,
const int  cell_deg,
const size_t  edge_deg,
const bool  use_threads = true,
std::ostream &  output = std::cout,
const bool  ortho = true 
)

Class constructor: initialises the data structure with the given mesh, and desired polynomial degrees of the basis functions.

The orthonormalisation comes at a cost in terms of manipulation of the basis functions. This should only be used when the polynomial degree is large and/or the cell is distorted. However, in these cases, it can make a huge difference on the observed convergence rate.

Parameters
mesh_ptrA pointer to the loaded mesh
cell_degThe degree of the cell polynomials
edge_degThe degree of the edge polynomials
use_threadsOptional argument to indicate if threads should be used
outputOptional argument for specifying outputs of messages
orthoOptional argument for choosing to orthonormalise or not the basis functions

◆ interpolate()

template<typename ContinuousFunction >
UVector HArDCore2D::HybridCore::interpolate ( const ContinuousFunction &  f,
const int  deg_cell,
const size_t  deg_edge,
size_t  doe 
) const

Compute the interpolant in the discrete space of a continuous function.

Returns
vector Xh of coefficients on the basis functions, as described in class UVector.
Parameters
ffunction to interpolate
deg_celldegree of cell polynomials for interpolation
deg_edgedegree of edge polynomials for interpolation
doedegree of exactness of the quadrature rules to compute the interpolate

◆ L2norm()

double HybridCore::L2norm ( const UVector Xh) const

Compute L2 norm of a discrete function (using cell values)

Parameters
XhVector of unknowns

◆ n_cell_dofs()

const size_t HArDCore2D::UVector::n_cell_dofs ( ) const
inline

Number of dofs in each cell.

◆ n_edge_dofs()

const size_t HArDCore2D::UVector::n_edge_dofs ( ) const
inline

Number of dofs on each edge.

◆ n_total_cell_dofs()

const size_t HArDCore2D::UVector::n_total_cell_dofs ( ) const
inline

Total number of cell dofs (in the vector, this is where the edge dofs start)

◆ operator()()

double HArDCore2D::UVector::operator() ( size_t  index) const
inline

Overloads the (): returns the corresponding coefficient.

◆ operator+()

UVector HArDCore2D::UVector::operator+ ( const UVector b)
inline

Overloads the addition: adds the coefficients.

◆ operator+=()

void HArDCore2D::UVector::operator+= ( const UVector b)
inline

Overloads the increment: adds the coefficients.

◆ operator-()

UVector HArDCore2D::UVector::operator- ( const UVector b)
inline

Overloads the subtraction: subtracts the coefficients.

◆ operator/=()

void HArDCore2D::UVector::operator/= ( const double &  r)
inline

Overloads the increment: adds the coefficients.

◆ restr()

Eigen::VectorXd UVector::restr ( size_t  iT) const

Extract the restriction of the unknowns corresponding to cell iT and its edges.

◆ UVector()

UVector::UVector ( const Eigen::VectorXd &  values,
const Mesh mesh,
const int  cell_deg,
const size_t  edge_deg 
)
Parameters
valuesvalues of the vector
meshreference to the mesh
cell_degpolynomial degrees in cell
edge_degpolynomial degrees on edge

◆ VertexValues()

Eigen::VectorXd HybridCore::VertexValues ( const UVector  Xh,
const std::string  from_dofs 
)

From a hybrid function, computes a vector of values at the vertices of the mesh.

Parameters
Xhvector of discrete unknowns on cell and edge polynomials
from_dofsType of unknowns to use: "cell" or "edge"