18 typedef std::function<Eigen::Vector2d(
const Eigen::Vector2d &)>
FunctionType;
25 const Eigen::MatrixXd & _serendipity,
26 const Eigen::MatrixXd & _extension,
27 const Eigen::MatrixXd & _reduction
47 return m_ddr_core.
mesh();
53 return m_ddr_core.
degree();
58 const int deg_quad = -1
63 inline const Eigen::MatrixXd &
ScurlCell(
size_t iT)
const
65 return (*m_cell_transfer_operators[iT]).serendipity;
69 inline const Eigen::MatrixXd &
EcurlCell(
size_t iT)
const
71 return (*m_cell_transfer_operators[iT]).extension;
75 inline const Eigen::MatrixXd &
RcurlCell(
size_t iT)
const
77 return (*m_cell_transfer_operators[iT]).reduction;
81 inline const Eigen::MatrixXd &
ScurlCell(
const Cell &
T)
const
87 inline const Eigen::MatrixXd &
EcurlCell(
const Cell &
T)
const
93 inline const Eigen::MatrixXd &
RcurlCell(
const Cell &
T)
const
102 inline const Eigen::MatrixXd
cellCurl(
size_t iT)
const
108 inline const Eigen::MatrixXd
cellCurl(
const Cell &
T)
const
128 const double & penalty_factor = 1.,
129 const Eigen::MatrixXd & mass_Pk2_T = Eigen::MatrixXd::Zero(1,1),
143 const std::string & side,
144 const double & penalty_factor = 1.,
145 const Eigen::MatrixXd & mass_Pk2_T = Eigen::MatrixXd::Zero(1,1),
161 return m_ddr_core.
cellBases(
T.global_index());
173 return m_ddr_core.
edgeBases(E.global_index());
177 TransferOperators _compute_cell_transfer_operators(
size_t iT);
184 std::vector<std::unique_ptr<TransferOperators> > m_cell_transfer_operators;
187 std::ostream & m_output;
Construct all polynomial spaces for the DDR sequence.
Definition ddrcore.hpp:63
Discrete Serendipity Hcurl space: local operators, L2 product and global interpolator.
Definition sxcurl.hpp:16
const Eigen::MatrixXd & RcurlCell(size_t iT) const
Return the reduction for the cell of index iT.
Definition sxcurl.hpp:75
const Eigen::MatrixXd & ScurlCell(size_t iT) const
Return the serendipity reconstruction for the cell of index iT.
Definition sxcurl.hpp:63
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the face of index iT.
Definition sxcurl.hpp:153
const Eigen::MatrixXd cellPotential(const Cell &T) const
Return the potential operator on cell T.
Definition sxcurl.hpp:120
Eigen::VectorXd interpolate(const FunctionType &v, const int deg_quad=-1) const
Interpolator of a continuous function.
Definition sxcurl.cpp:51
const Eigen::MatrixXd cellCurl(size_t iT) const
Return the full curl operator on the cell of index iT.
Definition sxcurl.hpp:102
const Eigen::MatrixXd cellPotential(size_t iT) const
Return the potential operator on the cell of index iT.
Definition sxcurl.hpp:114
const Eigen::MatrixXd cellCurl(const Cell &T) const
Return the full curl operator on cell T.
Definition sxcurl.hpp:108
const Eigen::MatrixXd & EcurlCell(size_t iT) const
Return the extension for the cell of index iT.
Definition sxcurl.hpp:69
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition sxcurl.hpp:159
const Eigen::MatrixXd & EcurlCell(const Cell &T) const
Return the extension for cell T.
Definition sxcurl.hpp:87
Eigen::MatrixXd computeL2ProductGradient(const size_t iT, const SXGrad &sx_grad, const std::string &side, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk2_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:186
const size_t & degree() const
Return the polynomial degree.
Definition sxcurl.hpp:51
const DDRCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition sxcurl.hpp:171
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition sxcurl.hpp:165
const Mesh & mesh() const
Return the mesh.
Definition sxcurl.hpp:45
std::function< Eigen::Vector2d(const Eigen::Vector2d &)> FunctionType
Definition sxcurl.hpp:18
const Eigen::MatrixXd & ScurlCell(const Cell &T) const
Return the serendipity reconstruction for cell T.
Definition sxcurl.hpp:81
Eigen::MatrixXd computeL2Product(const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk2_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const
Compute the matrix of the (weighted) L2-product.
Definition sxcurl.hpp:126
const Eigen::MatrixXd & RcurlCell(const Cell &T) const
Return the reduction for cell T.
Definition sxcurl.hpp:93
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:20
Base class for global DOF spaces.
Definition variabledofspace.hpp:17
Discrete Hcurl space: local operators, L2 product and global interpolator.
Definition xcurl.hpp:18
const Mesh & mesh() const
Return a const reference to the mesh.
Definition ddrcore.hpp:116
const size_t & degree() const
Return the polynomial degree.
Definition ddrcore.hpp:122
Eigen::MatrixXd potential
Definition xcurl.hpp:36
const LocalOperators & cellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition xcurl.hpp:61
Eigen::MatrixXd computeL2Product(const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk2_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 xcurl.cpp:189
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition ddrcore.hpp:128
const EdgeBases & edgeBases(size_t iE) const
Return edge bases for edge iE.
Definition ddrcore.hpp:136
Eigen::MatrixXd curl
Definition xcurl.hpp:35
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
Structure to store element bases.
Definition ddrcore.hpp:81
Structure to store edge bases.
Definition ddrcore.hpp:103
Structure for weights (scalar, at the moment) in integral.
Definition integralweight.hpp:33
A structure to store the serendipity, extension and reduction operators.
Definition sxcurl.hpp:23
TransferOperators(const Eigen::MatrixXd &_serendipity, const Eigen::MatrixXd &_extension, const Eigen::MatrixXd &_reduction)
Definition sxcurl.hpp:24
Eigen::MatrixXd extension
Definition sxcurl.hpp:37
Eigen::MatrixXd reduction
Definition sxcurl.hpp:38
Eigen::MatrixXd serendipity
Definition sxcurl.hpp:36