HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
sxcurl.hpp
Go to the documentation of this file.
1 #ifndef SXCURL_HPP
2 #define SXCURL_HPP
3 
4 #include <variabledofspace.hpp>
5 #include <ddrcore.hpp>
6 #include <integralweight.hpp>
7 #include <xcurl.hpp>
8 #include <sxgrad.hpp>
10 
11 namespace HArDCore3D
12 {
19 
20  class SXCurl : public VariableDOFSpace
21  {
22  public:
23  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType;
24 
26 
28  {
30  const Eigen::MatrixXd & _serendipity,
31  const Eigen::MatrixXd & _extension,
32  const Eigen::MatrixXd & _reduction
33  )
34  : serendipity(_serendipity),
35  extension(_extension),
36  reduction(_reduction)
37  {
38  // Do nothing
39  }
40 
41  Eigen::MatrixXd serendipity;
42  Eigen::MatrixXd extension;
43  Eigen::MatrixXd reduction;
44  };
45 
47  SXCurl(const DDRCore & ddr_core, const SerendipityProblem & ser_pro, bool use_threads = true, std::ostream & output = std::cout);
48 
50  const Mesh & mesh() const
51  {
52  return m_ddr_core.mesh();
53  }
54 
56  const size_t & degree() const
57  {
58  return m_ddr_core.degree();
59  }
60 
61  const SerendipityProblem & serPro() const
62  {
63  return m_ser_pro;
64  }
65 
67  Eigen::VectorXd interpolate(
68  const FunctionType & v,
69  const int doe_cell = -1,
70  const int doe_face = -1,
71  const int doe_edge = -1
72  ) const;
73 
74  //---------------------------------//
75  //---- Face transfer operators-----//
76 
78  inline const Eigen::MatrixXd & ScurlFace(size_t iF) const
79  {
80  return (*m_face_transfer_operators[iF]).serendipity;
81  }
82 
84  inline const Eigen::MatrixXd & EcurlFace(size_t iF) const
85  {
86  return (*m_face_transfer_operators[iF]).extension;
87  }
88 
90  inline const Eigen::MatrixXd & RcurlFace(size_t iF) const
91  {
92  return (*m_face_transfer_operators[iF]).reduction;
93  }
94 
96  inline const Eigen::MatrixXd & ScurlFace(const Face & F) const
97  {
98  return ScurlFace(F.global_index());
99  }
100 
102  inline const Eigen::MatrixXd & EcurlFace(const Face & F) const
103  {
104  return EcurlFace(F.global_index());
105  }
106 
108  inline const Eigen::MatrixXd & RcurlFace(const Face & F) const
109  {
110  return RcurlFace(F.global_index());
111  }
112 
113  //----------------------------------//
114  //---- Cell transfer operators -----//
116  inline const Eigen::MatrixXd & ScurlCell(size_t iT) const
117  {
118  return (*m_cell_transfer_operators[iT]).serendipity;
119  }
120 
122  inline const Eigen::MatrixXd & EcurlCell(size_t iT) const
123  {
124  return (*m_cell_transfer_operators[iT]).extension;
125  }
126 
128  inline const Eigen::MatrixXd & RcurlCell(size_t iT) const
129  {
130  return (*m_cell_transfer_operators[iT]).reduction;
131  }
132 
134  inline const Eigen::MatrixXd & ScurlCell(const Cell & T) const
135  {
136  return ScurlCell(T.global_index());
137  }
138 
140  inline const Eigen::MatrixXd & EcurlCell(const Cell & T) const
141  {
142  return EcurlCell(T.global_index());
143  }
144 
146  inline const Eigen::MatrixXd & RcurlCell(const Cell & T) const
147  {
148  return RcurlCell(T.global_index());
149  }
150 
151  //-----------------------------------------------------------------//
152  //---- Full curl and potential reconstructions, and L2 product ----//
153 
155  inline const Eigen::MatrixXd faceCurl(size_t iF) const
156  {
157  return m_xcurl.faceOperators(iF).curl * EcurlFace(iF);
158  }
159 
161  inline const Eigen::MatrixXd faceCurl(const Face & F) const
162  {
163  return faceCurl(F.global_index());
164  }
165 
167  inline const Eigen::MatrixXd facePotential(size_t iF) const
168  {
169  return m_xcurl.faceOperators(iF).potential * EcurlFace(iF);
170  }
171 
173  inline const Eigen::MatrixXd facePotential(const Face & F) const
174  {
175  return facePotential(F.global_index());
176  }
177 
179  inline const Eigen::MatrixXd cellCurl(size_t iT) const
180  {
181  return m_xcurl.cellOperators(iT).curl * EcurlCell(iT);
182  }
183 
185  inline const Eigen::MatrixXd cellCurl(const Cell & T) const
186  {
187  return cellCurl(T.global_index());
188  }
189 
191  inline const Eigen::MatrixXd cellPotential(size_t iT) const
192  {
193  return m_xcurl.cellOperators(iT).potential * EcurlCell(iT);
194  }
195 
197  inline const Eigen::MatrixXd cellPotential(const Cell & T) const
198  {
199  return cellPotential(T.global_index());
200  }
201 
203  Eigen::MatrixXd computeL2Product(
204  const size_t iT,
205  const double & penalty_factor = 1.,
206  const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
207  const IntegralWeight & weight = IntegralWeight(1.)
208  ) const
209  {
210  return EcurlCell(iT).transpose()
211  * m_xcurl.computeL2Product(iT, penalty_factor, mass_Pk3_T, weight)
212  *EcurlCell(iT);
213  }
214 
215 
217  Eigen::MatrixXd computeL2ProductGradient(
218  const size_t iT,
219  const SXGrad & sx_grad,
220  const std::string & side,
221  const double & penalty_factor = 1.,
222  const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
223  const IntegralWeight & weight = IntegralWeight(1.)
224  ) const;
225 
226  //-----------------------//
227  //---- Getters ----------//
228 
230  inline const DDRCore::CellBases & cellBases(size_t iT) const
231  {
232  return m_ddr_core.cellBases(iT);
233  }
234 
236  inline const DDRCore::CellBases & cellBases(const Cell & T) const
237  {
238  return m_ddr_core.cellBases(T.global_index());
239  }
240 
242  inline const DDRCore::FaceBases & faceBases(size_t iF) const
243  {
244  return m_ddr_core.faceBases(iF);
245  }
246 
248  inline const DDRCore::FaceBases & faceBases(const Face & F) const
249  {
250  return m_ddr_core.faceBases(F.global_index());
251  }
252 
254  inline const DDRCore::EdgeBases & edgeBases(size_t iE) const
255  {
256  return m_ddr_core.edgeBases(iE);
257  }
258 
260  inline const DDRCore::EdgeBases & edgeBases(const Edge & E) const
261  {
262  return m_ddr_core.edgeBases(E.global_index());
263  }
264 
265  private:
266  TransferOperators _compute_face_transfer_operators(size_t iF);
267  TransferOperators _compute_cell_transfer_operators(size_t iT);
268 
269  const DDRCore & m_ddr_core;
270  const SerendipityProblem & m_ser_pro;
271  const XCurl m_xcurl;
272 
273  // Containers for serendipity, extension and reduction operators
274  std::vector<std::unique_ptr<TransferOperators> > m_face_transfer_operators;
275  std::vector<std::unique_ptr<TransferOperators> > m_cell_transfer_operators;
276 
277  bool m_use_threads;
278  std::ostream & m_output;
279 
280  };
281 
282 } // end of namespace HArDCore3D
283 #endif
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