19 typedef std::function<double(
const Eigen::Vector2d &)>
FunctionType;
25 const Eigen::MatrixXd & _gradient,
26 const Eigen::MatrixXd & _potential
44 return m_ddr_core.
mesh();
50 return m_ddr_core.
degree();
56 const int deg_quad = -1
62 return *m_edge_operators[iE];
68 return *m_edge_operators[E.global_index()];
74 return *m_cell_operators[iT];
80 return *m_cell_operators[
T.global_index()];
104 return m_ddr_core.
edgeBases(E.global_index());
112 const double & penalty_factor = 1.,
113 const Eigen::MatrixXd & mass_Pkpo_T = Eigen::MatrixXd::Zero(1,1),
120 const Eigen::VectorXd & vT,
128 LocalOperators _compute_edge_gradient_potential(
size_t iE);
129 LocalOperators _compute_cell_gradient_potential(
size_t iT);
130 double _compute_squared_l2_norm(
size_t iT,
const Eigen::VectorXd & vT)
const;
134 std::ostream & m_output;
137 std::vector<std::unique_ptr<LocalOperators> > m_edge_operators;
138 std::vector<std::unique_ptr<LocalOperators> > m_cell_operators;
Construct all polynomial spaces for the DDR sequence.
Definition ddrcore.hpp:63
Base class for global DOF spaces. Provides functions to manipulate global DOFs (the local version bei...
Definition globaldofspace.hpp:16
Discrete H1 space: local operators, L2 product and global interpolator.
Definition xgrad.hpp:17
Eigen::Vector2d VectorRd
Definition basis.hpp:55
std::function< double(const Eigen::Vector2d &)> FunctionType
Definition xgrad.hpp:19
const Mesh & mesh() const
Return a const reference to the mesh.
Definition ddrcore.hpp:115
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition xgrad.hpp:96
const size_t & degree() const
Return the polynomial degree.
Definition ddrcore.hpp:121
Eigen::MatrixXd potential
Definition xgrad.hpp:35
LocalOperators(const Eigen::MatrixXd &_gradient, const Eigen::MatrixXd &_potential)
Definition xgrad.hpp:24
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the cell of index iT.
Definition xgrad.hpp:84
const size_t & degree() const
Return the polynomial degree.
Definition xgrad.hpp:48
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition ddrcore.hpp:127
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition xgrad.hpp:90
const LocalOperators & cellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition xgrad.hpp:72
Eigen::VectorXd interpolate(const FunctionType &q, const int deg_quad=-1) const
Interpolator of a continuous function.
Definition xgrad.cpp:62
Eigen::MatrixXd 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.
Definition xgrad.cpp:273
double evaluatePotential(const size_t iT, const Eigen::VectorXd &vT, const VectorRd &x) const
Evaluate the value of the potential at a point x.
Definition xgrad.cpp:354
const DDRCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition xgrad.hpp:102
const LocalOperators & cellOperators(const Cell &T) const
Return cell operators for cell T.
Definition xgrad.hpp:78
const LocalOperators & edgeOperators(size_t iE) const
Return edge operators for the edge of index iE.
Definition xgrad.hpp:60
const LocalOperators & edgeOperators(const Edge &E) const
Return edge operators for edge E.
Definition xgrad.hpp:66
Eigen::MatrixXd gradient
Definition xgrad.hpp:34
const Mesh & mesh() const
Return the mesh.
Definition xgrad.hpp:42
double computeL2Norm(const Eigen::VectorXd &v) const
Compute the L2-norm of a vector of the space.
Definition xgrad.cpp:378
const EdgeBases & edgeBases(size_t iE) const
Return edge bases for edge iE.
Definition ddrcore.hpp:135
bool use_threads
Definition HHO_DiffAdvecReac.hpp:47
if(strcmp(field, 'real')) % real valued entries T
Definition mmread.m:93
Definition ddr-klplate.hpp:27
static auto v
Definition ddrcore-test.hpp:32
static auto q
Definition ddrcore-test.hpp:14
Structure to store element bases.
Definition ddrcore.hpp:81
Structure to store edge bases.
Definition ddrcore.hpp:102
Structure for weights (scalar, at the moment) in integral.
Definition integralweight.hpp:33
A structure to store local operators (gradient and potential)
Definition xgrad.hpp:23