HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
|
Classes defining the HHO method (scalar and vector-valued) More...
Classes | |
class | HArDCore2D::HHOSpace |
Class definition: polynomial bases and operators. More... | |
struct | HArDCore2D::HHOSpace::CellBases |
Structure to store element bases. More... | |
struct | HArDCore2D::HHOSpace::EdgeBases |
Structure to store edge bases. More... | |
struct | HArDCore2D::HHOSpace::LocalOperators |
A structure to store local operators (gradient, potential, stabilisation) More... | |
class | HArDCore2D::VHHOSpace |
Class definition: polynomial bases and operators. More... | |
struct | HArDCore2D::VHHOSpace::CellBases |
Structure to store element bases. More... | |
struct | HArDCore2D::VHHOSpace::EdgeBases |
Structure to store edge bases. More... | |
struct | HArDCore2D::VHHOSpace::LocalOperators |
A structure to store local operators (gradient, potential, stabilisation) More... | |
Functions | |
HArDCore2D::HHOSpace::LocalOperators::LocalOperators (const Eigen::MatrixXd &_gradient, const Eigen::MatrixXd &_potential, const Eigen::MatrixXd &_stabilisation) | |
HArDCore2D::HHOSpace::HHOSpace (const Mesh &mesh, size_t K, bool use_threads=true, std::ostream &output=std::cout) | |
Constructor. | |
const Mesh & | HArDCore2D::HHOSpace::mesh () const |
Return a const reference to the mesh. | |
const size_t & | HArDCore2D::HHOSpace::degree () const |
Return the polynomial degree (common edge and elements) | |
Eigen::VectorXd | HArDCore2D::HHOSpace::interpolate (const FunctionType &q, const int doe_cell=-1, const int doe_edge=-1) const |
Interpolator of a continuous function. | |
const CellBases & | HArDCore2D::HHOSpace::cellBases (size_t iT) const |
Return cell bases for element iT. | |
const CellBases & | HArDCore2D::HHOSpace::cellBases (const Cell &T) const |
Return cell bases for cell T. | |
const EdgeBases & | HArDCore2D::HHOSpace::edgeBases (size_t iE) const |
Return edge bases for edge iE. | |
const EdgeBases & | HArDCore2D::HHOSpace::edgeBases (const Edge &E) const |
Return cell bases for edge E. | |
const LocalOperators & | HArDCore2D::HHOSpace::operators (size_t iT) const |
Return operators for the cell of index iT. | |
const LocalOperators & | HArDCore2D::HHOSpace::operators (const Cell &T) const |
Return cell operators for cell T. | |
std::vector< std::pair< double, double > > | HArDCore2D::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. | |
Eigen::VectorXd | HArDCore2D::HHOSpace::computeVertexValues (const Eigen::VectorXd &u) const |
Computes the values of the potential reconstruction at the mesh vertices. | |
HArDCore2D::VHHOSpace::LocalOperators::LocalOperators (const Eigen::MatrixXd &_gradient, const Eigen::MatrixXd &_symmetric_gradient, const Eigen::MatrixXd &_potential, const Eigen::MatrixXd &_potential_sym, const Eigen::MatrixXd &_stabilisation, const Eigen::MatrixXd &_difference) | |
HArDCore2D::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) | |
HArDCore2D::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 & | HArDCore2D::VHHOSpace::mesh () const |
Return a const reference to the mesh. | |
const size_t & | HArDCore2D::VHHOSpace::degree () const |
Return the polynomial degree (common face and elements) | |
const CellSelection & | HArDCore2D::VHHOSpace::boundaryStab () const |
Return the function to select the cells with boundary stabilisation. | |
Eigen::VectorXd | HArDCore2D::VHHOSpace::interpolate (const FunctionType &q, const int doe_cell=-1, const int doe_face=-1) const |
Interpolator of a continuous function. | |
Eigen::VectorXd | HArDCore2D::VHHOSpace::local_interpolate (size_t iT, const FunctionType &q, const int doe_cell=-1, const int doe_face=-1) const |
Local Interpolator of a continuous function on the element iT. | |
Eigen::VectorXd | HArDCore2D::VHHOSpace::components (size_t iT, const FunctionType &q, const int doe=-1) const |
Discrete components of continuous function in [P^k(T)]^d TODO more general function, with as input the basis space. | |
Eigen::VectorXd | HArDCore2D::VHHOSpace::components (size_t iT, const GradientType &q, const int doe=-1) const |
Discrete components of continuous function in [P^k(T;Symm)]^dxd TODO more general function, with as input the basis space. | |
const CellBases & | HArDCore2D::VHHOSpace::cellBases (size_t iT) const |
Return cell bases for element iT. | |
const CellBases & | HArDCore2D::VHHOSpace::cellBases (const Cell &T) const |
Return cell bases for cell T. | |
const EdgeBases & | HArDCore2D::VHHOSpace::edgeBases (size_t iE) const |
Return edge bases for edge iE. | |
const EdgeBases & | HArDCore2D::VHHOSpace::edgeBases (const Edge &E) const |
Return cell bases for edge E. | |
const LocalOperators & | HArDCore2D::VHHOSpace::operators (size_t iT) const |
Return operators for the cell of index iT. | |
const LocalOperators & | HArDCore2D::VHHOSpace::operators (const Cell &T) const |
Return cell operators for cell T. | |
std::vector< std::pair< double, double > > | HArDCore2D::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< Eigen::VectorXd > | HArDCore2D::VHHOSpace::computeVertexValues (const Eigen::VectorXd &u) const |
Classes defining the HHO method (scalar and vector-valued)
typedef std::function<bool(const Cell &)> HArDCore2D::CellSelection |
typedef std::function<double(const VectorRd &)> HArDCore2D::HHOSpace::FunctionType |
typedef std::function<VectorRd(const VectorRd &)> HArDCore2D::VHHOSpace::FunctionType |
typedef Cell HArDCore2D::HHOSpace::CellBases::GeometricSupport |
Geometric support.
typedef Edge HArDCore2D::HHOSpace::EdgeBases::GeometricSupport |
Geometric support.
typedef Cell HArDCore2D::VHHOSpace::CellBases::GeometricSupport |
Geometric support.
typedef Edge HArDCore2D::VHHOSpace::EdgeBases::GeometricSupport |
Geometric support.
typedef std::function<MatrixRd(const VectorRd &)> HArDCore2D::VHHOSpace::GradientType |
typedef TensorizedVectorFamily<PolyBasisCellType, dimspace> HArDCore2D::HHOSpace::PolydBasisCellType |
typedef TensorizedVectorFamily<PolyBasisCellType, dimspace> HArDCore2D::VHHOSpace::PolydBasisCellType |
typedef Family<TensorizedVectorFamily<PolyBasisEdgeType, dimspace> > HArDCore2D::VHHOSpace::PolydBasisEdgeType |
|
inline |
Return the function to select the cells with boundary stabilisation.
|
inline |
Return cell bases for cell T.
|
inline |
Return cell bases for cell T.
|
inline |
Return cell bases for element iT.
|
inline |
Return cell bases for element iT.
Eigen::VectorXd VHHOSpace::components | ( | size_t | iT, |
const FunctionType & | q, | ||
const int | doe = -1 |
||
) | const |
Discrete components of continuous function in [P^k(T)]^d TODO more general function, with as input the basis space.
q | The function to interpolate |
doe | The optional degre of cell quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used. |
Eigen::VectorXd VHHOSpace::components | ( | size_t | iT, |
const GradientType & | q, | ||
const int | doe = -1 |
||
) | const |
Discrete components of continuous function in [P^k(T;Symm)]^dxd TODO more general function, with as input the basis space.
q | The function to interpolate |
doe | The optional degre of cell quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used. |
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 |
Eigen::VectorXd HHOSpace::computeVertexValues | ( | const Eigen::VectorXd & | u | ) | const |
Computes the values of the potential reconstruction at the mesh vertices.
u | DOFs in the discrete space |
std::vector< Eigen::VectorXd > VHHOSpace::computeVertexValues | ( | const Eigen::VectorXd & | u | ) | const |
Computes the values of the potential reconstruction at the mesh vertices and stores them in a way that is compatible with vtu_writer
u | DOFs in the discrete space |
|
inline |
Return the polynomial degree (common edge and elements)
|
inline |
Return the polynomial degree (common face and elements)
|
inline |
Return cell bases for edge E.
|
inline |
Return cell bases for edge E.
|
inline |
Return edge bases for edge iE.
|
inline |
Return edge bases for edge iE.
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_edge = -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_edge | The optional degre of edge 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. |
Eigen::VectorXd VHHOSpace::local_interpolate | ( | size_t | iT, |
const FunctionType & | q, | ||
const int | doe_cell = -1 , |
||
const int | doe_face = -1 |
||
) | const |
Local Interpolator of a continuous function on the element iT.
iT | The element of the local interpolator |
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 |
_symmetric_gradient | Symmetric Gradient |
_potential | HHO Potential operator |
_stabilisation | H1 Stabilisation bilinear form associated to the HHO potential |
|
inline |
Return a const reference to the mesh.
|
inline |
Return a const reference to the mesh.
|
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)
|
static |
Eigen::MatrixXd HArDCore2D::VHHOSpace::LocalOperators::difference |
Eigen::MatrixXd HArDCore2D::HHOSpace::LocalOperators::gradient |
Eigen::MatrixXd HArDCore2D::VHHOSpace::LocalOperators::gradient |
std::unique_ptr<PolyBasisCellType> HArDCore2D::HHOSpace::CellBases::Polyk |
std::unique_ptr<PolyBasisEdgeType> HArDCore2D::HHOSpace::EdgeBases::Polyk |
std::unique_ptr<PolyBasisCellType> HArDCore2D::VHHOSpace::CellBases::Polyk |
std::unique_ptr<PolyBasisEdgeType> HArDCore2D::VHHOSpace::EdgeBases::Polyk |
std::unique_ptr<PolydBasisCellType> HArDCore2D::HHOSpace::CellBases::Polykd |
std::unique_ptr<PolydBasisCellType> HArDCore2D::VHHOSpace::CellBases::Polykd |
std::unique_ptr<PolydBasisEdgeType> HArDCore2D::VHHOSpace::EdgeBases::Polykd |
std::unique_ptr<PolydxdBasisCellType> HArDCore2D::VHHOSpace::CellBases::Polykdxd |
std::unique_ptr<PolyBasisCellType> HArDCore2D::HHOSpace::CellBases::Polykpo |
std::unique_ptr<PolyBasisCellType> HArDCore2D::VHHOSpace::CellBases::Polykpo |
std::unique_ptr<PolydBasisCellType> HArDCore2D::VHHOSpace::CellBases::Polykpod |
std::unique_ptr<PolySymdxdBasisCellType> HArDCore2D::VHHOSpace::CellBases::PolySymkdxd |
Eigen::MatrixXd HArDCore2D::HHOSpace::LocalOperators::potential |
Eigen::MatrixXd HArDCore2D::VHHOSpace::LocalOperators::potential |
Eigen::MatrixXd HArDCore2D::VHHOSpace::LocalOperators::potential_sym |
Eigen::MatrixXd HArDCore2D::HHOSpace::LocalOperators::stabilisation |
Eigen::MatrixXd HArDCore2D::VHHOSpace::LocalOperators::stabilisation |
Eigen::MatrixXd HArDCore2D::VHHOSpace::LocalOperators::symmetric_gradient |