HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
|
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< MonomialScalarBasisCell > | HArDCore2D::HybridCore::PolyCellBasisType |
type for cell basis | |
typedef Family< MonomialScalarBasisEdge > | HArDCore2D::HybridCore::PolyEdgeBasisType |
type for edge basis | |
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. | |
const QuadratureRule & | HArDCore2D::ElementQuad::get_quadT () const |
Returns quadrature rules in cell. | |
const QuadratureRule & | HArDCore2D::ElementQuad::get_quadF (size_t ilF) const |
Returns quadrature rules on edge with local number ilF. | |
const boost::multi_array< double, 2 > & | HArDCore2D::ElementQuad::get_phiT_quadT () const |
Returns values of cell basis functions at cell quadrature nodes. | |
const boost::multi_array< VectorRd, 2 > & | HArDCore2D::ElementQuad::get_dphiT_quadT () const |
Returns values of gradients of cell basis functions at cell quadrature nodes. | |
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. | |
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. | |
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. | |
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)\). | |
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. | |
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. | |
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. | |
template<> | |
const size_t | HArDCore2D::DimPoly< Edge > (const int m) |
Compute the size of the basis of 1-variate polynomials up to degree m. | |
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. | |
const int | HArDCore2D::UVector::get_cell_deg () const |
Return the cell degree. | |
const size_t | HArDCore2D::UVector::get_edge_deg () const |
Return the edge degree. | |
const size_t | HArDCore2D::UVector::n_cell_dofs () const |
Number of dofs in each cell. | |
const size_t | HArDCore2D::UVector::n_edge_dofs () const |
Number of dofs on each edge. | |
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) | |
Eigen::VectorXd | HArDCore2D::UVector::restr (size_t iT) const |
Extract the restriction of the unknowns corresponding to cell iT and its edges. | |
UVector | HArDCore2D::UVector::operator+ (const UVector &b) |
Overloads the addition: adds the coefficients. | |
UVector | HArDCore2D::UVector::operator- (const UVector &b) |
Overloads the subtraction: subtracts the coefficients. | |
void | HArDCore2D::UVector::operator+= (const UVector &b) |
Overloads the increment: adds the coefficients. | |
void | HArDCore2D::UVector::operator/= (const double &r) |
Overloads the increment: adds the coefficients. | |
double | HArDCore2D::UVector::operator() (size_t index) const |
Overloads the (): returns the corresponding coefficient. | |
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. | |
const Mesh * | HArDCore2D::HybridCore::get_mesh () const |
Returns a pointer to the mesh. | |
const int | HArDCore2D::HybridCore::CellDegree () const |
Return the degree of cell polynomials. | |
const int | HArDCore2D::HybridCore::CellDegreePos () const |
const size_t | HArDCore2D::HybridCore::EdgeDegree () const |
Return the degree of edge polynomials. | |
const PolyCellBasisType & | HArDCore2D::HybridCore::CellBasis (size_t iT) const |
Return cell basis for element with global index iT. | |
const PolyEdgeBasisType & | HArDCore2D::HybridCore::EdgeBasis (size_t iE) const |
Return edge basis for edge with global index iE. | |
double | HArDCore2D::HybridCore::L2norm (const UVector &Xh) const |
Compute L2 norm of a discrete function (using cell values) | |
double | HArDCore2D::HybridCore::H1norm (const UVector &Xh) const |
Compute discrete H1 norm of a discrete function. | |
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. | |
Eigen::VectorXd | HArDCore2D::HybridCore::compute_weights (size_t iT) const |
Computes the weights to get cell values from edge values when l=-1. | |
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. | |
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. | |
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. | |
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.
type for cell basis
type for edge basis
|
inline |
Return the values as an Eigen vector.
|
inline |
Return cell basis for element with global index iT.
|
inline |
Return the degree of cell polynomials.
|
inline |
Eigen::VectorXd HybridCore::compute_weights | ( | size_t | iT | ) | const |
Computes the weights to get cell values from edge values when l=-1.
|
inline |
|
inline |
Compute the size of the basis of 2-variate polynomials up to degree m.
m | Polynomial degree |
|
inline |
Compute the size of the basis of 1-variate polynomials up to degree m.
m | Polynomial degree |
|
inline |
Return edge basis for edge with global index iE.
|
inline |
Return the degree of edge polynomials.
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.
hho | A reference to the hybridcore instance |
iT | Global index of cell |
doeT | The degree of exactness for cell quadratures |
doeF | The degree of exactness of edge quadratures |
Evaluates a discrete function in the cell iT at point x.
Evaluates a discrete function on the edge iE at point x.
|
inline |
Return the cell degree.
|
inline |
Returns values of gradients of cell basis functions at edge quadrature nodes, for edge with local number ilF.
|
inline |
Returns values of gradients of cell basis functions at cell quadrature nodes.
|
inline |
Return the edge degree.
|
inline |
Returns a pointer to the mesh.
|
inline |
Returns values of edge basis functions at edge quadrature nodes, for edge with local number ilF.
|
inline |
Returns values of cell basis functions at edge quadrature nodes, for edge with local number ilF.
|
inline |
Returns values of cell basis functions at cell quadrature nodes.
|
inline |
Returns quadrature rules on edge with local number ilF.
|
inline |
Returns quadrature rules in cell.
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.
degree | local number of edge required degree of basis function |
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.
degree | local number of edge maximum degree of basis functions required |
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)\).
degree | maximal degree of basis functions that is required |
double HybridCore::H1norm | ( | const UVector & | Xh | ) | const |
Compute discrete H1 norm of a discrete function.
Xh | Vector of unknowns |
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.
mesh_ptr | A pointer to the loaded mesh |
cell_deg | The degree of the cell polynomials |
edge_deg | The degree of the edge polynomials |
use_threads | Optional argument to indicate if threads should be used |
output | Optional argument for specifying outputs of messages |
ortho | Optional argument for choosing to orthonormalise or not the basis functions |
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.
f | function to interpolate |
deg_cell | degree of cell polynomials for interpolation |
deg_edge | degree of edge polynomials for interpolation |
doe | degree of exactness of the quadrature rules to compute the interpolate |
double HybridCore::L2norm | ( | const UVector & | Xh | ) | const |
Compute L2 norm of a discrete function (using cell values)
Xh | Vector of unknowns |
|
inline |
Number of dofs in each cell.
|
inline |
Number of dofs on each edge.
|
inline |
Total number of cell dofs (in the vector, this is where the edge dofs start)
|
inline |
Overloads the (): returns the corresponding coefficient.
Overloads the addition: adds the coefficients.
|
inline |
Overloads the increment: adds the coefficients.
Overloads the subtraction: subtracts the coefficients.
|
inline |
Overloads the increment: adds the coefficients.
Eigen::VectorXd UVector::restr | ( | size_t | iT | ) | const |
Extract the restriction of the unknowns corresponding to cell iT and its edges.
UVector::UVector | ( | const Eigen::VectorXd & | values, |
const Mesh & | mesh, | ||
const int | cell_deg, | ||
const size_t | edge_deg | ||
) |
values | values of the vector |
mesh | reference to the mesh |
cell_deg | polynomial degrees in cell |
edge_deg | polynomial degrees on edge |
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.
Xh | vector of discrete unknowns on cell and edge polynomials |
from_dofs | Type of unknowns to use: "cell" or "edge" |