| 
    HArD::Core2D
    
   Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns 
   | 
 
#include <lepnccore.hpp>


Public Types | |
| using | cell_basis_type = std::function< double(double, double)> | 
| type for cell basis   | |
| using | cell_gradient_type = std::function< VectorRd(double, double)> | 
| type for gradients of cell basis   | |
  Public Types inherited from HArDCore2D::HybridCore | |
| typedef Family< MonomialScalarBasisCell > | PolyCellBasisType | 
| type for cell basis   | |
| typedef Family< MonomialScalarBasisEdge > | PolyEdgeBasisType | 
| type for edge basis   | |
Public Member Functions | |
| LEPNCCore (const Mesh *mesh_ptr, const size_t K, const size_t L, std::ostream &output=std::cout) | |
| Class constructor: initialise the LEPNCCore class, and creates non-conforming basis functions (and gradients)   | |
| const cell_basis_type & | nc_basis (size_t iT, size_t i) const | 
| Return a reference to the i'th non-conforming basis function of the cell iT.   | |
| const cell_gradient_type & | nc_gradient (size_t iT, size_t i) const | 
| Return a reference to the gradient of the i'th non-conforming basis function of the cell iT.   | |
| const VectorRd | ml_node (size_t iT, size_t i) const | 
| Return the i'th nodal point in cell iT (for mass lumping)   | |
| const std::vector< Eigen::ArrayXd > | nc_basis_quad (const size_t iT, const QuadratureRule quad) const | 
| Computes non-conforming basis functions at the given quadrature nodes.   | |
| const std::vector< Eigen::ArrayXXd > | grad_nc_basis_quad (const size_t iT, const QuadratureRule quad) const | 
| Compute \((\nabla \phi_i^{nc})_{i\in I}\) at the quadrature nodes, where \((\phi_i^{nc})_{i\in I}\) are the cell basis functions.   | |
| Eigen::MatrixXd | gram_matrix (const std::vector< Eigen::ArrayXd > &f_quad, const std::vector< Eigen::ArrayXd > &g_quad, const size_t &nrows, const size_t &ncols, const QuadratureRule &quad, const bool &sym, std::vector< double > L2weight={}) const | 
| Eigen::MatrixXd | gram_matrix (const std::vector< Eigen::ArrayXXd > &F_quad, const std::vector< Eigen::ArrayXXd > &G_quad, const size_t &nrows, const size_t &ncols, const QuadratureRule &quad, const bool &sym, std::vector< Eigen::Matrix2d > L2Weight={}) const | 
| Overloaded version of the previous one for vector-valued functions: the functions (F_i) and (G_j) are vector-valued functions.   | |
| template<typename Function > | |
| Eigen::VectorXd | nc_interpolate_moments (const Function &f, size_t doe) const | 
| Interpolates a continuous function on the degrees of freedom, using the moments on the basis functions. The first ones are the cell DOFs (DimPoly<Cell>(1) for each cell), the last ones are the edge DOFs (one for each edge)   | |
| template<typename Function > | |
| Eigen::VectorXd | nc_interpolate_ml (const Function &f, size_t doe) const | 
| Interpolates a continuous function on the degrees of freedom, using the moments on the basis functions associated to the edges and the nodal values (corresponding to mass-lumping) on the cell basis functions. The first ones are the cell DOFs (DimPoly<Cell>(1) for each cell), the last ones are the edge DOFs (one for each edge)   | |
| Eigen::VectorXd | nc_restr (const Eigen::VectorXd &Xh, size_t iT) const | 
| Extract from a global vector Xh of unknowns the non-conforming unknowns corresponding to cell iT.   | |
| double | nc_L2norm (const Eigen::VectorXd &Xh) const | 
| Compute L2 norm of a discrete function (given by coefficients on the basis functions)   | |
| double | nc_L2norm_ml (const Eigen::VectorXd &Xh) const | 
| Compute L2 norm of the mass-lumped discrete function (given by coefficients on the basis functions)   | |
| double | nc_H1norm (const Eigen::VectorXd &Xh) const | 
| Compute broken H1 norm of a discrete function (given by coefficients on the basis functions)   | |
| double | nc_evaluate_in_cell (const Eigen::VectorXd XTF, size_t iT, double x, double y) const | 
| Evaluates a non-conforming discrete function in the cell iT at point (x,y)   | |
| Eigen::VectorXd | nc_VertexValues (const Eigen::VectorXd Xh, const double weight=0) | 
| From a non-conforming discrete function, computes a vector of values at the vertices of the mesh.   | |
  Public Member Functions inherited from 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 * | get_mesh () const | 
| Returns a pointer to the mesh.   | |
| const int | CellDegree () const | 
| Return the degree of cell polynomials.   | |
| const int | CellDegreePos () const | 
| const size_t | EdgeDegree () const | 
| Return the degree of edge polynomials.   | |
| const PolyCellBasisType & | CellBasis (size_t iT) const | 
| Return cell basis for element with global index iT.   | |
| const PolyEdgeBasisType & | EdgeBasis (size_t iE) const | 
| Return edge basis for edge with global index iE.   | |
| double | L2norm (const UVector &Xh) const | 
| Compute L2 norm of a discrete function (using cell values)   | |
| double | H1norm (const UVector &Xh) const | 
| Compute discrete H1 norm of a discrete function.   | |
| template<typename ContinuousFunction > | |
| UVector | 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 | compute_weights (size_t iT) const | 
| Computes the weights to get cell values from edge values when l=-1.   | |
| double | evaluate_in_cell (const UVector Xh, size_t iT, VectorRd x) const | 
| Evaluates a discrete function in the cell iT at point x.   | |
| double | 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 | VertexValues (const UVector Xh, const std::string from_dofs) | 
| From a hybrid function, computes a vector of values at the vertices of the mesh.   | |
The LEPNCCore class provides basis functions for non-conforming schemes on generic polygonal meshes