HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
|
#include <elementquad.hpp>
Public Member Functions | |
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 & | get_quadT () const |
Returns quadrature rules in cell. | |
const QuadratureRule & | get_quadF (size_t ilF) const |
Returns quadrature rules on edge with local number ilF. | |
const boost::multi_array< double, 2 > & | get_phiT_quadT () const |
Returns values of cell basis functions at cell quadrature nodes. | |
const boost::multi_array< VectorRd, 2 > & | get_dphiT_quadT () const |
Returns values of gradients of cell basis functions at cell quadrature nodes. | |
const boost::multi_array< double, 2 > & | 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 > | 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 > & | 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 > | 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 > | 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 > | 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. | |
The ElementQuad class creates cell and edge quadrature rules, and vectors of values of basis functions and gradients at these points