21 typedef std::function<Eigen::Vector3d(
const Eigen::Vector3d &)>
FunctionType;
27 const Eigen::MatrixXd & _curl,
28 const Eigen::MatrixXd & _potential
46 return m_ddr_core.
mesh();
52 return m_ddr_core.
degree();
58 const int doe_cell = -1,
59 const int doe_face = -1,
60 const int doe_edge = -1
66 return *m_cell_operators[iT];
72 return * m_cell_operators[T.global_index()];
78 return *m_face_operators[iF];
84 return * m_face_operators[F.global_index()];
96 return m_ddr_core.
cellBases(T.global_index());
108 return m_ddr_core.
faceBases(F.global_index());
120 return m_ddr_core.
edgeBases(E.global_index());
128 const double & penalty_factor = 1.,
129 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
136 const XGrad & x_grad,
137 const std::string & side,
138 const double & penalty_factor = 1.,
139 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
146 const std::vector<Eigen::MatrixXd> & leftOp,
147 const std::vector<Eigen::MatrixXd> & rightOp,
148 const double & penalty_factor,
149 const Eigen::MatrixXd & w_mass_Pk3_T,
157 const double & penalty_factor = 1.,
158 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
163 LocalOperators _compute_face_curl_potential(
size_t iF);
164 LocalOperators _compute_cell_curl_potential(
size_t iT);
168 std::ostream & m_output;
171 std::vector<std::unique_ptr<LocalOperators> > m_cell_operators;
172 std::vector<std::unique_ptr<LocalOperators> > m_face_operators;
Construct all polynomial spaces for the DDR sequence.
Definition: ddrcore.hpp:62
Base class for global DOF spaces. Provides functions to manipulate global DOFs (the local version bei...
Definition: globaldofspace.hpp:16
Discrete Hcurl space: local operators, L2 product and global interpolator.
Definition: xcurl.hpp:19
Discrete H1 space: local operators, L2 product and global interpolator.
Definition: xgrad.hpp:18
Class to describe a mesh.
Definition: MeshND.hpp:17
const size_t & degree() const
Return the polynomial degree.
Definition: ddrcore.hpp:140
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the face of index iT.
Definition: xcurl.hpp:88
const DDRCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition: xcurl.hpp:118
Eigen::MatrixXd 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_Pk3_T, const IntegralWeight &weight) const
Compute the matrix of the L2 product, applying leftOp and rightOp to the variables....
Definition: xcurl.cpp:470
Eigen::MatrixXd computeL2ProductDOFs(size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk3_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.))
Alternate version of L2 product, using only cell potentials and 'DOFs' (as in Di-Pietro & Droniou,...
Definition: xcurl.cpp:589
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition: ddrcore.hpp:146
const LocalOperators & faceOperators(size_t iF) const
Return face operators for the face of index iF.
Definition: xcurl.hpp:76
Eigen::VectorXd interpolate(const FunctionType &v, const int doe_cell=-1, const int doe_face=-1, const int doe_edge=-1) const
Interpolator of a continuous function.
Definition: xcurl.cpp:63
const DDRCore::FaceBases & faceBases(size_t iF) const
Return face bases for the face of index iF.
Definition: xcurl.hpp:100
Eigen::MatrixXd computeL2ProductGradient(const size_t iT, const XGrad &x_grad, const std::string &side, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk3_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 discre...
Definition: xcurl.cpp:395
const LocalOperators & cellOperators(const Cell &T) const
Return cell operators for cell T.
Definition: xcurl.hpp:70
const LocalOperators & cellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition: xcurl.hpp:64
const Mesh & mesh() const
Return a const reference to the mesh.
Definition: ddrcore.hpp:134
const EdgeBases & edgeBases(size_t iE) const
Return edge bases for edge iE.
Definition: ddrcore.hpp:162
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType
Definition: xcurl.hpp:21
const LocalOperators & faceOperators(const Face &F) const
Return face operators for face F.
Definition: xcurl.hpp:82
const DDRCore::FaceBases & faceBases(const Face &F) const
Return cell bases for face F.
Definition: xcurl.hpp:106
Eigen::MatrixXd potential
Definition: xcurl.hpp:37
const Mesh & mesh() const
Return the mesh.
Definition: xcurl.hpp:44
Eigen::MatrixXd computeL2Product(const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk3_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. The stabilisation here is b...
Definition: xcurl.cpp:346
Eigen::MatrixXd curl
Definition: xcurl.hpp:36
XCurl(const DDRCore &ddr_core, bool use_threads=true, std::ostream &output=std::cout)
Constructor.
Definition: xcurl.cpp:13
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition: xcurl.hpp:94
const FaceBases & faceBases(size_t iF) const
Return face bases for face iF.
Definition: ddrcore.hpp:154
const size_t & degree() const
Return the polynomial degree.
Definition: xcurl.hpp:50
LocalOperators(const Eigen::MatrixXd &_curl, const Eigen::MatrixXd &_potential)
Definition: xcurl.hpp:26
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition: xcurl.hpp:112
bool use_threads
Definition: HHO_DiffAdvecReac.hpp:47
Definition: ddr-magnetostatics.hpp:40
MeshND::Edge< 2 > Edge
Definition: Mesh2D.hpp:11
MeshND::Face< 2 > Face
Definition: Mesh2D.hpp:12
MeshND::Cell< 2 > Cell
Definition: Mesh2D.hpp:13
Structure to store element bases.
Definition: ddrcore.hpp:86
Structure to store edge bases.
Definition: ddrcore.hpp:121
Structure to store face bases.
Definition: ddrcore.hpp:105
Structure for weights (scalar, at the moment) in integral.
Definition: integralweight.hpp:36
A structure to store the local operators (curl and potential)
Definition: xcurl.hpp:25