HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
hho-fullgradientdiff.hpp
Go to the documentation of this file.
1 // Implementation of the HHO scheme, with full gradient, for -div(K nabla u)=f with K scalar piecewise function.
2 //
3 // Author: Jerome Droniou (jerome.droniou@monash.edu)
4 //
5 
6 /*
7  *
8  * This implementation of HHO was developped following the principles described in
9  * Appendix B of the book
10  *
11  * The Hybrid High-Order Method for Polytopal Meshes: Design, Analysis, and Applications.
12  * D. A. Di Pietro and J. Droniou. Modeling, Simulation and Applications, vol. 19.
13  * Springer International Publishing, 2020, xxxi + 525p. doi: 10.1007/978-3-030-37203-3.
14  * url: https://hal.archives-ouvertes.fr/hal-02151813.
15  *
16  * If you use this code or part of it for a scientific publication, please cite the book
17  * above as a reference for the implementation.
18  *
19  */
20 
21 
22 #ifndef HHOFULLGRADIENTDIFF_HPP
23 #define HHOFULLGRADIENTDIFF_HPP
24 
25 #include <iostream>
26 
27 #include <boost/math/constants/constants.hpp>
28 
29 #include <Eigen/Sparse>
30 #include <unsupported/Eigen/SparseExtra>
31 
32 #include <mesh.hpp>
33 #include <mesh_builder.hpp>
35 
36 #include <hhospace.hpp>
38 #include <integralweight.hpp>
39 
45 namespace HArDCore3D
46 {
47 
55  {
56  typedef Eigen::SparseMatrix<double> SystemMatrixType;
57 
58  typedef std::function<double(const VectorRd &)> ForcingTermType;
59  typedef std::function<double(const VectorRd &)> SolutionType;
60  typedef std::function<VectorRd(const VectorRd &)> SolutionGradientType;
62 
65  const HHOSpace & hho_space,
66  const BoundaryConditions & BC,
67  bool use_threads,
68  std::ostream & output = std::cout
69  );
70 
73  const ForcingTermType & f,
74  const PermeabilityType & kappa,
75  const SolutionType & u,
76  Eigen::VectorXd & UDir
77  );
78 
80  inline size_t numSCDOFs() const
81  {
82  return m_hhospace.mesh().n_cells() * m_nloc_sc;
83  }
84 
86  inline size_t numDirDOFs() const
87  {
88  return m_BC.n_dir_faces()*m_hhospace.numLocalDofsFace();
89  }
90 
92  inline size_t numSkeletalDOFs() const
93  {
94  return m_hhospace.dimension() - numSCDOFs();
95  }
96 
98  inline size_t sizeSystem() const
99  {
100  return m_hhospace.dimension() - numSCDOFs() - numDirDOFs();
101  }
102 
104  inline const HHOSpace & hhospace() const
105  {
106  return m_hhospace;
107  }
108 
110  inline const Mesh & mesh() const
111  {
112  return m_hhospace.mesh();
113  }
114 
116  inline const SystemMatrixType & systemMatrix() const {
117  return m_A;
118  }
119 
122  return m_A;
123  }
124 
126  inline const Eigen::VectorXd & systemVector() const {
127  return m_b;
128  }
129 
131  inline Eigen::VectorXd & systemVector() {
132  return m_b;
133  }
134 
136  inline const double & stabilizationParameter() const {
137  return m_stab_par;
138  }
139 
141  inline double & stabilizationParameter() {
142  return m_stab_par;
143  }
144 
146  inline const SystemMatrixType & scMatrix() const {
147  return m_sc_A;
148  }
149 
151  inline Eigen::VectorXd & scVector() {
152  return m_sc_b;
153  }
154 
156  std::vector<double> computeEnergyNorms(
157  const std::vector<Eigen::VectorXd> & list_dofs
158  ) const;
159 
160  private:
162  std::pair<Eigen::MatrixXd, Eigen::VectorXd>
163  _compute_local_contribution(
164  size_t iT,
165  const ForcingTermType & f,
166  const PermeabilityType & kappa
167  );
168 
170  LocalStaticCondensation _compute_static_condensation(const size_t & iT) const;
171 
173  void _assemble_local_contribution(
174  size_t iT,
175  const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT,
176  std::list<Eigen::Triplet<double> > & A1,
177  Eigen::VectorXd & b1,
178  std::list<Eigen::Triplet<double> > & A2,
179  Eigen::VectorXd & b2
180  );
181 
182  const HHOSpace & m_hhospace;
183  const BoundaryConditions & m_BC;
184  bool m_use_threads;
185  std::ostream & m_output;
186  const size_t m_nloc_sc; // Number of statically condensed DOFs in each cell
187  SystemMatrixType m_A; // Matrix and RHS of statically condensed system
188  Eigen::VectorXd m_b;
189  SystemMatrixType m_sc_A; // Static condensation operator and RHS (to recover statically condensed DOFs)
190  Eigen::VectorXd m_sc_b;
191  double m_stab_par;
192  };
193 
194  //------------------------------------------------------------------------------
195  // Exact solutions
196  //------------------------------------------------------------------------------
197 
198  static const double PI = boost::math::constants::pi<double>();
199  using std::sin;
200  using std::cos;
201 
202  //------------------------------------------------------------------------------
203 
204  static const VectorRd vec_a = VectorRd(1.,2.,3.);
205 
207  linear_u = [](const VectorRd & x) -> double {
208  return vec_a.dot(x) + 4.;
209  };
210 
212  linear_gradu = [](const VectorRd & x) -> VectorRd {
213  return vec_a;
214  };
215 
217  linear_f = [](const VectorRd & x) -> double {
218  return 0.;
219  };
220 
223 
224 
225  //------------------------------------------------------------------------------
226 
228  trigonometric_u = [](const VectorRd & x) -> double {
229  return sin(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2));
230  };
231 
233  trigonometric_gradu = [](const VectorRd & x) -> VectorRd {
234  return PI * VectorRd(
235  cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
236  sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
237  sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
238  );
239  };
240 
242  trigonometric_f = [](const VectorRd & x) -> double {
243  return 3. * std::pow(PI, 2) * trigonometric_u(x);
244  };
245 
248 
249 
250 } // end of namespace HArDCore3D
251 
252 #endif
The BoundaryConditions class provides definition of boundary conditions.
Definition: BoundaryConditions.hpp:43
const size_t n_dir_faces() const
Returns the number of Dirichlet faces.
Definition: BoundaryConditions.hpp:63
Class definition: polynomial bases and operators.
Definition: hhospace.hpp:53
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: localdofspace.hpp:80
size_t numLocalDofsFace() const
Returns the number of local face DOFs.
Definition: localdofspace.hpp:49
static const double PI
Definition: ddr-magnetostatics.hpp:186
static Magnetostatics::SolutionPotentialType linear_u
Definition: ddr-magnetostatics.hpp:213
static Magnetostatics::SolutionPotentialType trigonometric_u
Definition: ddr-magnetostatics.hpp:237
static Magnetostatics::ForcingTermType trigonometric_f
Definition: ddr-magnetostatics.hpp:255
static Magnetostatics::ForcingTermType linear_f
Definition: ddr-magnetostatics.hpp:227
const Mesh & mesh() const
Return a const reference to the mesh.
Definition: hhospace.hpp:114
bool use_threads
Definition: HHO_DiffAdvecReac.hpp:47
static Brinkman::VelocityGradientType linear_gradu
Definition: hho-brinkman.hpp:369
static Brinkman::VelocityGradientType trigonometric_gradu
Definition: hho-brinkman.hpp:514
std::function< double(const VectorRd &)> SolutionType
Definition: hho-fullgradientdiff.hpp:59
Eigen::VectorXd & systemVector()
Returns the linear system right-hand side vector.
Definition: hho-fullgradientdiff.hpp:131
std::function< double(const VectorRd &)> ForcingTermType
Definition: hho-fullgradientdiff.hpp:58
const Mesh & mesh() const
Returns the mesh.
Definition: hho-fullgradientdiff.hpp:110
std::vector< double > computeEnergyNorms(const std::vector< Eigen::VectorXd > &list_dofs) const
Compute the discrete energy norm of a family of vectors representing the dofs.
Definition: hho-fullgradientdiff.cpp:452
Eigen::SparseMatrix< double > SystemMatrixType
Definition: hho-fullgradientdiff.hpp:56
const Eigen::VectorXd & systemVector() const
Returns the linear system right-hand side vector.
Definition: hho-fullgradientdiff.hpp:126
const SystemMatrixType & systemMatrix() const
Returns the linear system matrix.
Definition: hho-fullgradientdiff.hpp:116
size_t numSkeletalDOFs() const
Returns the number of DOFs after SC but before eliminating Dirichlet DOFs.
Definition: hho-fullgradientdiff.hpp:92
static const VectorRd vec_a
Definition: hho-fullgradientdiff.hpp:204
const HHOSpace & hhospace() const
Returns the HHO space.
Definition: hho-fullgradientdiff.hpp:104
static FullGradientDiffusion::PermeabilityType trigonometric_kappa
Definition: hho-fullgradientdiff.hpp:247
size_t sizeSystem() const
Returns the size of the final system, after application of SC and removal of Dirichlet BCs.
Definition: hho-fullgradientdiff.hpp:98
FullGradientDiffusion(const HHOSpace &hho_space, const BoundaryConditions &BC, bool use_threads, std::ostream &output=std::cout)
Constructor.
Definition: hho-fullgradientdiff.cpp:264
size_t numDirDOFs() const
Returns the number of Dirichlet DOFs.
Definition: hho-fullgradientdiff.hpp:86
size_t numSCDOFs() const
Returns the number of statically condensed DOFs (the cell DOFs)
Definition: hho-fullgradientdiff.hpp:80
Eigen::VectorXd & scVector()
Returns the static condensation rhs.
Definition: hho-fullgradientdiff.hpp:151
IntegralWeight PermeabilityType
Definition: hho-fullgradientdiff.hpp:61
double & stabilizationParameter()
Returns the stabilization parameter.
Definition: hho-fullgradientdiff.hpp:141
const SystemMatrixType & scMatrix() const
Returns the static condensation recovery operator.
Definition: hho-fullgradientdiff.hpp:146
static FullGradientDiffusion::PermeabilityType linear_kappa
Definition: hho-fullgradientdiff.hpp:222
SystemMatrixType & systemMatrix()
Returns the linear system matrix.
Definition: hho-fullgradientdiff.hpp:121
std::function< VectorRd(const VectorRd &)> SolutionGradientType
Definition: hho-fullgradientdiff.hpp:60
const double & stabilizationParameter() const
Returns the stabilization parameter (scaling)
Definition: hho-fullgradientdiff.hpp:136
void assembleLinearSystem(const ForcingTermType &f, const PermeabilityType &kappa, const SolutionType &u, Eigen::VectorXd &UDir)
Assemble the global system
Definition: hho-fullgradientdiff.cpp:289
std::size_t n_cells() const
number of cells in the mesh.
Definition: MeshND.hpp:60
Definition: ddr-magnetostatics.hpp:40
MeshND::VectorRd< 2 > VectorRd
Definition: Mesh2D.hpp:14
Assemble a diffusion problem.
Definition: hho-fullgradientdiff.hpp:55
Structure for weights (scalar, at the moment) in integral.
Definition: integralweight.hpp:36
Structure to store information for, and perform, local static condensation.
Definition: local_static_condensation.hpp:25