23 typedef std::function<Eigen::Vector3d(
const Eigen::Vector3d &)>
FunctionType;
52 return m_ddr_core.
mesh();
58 return m_ddr_core.
degree();
80 return (*m_face_transfer_operators[
iF]).serendipity;
86 return (*m_face_transfer_operators[
iF]).extension;
92 return (*m_face_transfer_operators[
iF]).reduction;
96 inline const Eigen::MatrixXd &
ScurlFace(
const Face &
F)
const
102 inline const Eigen::MatrixXd &
EcurlFace(
const Face &
F)
const
108 inline const Eigen::MatrixXd &
RcurlFace(
const Face &
F)
const
118 return (*m_cell_transfer_operators[
iT]).serendipity;
124 return (*m_cell_transfer_operators[
iT]).extension;
130 return (*m_cell_transfer_operators[
iT]).reduction;
134 inline const Eigen::MatrixXd &
ScurlCell(
const Cell & T)
const
140 inline const Eigen::MatrixXd &
EcurlCell(
const Cell & T)
const
146 inline const Eigen::MatrixXd &
RcurlCell(
const Cell & T)
const
161 inline const Eigen::MatrixXd
faceCurl(
const Face &
F)
const
185 inline const Eigen::MatrixXd
cellCurl(
const Cell & T)
const
206 const Eigen::MatrixXd &
mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
220 const std::string &
side,
222 const Eigen::MatrixXd &
mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
228 const Eigen::VectorXd & u
243 return m_ddr_core.
cellBases(T.global_index());
255 return m_ddr_core.
faceBases(
F.global_index());
267 return m_ddr_core.
edgeBases(E.global_index());
271 TransferOperators _compute_face_transfer_operators(
size_t iF);
272 TransferOperators _compute_cell_transfer_operators(
size_t iT);
279 std::vector<std::unique_ptr<TransferOperators> > m_face_transfer_operators;
280 std::vector<std::unique_ptr<TransferOperators> > m_cell_transfer_operators;
283 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
@ Matrix
Definition basis.hpp:67
const LocalOperators & faceOperators(size_t iF) const
Return face operators for the face of index iF.
Definition xcurl.hpp:76
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 Mesh & mesh() const
Return the mesh.
Definition sxcurl.hpp:50
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 cellPotential(const Cell &T) const
Return the potential operator on cell T.
Definition sxcurl.hpp:197
const Eigen::MatrixXd faceCurl(size_t iF) const
Return the full curl operator on the face of index iF.
Definition sxcurl.hpp:155
const Eigen::MatrixXd & RcurlCell(const Cell &T) const
Return the reduction for cell T.
Definition sxcurl.hpp:146
const Eigen::MatrixXd & EcurlCell(const Cell &T) const
Return the extension for cell T.
Definition sxcurl.hpp:140
const Eigen::MatrixXd & RcurlFace(size_t iF) const
Return the reduction for the face of index iF.
Definition sxcurl.hpp:90
const Eigen::MatrixXd & EcurlFace(const Face &F) const
Return the extension for face F.
Definition sxcurl.hpp:102
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 Mesh & mesh() const
Return a const reference to the mesh.
Definition ddrcore.hpp:134
const Eigen::MatrixXd faceCurl(const Face &F) const
Return the full curl operator on face F.
Definition sxcurl.hpp:161
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 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 Eigen::MatrixXd & ScurlCell(size_t iT) const
Return the serendipity reconstruction for the cell of index iT.
Definition sxcurl.hpp:116
std::vector< VectorRd > computeVertexValues(const Eigen::VectorXd &u) const
Computes the values of the potential reconstruction at the mesh vertices.
Definition sxcurl.cpp:405
const Eigen::MatrixXd & ScurlCell(const Cell &T) const
Return the serendipity reconstruction for cell T.
Definition sxcurl.hpp:134
const DDRCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition sxcurl.hpp:265
const Eigen::MatrixXd & RcurlCell(size_t iT) const
Return the reduction for the cell of index iT.
Definition sxcurl.hpp:128
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the face of index iT.
Definition sxcurl.hpp:235
const DDRCore::FaceBases & faceBases(const Face &F) const
Return cell bases for face F.
Definition sxcurl.hpp:253
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition sxcurl.hpp:259
const DDRCore::FaceBases & faceBases(size_t iF) const
Return face bases for the face of index iF.
Definition sxcurl.hpp:247
const LocalOperators & cellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition xcurl.hpp:64
Eigen::MatrixXd potential
Definition xcurl.hpp:37
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition sxcurl.hpp:241
Eigen::MatrixXd extension
Definition sxcurl.hpp:42
const size_t & degree() const
Return the polynomial degree.
Definition ddrcore.hpp:140
const SerendipityProblem & serPro() const
Definition sxcurl.hpp:61
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
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 & ScurlFace(const Face &F) const
Return the serendipity reconstruction for face F.
Definition sxcurl.hpp:96
const Eigen::MatrixXd cellCurl(const Cell &T) const
Return the full curl operator on cell T.
Definition sxcurl.hpp:185
const size_t & degree() const
Return the polynomial degree.
Definition sxcurl.hpp:56
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 & EcurlCell(size_t iT) const
Return the extension for the cell of index iT.
Definition sxcurl.hpp:122
const Eigen::MatrixXd & EcurlFace(size_t iF) const
Return the extension for the face of index iF.
Definition sxcurl.hpp:84
const FaceBases & faceBases(size_t iF) const
Return face bases for face iF.
Definition ddrcore.hpp:154
const Eigen::MatrixXd & ScurlFace(size_t iF) const
Return the serendipity reconstruction for the face of index iF.
Definition sxcurl.hpp:78
const Eigen::MatrixXd & RcurlFace(const Face &F) const
Return cell reduction for cell T.
Definition sxcurl.hpp:108
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition ddrcore.hpp:146
Eigen::MatrixXd reduction
Definition sxcurl.hpp:43
const Eigen::MatrixXd facePotential(const Face &F) const
Return the potential operator on face F.
Definition sxcurl.hpp:173
bool use_threads
Definition HHO_DiffAdvecReac.hpp:47
Definition ddr-magnetostatics.hpp:41
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