HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
vdiv.hpp
Go to the documentation of this file.
1 #ifndef VDIV_HPP
2 #define VDIV_HPP
3 
4 #include <globaldofspace.hpp>
5 #include <vemcore.hpp>
6 #include <integralweight.hpp>
7 
8 namespace HArDCore3D
9 {
16 
25  class VDiv : public GlobalDOFSpace
26  {
27  public:
28  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType;
29  typedef std::function<double(const Eigen::Vector3d &)> DivergenceType;
30 
33  {
35  const Eigen::MatrixXd & _proj_divergence,
36  const Eigen::MatrixXd & _proj_function
37  )
38  : proj_divergence(_proj_divergence),
39  proj_function(_proj_function)
40  {
41  // Do nothing
42  }
43 
44  Eigen::MatrixXd proj_divergence;
45  Eigen::MatrixXd proj_function;
46  };
47 
49  VDiv(const VEMCore & vem_core, bool use_threads = true, std::ostream & output = std::cout);
50 
52  const Mesh & mesh() const
53  {
54  return m_vem_core.mesh();
55  }
56 
58  const size_t & degree() const
59  {
60  return m_vem_core.degree();
61  }
62 
64  Eigen::VectorXd interpolate(
65  const FunctionType & v,
66  const DivergenceType & div_v,
67  const int doe_cell = -1,
68  const int doe_face = -1
69  ) const;
70 
72  inline const LocalOperators & cellOperators(size_t iT) const
73  {
74  return *m_cell_operators[iT];
75  }
76 
78  inline const LocalOperators & cellOperators(const Cell & T) const
79  {
80  return * m_cell_operators[T.global_index()];
81  }
82 
84  inline const LocalOperators & faceOperators(size_t iF) const
85  {
86  return *m_face_operators[iF];
87  }
88 
90  inline const LocalOperators & faceOperators(const Face & F) const
91  {
92  return * m_face_operators[F.global_index()];
93  }
94 
96  inline const VEMCore::CellBases & cellBases(size_t iT) const
97  {
98  return m_vem_core.cellBases(iT);
99  }
100 
102  inline const VEMCore::CellBases & cellBases(const Cell & T) const
103  {
104  return m_vem_core.cellBases(T.global_index());
105  }
106 
108  inline const VEMCore::FaceBases & faceBases(size_t iF) const
109  {
110  return m_vem_core.faceBases(iF);
111  }
112 
114  inline const VEMCore::FaceBases & faceBases(const Face & F) const
115  {
116  return m_vem_core.faceBases(F.global_index());
117  }
118 
120  inline const VEMCore::EdgeBases & edgeBases(size_t iE) const
121  {
122  return m_vem_core.edgeBases(iE);
123  }
124 
126  inline const VEMCore::EdgeBases & edgeBases(const Edge & E) const
127  {
128  return m_vem_core.edgeBases(E.global_index());
129  }
130 
132  // The mass matrix of P^k(T)^3 is the most expensive mass matrix in the calculation of this norm, which
133  // is why there's the option of passing it as parameter if it's been already pre-computed when the norm is called.
134  Eigen::MatrixXd computeL2Product(
135  const size_t iT,
136  const double & penalty_factor = 1.,
137  const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
138  const IntegralWeight & weight = IntegralWeight(1.)
139  ) const;
140 
141  private:
142  LocalOperators _compute_cell_operators(size_t iT);
143 
144  const VEMCore & m_vem_core;
145  bool m_use_threads;
146  std::ostream & m_output;
147 
148  // Containers for local operators
149  std::vector<std::unique_ptr<LocalOperators> > m_cell_operators;
150  std::vector<std::unique_ptr<LocalOperators> > m_face_operators;
151 
152  };
153 
154 } // end of namespace HArDCore3D
155 #endif
Base class for global DOF spaces. Provides functions to manipulate global DOFs (the local version bei...
Definition: globaldofspace.hpp:16
Virtual Hdiv space: local operators, L2 product and global interpolator.
Definition: vdiv.hpp:26
Construct all polynomial spaces for the VEM sequence.
Definition: vemcore.hpp:39
Class to describe a mesh.
Definition: MeshND.hpp:17
bool use_threads
Definition: HHO_DiffAdvecReac.hpp:47
const Mesh & mesh() const
Return the mesh.
Definition: vdiv.hpp:52
Eigen::VectorXd interpolate(const FunctionType &v, const DivergenceType &div_v, const int doe_cell=-1, const int doe_face=-1) const
Interpolator of a continuous function.
Definition: vdiv.cpp:52
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition: vemcore.hpp:129
const VEMCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition: vdiv.hpp:102
const LocalOperators & faceOperators(const Face &F) const
Return face operators for face F.
Definition: vdiv.hpp:90
Eigen::MatrixXd proj_function
Definition: vdiv.hpp:45
const Mesh & mesh() const
Return a const reference to the mesh.
Definition: vemcore.hpp:117
const VEMCore::CellBases & cellBases(size_t iT) const
Return cell bases for the face of index iT.
Definition: vdiv.hpp:96
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: vdiv.cpp:207
const FaceBases & faceBases(size_t iF) const
Return face bases for face iF.
Definition: vemcore.hpp:137
VDiv(const VEMCore &vem_core, bool use_threads=true, std::ostream &output=std::cout)
Constructor.
Definition: vdiv.cpp:13
const size_t & degree() const
Return the polynomial degree.
Definition: vemcore.hpp:123
Eigen::MatrixXd proj_divergence
Definition: vdiv.hpp:44
const LocalOperators & cellOperators(const Cell &T) const
Return cell operators for cell T.
Definition: vdiv.hpp:78
const LocalOperators & faceOperators(size_t iF) const
Return face operators for the face of index iF.
Definition: vdiv.hpp:84
std::function< double(const Eigen::Vector3d &)> DivergenceType
Definition: vdiv.hpp:29
const size_t & degree() const
Return the polynomial degree.
Definition: vdiv.hpp:58
LocalOperators(const Eigen::MatrixXd &_proj_divergence, const Eigen::MatrixXd &_proj_function)
Definition: vdiv.hpp:34
const VEMCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition: vdiv.hpp:126
const VEMCore::FaceBases & faceBases(const Face &F) const
Return cell bases for face F.
Definition: vdiv.hpp:114
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType
Definition: vdiv.hpp:28
const LocalOperators & cellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition: vdiv.hpp:72
const VEMCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition: vdiv.hpp:120
const EdgeBases & edgeBases(size_t iE) const
Return edge bases for edge iE.
Definition: vemcore.hpp:145
const VEMCore::FaceBases & faceBases(size_t iF) const
Return face bases for the face of index iF.
Definition: vdiv.hpp:108
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 for weights (scalar, at the moment) in integral.
Definition: integralweight.hpp:36
A structure to store the local operators (projections of div and function)
Definition: vdiv.hpp:33
Structure to store element bases.
Definition: vemcore.hpp:63
Structure to store edge bases.
Definition: vemcore.hpp:104
Structure to store face bases.
Definition: vemcore.hpp:84