6 #include <unsupported/Eigen/KroneckerProduct>
20 typedef std::function<Eigen::Vector3d(
const Eigen::Vector3d &)>
FunctionType;
28 const Eigen::MatrixXd & _serendipity,
29 const Eigen::MatrixXd & _extension,
30 const Eigen::MatrixXd & _reduction
50 return m_sxcurl.
mesh();
74 const int doe_cell = -1,
75 const int doe_face = -1,
76 const int doe_edge = -1
82 inline const Eigen::VectorXd
LaxcurlI(
size_t iG, Eigen::VectorXd & v)
const
84 return Eigen::Map<Eigen::VectorXd,0,Eigen::InnerStride<>>(v.data()+iG, m_sxcurl.
dimension(),Eigen::InnerStride<>(m_lie_algebra.
dimension()));
91 inline const Eigen::MatrixXd &
ScurlFace(
size_t iF)
const
93 return (*m_face_transfer_operators[iF]).serendipity;
97 inline const Eigen::MatrixXd &
EcurlFace(
size_t iF)
const
99 return (*m_face_transfer_operators[iF]).extension;
103 inline const Eigen::MatrixXd &
RcurlFace(
size_t iF)
const
105 return (*m_face_transfer_operators[iF]).reduction;
129 inline const Eigen::MatrixXd &
ScurlCell(
size_t iT)
const
131 return (*m_cell_transfer_operators[iT]).serendipity;
135 inline const Eigen::MatrixXd &
EcurlCell(
size_t iT)
const
137 return (*m_cell_transfer_operators[iT]).extension;
141 inline const Eigen::MatrixXd &
RcurlCell(
size_t iT)
const
143 return (*m_cell_transfer_operators[iT]).reduction;
168 inline const Eigen::MatrixXd
faceCurl(
size_t iF)
const
170 Eigen::MatrixXd
id = Eigen::MatrixXd::Identity(m_lie_algebra.
dimension(), m_lie_algebra.
dimension());
171 return Eigen::KroneckerProduct(m_sxcurl.
faceCurl(iF),
id);
183 Eigen::MatrixXd
id = Eigen::MatrixXd::Identity(m_lie_algebra.
dimension(), m_lie_algebra.
dimension());
184 return Eigen::KroneckerProduct(m_sxcurl.
facePotential(iF),
id);
194 inline const Eigen::MatrixXd
cellCurl(
size_t iT)
const
196 Eigen::MatrixXd
id = Eigen::MatrixXd::Identity(m_lie_algebra.
dimension(), m_lie_algebra.
dimension());
197 return Eigen::KroneckerProduct(m_sxcurl.
cellCurl(iT),
id);
209 Eigen::MatrixXd
id = Eigen::MatrixXd::Identity(m_lie_algebra.
dimension(), m_lie_algebra.
dimension());
210 return Eigen::KroneckerProduct(m_sxcurl.
cellPotential(iT),
id);
222 const double & penalty_factor = 1.,
223 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
227 return Eigen::KroneckerProduct(m_sxcurl.
computeL2Product(iT, penalty_factor, mass_Pk3_T, weight), m_lie_algebra.
massMatrix());
235 const std::string & side,
236 const double & penalty_factor = 1.,
237 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
253 return m_sxcurl.
cellBases(T.global_index());
265 return m_sxcurl.
faceBases(F.global_index());
277 return m_sxcurl.
edgeBases(E.global_index());
281 TransferOperators _compute_face_transfer_operators(
size_t iF);
282 TransferOperators _compute_cell_transfer_operators(
size_t iT);
288 std::vector<std::unique_ptr<TransferOperators> > m_face_transfer_operators;
289 std::vector<std::unique_ptr<TransferOperators> > m_cell_transfer_operators;
292 std::ostream & m_output;
Discrete Lie algebra valued Serendipity Hcurl space: local operators, L2 product and global interpola...
Definition: lasxcurl.hpp:18
Lie algebra class: mass matrix, structure constants and Lie bracket.
Definition: liealgebra.hpp:17
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
Class to describe a mesh.
Definition: MeshND.hpp:17
size_t dimension() const
Returns the dimension of the global space (all DOFs for all geometric entities)
Definition: variabledofspace.hpp:158
const Eigen::MatrixXd cellCurl(size_t iT) const
Return the full curl operator on the cell of index iT.
Definition: sxcurl.hpp:179
const Eigen::MatrixXd cellPotential(size_t iT) const
Return the potential operator on the cell of index iT.
Definition: sxcurl.hpp:191
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
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 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 Mesh & mesh() const
Return the mesh.
Definition: sxcurl.hpp:50
const SerendipityProblem & serPro() const
Definition: sxcurl.hpp:61
const Eigen::MatrixXd facePotential(size_t iF) const
Return the potential operator on the face of index iF.
Definition: sxcurl.hpp:167
bool use_threads
Definition: HHO_DiffAdvecReac.hpp:47
Eigen::MatrixXd extension
Definition: lasxcurl.hpp:40
const DDRCore::FaceBases & faceBases(const Face &F) const
Return cell bases for face F.
Definition: lasxcurl.hpp:263
const Mesh & mesh() const
Return the mesh.
Definition: lasxcurl.hpp:48
const Eigen::MatrixXd cellCurl(size_t iT) const
Return the full curl operator on the cell of index iT.
Definition: lasxcurl.hpp:194
const Eigen::MatrixXd facePotential(const Face &F) const
Return the potential operator on face F.
Definition: lasxcurl.hpp:188
const Eigen::MatrixXd & ScurlFace(size_t iF) const
Return the serendipity reconstruction for the face of index iF.
Definition: lasxcurl.hpp:91
const DDRCore::FaceBases & faceBases(size_t iF) const
Return face bases for the face of index iF.
Definition: lasxcurl.hpp:257
const Eigen::MatrixXd cellPotential(const Cell &T) const
Return the potential operator on cell T.
Definition: lasxcurl.hpp:214
const Eigen::MatrixXd & ScurlFace(const Face &F) const
Return the serendipity reconstruction for face F.
Definition: lasxcurl.hpp:109
std::vector< FunctionType > LAFunctionType
Definition: lasxcurl.hpp:21
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: lasxcurl.hpp:220
const size_t & degree() const
Return the polynomial degree.
Definition: lasxcurl.hpp:54
const Eigen::MatrixXd & massMatrix() const
Returns the Gram matrix of the Lie algebra.
Definition: liealgebra.hpp:40
const Eigen::MatrixXd & RcurlFace(const Face &F) const
Return cell reduction for cell T.
Definition: lasxcurl.hpp:121
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition: lasxcurl.hpp:251
const Eigen::MatrixXd & EcurlFace(const Face &F) const
Return the extension for face F.
Definition: lasxcurl.hpp:115
const Eigen::MatrixXd & EcurlFace(size_t iF) const
Return the extension for the face of index iF.
Definition: lasxcurl.hpp:97
const Eigen::MatrixXd cellCurl(const Cell &T) const
Return the full curl operator on cell T.
Definition: lasxcurl.hpp:201
const Eigen::MatrixXd facePotential(size_t iF) const
Return the potential operator on the face of index iF.
Definition: lasxcurl.hpp:181
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the face of index iT.
Definition: lasxcurl.hpp:245
Eigen::MatrixXd serendipity
Definition: lasxcurl.hpp:39
const LieAlgebra & lieAlg() const
Return the Lie algebra.
Definition: lasxcurl.hpp:60
const Eigen::VectorXd LaxcurlI(size_t iG, Eigen::VectorXd &v) const
Definition: lasxcurl.hpp:82
const DDRCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition: lasxcurl.hpp:275
const Eigen::MatrixXd & ScurlCell(size_t iT) const
Return the serendipity reconstruction for the cell of index iT.
Definition: lasxcurl.hpp:129
const Eigen::MatrixXd faceCurl(size_t iF) const
Return the full curl operator on the face of index iF.
Definition: lasxcurl.hpp:168
const size_t dimension() const
Assuming orthonormal basis.
Definition: liealgebra.hpp:28
const Eigen::MatrixXd & EcurlCell(const Cell &T) const
Return the extension for cell T.
Definition: lasxcurl.hpp:153
TransferOperators(const Eigen::MatrixXd &_serendipity, const Eigen::MatrixXd &_extension, const Eigen::MatrixXd &_reduction)
Definition: lasxcurl.hpp:27
Eigen::MatrixXd reduction
Definition: lasxcurl.hpp:41
const Eigen::MatrixXd & EcurlCell(size_t iT) const
Return the extension for the cell of index iT.
Definition: lasxcurl.hpp:135
const Eigen::MatrixXd & RcurlFace(size_t iF) const
Return the reduction for the face of index iF.
Definition: lasxcurl.hpp:103
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition: lasxcurl.hpp:269
const Eigen::MatrixXd cellPotential(size_t iT) const
Return the potential operator on the cell of index iT.
Definition: lasxcurl.hpp:207
const Eigen::MatrixXd & RcurlCell(size_t iT) const
Return the reduction for the cell of index iT.
Definition: lasxcurl.hpp:141
const SerendipityProblem & serPro() const
Return the serendipity operators.
Definition: lasxcurl.hpp:66
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: lasxcurl.cpp:100
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType
Definition: lasxcurl.hpp:20
const Eigen::MatrixXd faceCurl(const Face &F) const
Return the full curl operator on face F.
Definition: lasxcurl.hpp:175
Eigen::VectorXd interpolate(const LAFunctionType &v, const int doe_cell=-1, const int doe_face=-1, const int doe_edge=-1) const
Interpolator of a continuous function.
Definition: lasxcurl.cpp:61
const Eigen::MatrixXd & RcurlCell(const Cell &T) const
Return the reduction for cell T.
Definition: lasxcurl.hpp:159
LASXCurl(const LieAlgebra &lie_algebra, const SXCurl &sx_curl, bool use_threads=true, std::ostream &output=std::cout)
Constructor.
Definition: lasxcurl.cpp:10
const Eigen::MatrixXd & ScurlCell(const Cell &T) const
Return the serendipity reconstruction for cell T.
Definition: lasxcurl.hpp:147
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: lasxcurl.hpp:26