HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
|
#include <basis.hpp>
Public Types | |
typedef BasisType::FunctionValue | FunctionValue |
typedef BasisType::GradientValue | GradientValue |
typedef BasisType::CurlValue | CurlValue |
typedef BasisType::DivergenceValue | DivergenceValue |
typedef BasisType::RotorValue | RotorValue |
typedef BasisType::HessianValue | HessianValue |
typedef BasisType::GeometricSupport | GeometricSupport |
typedef BasisType | AncestorType |
Public Member Functions | |
Family (const BasisType &basis, const Eigen::MatrixXd &matrix) | |
Constructor. | |
size_t | dimension () const |
Dimension of the family. This is actually the number of functions in the family, not necessarily linearly independent. | |
FunctionValue | function (size_t i, const VectorRd &x) const |
Evaluate the i-th function at point x. | |
FunctionValue | function (size_t i, size_t iqn, const boost::multi_array< FunctionValue, 2 > &ancestor_value_quad) const |
Evaluate the i-th function at a quadrature point iqn, knowing all the values of ancestor basis functions at the quadrature nodes (provided by eval_quad) | |
GradientValue | gradient (size_t i, const VectorRd &x) const |
Evaluate the gradient of the i-th function at point x. | |
GradientValue | gradient (size_t i, size_t iqn, const boost::multi_array< GradientValue, 2 > &ancestor_gradient_quad) const |
Evaluate the gradient of the i-th function at a quadrature point iqn, knowing all the gradients of ancestor basis functions at the quadrature nodes (provided by eval_quad) | |
CurlValue | curl (size_t i, const VectorRd &x) const |
Evaluate the curl of the i-th function at point x. | |
CurlValue | curl (size_t i, size_t iqn, const boost::multi_array< CurlValue, 2 > &ancestor_curl_quad) const |
Evaluate the curl of the i-th function at a quadrature point iqn, knowing all the curls of ancestor basis functions at the quadrature nodes (provided by eval_quad) | |
DivergenceValue | divergence (size_t i, const VectorRd &x) const |
Evaluate the divergence of the i-th function at point x. | |
DivergenceValue | divergence (size_t i, size_t iqn, const boost::multi_array< DivergenceValue, 2 > &ancestor_divergence_quad) const |
Evaluate the divergence of the i-th function at a quadrature point iqn, knowing all the divergences of ancestor basis functions at the quadrature nodes (provided by eval_quad) | |
DivergenceValue | rotor (size_t i, const VectorRd &x) const |
Evaluate the scalar rotor of the i-th function at point x. | |
RotorValue | rotor (size_t i, size_t iqn, const boost::multi_array< RotorValue, 2 > &ancestor_rotor_quad) const |
Evaluate the scalar rotor of the i-th function at a quadrature point iqn, knowing all the scalar rotors of ancestor basis functions at the quadrature nodes (provided by eval_quad) | |
HessianValue | hessian (size_t i, const VectorRd &x) const |
Evaluate the Hessian of the i-th basis function at point x. | |
HessianValue | hessian (size_t i, size_t iqn, const boost::multi_array< HessianValue, 2 > &ancestor_hessian_quad) const |
Evaluate the hessian of the i-th function at a quadrature point iqn, knowing all the hessian of ancestor basis functions at the quadrature nodes (provided by eval_quad) | |
const Eigen::MatrixXd & | matrix () const |
Return the coefficient matrix. | |
constexpr const BasisType & | ancestor () const |
Return the ancestor. | |
size_t | max_degree () const |
Returns the maximum degree of the basis functions. | |
Static Public Attributes | |
static constexpr const TensorRankE | tensorRank = BasisType::tensorRank |
static constexpr const bool | hasAncestor = true |
static const bool | hasFunction = BasisType::hasFunction |
static const bool | hasGradient = BasisType::hasGradient |
static const bool | hasCurl = BasisType::hasCurl |
static const bool | hasDivergence = BasisType::hasDivergence |
static const bool | hasRotor = BasisType::hasRotor |
static const bool | hasHessian = BasisType::hasHessian |
Family of functions expressed as linear combination of the functions of a given basis
If \((f_1,...,f_r)\) is the basis, then the family is \((\phi_1,...,\phi_l)\) where \(\phi_i = \sum_j M_ij f_j\). The matrix \(M\) is the member m_matrix.