HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
|
Classes providing tools to the Discrete De Rham sequence. More...
Classes | |
class | HArDCore2D::DDRCore |
Construct all polynomial spaces for the DDR sequence. More... | |
struct | HArDCore2D::DDRCore::CellBases |
Structure to store element bases. More... | |
struct | HArDCore2D::DDRCore::EdgeBases |
Structure to store edge bases. More... | |
class | HArDCore2D::EXCurl |
Extended XCurl space, with vector-valued polynomials on the edges. More... | |
struct | HArDCore2D::EXCurl::hhoLocalOperators |
A structure to store the local HHO operators. More... | |
class | HArDCore2D::SerendipityProblem |
Construct all polynomial spaces for the DDR sequence. More... | |
class | HArDCore2D::SXGrad |
Discrete Serendipity Hgrad space: local operators, L2 product and global interpolator. More... | |
struct | HArDCore2D::SXGrad::TransferOperators |
A structure to store the serendipity, extension and reduction operators. More... | |
class | HArDCore2D::VSXGrad |
Vector version of sXgrad, the arbitrary order space with nodal primal unknowns. More... | |
struct | HArDCore2D::VSXGrad::LocalOperators |
struct | HArDCore2D::VSXGrad::CellBases |
struct | HArDCore2D::VSXGrad::EdgeBases |
class | HArDCore2D::XCurl |
Discrete Hcurl space: local operators, L2 product and global interpolator. More... | |
struct | HArDCore2D::XCurl::LocalOperators |
A structure to store the local operators (curl and potential) More... | |
class | HArDCore2D::XGrad |
Discrete H1 space: local operators, L2 product and global interpolator. More... | |
struct | HArDCore2D::XGrad::LocalOperators |
A structure to store local operators (gradient and potential) More... | |
class | HArDCore2D::XHess |
Discrete H2 space: local operators, L2 product and global interpolator. More... | |
struct | HArDCore2D::XHess::LocalOperators |
A structure to store local operators (gradient and potential) More... | |
struct | HArDCore2D::XHess::TransferOperators |
class | HArDCore2D::XRot |
Discrete H1 space: local operators, L2 product and global interpolator. More... | |
struct | HArDCore2D::XRot::LocalOperators |
A structure to store local operators (vector rotor and potential) More... | |
class | HArDCore2D::XRotRot |
Discrete HRotRot space: local operators, L2 product and global interpolator. More... | |
struct | HArDCore2D::XRotRot::LocalOperators |
A structure to store the local operators (scalar rotor and potential) More... | |
Functions | |
HArDCore2D::DDRCore::DDRCore (const Mesh &mesh, size_t K, bool use_threads=true, std::ostream &output=std::cout) | |
Constructor. | |
const Mesh & | HArDCore2D::DDRCore::mesh () const |
Return a const reference to the mesh. | |
const size_t & | HArDCore2D::DDRCore::degree () const |
Return the polynomial degree. | |
const CellBases & | HArDCore2D::DDRCore::cellBases (size_t iT) const |
Return cell bases for element iT. | |
const EdgeBases & | HArDCore2D::DDRCore::edgeBases (size_t iE) const |
Return edge bases for edge iE. | |
HArDCore2D::EXCurl::hhoLocalOperators::hhoLocalOperators (const Eigen::MatrixXd &_gradient, const Eigen::MatrixXd &_symmetric_gradient, const Eigen::MatrixXd &_potential, const Eigen::MatrixXd &_stabilisation) | |
HArDCore2D::EXCurl::EXCurl (const DDRCore &ddr_core, bool use_threads=true, std::ostream &output=std::cout) | |
Constructor. | |
const Mesh & | HArDCore2D::EXCurl::mesh () const |
Return the mesh. | |
const size_t & | HArDCore2D::EXCurl::degree () const |
Return the polynomial degree. | |
Eigen::VectorXd | HArDCore2D::EXCurl::interpolate (const FunctionType &v, const int deg_quad=-1) const |
Interpolator of a continuous function. | |
const Eigen::MatrixXd | HArDCore2D::EXCurl::ddrPotential (size_t iT) const |
Return cell potential for the cell of index iT. | |
const Eigen::MatrixXd | HArDCore2D::EXCurl::ddrCurl (size_t iT) const |
Return cell (scalar) curl for the cell of index iT. | |
const Eigen::MatrixXd | HArDCore2D::EXCurl::ddrPotential (const Cell &T) const |
Return cell potential for cell T. | |
const Eigen::MatrixXd | HArDCore2D::EXCurl::ddrCurl (const Cell &T) const |
Return cell (scalar) curl for cell T. | |
const hhoLocalOperators & | HArDCore2D::EXCurl::hhoOperators (size_t iT) const |
Return hho operators for the cell of index iT. | |
const hhoLocalOperators & | HArDCore2D::EXCurl::hhoOperators (const Cell &T) const |
Return cell operators for cell T. | |
const DDRCore::CellBases & | HArDCore2D::EXCurl::cellBases (size_t iT) const |
Return ddrcore cell bases for the cell of index iT. | |
const DDRCore::CellBases & | HArDCore2D::EXCurl::cellBases (const Cell &T) const |
Return ddrcore cell bases for cell T. | |
const Polyk2x2Type & | HArDCore2D::EXCurl::Polyk2x2 (size_t iT) const |
Return basis for Matricial P^k space. | |
const Polykpo2Type & | HArDCore2D::EXCurl::Polykpo2 (size_t iT) const |
Return basis for the (P^{k+1})^2 space. | |
const DDRCore::EdgeBases & | HArDCore2D::EXCurl::edgeBases (size_t iE) const |
Return ddrcore edge bases for the edge of index iE. | |
const DDRCore::EdgeBases & | HArDCore2D::EXCurl::edgeBases (const Edge &E) const |
Return ddrcore edge bases for edge E. | |
Eigen::MatrixXd | HArDCore2D::EXCurl::computeL2Product (const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk2_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const |
Compute the matrix of the (weighted) L2-product for the cell of index iT. | |
Eigen::MatrixXd | HArDCore2D::EXCurl::computeL2ProductGradient (const size_t iT, const XGrad &x_grad, const std::string &side, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk2_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const |
Compute the matrix of the (weighted) L2-product as 'computeL2Product', with application of the discrete gradient on the left/right/both sides (depending on argument "side"). | |
Eigen::MatrixXd | HArDCore2D::EXCurl::computeL2Product_with_Ops (const size_t iT, const std::vector< Eigen::MatrixXd > &leftOp, const std::vector< Eigen::MatrixXd > &rightOp, const double &penalty_factor, const Eigen::MatrixXd &w_mass_Pk2_T, const IntegralWeight &weight) const |
Compute the matrix of the L2 product, applying leftOp and rightOp to the variables. Probably not directly called, mostly invoked through the wrappers computeL2Product and computeL2ProductGradient. | |
HArDCore2D::SerendipityProblem::SerendipityProblem (const DDRCore &ddrcore, bool use_threads=true, std::ostream &output=std::cout) | |
Constructor. | |
const std::vector< size_t > & | HArDCore2D::SerendipityProblem::serendipityEdges (size_t iT) const |
Return the list of serendipity edges in a cell. | |
const int | HArDCore2D::SerendipityProblem::n_serendipityEdges (size_t iT) const |
Return the number of serendipity edges in a cell. | |
const int | HArDCore2D::SerendipityProblem::serDegreeCell (size_t iT) const |
Return the serendipity degree ell_T in a cell. | |
size_t | HArDCore2D::SerendipityProblem::dimCellPolyl (size_t iT) const |
Return the dimension of P^l on cell of index iT. | |
size_t | HArDCore2D::SerendipityProblem::dimCellPolyl (const Cell &T) const |
Return the dimension of P^{l+1} on cell T. | |
const PolylBasisCellType & | HArDCore2D::SerendipityProblem::cellBasisPolyl (size_t iT) const |
Return the basis of P^l on cell of index iT. | |
const PolylBasisCellType & | HArDCore2D::SerendipityProblem::cellBasisPolyl (const Cell &T) const |
Return the basis of P^l on cell T. | |
Eigen::VectorXd | HArDCore2D::SerendipityProblem::nDOFs_cells_SXGrad () const |
Number of DOFs on cells for serendipity XGrad space. | |
size_t | HArDCore2D::SerendipityProblem::dimCellRolyCompllpo (size_t iT) const |
Return the dimension of R^{c,l+1} on cell of index iT. | |
size_t | HArDCore2D::SerendipityProblem::dimCellRolyCompllpo (const Cell &T) const |
Return the dimension of R^{c,l+1} on cell T. | |
const RolyCompllpoBasisCellType & | HArDCore2D::SerendipityProblem::cellBasisRolyCompllpo (size_t iT) const |
Return the basis of R^{c,l+1} on cell of index iT. | |
const RolyCompllpoBasisCellType & | HArDCore2D::SerendipityProblem::cellBasisRolyCompllpo (const Cell &T) const |
Return the basis of R^{c,l+1} on cell T. | |
Eigen::VectorXd | HArDCore2D::SerendipityProblem::nDOFs_cells_SXCurl () const |
Number of DOFs on cells for serendipity XCurl space. | |
Eigen::VectorXd | HArDCore2D::SerendipityProblem::nDOFs_cells_SXRotRot () const |
Number of DOFs on cells for serendipity XRotRot space. | |
const Eigen::MatrixXd | HArDCore2D::SerendipityProblem::SerendipityOperatorCell (const size_t iT, const Eigen::MatrixXd <) const |
Compute the serendipity operator on the cell of index iT. | |
const Eigen::MatrixXd | HArDCore2D::SerendipityProblem::SerendipityOperatorCell (const Cell &T, const Eigen::MatrixXd <) const |
Compute the serendipity operator on the Cell T. | |
const Mesh & | HArDCore2D::SerendipityProblem::mesh () const |
Return a const reference to the mesh. | |
const DDRCore & | HArDCore2D::SerendipityProblem::ddrCore () const |
Return a const reference to the underlying DDR core. | |
HArDCore2D::SXGrad::TransferOperators::TransferOperators (const Eigen::MatrixXd &_serendipity, const Eigen::MatrixXd &_extension, const Eigen::MatrixXd &_reduction) | |
HArDCore2D::SXGrad::SXGrad (const DDRCore &ddr_core, const SerendipityProblem &ser_pro, bool use_threads=true, std::ostream &output=std::cout) | |
Constructor. | |
const Mesh & | HArDCore2D::SXGrad::mesh () const |
Return the mesh. | |
const size_t & | HArDCore2D::SXGrad::degree () const |
Return the polynomial degree. | |
Eigen::VectorXd | HArDCore2D::SXGrad::interpolate (const FunctionType &q, const int deg_quad=-1) const |
const Eigen::MatrixXd & | HArDCore2D::SXGrad::SgradCell (size_t iT) const |
Return the serendipity reconstruction for the cell of index iT. | |
const Eigen::MatrixXd & | HArDCore2D::SXGrad::EgradCell (size_t iT) const |
Return the extension for the cell of index iT. | |
const Eigen::MatrixXd & | HArDCore2D::SXGrad::RgradCell (size_t iT) const |
Return the reduction for the cell of index iT. | |
const Eigen::MatrixXd & | HArDCore2D::SXGrad::SgradCell (const Cell &T) const |
Return the serendipity reconstruction for cell T. | |
const Eigen::MatrixXd & | HArDCore2D::SXGrad::EgradCell (const Cell &T) const |
Return the extension for cell T. | |
const Eigen::MatrixXd & | HArDCore2D::SXGrad::RgradCell (const Cell &T) const |
Return the reduction for cell T. | |
const TransferOperators & | HArDCore2D::SXGrad::cellOperators (size_t iT) const |
Return cell operators for the cell of index iT. | |
const TransferOperators & | HArDCore2D::SXGrad::cellOperators (const Cell &T) const |
Return cell operators for cell T. | |
const Eigen::MatrixXd | HArDCore2D::SXGrad::edgeGradient (size_t iE) const |
Return the full gradient operator on the edge of index iE. | |
const Eigen::MatrixXd | HArDCore2D::SXGrad::edgeGradient (const Edge &E) const |
Return the full gradient operator on edge E. | |
const Eigen::MatrixXd | HArDCore2D::SXGrad::edgePotential (size_t iE) const |
Return the potential operator on the edge of index iE. | |
const Eigen::MatrixXd | HArDCore2D::SXGrad::edgePotential (const Edge &E) const |
Return the potential operator on edge E. | |
const Eigen::MatrixXd | HArDCore2D::SXGrad::cellGradient (size_t iT) const |
Return the full gradient operator on the cell of index iT. | |
const Eigen::MatrixXd | HArDCore2D::SXGrad::cellGradient (const Cell &T) const |
Return the full gradient operator on cell T. | |
const Eigen::MatrixXd | HArDCore2D::SXGrad::cellPotential (size_t iT) const |
Return the potential operator on the cell of index iT. | |
const Eigen::MatrixXd | HArDCore2D::SXGrad::cellPotential (const Cell &T) const |
Return the potential operator on cell T. | |
Eigen::MatrixXd | HArDCore2D::SXGrad::computeL2Product (const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk2_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const |
Compute the matrix of the (weighted) L2-product. | |
double | HArDCore2D::SXGrad::computeL2Norm (const Eigen::VectorXd &v) const |
Compute the L2-norm of a vector of the space. | |
Eigen::MatrixXd | HArDCore2D::SXGrad::computeStabilisation (const size_t iT, const IntegralWeight &weight=IntegralWeight(1.)) const |
Computes only the stabilisation matrix of the (weighted) L2-product for the cell of index iT. | |
std::vector< double > | HArDCore2D::SXGrad::computeVertexValues (const Eigen::VectorXd &u) const |
Computes the values of the potential reconstruction at the mesh vertices. | |
const DDRCore::CellBases & | HArDCore2D::SXGrad::cellBases (size_t iT) const |
Return cell bases for the cell of index iT. | |
const DDRCore::CellBases & | HArDCore2D::SXGrad::cellBases (const Cell &T) const |
Return cell bases for cell T. | |
const DDRCore::EdgeBases & | HArDCore2D::SXGrad::edgeBases (size_t iE) const |
Return edge bases for the edge of index iE. | |
const DDRCore::EdgeBases & | HArDCore2D::SXGrad::edgeBases (const Edge &E) const |
Return edge bases for edge E. | |
HArDCore2D::VSXGrad::LocalOperators::LocalOperators (const Eigen::MatrixXd &_rot, const Eigen::MatrixXd &_potential) | |
HArDCore2D::VSXGrad::VSXGrad (const DDRCore &ddr_core, const SerendipityProblem &sp, bool use_threads=true, std::ostream &output=std::cout) | |
Constructor. | |
const Mesh & | HArDCore2D::VSXGrad::mesh () const |
Return the mesh. | |
const SXGrad & | HArDCore2D::VSXGrad::sXgrad () const |
Return the underlying sXgrad space. | |
const size_t & | HArDCore2D::VSXGrad::degree () const |
Return the polynomial degree. | |
Eigen::VectorXd | HArDCore2D::VSXGrad::interpolate (const FunctionType &q, const int doe_cell=-1) const |
Interpolator of a continuous function. | |
auto | HArDCore2D::VSXGrad::PolykE (size_t iE) const -> const PolyBasisEdgeType & |
auto | HArDCore2D::VSXGrad::PolykpoE (size_t iE) const -> const PolyBasisEdgeType & |
auto | HArDCore2D::VSXGrad::PolykptwoE (size_t iE) const -> const PolyBasisEdgeType & |
auto | HArDCore2D::VSXGrad::Polykpo (size_t iT) const -> const PolyBasisCellType & |
auto | HArDCore2D::VSXGrad::Polykmo (size_t iT) const -> const PolyBasisCellType & |
auto | HArDCore2D::VSXGrad::Polyk (size_t iT) const -> const PolyBasisCellType & |
auto | HArDCore2D::VSXGrad::Polykptwo (size_t iT) const -> const PolyBasisCellType & |
auto | HArDCore2D::VSXGrad::Polyk2 (size_t iT) const -> const Poly2BasisCellType & |
const Eigen::MatrixXd | HArDCore2D::VSXGrad::edgeGradient (size_t iE) const |
Return the gradient operator (expressed on Polyk3) on the edge of index iE. | |
const Eigen::MatrixXd | HArDCore2D::VSXGrad::edgeGradient (const Edge &E) const |
Return the gradient operator (expressed on Polyk3) on edge E. | |
const Eigen::MatrixXd | HArDCore2D::VSXGrad::edgePotential (size_t iE) const |
Return the potential operator (expressed on Polykpo3) on the edge of index iE. | |
const Eigen::MatrixXd | HArDCore2D::VSXGrad::edgePotential (const Edge &E) const |
Return the potential operator (expressed on Polykpo3) on edge E. | |
const Eigen::MatrixXd | HArDCore2D::VSXGrad::cellGradient (size_t iT) const |
Return the gradient operator (expressed on Polyk3x3) on the cell of index iT. | |
const Eigen::MatrixXd | HArDCore2D::VSXGrad::cellGradient (const Cell &T) const |
Return the gradient operator (expressed on Polyk3x3) on cell T. | |
const Eigen::MatrixXd | HArDCore2D::VSXGrad::cellSymGradient (size_t iT) const |
Return the symmetric gradient operator (expressed on Polyk3x3) on the cell of index iT. | |
const Eigen::MatrixXd | HArDCore2D::VSXGrad::cellSymGradient (const Cell &T) const |
Return the symmetric gradient operator (expressed on Polyk3x3) on cell T. | |
const Eigen::MatrixXd | HArDCore2D::VSXGrad::cellDivergence (size_t iT) const |
Return the divergence operator (expressed on the ancester of Polyk3x3, which should be m_sxgrad.cellBases(iT).Polyk) on the cell of index iT. | |
const Eigen::MatrixXd | HArDCore2D::VSXGrad::cellDivergence (const Cell &T) const |
Return the divergence operator (expressed on the ancester of Polyk3x3, which should be m_sxgrad.cellBases(iT).Polyk) on the cell of index iT. | |
const Eigen::MatrixXd | HArDCore2D::VSXGrad::cellPotential (size_t iT) const |
Return the potential operator (expressed on Polykpo3) on the cell of index iT. | |
const Eigen::MatrixXd | HArDCore2D::VSXGrad::cellPotential (const Cell &T) const |
Return the potential operator on cell T. | |
Eigen::MatrixXd | HArDCore2D::VSXGrad::_permute_bases (size_t iT) |
Eigen::MatrixXd | HArDCore2D::VSXGrad::computeStabilisation (const size_t iT, const IntegralWeight &weight=IntegralWeight(1.)) const |
Compute the H^1-stabilisation (based on the L^2 stabilisation implemented in SXgrad) | |
std::vector< double > | HArDCore2D::VSXGrad::computeH1norms (const std::vector< Eigen::VectorXd > &list_dofs, const std::string typegrad="full") const |
Compute the H^1 discrete norm, with full or symmetric gradient. | |
std::vector< VectorRd > | HArDCore2D::VSXGrad::computeVertexValues (const Eigen::VectorXd &u) const |
Computes the values of the potential reconstruction at the mesh vertices. | |
const CellBases & | HArDCore2D::VSXGrad::cellBases (size_t iT) const |
Return vector/matrix cell bases for the face of index iT. | |
const CellBases & | HArDCore2D::VSXGrad::cellBases (const Cell &T) const |
Return vector/matrix cell bases for cell T. | |
const EdgeBases & | HArDCore2D::VSXGrad::edgeBases (size_t iE) const |
Return edge bases for the edge of index iE. | |
const EdgeBases & | HArDCore2D::VSXGrad::edgeBases (const Edge &E) const |
Return edge bases for edge E. | |
const LocalOperators & | HArDCore2D::VSXGrad::cellOperators (size_t iT) const |
Return cell operators for the cell of index iT. | |
const LocalOperators & | HArDCore2D::VSXGrad::cellOperators (const Cell &T) const |
Return cell operators for cell T. | |
HArDCore2D::XCurl::LocalOperators::LocalOperators (const Eigen::MatrixXd &_curl, const Eigen::MatrixXd &_potential) | |
HArDCore2D::XCurl::XCurl (const DDRCore &ddr_core, bool use_threads=true, std::ostream &output=std::cout) | |
Constructor. | |
const Mesh & | HArDCore2D::XCurl::mesh () const |
Return the mesh. | |
const size_t & | HArDCore2D::XCurl::degree () const |
Return the polynomial degree. | |
Eigen::VectorXd | HArDCore2D::XCurl::interpolate (const FunctionType &v, const int deg_quad=-1) const |
Interpolator of a continuous function. | |
const LocalOperators & | HArDCore2D::XCurl::cellOperators (size_t iT) const |
Return cell operators for the cell of index iT. | |
const LocalOperators & | HArDCore2D::XCurl::cellOperators (const Cell &T) const |
Return cell operators for cell T. | |
const DDRCore::CellBases & | HArDCore2D::XCurl::cellBases (size_t iT) const |
Return cell bases for the cell of index iT. | |
const DDRCore::CellBases & | HArDCore2D::XCurl::cellBases (const Cell &T) const |
Return cell bases for cell T. | |
const DDRCore::EdgeBases & | HArDCore2D::XCurl::edgeBases (size_t iE) const |
Return edge bases for the edge of index iE. | |
const DDRCore::EdgeBases & | HArDCore2D::XCurl::edgeBases (const Edge &E) const |
Return edge bases for edge E. | |
Eigen::MatrixXd | HArDCore2D::XCurl::computeL2Product (const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk2_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const |
Compute the matrix of the (weighted) L2-product for the cell of index iT. | |
Eigen::MatrixXd | HArDCore2D::XCurl::computeL2ProductGradient (const size_t iT, const XGrad &x_grad, const std::string &side, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk2_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const |
Compute the matrix of the (weighted) L2-product as 'computeL2Product', with application of the discrete gradient on the left/right/both sides (depending on argument "side"). | |
Eigen::MatrixXd | HArDCore2D::XCurl::computeL2Product_with_Ops (const size_t iT, const std::vector< Eigen::MatrixXd > &leftOp, const std::vector< Eigen::MatrixXd > &rightOp, const double &penalty_factor, const Eigen::MatrixXd &w_mass_Pk2_T, const IntegralWeight &weight) const |
Compute the matrix of the L2 product, applying leftOp and rightOp to the variables. Probably not directly called, mostly invoked through the wrappers computeL2Product and computeL2ProductGradient. | |
HArDCore2D::XGrad::LocalOperators::LocalOperators (const Eigen::MatrixXd &_gradient, const Eigen::MatrixXd &_potential) | |
HArDCore2D::XGrad::XGrad (const DDRCore &ddr_core, bool use_threads=true, std::ostream &output=std::cout) | |
Constructor. | |
const Mesh & | HArDCore2D::XGrad::mesh () const |
Return the mesh. | |
const size_t & | HArDCore2D::XGrad::degree () const |
Return the polynomial degree. | |
Eigen::VectorXd | HArDCore2D::XGrad::interpolate (const FunctionType &q, const int deg_quad=-1) const |
Interpolator of a continuous function. | |
const LocalOperators & | HArDCore2D::XGrad::edgeOperators (size_t iE) const |
Return edge operators for the edge of index iE. | |
const LocalOperators & | HArDCore2D::XGrad::edgeOperators (const Edge &E) const |
Return edge operators for edge E. | |
const LocalOperators & | HArDCore2D::XGrad::cellOperators (size_t iT) const |
Return cell operators for the cell of index iT. | |
const LocalOperators & | HArDCore2D::XGrad::cellOperators (const Cell &T) const |
Return cell operators for cell T. | |
const DDRCore::CellBases & | HArDCore2D::XGrad::cellBases (size_t iT) const |
Return cell bases for the cell of index iT. | |
const DDRCore::CellBases & | HArDCore2D::XGrad::cellBases (const Cell &T) const |
Return cell bases for cell T. | |
const DDRCore::EdgeBases & | HArDCore2D::XGrad::edgeBases (size_t iE) const |
Return edge bases for the edge of index iE. | |
const DDRCore::EdgeBases & | HArDCore2D::XGrad::edgeBases (const Edge &E) const |
Return edge bases for edge E. | |
Eigen::MatrixXd | HArDCore2D::XGrad::computeStabilisation (const size_t iT, const IntegralWeight &weight=IntegralWeight(1.)) const |
Computes only the stabilisation matrix of the (weighted) L2-product for the cell of index iT. This stabilisation is based on cell and face potentials. | |
Eigen::MatrixXd | HArDCore2D::XGrad::computeL2Product (const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pkpo_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const |
Compute the matrix of the (weighted) L2-product for the cell of index iT. | |
double | HArDCore2D::XGrad::evaluatePotential (const size_t iT, const Eigen::VectorXd &vT, const VectorRd &x) const |
Evaluate the value of the potential at a point x. | |
double | HArDCore2D::XGrad::computeL2Norm (const Eigen::VectorXd &v) const |
Compute the L2-norm of a vector of the space. | |
HArDCore2D::XHess::LocalOperators::LocalOperators (const Eigen::MatrixXd &_gradient, const Eigen::MatrixXd &_potential) | |
HArDCore2D::XHess::TransferOperators::TransferOperators (const Eigen::MatrixXd &_serendipity) | |
HArDCore2D::XHess::XHess (const DDRCore &ddr_core, const SerendipityProblem &ser_pro, bool use_threads=true, std::ostream &output=std::cout) | |
Constructor. | |
const Mesh & | HArDCore2D::XHess::mesh () const |
Return the mesh. | |
const size_t & | HArDCore2D::XHess::degree () const |
Return the polynomial degree. | |
Eigen::VectorXd | HArDCore2D::XHess::interpolate (const FunctionType &q, const GradFunctionType &Dq, const int deg_quad=-1) const |
Interpolator of a continuous function. | |
const Eigen::MatrixXd & | HArDCore2D::XHess::SgradCell (size_t iT) const |
Return the serendipity reconstruction for the cell of index iT. | |
const Eigen::MatrixXd & | HArDCore2D::XHess::SgradCell (const Cell &T) const |
Return the serendipity reconstruction for cell T. | |
const TransferOperators & | HArDCore2D::XHess::TcellOperators (size_t iT) const |
Return cell operators for the cell of index iT. | |
const TransferOperators & | HArDCore2D::XHess::TcellOperators (const Cell &T) const |
Return cell operators for cell T. | |
const LocalOperators & | HArDCore2D::XHess::edgeOperators (size_t iE) const |
Return edge operators for the edge of index iE. | |
const LocalOperators & | HArDCore2D::XHess::edgeOperators (const Edge &E) const |
Return edge operators for edge E. | |
const LocalOperators & | HArDCore2D::XHess::cellOperators (size_t iT) const |
Return cell operators for the cell of index iT. | |
const LocalOperators & | HArDCore2D::XHess::cellOperators (const Cell &T) const |
Return cell operators for cell T. | |
const DDRCore::CellBases & | HArDCore2D::XHess::cellBases (size_t iT) const |
Return cell bases for the cell of index iT. | |
const DDRCore::CellBases & | HArDCore2D::XHess::cellBases (const Cell &T) const |
Return cell bases for cell T. | |
const DDRCore::EdgeBases & | HArDCore2D::XHess::edgeBases (size_t iE) const |
Return edge bases for the edge of index iE. | |
const DDRCore::EdgeBases & | HArDCore2D::XHess::edgeBases (const Edge &E) const |
Return edge bases for edge E. | |
Eigen::MatrixXd | HArDCore2D::XHess::computeL2Product (const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pkpo_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const |
Compute the matrix of the (weighted) L2-product for the cell of index iT. | |
double | HArDCore2D::XHess::evaluatePotential (const size_t iT, const Eigen::VectorXd &vT, const VectorRd &x) const |
Evaluate the value of the potential at a point x. | |
double | HArDCore2D::XHess::computeL2Norm (const Eigen::VectorXd &v) const |
Compute the L2-norm of a vector of the space. | |
HArDCore2D::XRot::LocalOperators::LocalOperators (const Eigen::MatrixXd &_rotor, const Eigen::MatrixXd &_rotor_rhs, const Eigen::MatrixXd &_potential) | |
HArDCore2D::XRot::XRot (const DDRCore &ddr_core, bool use_threads=true, std::ostream &output=std::cout) | |
Constructor. | |
const Mesh & | HArDCore2D::XRot::mesh () const |
Return the mesh. | |
const bool | HArDCore2D::XRot::useThreads () const |
Return true if we use thread-based parallelism. | |
const size_t & | HArDCore2D::XRot::degree () const |
Return the polynomial degree. | |
Eigen::VectorXd | HArDCore2D::XRot::interpolate (const FunctionType &q, const int deg_quad=-1) const |
Interpolator of a continuous function. | |
const Eigen::MatrixXd & | HArDCore2D::XRot::edgePotential (size_t iE) const |
Return edge potential for the edge of index iE. | |
const Eigen::MatrixXd & | HArDCore2D::XRot::edgePotential (const Edge &E) const |
Return edge potential for the edge E. | |
const LocalOperators & | HArDCore2D::XRot::cellOperators (size_t iT) const |
Return cell operators for the cell of index iT. | |
const LocalOperators & | HArDCore2D::XRot::cellOperators (const Cell &T) const |
Return cell operators for cell T. | |
const DDRCore::CellBases & | HArDCore2D::XRot::cellBases (size_t iT) const |
Return cell bases for the cell of index iT. | |
const DDRCore::CellBases & | HArDCore2D::XRot::cellBases (const Cell &T) const |
Return cell bases for cell T. | |
const DDRCore::EdgeBases & | HArDCore2D::XRot::edgeBases (size_t iE) const |
Return edge bases for the edge of index iE. | |
const DDRCore::EdgeBases & | HArDCore2D::XRot::edgeBases (const Edge &E) const |
Return edge bases for edge E. | |
Eigen::MatrixXd | HArDCore2D::XRot::computeRotorL2Product (const size_t iT, const double &penalty_factor_cell=1., const double &penalty_factor_edge=1.) const |
Compute rotor L2-product. | |
double | HArDCore2D::XRot::computeRotorL2Norm (const Eigen::VectorXd &v) const |
Compute rotor L2-norm. | |
HArDCore2D::XRotRot::LocalOperators::LocalOperators (const Eigen::MatrixXd &_rotor, const Eigen::MatrixXd &_potential) | |
HArDCore2D::XRotRot::XRotRot (const DDRCore &ddr_core, bool use_threads=true, std::ostream &output=std::cout) | |
Constructor. | |
const Mesh & | HArDCore2D::XRotRot::mesh () const |
Return the mesh. | |
const size_t & | HArDCore2D::XRotRot::degree () const |
Return the polynomial degree. | |
Eigen::VectorXd | HArDCore2D::XRotRot::interpolate (const FunctionType &v, const RotType &rot_v, const int deg_quad=-1) const |
Interpolator of a continuous function. | |
const LocalOperators & | HArDCore2D::XRotRot::cellOperators (size_t iT) const |
Return cell operators for the cell of index iT. | |
const LocalOperators & | HArDCore2D::XRotRot::cellOperators (const Cell &T) const |
Return cell operators for cell T. | |
Eigen::MatrixXd | HArDCore2D::XRotRot::cellRotor (size_t iT) const |
Return cell rotor for cell of index iT. | |
Eigen::MatrixXd | HArDCore2D::XRotRot::cellRotor (const Cell &T) const |
Return cell rotor for cell T. | |
const DDRCore::CellBases & | HArDCore2D::XRotRot::cellBases (size_t iT) const |
Return cell bases for the cell of index iT. | |
const DDRCore::CellBases & | HArDCore2D::XRotRot::cellBases (const Cell &T) const |
Return cell bases for cell T. | |
const DDRCore::EdgeBases & | HArDCore2D::XRotRot::edgeBases (size_t iE) const |
Return edge bases for the edge of index iE. | |
const DDRCore::EdgeBases & | HArDCore2D::XRotRot::edgeBases (const Edge &E) const |
Return edge bases for edge E. | |
Eigen::MatrixXd | HArDCore2D::XRotRot::computeL2Product (size_t iT, const double &penalty_factor=1., const IntegralWeight &weight=IntegralWeight(1.)) const |
Compute the matrix of the (weighted) L2-product for the cell of index iT. | |
Eigen::MatrixXd | HArDCore2D::XRotRot::computeGradientPotentialL2Product (size_t iT, const XGrad *x_grad, const double &penalty_factor=1., const IntegralWeight &weight=IntegralWeight(1.)) const |
Compute the matrix of the (weighted) L2-product as 'computeL2Product', with application of the discrete gradient on the left side. | |
Eigen::MatrixXd | HArDCore2D::XRotRot::computeGradientL2Product (size_t iT, const XGrad *x_grad, const double &penalty_factor=1., const IntegralWeight &weight=IntegralWeight(1.)) const |
Compute the matrix of the (weighted) gradient-gradient L2-product. | |
template<typename LeftOperatorFillerType , typename RightOperatorFillerType > | |
Eigen::MatrixXd | HArDCore2D::XRotRot::computeL2ProductWithOperatorFillers (size_t iT, const double &penalty_factor, const IntegralWeight &weight, LeftOperatorFillerType fillLeftOp, RightOperatorFillerType fillRightOp) const |
Compute the matrix of the L2 product given operators filling functions. | |
double | HArDCore2D::XRotRot::computeL2Norm (const Eigen::VectorXd &v) const |
Compute the L2-norm of a vector of the space. | |
double | HArDCore2D::XRotRot::computeGradientL2Norm (const Eigen::VectorXd &v, const XGrad *x_grad) const |
Compute the L2-norm of the discrete gradient of a vector in XGrad. | |
Classes providing tools to the Discrete De Rham sequence.
typedef std::function<Eigen::Vector2d(const Eigen::Vector2d &)> HArDCore2D::EXCurl::FunctionType |
typedef std::function<double(const Eigen::Vector2d &)> HArDCore2D::SXGrad::FunctionType |
typedef std::function<Eigen::Vector2d(const Eigen::Vector2d &)> HArDCore2D::VSXGrad::FunctionType |
typedef std::function<Eigen::Vector2d(const Eigen::Vector2d &)> HArDCore2D::XCurl::FunctionType |
typedef std::function<double(const Eigen::Vector2d &)> HArDCore2D::XGrad::FunctionType |
typedef std::function<double(const Eigen::Vector2d &)> HArDCore2D::XHess::FunctionType |
typedef std::function<double(const Eigen::Vector2d &)> HArDCore2D::XRot::FunctionType |
typedef std::function<Eigen::Vector2d(const Eigen::Vector2d &)> HArDCore2D::XRotRot::FunctionType |
typedef Cell HArDCore2D::DDRCore::CellBases::GeometricSupport |
Geometric support.
typedef Edge HArDCore2D::DDRCore::EdgeBases::GeometricSupport |
Geometric support.
typedef Family<GradientBasis<ShiftedBasis<MonomialScalarBasisCell> > > HArDCore2D::DDRCore::GolyBasisCellType |
typedef std::function<Eigen::Vector2d(const Eigen::Vector2d &)> HArDCore2D::XHess::GradFunctionType |
typedef Eigen::FullPivLU<Eigen::MatrixXd> HArDCore2D::SerendipityProblem::InverseProblem |
Type for inverses of matrix for serendipity problem.
typedef TensorizedVectorFamily<DDRCore::PolyBasisCellType, dimspace> HArDCore2D::VSXGrad::Poly2BasisCellType |
typedef TensorizedVectorFamily<DDRCore::PolyBasisEdgeType, dimspace> HArDCore2D::VSXGrad::Poly2BasisEdgeType |
typedef MatrixFamily<DDRCore::PolyBasisCellType, dimspace> HArDCore2D::VSXGrad::Poly2x2BasisCellType |
typedef TensorizedVectorFamily<DDRCore::PolyBasisCellType, dimspace> HArDCore2D::EXCurl::Polykpo2Type |
typedef RestrictedBasis<DDRCore::PolyBasisCellType > HArDCore2D::SerendipityProblem::PolylBasisCellType |
typedef Family<CurlBasis<ShiftedBasis<MonomialScalarBasisCell> > > HArDCore2D::DDRCore::RolyBasisCellType |
typedef RestrictedBasis<DDRCore::RolyComplBasisCellType > HArDCore2D::SerendipityProblem::RolyCompllpoBasisCellType |
typedef std::function<double(const Eigen::Vector2d &)> HArDCore2D::XRotRot::RotType |
Eigen::MatrixXd VSXGrad::_permute_bases | ( | size_t | iT | ) |
|
inline |
Return ddrcore cell bases for cell T.
|
inline |
Return cell bases for cell T.
|
inline |
Return vector/matrix cell bases for cell T.
|
inline |
Return cell bases for cell T.
|
inline |
Return cell bases for cell T.
|
inline |
Return cell bases for cell T.
|
inline |
Return cell bases for cell T.
|
inline |
Return cell bases for cell T.
|
inline |
Return cell bases for element iT.
|
inline |
Return ddrcore cell bases for the cell of index iT.
|
inline |
Return cell bases for the cell of index iT.
|
inline |
Return vector/matrix cell bases for the face of index iT.
|
inline |
Return cell bases for the cell of index iT.
|
inline |
Return cell bases for the cell of index iT.
|
inline |
Return cell bases for the cell of index iT.
|
inline |
Return cell bases for the cell of index iT.
|
inline |
Return cell bases for the cell of index iT.
|
inline |
Return the basis of P^l on cell T.
|
inline |
Return the basis of P^l on cell of index iT.
|
inline |
Return the basis of R^{c,l+1} on cell T.
|
inline |
Return the basis of R^{c,l+1} on cell of index iT.
|
inline |
Return the divergence operator (expressed on the ancester of Polyk3x3, which should be m_sxgrad.cellBases(iT).Polyk) on the cell of index iT.
|
inline |
Return the divergence operator (expressed on the ancester of Polyk3x3, which should be m_sxgrad.cellBases(iT).Polyk) on the cell of index iT.
|
inline |
Return the full gradient operator on cell T.
|
inline |
Return the gradient operator (expressed on Polyk3x3) on cell T.
|
inline |
Return the full gradient operator on the cell of index iT.
|
inline |
Return the gradient operator (expressed on Polyk3x3) on the cell of index iT.
|
inline |
Return cell operators for cell T.
|
inline |
Return cell operators for cell T.
|
inline |
Return cell operators for cell T.
|
inline |
Return cell operators for cell T.
|
inline |
Return cell operators for cell T.
|
inline |
Return cell operators for cell T.
|
inline |
Return cell operators for cell T.
|
inline |
Return cell operators for the cell of index iT.
|
inline |
Return cell operators for the cell of index iT.
|
inline |
Return cell operators for the cell of index iT.
|
inline |
Return cell operators for the cell of index iT.
|
inline |
Return cell operators for the cell of index iT.
|
inline |
Return cell operators for the cell of index iT.
|
inline |
Return cell operators for the cell of index iT.
|
inline |
Return the potential operator on cell T.
|
inline |
Return the potential operator on cell T.
|
inline |
Return the potential operator on the cell of index iT.
|
inline |
Return the potential operator (expressed on Polykpo3) on the cell of index iT.
|
inline |
Return cell rotor for cell T.
|
inline |
Return cell rotor for cell of index iT.
|
inline |
Return the symmetric gradient operator (expressed on Polyk3x3) on cell T.
|
inline |
Return the symmetric gradient operator (expressed on Polyk3x3) on the cell of index iT.
double XRotRot::computeGradientL2Norm | ( | const Eigen::VectorXd & | v, |
const XGrad * | x_grad | ||
) | const |
Compute the L2-norm of the discrete gradient of a vector in XGrad.
Eigen::MatrixXd XRotRot::computeGradientL2Product | ( | size_t | iT, |
const XGrad * | x_grad, | ||
const double & | penalty_factor = 1. , |
||
const IntegralWeight & | weight = IntegralWeight(1.) |
||
) | const |
Compute the matrix of the (weighted) gradient-gradient L2-product.
iT | Index of the cell |
x_grad | Instance of XGrad to access the full gradients |
penalty_factor | Pre-factor for stabilisation term |
weight | Weight function in the L2 product, defaults to constant 1. |
Eigen::MatrixXd XRotRot::computeGradientPotentialL2Product | ( | size_t | iT, |
const XGrad * | x_grad, | ||
const double & | penalty_factor = 1. , |
||
const IntegralWeight & | weight = IntegralWeight(1.) |
||
) | const |
Compute the matrix of the (weighted) L2-product as 'computeL2Product', with application of the discrete gradient on the left side.
iT | Index of the cell |
x_grad | Instance of XGrad to access the full gradients |
penalty_factor | Pre-factor for stabilisation term |
weight | Weight function in the L2 product, defaults to constant 1. |
std::vector< double > VSXGrad::computeH1norms | ( | const std::vector< Eigen::VectorXd > & | list_dofs, |
const std::string | typegrad = "full" |
||
) | const |
Compute the H^1 discrete norm, with full or symmetric gradient.
list_dofs | The list of vectors representing the dofs |
typegrad | computes the norm using the "full" or "symmetric" gradient |
double SXGrad::computeL2Norm | ( | const Eigen::VectorXd & | v | ) | const |
Compute the L2-norm of a vector of the space.
double XGrad::computeL2Norm | ( | const Eigen::VectorXd & | v | ) | const |
Compute the L2-norm of a vector of the space.
double XHess::computeL2Norm | ( | const Eigen::VectorXd & | v | ) | const |
Compute the L2-norm of a vector of the space.
double XRotRot::computeL2Norm | ( | const Eigen::VectorXd & | v | ) | const |
Compute the L2-norm of a vector of the space.
Eigen::MatrixXd EXCurl::computeL2Product | ( | const size_t | iT, |
const double & | penalty_factor = 1. , |
||
const Eigen::MatrixXd & | mass_Pk2_T = Eigen::MatrixXd::Zero(1,1) , |
||
const IntegralWeight & | weight = IntegralWeight(1.) |
||
) | const |
Compute the matrix of the (weighted) L2-product for the cell of index iT.
iT | index of the cell |
penalty_factor | pre-factor for stabilisation term |
mass_Pk2_T | if pre-computed, the mass matrix of (P^k(T))^3; if none is pre-computed, passing Eigen::MatrixXd::Zero(1,1) will force the calculation |
weight | weight function in the L2 product, defaults to 1 |
|
inline |
Compute the matrix of the (weighted) L2-product.
iT | index of the cell |
penalty_factor | pre-factor for stabilisation term |
mass_Pk2_T | if pre-computed, the mass matrix of (P^k(T))^3; if none is pre-computed, passing Eigen::MatrixXd::Zero(1,1) will force the calculation |
weight | weight function in the L2 product, defaults to 1 |
Eigen::MatrixXd XCurl::computeL2Product | ( | const size_t | iT, |
const double & | penalty_factor = 1. , |
||
const Eigen::MatrixXd & | mass_Pk2_T = Eigen::MatrixXd::Zero(1,1) , |
||
const IntegralWeight & | weight = IntegralWeight(1.) |
||
) | const |
Compute the matrix of the (weighted) L2-product for the cell of index iT.
iT | index of the cell |
penalty_factor | pre-factor for stabilisation term |
mass_Pk2_T | if pre-computed, the mass matrix of (P^k(T))^3; if none is pre-computed, passing Eigen::MatrixXd::Zero(1,1) will force the calculation |
weight | weight function in the L2 product, defaults to 1 |
Eigen::MatrixXd XGrad::computeL2Product | ( | const size_t | iT, |
const double & | penalty_factor = 1. , |
||
const Eigen::MatrixXd & | mass_Pkpo_T = Eigen::MatrixXd::Zero(1,1) , |
||
const IntegralWeight & | weight = IntegralWeight(1.) |
||
) | const |
Compute the matrix of the (weighted) L2-product for the cell of index iT.
iT | index of the cell |
penalty_factor | pre-factor for stabilisation term |
mass_Pkpo_T | if pre-computed, the mass matrix of P^{k+1}(T); if none is pre-computed, passing Eigen::MatrixXd::Zero(1,1) will force the calculation |
weight | weight function in the L2 product, defaults to 1 |
Eigen::MatrixXd XHess::computeL2Product | ( | const size_t | iT, |
const double & | penalty_factor = 1. , |
||
const Eigen::MatrixXd & | mass_Pkpo_T = Eigen::MatrixXd::Zero(1,1) , |
||
const IntegralWeight & | weight = IntegralWeight(1.) |
||
) | const |
Compute the matrix of the (weighted) L2-product for the cell of index iT.
iT | index of the cell |
penalty_factor | pre-factor for stabilisation term |
mass_Pkpo_T | if pre-computed, the mass matrix of P^{k+1}(T); if none is pre-computed, passing Eigen::MatrixXd::Zero(1,1) will force the calculation |
weight | weight function in the L2 product, defaults to 1 |
Eigen::MatrixXd XRotRot::computeL2Product | ( | size_t | iT, |
const double & | penalty_factor = 1. , |
||
const IntegralWeight & | weight = IntegralWeight(1.) |
||
) | const |
Compute the matrix of the (weighted) L2-product for the cell of index iT.
iT | Index of the cell |
penalty_factor | Pre-factor for stabilisation term |
weight | Weight function in the L2 product, defaults to 1 |
Eigen::MatrixXd HArDCore2D::EXCurl::computeL2Product_with_Ops | ( | const size_t | iT, |
const std::vector< Eigen::MatrixXd > & | leftOp, | ||
const std::vector< Eigen::MatrixXd > & | rightOp, | ||
const double & | penalty_factor, | ||
const Eigen::MatrixXd & | w_mass_Pk2_T, | ||
const IntegralWeight & | weight | ||
) | const |
Compute the matrix of the L2 product, applying leftOp and rightOp to the variables. Probably not directly called, mostly invoked through the wrappers computeL2Product and computeL2ProductGradient.
iT | index of the cell |
leftOp | edge and element operators to apply on the left |
rightOp | edge and element operators to apply on the right |
penalty_factor | pre-factor for stabilisation term |
w_mass_Pk2_T | mass matrix of (P^k(T))^3 weighted by weight |
weight | weight function in the L2 product |
Eigen::MatrixXd XCurl::computeL2Product_with_Ops | ( | const size_t | iT, |
const std::vector< Eigen::MatrixXd > & | leftOp, | ||
const std::vector< Eigen::MatrixXd > & | rightOp, | ||
const double & | penalty_factor, | ||
const Eigen::MatrixXd & | w_mass_Pk2_T, | ||
const IntegralWeight & | weight | ||
) | const |
Compute the matrix of the L2 product, applying leftOp and rightOp to the variables. Probably not directly called, mostly invoked through the wrappers computeL2Product and computeL2ProductGradient.
iT | index of the cell |
leftOp | edge and element operators to apply on the left |
rightOp | edge and element operators to apply on the right |
penalty_factor | pre-factor for stabilisation term |
w_mass_Pk2_T | mass matrix of (P^k(T))^3 weighted by weight |
weight | weight function in the L2 product |
Eigen::MatrixXd EXCurl::computeL2ProductGradient | ( | const size_t | iT, |
const XGrad & | x_grad, | ||
const std::string & | side, | ||
const double & | penalty_factor = 1. , |
||
const Eigen::MatrixXd & | mass_Pk2_T = Eigen::MatrixXd::Zero(1,1) , |
||
const IntegralWeight & | weight = IntegralWeight(1.) |
||
) | const |
Compute the matrix of the (weighted) L2-product as 'computeL2Product', with application of the discrete gradient on the left/right/both sides (depending on argument "side").
iT | index of the cell |
x_grad | instance of XGrad to access the full gradients |
side | which side (left, right, both) we apply the gradient to |
penalty_factor | pre-factor for stabilisation term |
mass_Pk2_T | if pre-computed, the mass matrix of (P^k(T))^3; if none is pre-computed, passing Eigen::MatrixXd::Zero(1,1) will force the calculation |
weight | weight function in the L2 product, defaults to constant 1. |
Eigen::MatrixXd XCurl::computeL2ProductGradient | ( | const size_t | iT, |
const XGrad & | x_grad, | ||
const std::string & | side, | ||
const double & | penalty_factor = 1. , |
||
const Eigen::MatrixXd & | mass_Pk2_T = Eigen::MatrixXd::Zero(1,1) , |
||
const IntegralWeight & | weight = IntegralWeight(1.) |
||
) | const |
Compute the matrix of the (weighted) L2-product as 'computeL2Product', with application of the discrete gradient on the left/right/both sides (depending on argument "side").
iT | index of the cell |
x_grad | instance of XGrad to access the full gradients |
side | which side (left, right, both) we apply the gradient to |
penalty_factor | pre-factor for stabilisation term |
mass_Pk2_T | if pre-computed, the mass matrix of (P^k(T))^3; if none is pre-computed, passing Eigen::MatrixXd::Zero(1,1) will force the calculation |
weight | weight function in the L2 product, defaults to constant 1. |
Eigen::MatrixXd HArDCore2D::XRotRot::computeL2ProductWithOperatorFillers | ( | size_t | iT, |
const double & | penalty_factor, | ||
const IntegralWeight & | weight, | ||
LeftOperatorFillerType | fillLeftOp, | ||
RightOperatorFillerType | fillRightOp | ||
) | const |
Compute the matrix of the L2 product given operators filling functions.
iT | Index of the cell |
penalty_factor | Penalty factor for stabilisation term |
weight | Weight function in the L2 product |
fillLeftOp | Function to fill the left-hand side operator |
fillRightOp | Function to fill the right-hand side operator |
double XRot::computeRotorL2Norm | ( | const Eigen::VectorXd & | v | ) | const |
Compute rotor L2-norm.
Eigen::MatrixXd XRot::computeRotorL2Product | ( | const size_t | iT, |
const double & | penalty_factor_cell = 1. , |
||
const double & | penalty_factor_edge = 1. |
||
) | const |
Compute rotor L2-product.
iT | Index of the cell |
penalty_factor_cell | Coefficient for the cell-based term |
penalty_factor_edge | Coefficient for the edge-based term |
|
inline |
Computes only the stabilisation matrix of the (weighted) L2-product for the cell of index iT.
iT | index of the cell |
weight | weight function in the stabilisation, defaults to 1 |
|
inline |
Compute the H^1-stabilisation (based on the L^2 stabilisation implemented in SXgrad)
iT | index of the cell |
weight | weight function in the L2 product, defaults to 1 |
Eigen::MatrixXd XGrad::computeStabilisation | ( | const size_t | iT, |
const IntegralWeight & | weight = IntegralWeight(1.) |
||
) | const |
Computes only the stabilisation matrix of the (weighted) L2-product for the cell of index iT. This stabilisation is based on cell and face potentials.
iT | index of the cell |
weight | weight function in the stabilisation, defaults to 1 |
std::vector< double > HArDCore2D::SXGrad::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< VectorRd > VSXGrad::computeVertexValues | ( | const Eigen::VectorXd & | u | ) | const |
Computes the values of the potential reconstruction at the mesh vertices.
u | DOFs in the discrete space |
|
inline |
Return a const reference to the underlying DDR core.
DDRCore::DDRCore | ( | const Mesh & | mesh, |
size_t | K, | ||
bool | use_threads = true , |
||
std::ostream & | output = std::cout |
||
) |
Constructor.
|
inline |
Return cell (scalar) curl for cell T.
|
inline |
Return cell (scalar) curl for the cell of index iT.
|
inline |
Return cell potential for cell T.
|
inline |
Return cell potential for the cell of index iT.
|
inline |
Return the polynomial degree.
|
inline |
Return the polynomial degree.
|
inline |
Return the polynomial degree.
|
inline |
Return the polynomial degree.
|
inline |
Return the polynomial degree.
|
inline |
Return the polynomial degree.
|
inline |
Return the polynomial degree.
|
inline |
Return the polynomial degree.
|
inline |
Return the polynomial degree.
|
inline |
Return the dimension of P^{l+1} on cell T.
|
inline |
Return the dimension of P^l on cell of index iT.
|
inline |
Return the dimension of R^{c,l+1} on cell T.
|
inline |
Return the dimension of R^{c,l+1} on cell of index iT.
|
inline |
Return ddrcore edge bases for edge E.
|
inline |
Return edge bases for edge E.
|
inline |
Return edge bases for edge E.
|
inline |
Return edge bases for edge E.
|
inline |
Return edge bases for edge E.
|
inline |
Return edge bases for edge E.
|
inline |
Return edge bases for edge E.
|
inline |
Return edge bases for edge E.
|
inline |
Return edge bases for edge iE.
|
inline |
Return ddrcore edge bases for the edge of index iE.
|
inline |
Return edge bases for the edge of index iE.
|
inline |
Return edge bases for the edge of index iE.
|
inline |
Return edge bases for the edge of index iE.
|
inline |
Return edge bases for the edge of index iE.
|
inline |
Return edge bases for the edge of index iE.
|
inline |
Return edge bases for the edge of index iE.
|
inline |
Return edge bases for the edge of index iE.
|
inline |
Return the full gradient operator on edge E.
|
inline |
Return the gradient operator (expressed on Polyk3) on edge E.
|
inline |
Return the full gradient operator on the edge of index iE.
|
inline |
Return the gradient operator (expressed on Polyk3) on the edge of index iE.
|
inline |
Return edge operators for edge E.
|
inline |
Return edge operators for edge E.
|
inline |
Return edge operators for the edge of index iE.
|
inline |
Return edge operators for the edge of index iE.
|
inline |
Return the potential operator on edge E.
|
inline |
Return the potential operator (expressed on Polykpo3) on edge E.
|
inline |
Return edge potential for the edge E.
|
inline |
Return the potential operator on the edge of index iE.
|
inline |
Return the potential operator (expressed on Polykpo3) on the edge of index iE.
|
inline |
Return edge potential for the edge of index iE.
|
inline |
Return the extension for cell T.
|
inline |
Return the extension for the cell of index iT.
double XGrad::evaluatePotential | ( | const size_t | iT, |
const Eigen::VectorXd & | vT, | ||
const VectorRd & | x | ||
) | const |
Evaluate the value of the potential at a point x.
iT | index of the cell in which to take the potential |
vT | vector of local DOFs |
x | point at which to evaluate the potential |
double XHess::evaluatePotential | ( | const size_t | iT, |
const Eigen::VectorXd & | vT, | ||
const VectorRd & | x | ||
) | const |
Evaluate the value of the potential at a point x.
iT | index of the cell in which to take the potential |
vT | vector of local DOFs |
x | point at which to evaluate the potential |
EXCurl::EXCurl | ( | const DDRCore & | ddr_core, |
bool | use_threads = true , |
||
std::ostream & | output = std::cout |
||
) |
Constructor.
|
inline |
_gradient | gradient operator |
_symmetric_gradient | symmetric gradient operator |
_potential | Potential operator |
_stabilisation | Stabilisation |
|
inline |
Return cell operators for cell T.
|
inline |
Return hho operators for the cell of index iT.
Eigen::VectorXd XHess::interpolate | ( | const FunctionType & | q, |
const GradFunctionType & | Dq, | ||
const int | deg_quad = -1 |
||
) | const |
Interpolator of a continuous function.
q | The function to interpolate |
Dq | The gradient of the function |
deg_quad | The optional degre of quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used. |
Eigen::VectorXd SXGrad::interpolate | ( | const FunctionType & | q, |
const int | deg_quad = -1 |
||
) | const |
q | The function to interpolate |
deg_quad | The optional degre of quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used. |
Eigen::VectorXd XGrad::interpolate | ( | const FunctionType & | q, |
const int | deg_quad = -1 |
||
) | const |
Interpolator of a continuous function.
q | The function to interpolate |
deg_quad | The optional degre of quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used. |
Eigen::VectorXd XRot::interpolate | ( | const FunctionType & | q, |
const int | deg_quad = -1 |
||
) | const |
Interpolator of a continuous function.
q | The function to interpolate |
deg_quad | The optional degre of quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used. |
Eigen::VectorXd VSXGrad::interpolate | ( | const FunctionType & | q, |
const int | doe_cell = -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. |
Eigen::VectorXd EXCurl::interpolate | ( | const FunctionType & | v, |
const int | deg_quad = -1 |
||
) | const |
Interpolator of a continuous function.
v | The function to interpolate |
deg_quad | The optional degre of quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used. |
Eigen::VectorXd XCurl::interpolate | ( | const FunctionType & | v, |
const int | deg_quad = -1 |
||
) | const |
Interpolator of a continuous function.
v | The function to interpolate |
deg_quad | The optional degre of quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used. |
Eigen::VectorXd XRotRot::interpolate | ( | const FunctionType & | v, |
const RotType & | rot_v, | ||
const int | deg_quad = -1 |
||
) | const |
Interpolator of a continuous function.
v | The function to interpolate |
rot_v | Rotor of the function to interpolate |
deg_quad | The optional degre of quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used. |
|
inline |
_curl | Curl operator |
_potential | Potential operator |
|
inline |
_gradient | Gradient operator |
_potential | Potential operator |
|
inline |
_gradient | Gradient operator |
_potential | Potential operator |
|
inline |
_rot | Curl operator |
_potential | Potential operator |
|
inline |
_rotor | Rot operator |
_potential | Potential operator |
|
inline |
_rotor | Vector rotor operator |
_rotor_rhs | Vector rotor right-hand side |
_potential | Potential operator |
|
inline |
Return a const reference to the mesh.
|
inline |
Return the mesh.
|
inline |
Return a const reference to the mesh.
|
inline |
Return the mesh.
|
inline |
Return the mesh.
|
inline |
Return the mesh.
|
inline |
Return the mesh.
|
inline |
Return the mesh.
|
inline |
Return the mesh.
|
inline |
Return the mesh.
|
inline |
Return the number of serendipity edges in a cell.
Eigen::VectorXd SerendipityProblem::nDOFs_cells_SXCurl | ( | ) | const |
Number of DOFs on cells for serendipity XCurl space.
Eigen::VectorXd SerendipityProblem::nDOFs_cells_SXGrad | ( | ) | const |
Number of DOFs on cells for serendipity XGrad space.
Eigen::VectorXd SerendipityProblem::nDOFs_cells_SXRotRot | ( | ) | const |
Number of DOFs on cells for serendipity XRotRot space.
|
inline |
|
inline |
|
inline |
Return basis for Matricial P^k space.
|
inline |
|
inline |
|
inline |
|
inline |
Return basis for the (P^{k+1})^2 space.
|
inline |
|
inline |
|
inline |
|
inline |
Return the reduction for cell T.
|
inline |
Return the reduction for the cell of index iT.
|
inline |
Return the serendipity degree ell_T in a cell.
|
inline |
Return the list of serendipity edges in a cell.
|
inline |
Compute the serendipity operator on the Cell T.
const Eigen::MatrixXd SerendipityProblem::SerendipityOperatorCell | ( | const size_t | iT, |
const Eigen::MatrixXd & | LT | ||
) | const |
Compute the serendipity operator on the cell of index iT.
SerendipityProblem::SerendipityProblem | ( | const DDRCore & | ddrcore, |
bool | use_threads = true , |
||
std::ostream & | output = std::cout |
||
) |
Constructor.
|
inline |
Return the serendipity reconstruction for cell T.
|
inline |
Return the serendipity reconstruction for cell T.
|
inline |
Return the serendipity reconstruction for the cell of index iT.
|
inline |
Return the serendipity reconstruction for the cell of index iT.
|
inline |
Return the underlying sXgrad space.
SXGrad::SXGrad | ( | const DDRCore & | ddr_core, |
const SerendipityProblem & | ser_pro, | ||
bool | use_threads = true , |
||
std::ostream & | output = std::cout |
||
) |
Constructor.
|
inline |
Return cell operators for cell T.
|
inline |
Return cell operators for the cell of index iT.
|
inline |
_serendipity | Serendipity reconstruction |
|
inline |
|
inline |
Return true if we use thread-based parallelism.
VSXGrad::VSXGrad | ( | const DDRCore & | ddr_core, |
const SerendipityProblem & | sp, | ||
bool | use_threads = true , |
||
std::ostream & | output = std::cout |
||
) |
Constructor.
XCurl::XCurl | ( | const DDRCore & | ddr_core, |
bool | use_threads = true , |
||
std::ostream & | output = std::cout |
||
) |
Constructor.
XGrad::XGrad | ( | const DDRCore & | ddr_core, |
bool | use_threads = true , |
||
std::ostream & | output = std::cout |
||
) |
Constructor.
XHess::XHess | ( | const DDRCore & | ddr_core, |
const SerendipityProblem & | ser_pro, | ||
bool | use_threads = true , |
||
std::ostream & | output = std::cout |
||
) |
Constructor.
XRot::XRot | ( | const DDRCore & | ddr_core, |
bool | use_threads = true , |
||
std::ostream & | output = std::cout |
||
) |
Constructor.
XRotRot::XRotRot | ( | const DDRCore & | ddr_core, |
bool | use_threads = true , |
||
std::ostream & | output = std::cout |
||
) |
Constructor.
Eigen::MatrixXd HArDCore2D::XCurl::LocalOperators::curl |
Eigen::MatrixXd HArDCore2D::SXGrad::TransferOperators::extension |
std::unique_ptr<GolyComplBasisCellType> HArDCore2D::DDRCore::CellBases::GolyComplkp2 |
Eigen::MatrixXd HArDCore2D::EXCurl::hhoLocalOperators::gradient |
Eigen::MatrixXd HArDCore2D::XGrad::LocalOperators::gradient |
Eigen::MatrixXd HArDCore2D::XHess::LocalOperators::gradient |
std::unique_ptr<PolyBasisCellType> HArDCore2D::DDRCore::CellBases::Polyk |
std::unique_ptr<PolyBasisEdgeType> HArDCore2D::DDRCore::EdgeBases::Polyk |
std::unique_ptr<Poly2BasisCellType> HArDCore2D::DDRCore::CellBases::Polyk2 |
std::unique_ptr<Poly2BasisEdgeType> HArDCore2D::VSXGrad::EdgeBases::Polyk2 |
std::unique_ptr<PolyBasisCellType> HArDCore2D::DDRCore::CellBases::Polykmo |
std::unique_ptr<PolyBasisEdgeType> HArDCore2D::DDRCore::EdgeBases::Polykmo |
std::unique_ptr<Poly2BasisCellType> HArDCore2D::DDRCore::CellBases::Polykmo2 |
std::unique_ptr<Poly2BasisCellType> HArDCore2D::VSXGrad::CellBases::Polykmo2 |
std::unique_ptr<Poly2BasisEdgeType> HArDCore2D::VSXGrad::EdgeBases::Polykmo2 |
std::unique_ptr<PolyBasisCellType> HArDCore2D::DDRCore::CellBases::Polykmtwo |
std::unique_ptr<PolyBasisCellType> HArDCore2D::DDRCore::CellBases::Polykpo |
std::unique_ptr<PolyBasisEdgeType> HArDCore2D::DDRCore::EdgeBases::Polykpo |
std::unique_ptr<Poly2BasisEdgeType> HArDCore2D::VSXGrad::EdgeBases::Polykpo2 |
std::unique_ptr<Poly2x2BasisCellType> HArDCore2D::VSXGrad::CellBases::Polykpo2x2 |
std::unique_ptr<Poly2BasisCellType> HArDCore2D::VSXGrad::CellBases::Polykptwo2 |
std::unique_ptr<Poly2BasisEdgeType> HArDCore2D::VSXGrad::EdgeBases::Polykptwo2 |
Eigen::MatrixXd HArDCore2D::EXCurl::hhoLocalOperators::potential |
Eigen::MatrixXd HArDCore2D::VSXGrad::LocalOperators::potential |
Eigen::MatrixXd HArDCore2D::XCurl::LocalOperators::potential |
Eigen::MatrixXd HArDCore2D::XGrad::LocalOperators::potential |
Eigen::MatrixXd HArDCore2D::XHess::LocalOperators::potential |
Eigen::MatrixXd HArDCore2D::XRot::LocalOperators::potential |
Eigen::MatrixXd HArDCore2D::XRotRot::LocalOperators::potential |
Eigen::MatrixXd HArDCore2D::SXGrad::TransferOperators::reduction |
std::unique_ptr<RolyComplBasisCellType> HArDCore2D::DDRCore::CellBases::RolyComplk |
std::unique_ptr<RolyComplBasisCellType> HArDCore2D::DDRCore::CellBases::RolyComplkmo |
std::unique_ptr<RolyComplBasisCellType> HArDCore2D::DDRCore::CellBases::RolyComplkp2 |
std::unique_ptr<RolyBasisCellType> HArDCore2D::DDRCore::CellBases::Rolykmo |
Eigen::MatrixXd HArDCore2D::VSXGrad::LocalOperators::rot |
Eigen::MatrixXd HArDCore2D::XRot::LocalOperators::rotor |
Eigen::MatrixXd HArDCore2D::XRotRot::LocalOperators::rotor |
Eigen::MatrixXd HArDCore2D::XRot::LocalOperators::rotor_rhs |
Eigen::MatrixXd HArDCore2D::SXGrad::TransferOperators::serendipity |
Eigen::MatrixXd HArDCore2D::XHess::TransferOperators::serendipity |
Eigen::MatrixXd HArDCore2D::EXCurl::hhoLocalOperators::stabilisation |
Eigen::MatrixXd HArDCore2D::EXCurl::hhoLocalOperators::symmetric_gradient |