|
HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
|
Classes defining the HHO method (scalar and vector-valued) More...
Classes | |
| class | HArDCore3D::HHOSpace |
| Class definition: polynomial bases and operators. More... | |
| struct | HArDCore3D::HHOSpace::CellBases |
| Structure to store element bases. More... | |
| struct | HArDCore3D::HHOSpace::FaceBases |
| Structure to store face bases. More... | |
| struct | HArDCore3D::HHOSpace::LocalOperators |
| A structure to store local operators (gradient, potential, stabilisation) More... | |
| class | HArDCore3D::VHHOSpace |
| Class definition: polynomial bases and operators. More... | |
| struct | HArDCore3D::VHHOSpace::CellBases |
| Structure to store element bases. More... | |
| struct | HArDCore3D::VHHOSpace::FaceBases |
| Structure to store face bases. More... | |
| struct | HArDCore3D::VHHOSpace::LocalOperators |
| A structure to store local operators (gradient, potential, stabilisation) More... | |
Functions | |
| HArDCore3D::HHOSpace::LocalOperators::LocalOperators (const Eigen::MatrixXd &_gradient, const Eigen::MatrixXd &_potential, const Eigen::MatrixXd &_stabilisation) | |
| HArDCore3D::HHOSpace::HHOSpace (const Mesh &mesh, size_t K, bool use_threads=true, std::ostream &output=std::cout) | |
| Constructor. | |
| const Mesh & | HArDCore3D::HHOSpace::mesh () const |
| Return a const reference to the mesh. | |
| const size_t & | HArDCore3D::HHOSpace::degree () const |
| Return the polynomial degree (common face and elements) | |
| Eigen::VectorXd | HArDCore3D::HHOSpace::interpolate (const FunctionType &q, const int doe_cell=-1, const int doe_face=-1) const |
| Interpolator of a continuous function. | |
| const CellBases & | HArDCore3D::HHOSpace::cellBases (size_t iT) const |
| Return cell bases for element iT. | |
| const CellBases & | HArDCore3D::HHOSpace::cellBases (const Cell &T) const |
| Return cell bases for cell T. | |
| const FaceBases & | HArDCore3D::HHOSpace::faceBases (size_t iF) const |
| Return face bases for face iF. | |
| const FaceBases & | HArDCore3D::HHOSpace::faceBases (const Face &F) const |
| Return cell bases for face F. | |
| const LocalOperators & | HArDCore3D::HHOSpace::operators (size_t iT) const |
| Return operators for the cell of index iT. | |
| const LocalOperators & | HArDCore3D::HHOSpace::operators (const Cell &T) const |
| Return cell operators for cell T. | |
| std::vector< std::pair< double, double > > | HArDCore3D::HHOSpace::computeNorms (const std::vector< Eigen::VectorXd > &list_dofs) const |
| Computes the discrete L2 (cell unknowns only) and H1 norms of a list of vectors. | |
| std::vector< double > | HArDCore3D::HHOSpace::computeVertexValues (const Eigen::VectorXd &u) const |
| Computes the values of the potential reconstruction at the mesh vertices. | |
| HArDCore3D::VHHOSpace::LocalOperators::LocalOperators (const Eigen::MatrixXd &_gradient, const Eigen::MatrixXd &_potential, const Eigen::MatrixXd &_stabilisation, const Eigen::MatrixXd &_potential_div, const Eigen::MatrixXd &_stabilisation_div) | |
| HArDCore3D::VHHOSpace::VHHOSpace (const Mesh &mesh, size_t K, const CellSelection &BoundaryStab, bool use_threads=true, std::ostream &output=std::cout) | |
| Constructor (with function to select cells in which boundary faces are used in the stabilisation) | |
| HArDCore3D::VHHOSpace::VHHOSpace (const Mesh &mesh, size_t K, bool use_threads=true, std::ostream &output=std::cout) | |
| Overloaded constructor when the selection of boundary stabilisation is not entered (all boundary faces are then used) | |
| const Mesh & | HArDCore3D::VHHOSpace::mesh () const |
| Return a const reference to the mesh. | |
| const size_t & | HArDCore3D::VHHOSpace::degree () const |
| Return the polynomial degree (common face and elements) | |
| const CellSelection & | HArDCore3D::VHHOSpace::boundaryStab () const |
| Return the function to select the cells with boundary stabilisation. | |
| Eigen::VectorXd | HArDCore3D::VHHOSpace::interpolate (const FunctionType &q, const int doe_cell=-1, const int doe_face=-1) const |
| Interpolator of a continuous function. | |
| const CellBases & | HArDCore3D::VHHOSpace::cellBases (size_t iT) const |
| Return cell bases for element iT. | |
| const CellBases & | HArDCore3D::VHHOSpace::cellBases (const Cell &T) const |
| Return cell bases for cell T. | |
| const FaceBases & | HArDCore3D::VHHOSpace::faceBases (size_t iF) const |
| Return face bases for face iF. | |
| const FaceBases & | HArDCore3D::VHHOSpace::faceBases (const Face &F) const |
| Return cell bases for face F. | |
| const LocalOperators & | HArDCore3D::VHHOSpace::operators (size_t iT) const |
| Return operators for the cell of index iT. | |
| const LocalOperators & | HArDCore3D::VHHOSpace::operators (const Cell &T) const |
| Return cell operators for cell T. | |
| std::vector< std::pair< double, double > > | HArDCore3D::VHHOSpace::computeNorms (const std::vector< Eigen::VectorXd > &list_dofs) const |
| Computes the discrete L2 (cell unknowns only) and H1 norms of a list of vectors. | |
| std::vector< VectorRd > | HArDCore3D::VHHOSpace::computeVertexValues (const Eigen::VectorXd &u) const |
| Computes the values of the potential reconstruction at the mesh vertices. | |
Classes defining the HHO method (scalar and vector-valued)
| typedef std::function<bool(const Cell &)> HArDCore3D::CellSelection |
| typedef std::function<double(const VectorRd &)> HArDCore3D::HHOSpace::FunctionType |
| typedef std::function<VectorRd(const VectorRd &)> HArDCore3D::VHHOSpace::FunctionType |
Geometric support.
Geometric support.
Geometric support.
Geometric support.
| typedef TensorizedVectorFamily<PolyBasisCellType, dimspace> HArDCore3D::HHOSpace::PolydBasisCellType |
| typedef TensorizedVectorFamily<PolyBasisCellType, dimspace> HArDCore3D::VHHOSpace::PolydBasisCellType |
| typedef Family<TensorizedVectorFamily<PolyBasisFaceType, dimspace> > HArDCore3D::VHHOSpace::PolydBasisFaceType |
|
inline |
Return the function to select the cells with boundary stabilisation.
Return cell bases for cell T.
Return cell bases for cell T.
Return cell bases for element iT.
Return cell bases for element iT.
| std::vector< std::pair< double, double > > HHOSpace::computeNorms | ( | const std::vector< Eigen::VectorXd > & | list_dofs | ) | const |
Computes the discrete L2 (cell unknowns only) and H1 norms of a list of vectors.
| list_dofs | The list of vectors representing the dofs |
| std::vector< std::pair< double, double > > VHHOSpace::computeNorms | ( | const std::vector< Eigen::VectorXd > & | list_dofs | ) | const |
Computes the discrete L2 (cell unknowns only) and H1 norms of a list of vectors.
| list_dofs | The list of vectors representing the dofs |
Computes the values of the potential reconstruction at the mesh vertices.
| u | DOFs in the discrete space |
| std::vector< VectorRd > VHHOSpace::computeVertexValues | ( | const Eigen::VectorXd & | u | ) | const |
Computes the values of the potential reconstruction at the mesh vertices.
| u | DOFs in the discrete space |
Return the polynomial degree (common face and elements)
Return the polynomial degree (common face and elements)
Return cell bases for face F.
Return cell bases for face F.
Return face bases for face iF.
Return face bases for face iF.
| HHOSpace::HHOSpace | ( | const Mesh & | mesh, |
| size_t | K, | ||
| bool | use_threads = true, |
||
| std::ostream & | output = std::cout |
||
| ) |
Constructor.
| Eigen::VectorXd HHOSpace::interpolate | ( | const FunctionType & | q, |
| const int | doe_cell = -1, |
||
| const int | doe_face = -1 |
||
| ) | const |
Interpolator of a continuous function.
| q | The function to interpolate |
| doe_cell | The optional degre of cell quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used. |
| doe_face | The optional degre of face quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used. |
| Eigen::VectorXd VHHOSpace::interpolate | ( | const FunctionType & | q, |
| const int | doe_cell = -1, |
||
| const int | doe_face = -1 |
||
| ) | const |
Interpolator of a continuous function.
| q | The function to interpolate |
| doe_cell | The optional degre of cell quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used. |
| doe_face | The optional degre of face quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used. |
|
inline |
| _gradient | Gradient operator |
| _potential | Potential operator |
| _stabilisation | Stabilisation bilinear form |
|
inline |
| _gradient | Gradient operator |
| _potential | HHO Potential operator |
| _stabilisation | H1 Stabilisation bilinear form associated to the HHO potential |
| _potential_div | Potential reconstructed from divergence |
| _stabilisation_div | L2 Stabilisation associated to _potential_div |
|
inline |
Return cell operators for cell T.
|
inline |
Return cell operators for cell T.
|
inline |
Return operators for the cell of index iT.
|
inline |
Return operators for the cell of index iT.
|
inline |
Overloaded constructor when the selection of boundary stabilisation is not entered (all boundary faces are then used)
| VHHOSpace::VHHOSpace | ( | const Mesh & | mesh, |
| size_t | K, | ||
| const CellSelection & | BoundaryStab, | ||
| bool | use_threads = true, |
||
| std::ostream & | output = std::cout |
||
| ) |
Constructor (with function to select cells in which boundary faces are used in the stabilisation)
| Eigen::MatrixXd HArDCore3D::HHOSpace::LocalOperators::gradient |
| Eigen::MatrixXd HArDCore3D::VHHOSpace::LocalOperators::gradient |
| std::unique_ptr<PolyBasisCellType> HArDCore3D::HHOSpace::CellBases::Polyk |
| std::unique_ptr<PolyBasisFaceType> HArDCore3D::HHOSpace::FaceBases::Polyk |
| std::unique_ptr<PolyBasisCellType> HArDCore3D::VHHOSpace::CellBases::Polyk |
| std::unique_ptr<PolyBasisFaceType> HArDCore3D::VHHOSpace::FaceBases::Polyk |
| std::unique_ptr<PolydBasisCellType> HArDCore3D::HHOSpace::CellBases::Polykd |
| std::unique_ptr<PolydBasisCellType> HArDCore3D::VHHOSpace::CellBases::Polykd |
| std::unique_ptr<PolydBasisFaceType> HArDCore3D::VHHOSpace::FaceBases::Polykd |
| std::unique_ptr<PolydxdBasisCellType> HArDCore3D::VHHOSpace::CellBases::Polykdxd |
| std::unique_ptr<PolyBasisCellType> HArDCore3D::HHOSpace::CellBases::Polykpo |
| std::unique_ptr<PolyBasisCellType> HArDCore3D::VHHOSpace::CellBases::Polykpo |
| std::unique_ptr<PolydBasisCellType> HArDCore3D::VHHOSpace::CellBases::Polykpod |
| Eigen::MatrixXd HArDCore3D::HHOSpace::LocalOperators::potential |
| Eigen::MatrixXd HArDCore3D::VHHOSpace::LocalOperators::potential |
| Eigen::MatrixXd HArDCore3D::VHHOSpace::LocalOperators::potential_div |
| Eigen::MatrixXd HArDCore3D::HHOSpace::LocalOperators::stabilisation |
| Eigen::MatrixXd HArDCore3D::VHHOSpace::LocalOperators::stabilisation |
| Eigen::MatrixXd HArDCore3D::VHHOSpace::LocalOperators::stabilisation_div |