Classes providing cell and face quadrature rules, and values of basis functions at the nodes.
More...
|
| | 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 QuadratureRule & | HArDCore3D::ElementQuad::get_quadT () const |
| | Returns quadrature rules in cell.
|
| |
| const QuadratureRule & | HArDCore3D::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 Mesh * | HArDCore3D::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 PolyCellBasisType & | HArDCore3D::HybridCore::CellBasis (size_t iT) const |
| | Return cell basis for element with global index iT.
|
| |
| const PolyFaceBasisType & | HArDCore3D::HybridCore::FaceBasis (size_t iF) const |
| | Return face basis for face with global index iF.
|
| |
| const PolyEdgeBasisType & | HArDCore3D::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.
|
| |
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.
| 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
-
| degree | maximal degree of basis functions that is required |