23 typedef std::function<Eigen::Vector3d(
const Eigen::Vector3d &)>
FunctionType;
30 const Eigen::MatrixXd & _serendipity,
31 const Eigen::MatrixXd & _extension,
32 const Eigen::MatrixXd & _reduction
52 return m_ddr_core.
mesh();
58 return m_ddr_core.
degree();
69 const int doe_cell = -1,
70 const int doe_face = -1,
71 const int doe_edge = -1
78 inline const Eigen::MatrixXd &
ScurlFace(
size_t iF)
const
80 return (*m_face_transfer_operators[iF]).serendipity;
84 inline const Eigen::MatrixXd &
EcurlFace(
size_t iF)
const
86 return (*m_face_transfer_operators[iF]).extension;
90 inline const Eigen::MatrixXd &
RcurlFace(
size_t iF)
const
92 return (*m_face_transfer_operators[iF]).reduction;
116 inline const Eigen::MatrixXd &
ScurlCell(
size_t iT)
const
118 return (*m_cell_transfer_operators[iT]).serendipity;
122 inline const Eigen::MatrixXd &
EcurlCell(
size_t iT)
const
124 return (*m_cell_transfer_operators[iT]).extension;
128 inline const Eigen::MatrixXd &
RcurlCell(
size_t iT)
const
130 return (*m_cell_transfer_operators[iT]).reduction;
155 inline const Eigen::MatrixXd
faceCurl(
size_t iF)
const
179 inline const Eigen::MatrixXd
cellCurl(
size_t iT)
const
205 const double & penalty_factor = 1.,
206 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
220 const std::string & side,
221 const double & penalty_factor = 1.,
222 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
238 return m_ddr_core.
cellBases(T.global_index());
250 return m_ddr_core.
faceBases(F.global_index());
262 return m_ddr_core.
edgeBases(E.global_index());
266 TransferOperators _compute_face_transfer_operators(
size_t iF);
267 TransferOperators _compute_cell_transfer_operators(
size_t iT);
274 std::vector<std::unique_ptr<TransferOperators> > m_face_transfer_operators;
275 std::vector<std::unique_ptr<TransferOperators> > m_cell_transfer_operators;
278 std::ostream & m_output;
Construct all polynomial spaces for the DDR sequence.
Definition: ddrcore.hpp:62
Discrete Serendipity Hcurl space: local operators, L2 product and global interpolator.
Definition: sxcurl.hpp:21
Discrete Serendipity Hgrad space: local operators, L2 product and global interpolator.
Definition: sxgrad.hpp:20
Construct all polynomial spaces for the DDR sequence.
Definition: serendipity_problem.hpp:40
Base class for global DOF spaces.
Definition: variabledofspace.hpp:17
Discrete Hcurl space: local operators, L2 product and global interpolator.
Definition: xcurl.hpp:19
Class to describe a mesh.
Definition: MeshND.hpp:17
const size_t & degree() const
Return the polynomial degree.
Definition: ddrcore.hpp:140
const Eigen::MatrixXd cellCurl(size_t iT) const
Return the full curl operator on the cell of index iT.
Definition: sxcurl.hpp:179
Eigen::MatrixXd serendipity
Definition: sxcurl.hpp:41
const Eigen::MatrixXd cellPotential(size_t iT) const
Return the potential operator on the cell of index iT.
Definition: sxcurl.hpp:191
const Eigen::MatrixXd & EcurlCell(size_t iT) const
Return the extension for the cell of index iT.
Definition: sxcurl.hpp:122
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition: ddrcore.hpp:146
const Eigen::MatrixXd & EcurlFace(size_t iF) const
Return the extension for the face of index iF.
Definition: sxcurl.hpp:84
const LocalOperators & faceOperators(size_t iF) const
Return face operators for the face of index iF.
Definition: xcurl.hpp:76
const Eigen::MatrixXd cellPotential(const Cell &T) const
Return the potential operator on cell T.
Definition: sxcurl.hpp:197
const DDRCore::FaceBases & faceBases(size_t iF) const
Return face bases for the face of index iF.
Definition: sxcurl.hpp:242
const Eigen::MatrixXd faceCurl(size_t iF) const
Return the full curl operator on the face of index iF.
Definition: sxcurl.hpp:155
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition: sxcurl.hpp:254
const Eigen::MatrixXd & RcurlFace(size_t iF) const
Return the reduction for the face of index iF.
Definition: sxcurl.hpp:90
const DDRCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition: sxcurl.hpp:260
const DDRCore::FaceBases & faceBases(const Face &F) const
Return cell bases for face F.
Definition: sxcurl.hpp:248
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: sxcurl.cpp:65
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 DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition: sxcurl.hpp:236
const Eigen::MatrixXd faceCurl(const Face &F) const
Return the full curl operator on face F.
Definition: sxcurl.hpp:161
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType
Definition: sxcurl.hpp:23
TransferOperators(const Eigen::MatrixXd &_serendipity, const Eigen::MatrixXd &_extension, const Eigen::MatrixXd &_reduction)
Definition: sxcurl.hpp:29
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.
Definition: sxcurl.hpp:203
const EdgeBases & edgeBases(size_t iE) const
Return edge bases for edge iE.
Definition: ddrcore.hpp:162
const size_t & degree() const
Return the polynomial degree.
Definition: sxcurl.hpp:56
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the face of index iT.
Definition: sxcurl.hpp:230
const Eigen::MatrixXd & ScurlFace(const Face &F) const
Return the serendipity reconstruction for face F.
Definition: sxcurl.hpp:96
SXCurl(const DDRCore &ddr_core, const SerendipityProblem &ser_pro, bool use_threads=true, std::ostream &output=std::cout)
Constructor.
Definition: sxcurl.cpp:13
const Mesh & mesh() const
Return the mesh.
Definition: sxcurl.hpp:50
const Eigen::MatrixXd & RcurlCell(const Cell &T) const
Return the reduction for cell T.
Definition: sxcurl.hpp:146
const Eigen::MatrixXd & RcurlFace(const Face &F) const
Return cell reduction for cell T.
Definition: sxcurl.hpp:108
Eigen::MatrixXd potential
Definition: xcurl.hpp:37
Eigen::MatrixXd extension
Definition: sxcurl.hpp:42
Eigen::MatrixXd computeL2ProductGradient(const size_t iT, const SXGrad &sx_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: sxcurl.cpp:327
const Eigen::MatrixXd & ScurlCell(const Cell &T) const
Return the serendipity reconstruction for cell T.
Definition: sxcurl.hpp:134
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
const Eigen::MatrixXd cellCurl(const Cell &T) const
Return the full curl operator on cell T.
Definition: sxcurl.hpp:185
const Eigen::MatrixXd & EcurlCell(const Cell &T) const
Return the extension for cell T.
Definition: sxcurl.hpp:140
const Eigen::MatrixXd & EcurlFace(const Face &F) const
Return the extension for face F.
Definition: sxcurl.hpp:102
const SerendipityProblem & serPro() const
Definition: sxcurl.hpp:61
const FaceBases & faceBases(size_t iF) const
Return face bases for face iF.
Definition: ddrcore.hpp:154
const Eigen::MatrixXd facePotential(size_t iF) const
Return the potential operator on the face of index iF.
Definition: sxcurl.hpp:167
const Eigen::MatrixXd & ScurlFace(size_t iF) const
Return the serendipity reconstruction for the face of index iF.
Definition: sxcurl.hpp:78
Eigen::MatrixXd reduction
Definition: sxcurl.hpp:43
const Eigen::MatrixXd & ScurlCell(size_t iT) const
Return the serendipity reconstruction for the cell of index iT.
Definition: sxcurl.hpp:116
const Eigen::MatrixXd facePotential(const Face &F) const
Return the potential operator on face F.
Definition: sxcurl.hpp:173
const Eigen::MatrixXd & RcurlCell(size_t iT) const
Return the reduction for the cell of index iT.
Definition: sxcurl.hpp:128
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 serendipity, extension and reduction operators.
Definition: sxcurl.hpp:28