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
HybridCore

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

Classes

class  HArDCore3D::ElementQuad
 
class  HArDCore3D::UVector
 
class  HArDCore3D::HybridCore
 

Typedefs

typedef Family< MonomialScalarBasisCellHArDCore3D::HybridCore::PolyCellBasisType
 type for cell basis
 
typedef Family< MonomialScalarBasisFaceHArDCore3D::HybridCore::PolyFaceBasisType
 type for face basis
 
typedef Family< MonomialScalarBasisEdgeHArDCore3D::HybridCore::PolyEdgeBasisType
 type for edge basis
 

Functions

 HArDCore3D::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.
 
const QuadratureRuleHArDCore3D::ElementQuad::get_quadT () const
 Returns quadrature rules in cell.
 
const QuadratureRuleHArDCore3D::ElementQuad::get_quadF (size_t ilF) const
 Returns quadrature rules on face with local number ilF.
 
const boost::multi_array< double, 2 > & HArDCore3D::ElementQuad::get_phiT_quadT () const
 Returns values of cell basis functions at cell quadrature nodes.
 
const boost::multi_array< VectorRd, 2 > & HArDCore3D::ElementQuad::get_dphiT_quadT () const
 Returns values of gradients of cell basis functions at cell quadrature nodes.
 
const boost::multi_array< double, 2 > & HArDCore3D::ElementQuad::get_phiT_quadF (size_t ilF) const
 Returns values of cell basis functions at face quadrature nodes, for face with local number ilF.
 
const boost::multi_array< double, 2 > HArDCore3D::ElementQuad::get_phiF_quadF (size_t ilF) const
 Returns values of face basis functions at face quadrature nodes, for face with local number ilF.
 
const boost::multi_array< VectorRd, 2 > & HArDCore3D::ElementQuad::get_dphiT_quadF (size_t ilF) const
 Returns values of gradients of cell basis functions at face quadrature nodes, for face with local number ilF.
 
boost::multi_array< VectorRd, 2 > HArDCore3D::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,0), (\phi_2,0,0), ..., (\phi_N,0,0), (0,\phi_1,0), (0,\phi_2,0) ... (0,\phi_N,0), (0,0,\phi_1) ... (0,0,\phi_N)\).
 
boost::multi_array< VectorRd, 2 > HArDCore3D::ElementQuad::get_vec_phiT_quadF (size_t ilF, size_t degree) const
 Builds on the fly the values of vector cell basis functions at face quadrature nodes. The vector basis is obtained by tensorization of the scalar one as in get_vec_phiT_quadT.
 
boost::multi_array< VectorRd, 2 > HArDCore3D::ElementQuad::get_vec_phiF_quadF (size_t ilF, size_t degree) const
 Builds on the fly the values of vector face basis functions at face quadrature nodes. The vector basis is obtained by tensorization of the scalar one as in get_vec_phiT_quadT.
 
template<typename GeometricSupport >
size_t const HArDCore3D::DimPoly (int m)
 
template<>
const size_t HArDCore3D::DimPoly< Cell > (const int m)
 Compute the size of the basis of 3-variate polynomials up to degree m.
 
template<>
const size_t HArDCore3D::DimPoly< Face > (const int m)
 Compute the size of the basis of 2-variate polynomials up to degree m.
 
template<>
const size_t HArDCore3D::DimPoly< Edge > (const int m)
 Compute the size of the basis of 1-variate polynomials up to degree m.
 
 HArDCore3D::UVector::UVector (const Eigen::VectorXd values, const Mesh &mesh, const int cell_deg, const size_t face_deg)
 
Eigen::VectorXd & HArDCore3D::UVector::asVectorXd () const
 Return the values as an Eigen vector.
 
const int HArDCore3D::UVector::get_cell_deg () const
 Return the cell degree.
 
const size_t HArDCore3D::UVector::get_face_deg () const
 Return the face degree.
 
const size_t HArDCore3D::UVector::n_cell_dofs () const
 Number of dofs in each cell.
 
const size_t HArDCore3D::UVector::n_face_dofs () const
 Number of dofs on each face.
 
const size_t HArDCore3D::UVector::n_total_cell_dofs () const
 Total number of cell dofs (in the vector, this is where the face dofs start)
 
Eigen::VectorXd HArDCore3D::UVector::restr (size_t iT) const
 Extract the restriction of the unknowns corresponding to cell iT and its faces.
 
UVector HArDCore3D::UVector::operator+ (const UVector &b)
 Overloads the addition: adds the coefficients.
 
UVector HArDCore3D::UVector::operator- (const UVector &b)
 Overloads the subtraction: subtracts the coefficients.
 
