20 typedef std::function<double(
const Eigen::Vector2d &)>
FunctionType;
26 const Eigen::MatrixXd & _rotor,
27 const Eigen::MatrixXd & _rotor_rhs,
28 const Eigen::MatrixXd & _potential
48 return m_ddr_core.
mesh();
60 return m_ddr_core.
degree();
66 const int deg_quad = -1
72 return *m_edge_potentials[iE];
78 return *m_edge_potentials[E.global_index()];
84 return *m_cell_operators[iT];
90 return *m_cell_operators[
T.global_index()];
102 return m_ddr_core.
cellBases(
T.global_index());
114 return m_ddr_core.
edgeBases(E.global_index());
120 const double & penalty_factor_cell = 1.,
121 const double & penalty_factor_edge = 1.
129 Eigen::MatrixXd _compute_edge_potential(
size_t iE);
130 LocalOperators _compute_cell_operators(
size_t iT);
134 std::ostream & m_output;
137 std::vector<std::unique_ptr<Eigen::MatrixXd> > m_edge_potentials;
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 xrot.hpp:18
const Mesh & mesh() const
Return a const reference to the mesh.
Definition ddrcore.hpp:115
const bool useThreads() const
Return true if we use thread-based parallelism.
Definition xrot.hpp:52
Eigen::MatrixXd computeRotorL2Product(const size_t iT, const double &penalty_factor_cell=1., const double &penalty_factor_edge=1.) const
Compute rotor L2-product.
Definition xrot.cpp:287
const size_t & degree() const
Return the polynomial degree.
Definition ddrcore.hpp:121
Eigen::MatrixXd rotor_rhs
Definition xrot.hpp:38
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition xrot.hpp:106
const size_t & degree() const
Return the polynomial degree.
Definition xrot.hpp:58
double computeRotorL2Norm(const Eigen::VectorXd &v) const
Compute rotor L2-norm.
Definition xrot.cpp:348
Eigen::MatrixXd rotor
Definition xrot.hpp:37
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition xrot.hpp:100
std::function< double(const Eigen::Vector2d &)> FunctionType
Definition xrot.hpp:20
Eigen::MatrixXd potential
Definition xrot.hpp:39
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the cell of index iT.
Definition xrot.hpp:94
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition ddrcore.hpp:127
const Eigen::MatrixXd & edgePotential(const Edge &E) const
Return edge potential for the edge E.
Definition xrot.hpp:76
Eigen::VectorXd interpolate(const FunctionType &q, const int deg_quad=-1) const
Interpolator of a continuous function.
Definition xrot.cpp:61
const Mesh & mesh() const
Return the mesh.
Definition xrot.hpp:46
const DDRCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition xrot.hpp:112
const Eigen::MatrixXd & edgePotential(size_t iE) const
Return edge potential for the edge of index iE.
Definition xrot.hpp:70
const EdgeBases & edgeBases(size_t iE) const
Return edge bases for edge iE.
Definition ddrcore.hpp:135
const LocalOperators & cellOperators(const Cell &T) const
Return cell operators for cell T.
Definition xrot.hpp:88
LocalOperators(const Eigen::MatrixXd &_rotor, const Eigen::MatrixXd &_rotor_rhs, const Eigen::MatrixXd &_potential)
Definition xrot.hpp:25
const LocalOperators & cellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition xrot.hpp:82
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
A structure to store local operators (vector rotor and potential)
Definition xrot.hpp:24