double HArDCore3D::UVector::operator() (size_t index) const
 Overloads the (): returns the corresponding coefficient.
 
 HArDCore3D::HybridCore::HybridCore (const Mesh *mesh_ptr, const int cell_deg, const size_t face_deg, const int 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.
 
const MeshHArDCore3D::HybridCore::get_mesh () const
 Returns a pointer to the mesh.
 
const int HArDCore3D::HybridCore::CellDegree () const
 Return the degree of cell polynomials.
 
const int HArDCore3D::HybridCore::CellDegreePos () const
 
const size_t HArDCore3D::HybridCore::FaceDegree () const
 Return the degree of face polynomials.
 
const PolyCellBasisTypeHArDCore3D::HybridCore::CellBasis (size_t iT) const
 Return cell basis for element with global index iT.
 
const PolyFaceBasisTypeHArDCore3D::HybridCore::FaceBasis (size_t iF) const
 Return face basis for face with global index iF.
 
const PolyEdgeBasisTypeHArDCore3D::HybridCore::EdgeBasis (size_t iE) const
 Return edge basis for edge with global index iE.
 
double HArDCore3D::HybridCore::L2norm (const UVector &Xh) const
 Compute L2 norm of a discrete function (using cell values)
 
double HArDCore3D::HybridCore::H1norm (const UVector &Xh) const
 Compute discrete H1 norm of a discrete function.
 
template<typename ContinuousFunction >
UVector HArDCore3D::HybridCore::interpolate (const ContinuousFunction &f, const int deg_cell, const size_t deg_face, size_t doe) const
 Compute the interpolant in the discrete space of a continuous function.
 
Eigen::VectorXd HArDCore3D::HybridCore::compute_weights (size_t iT) const
 Computes the weights to get cell values from face values when l=-1.
 
double HArDCore3D::HybridCore::evaluate_in_cell (const UVector Xh, size_t iT, VectorRd x) const
 Evaluates a discrete function in the cell iT at point x.
 
double HArDCore3D::HybridCore::evaluate_in_face (const UVector Xh, size_t iF, VectorRd x) const
 Evaluates a discrete function on the face iF at point x.
 
Eigen::VectorXd HArDCore3D::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.
 

Detailed Description

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

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

Typedef Documentation

◆ PolyCellBasisType

type for cell basis

◆ PolyEdgeBasisType

type for edge basis

◆ PolyFaceBasisType

type for face basis

Function Documentation

◆ asVectorXd()

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

Return the values as an Eigen vector.

◆ CellBasis()

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

Return cell basis for element with global index iT.

◆ CellDegree()

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

Return the degree of cell polynomials.

◆ CellDegreePos()

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

◆ compute_weights()

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

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

◆ DimPoly()

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

◆ DimPoly< Cell >()

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

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

Parameters
mPolynomial degree

◆ DimPoly< Edge >()

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

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

Parameters
mPolynomial degree

◆ DimPoly< Face >()

template<>
const size_t HArDCore3D::DimPoly< Face > ( const int  m)
inline

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

Parameters
mPolynomial degree

◆ EdgeBasis()

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

Return edge basis for edge with global index iE.

◆ 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 face 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_face()

double HybridCore::evaluate_in_face ( const UVector  Xh,
size_t  iF,
VectorRd  x 
) const

Evaluates a discrete function on the face iF at point x.

◆ FaceBasis()

const PolyFaceBasisType & HArDCore3D::HybridCore::FaceBasis ( size_t  iF) const
inline

Return face basis for face with global index iF.

◆ FaceDegree()

const size_t HArDCore3D::HybridCore::FaceDegree ( ) const
inline

Return the degree of face polynomials.

◆ get_cell_deg()

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

Return the cell degree.

◆ get_dphiT_quadF()

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

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

◆ get_dphiT_quadT()

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

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

◆ get_face_deg()

const size_t HArDCore3D::UVector::get_face_deg ( ) const
inline

Return the face degree.

◆ get_mesh()

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

Returns a pointer to the mesh.

◆ get_phiF_quadF()

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

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

◆ get_phiT_quadF()

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

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

◆ get_phiT_quadT()

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

Returns values of cell basis functions at cell quadrature nodes.

◆ get_quadF()

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

Returns quadrature rules on face with local number ilF.

◆ get_quadT()

const QuadratureRule & HArDCore3D::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 face basis functions at face 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 face 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 face 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 face 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,0), (\phi_2,0,0), ..., (\phi_N,0,0), (0,\phi_1,0), (0,\phi_2,0) ... (0,\phi_N,0), (0,0,\phi_1) ... (0,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  face_deg,
const int  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
face_degThe degree of the face 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()

UVector HArDCore3D::HybridCore::interpolate ( const ContinuousFunction f,
const int  deg_cell,
const size_t  deg_face,
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_facedegree of face 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 HArDCore3D::UVector::n_cell_dofs ( ) const
inline

Number of dofs in each cell.

◆ n_face_dofs()

const size_t HArDCore3D::UVector::n_face_dofs ( ) const
inline

Number of dofs on each face.

◆ n_total_cell_dofs()

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

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

◆ operator()()

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

Overloads the (): returns the corresponding coefficient.

◆ operator+()

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

Overloads the addition: adds the coefficients.

◆ operator-()

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

Overloads the subtraction: subtracts the coefficients.

◆ restr()

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

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

◆ UVector()

UVector::UVector ( const Eigen::VectorXd  values,
const Mesh mesh,
const int  cell_deg,
const size_t  face_deg 
)
Parameters
valuesvalues of the vector
meshreference to the mesh
cell_degpolynomial degrees in cell
face_degpolynomial degrees on faces

◆ 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 face polynomials
from_dofsType of unknowns to use: "cell" or "face"