HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
hho-brinkman.hpp
Go to the documentation of this file.
1 // Implementation of the HHO scheme, with full gradient, for the Brinkman equations -div(mu grad u) + nu u + grad p=f, div u = g.
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 HHOBRINKMAN_HPP
23 #define HHOBRINKMAN_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 <vhhospace.hpp>
38 #include <integralweight.hpp>
39 
45 namespace HArDCore3D
46 {
47 
53  //------------------------------------------------------------------------------
54  // Brinkman model
55  //------------------------------------------------------------------------------
56 
59  {
62 
65  const ViscosityType & _mu,
66  const PermeabilityInvType & _kappainv
67  )
68  : mu(_mu),
69  kappainv(_kappainv)
70  {
71  // do nothing
72  }
73 
75  double Cf(const Cell & T) const
76  {
77  double muT = mu.value(T, T.center_mass());
78  double kappainvT = kappainv.value(T, T.center_mass());
79 
80  // Thres is the level at which the parameter is put at 0; note that Cf is also bounded below by thres
81  // to avoid singularities when taking Cf^{-1}.
82  double thres = 1e-14;
83  return (muT > thres ? std::max(thres, kappainvT * std::pow(T.diam(), 2) / muT) : 1./thres);
84  }
85 
88  };
89 
92  {
94  BrinkmanNorms(double norm_u, double norm_p):
95  u(norm_u),
96  p(norm_p),
97  energy( std::sqrt( std::pow(norm_u, 2) + std::pow(norm_p, 2) ) )
98  {
99  // do nothing
100  };
101 
102  double u;
103  double p;
104  double energy;
105  };
106 
107 
109 
117  struct Brinkman
118  {
119  typedef Eigen::SparseMatrix<double> SystemMatrixType;
120 
121  typedef std::function<VectorRd(const VectorRd &)> MomentumForcingTermType;
122  typedef std::function<double(const VectorRd &)> CompressibilityForcingTermType;
123  typedef std::function<VectorRd(const VectorRd &)> VelocityType;
124  typedef std::function<MatrixRd(const VectorRd &)> VelocityGradientType;
125  typedef std::function<double(const VectorRd &)> PressureType;
126  typedef std::function<VectorRd(const VectorRd &)> PressureGradientType;
127 
129  Brinkman(
130  const VHHOSpace & vhho_space,
131  const GlobalDOFSpace & p_space,
132  const BoundaryConditions & BC,
133  bool use_threads,
134  std::ostream & output = std::cout
135  );
136 
139  const MomentumForcingTermType & f,
141  const BrinkmanParameters & para,
142  const VelocityType & u,
143  Eigen::VectorXd & UDir
144  );
145 
147  inline size_t nloc_sc_u() const
148  {
149  return m_nloc_sc_u;
150  }
151 
153  inline size_t nloc_sc_p() const
154  {
155  return m_nloc_sc_p;
156  }
157 
159  inline size_t numSCDOFs_u() const
160  {
161  return mesh().n_cells() * m_nloc_sc_u;
162  }
163 
165  inline size_t numSCDOFs_p() const
166  {
167  return mesh().n_cells() * m_nloc_sc_p;
168  }
169 
171  inline size_t numSCDOFs() const
172  {
173  return numSCDOFs_u() + numSCDOFs_p();
174  }
175 
177  inline size_t numDirDOFs() const
178  {
179  return m_BC.n_dir_faces()*m_vhhospace.numLocalDofsFace();
180  }
181 
183  inline size_t dimVelocity() const
184  {
185  return m_vhhospace.dimension();
186  }
187 
189  inline size_t dimPressure() const
190  {
191  return m_pspace.dimension();
192  }
193 
195  inline size_t numNonSCDOFs() const
196  {
197  return dimVelocity() + dimPressure() - numSCDOFs() + 1;
198  }
199 
201  inline size_t sizeSystem() const
202  {
203  return numNonSCDOFs() - numDirDOFs();
204  }
205 
207  inline const VHHOSpace & vhhospace() const
208  {
209  return m_vhhospace;
210  }
211 
213  inline const GlobalDOFSpace & pspace() const
214  {
215  return m_pspace;
216  }
217 
219  inline const Mesh & mesh() const
220  {
221  return m_vhhospace.mesh();
222  }
223 
225  inline const SystemMatrixType & systemMatrix() const {
226  return m_A;
227  }
228 
231  return m_A;
232  }
233 
235  inline const Eigen::VectorXd & systemVector() const {
236  return m_b;
237  }
238 
240  inline Eigen::VectorXd & systemVector() {
241  return m_b;
242  }
243 
245  inline const std::pair<double,double> & stabilizationParameter() const {
246  return m_stab_par;
247  }
248 
250  inline std::pair<double,double> & stabilizationParameter() {
251  return m_stab_par;
252  }
253 
255  inline const SystemMatrixType & scMatrix() const {
256  return m_sc_A;
257  }
258 
260  inline Eigen::VectorXd & scVector() {
261  return m_sc_b;
262  }
263 
265  Eigen::VectorXd interpolate(
266  const VelocityType & u,
267  const PressureType & p,
268  const int doe_cell = -1,
269  const int doe_face = -1
270  ) const;
271 
273  std::vector<BrinkmanNorms> computeBrinkmanNorms(
274  const std::vector<Eigen::VectorXd> & list_dofs,
275  const BrinkmanParameters & para
276  ) const;
277 
279  std::vector<double> pressureVertexValues(
280  const Eigen::VectorXd & p
281  ) const;
282 
284  std::pair<double, double> computeConditionNum() const;
285 
287  std::pair<double,double> computeFlux(
288  const Eigen::VectorXd & u,
289  const std::pair<std::vector<VectorRd>, VectorRd> & surf
290  ) const;
291 
292  private:
294  std::pair<Eigen::MatrixXd, Eigen::VectorXd>
295  _compute_local_contribution(
296  size_t iT,
297  const MomentumForcingTermType & f,
299  const BrinkmanParameters & para
300  );
301 
303  LocalStaticCondensation _compute_static_condensation(const size_t & iT) const;
304 
306  void _assemble_local_contribution(
307  size_t iT,
308  const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT,
309  std::list<Eigen::Triplet<double> > & A1,
310  Eigen::VectorXd & b1,
311  std::list<Eigen::Triplet<double> > & A2,
312  Eigen::VectorXd & b2
313  );
314 
316  const VHHOSpace & m_vhhospace;
317  const GlobalDOFSpace & m_pspace;
318 
319  const BoundaryConditions & m_BC;
320  bool m_use_threads;
321  std::ostream & m_output;
322  const size_t m_nloc_sc_u; // Number of statically condensed velocity DOFs in each cell
323  const size_t m_nloc_sc_p; // Number of statically condensed pressure DOFs in each cell
324  SystemMatrixType m_A; // Matrix and RHS of statically condensed system
325  Eigen::VectorXd m_b;
326  SystemMatrixType m_sc_A; // Static condensation operator and RHS (to recover statically condensed DOFs)
327  Eigen::VectorXd m_sc_b;
328  std::pair<double,double> m_stab_par; // Stabilisation parameters: first for Stokes, second for Darcy
329  };
330 
331 
332  //------------------------------------------------------------------------------
333  // Exact solutions
334  //------------------------------------------------------------------------------
335 
336  static const double PI = boost::math::constants::pi<double>();
337  using std::sin;
338  using std::cos;
339 
340  double scaling_mu = 1.;
341  double scaling_kappainv = 1.;
342 
343  // Threshold for checking that quantities are zero or not
344  double constexpr eps=1e-14;
345 
346  //---------------- Linear -----------------------------------------------------------
347 
348 
351 
354 
355  static const MatrixRd mat_u = (MatrixRd() << 1.,2.,3., 0.,-1.,1. , 0.,1.,1.).finished();
356  static const VectorRd vec_p = VectorRd(-1., 2., -5.);
357 
359  linear_u = [](const VectorRd & x) -> VectorRd {
360  return mat_u * x;
361  };
362 
364  linear_p = [](const VectorRd & x) -> double {
365  return vec_p.dot(x - VectorRd(.5,.5,.5)); // Assuming domain (0,1)^3 to have zero average
366  };
367 
369  linear_gradu = [](const VectorRd & x) -> MatrixRd {
370  return mat_u;
371  };
372 
374  linear_gradp = [](const VectorRd & x) -> VectorRd {
375  return vec_p;
376  };
377 
379  linear_f = [](const VectorRd & x) -> VectorRd {
380  return VectorRd::Zero() + scaling_kappainv*linear_u(x) + linear_gradp(x);
381  };
382 
384  linear_g = [](const VectorRd & x)->double { return linear_gradu(x).trace();};
385 
386 
387  //---------------- Quadratic ------------------------------------------------
390 
393 
395  quadratic_u = [](const VectorRd & x) -> VectorRd {
396  return VectorRd(
397  std::pow(x(0), 2),
398  std::pow(x(1), 2),
399  std::pow(x(2), 2)
400  );
401  };
402 
404  quadratic_gradu = [](const VectorRd & x) -> MatrixRd {
405  MatrixRd M = MatrixRd::Zero();
406  M.row(0) << 2. * x(0),
407  0.,
408  0.;
409  M.row(1) << 0.,
410  2 * x(1),
411  0.;
412  M.row(2) << 0.,
413  0.,
414  2 * x(2);
415 
416  return M;
417  };
418 
420  quadratic_p = [](const VectorRd & x) -> double {
421  return std::pow(x(0), 2) + x(1)*x(2) - 1./3. - 1./4.;
422  };
423 
425  quadratic_gradp = [](const VectorRd & x) -> VectorRd {
426  return VectorRd( 2 * x(0),
427  x(2),
428  x(1)
429  );
430  };
431 
433  quadratic_f = [](const VectorRd & x) -> VectorRd {
434  return - scaling_mu * VectorRd( 2., 2., 2.) + scaling_kappainv * quadratic_u(x) + quadratic_gradp(x);
435  };
436 
438  quadratic_g = [](const VectorRd & x)->double { return quadratic_gradu(x).trace();};
439 
440 
441  //---------------- Cubic ------------------------------------------------
444 
447 
449  cubic_u = [](const VectorRd & x) -> VectorRd {
450  return VectorRd(
451  std::pow(x(0), 3) + std::pow(x(0), 2)*x(1),
452  std::pow(x(1), 2) + x(0) * x(1) * x(2),
453  std::pow(x(2), 3)
454  );
455  };
456 
458  cubic_gradu = [](const VectorRd & x) -> MatrixRd {
459  MatrixRd M = MatrixRd::Zero();
460  M.row(0) << 3. * std::pow(x(0), 2) + 2 * x(0) * x(1),
461  std::pow(x(0), 2),
462  0.;
463  M.row(1) << x(1) * x(2),
464  2 * x(1) + x(0) * x(2),
465  x(0) * x(1);
466  M.row(2) << 0.,
467  0.,
468  3 * std::pow(x(2), 2);
469 
470  return M;
471  };
472 
474  cubic_p = [](const VectorRd & x) -> double {
475  return std::pow(x(0), 3) + x(1)*x(2) - 1./4. - 1./4.;
476  };
477 
479  cubic_gradp = [](const VectorRd & x) -> VectorRd {
480  return VectorRd( 3 * std::pow(x(0), 2),
481  x(2),
482  x(1)
483  );
484  };
485 
487  cubic_f = [](const VectorRd & x) -> VectorRd {
488  return - scaling_mu * VectorRd( 6 * x(0) + 2 * x(1), 2., 6 * x(2)) + scaling_kappainv * cubic_u(x) + cubic_gradp(x);
489  };
490 
492  cubic_g = [](const VectorRd & x)->double { return cubic_gradu(x).trace();};
493 
494 
495  //---------------- Trigonometric ------------------------------------------------
498 
501 
502  double constexpr pressure_scaling = 1.;
503 
505  trigonometric_u = [](const VectorRd & x) -> VectorRd {
506  return VectorRd(
507  0.5 * sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2)),
508  0.5 * cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
509  -cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2))
510  );
511  };
512 
514  trigonometric_gradu = [](const VectorRd & x) -> MatrixRd {
515  MatrixRd M = MatrixRd::Zero();
516  M.row(0) << PI * cos(2 * PI * x(0)) * cos(2 * PI * x(1)) * cos(2 * PI * x(2)),
517  -PI * sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
518  -PI * sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2));
519  M.row(1) << - PI * sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
520  PI * cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2)),
521  -PI * cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2));
522  M.row(2) << 2 * PI * sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
523  2 * PI * cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
524  - 2 * PI * cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2));
525 
526  return M;
527  };
528 
530  trigonometric_p = [](const VectorRd & x) -> double {
531  return pressure_scaling * sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2));
532  };
533 
535  trigonometric_gradp = [](const VectorRd & x) -> VectorRd {
536  return pressure_scaling * 2. * PI * VectorRd(
537  cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
538  sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
539  sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2))
540  );
541  };
542 
544  trigonometric_f = [](const VectorRd & x) -> VectorRd {
545  return scaling_mu * 6. * std::pow(PI, 2) * VectorRd(
546  sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2)),
547  cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
548  -2. * cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2))
549  )
551  };
552 
554  trigonometric_g = [](const VectorRd & x)->double { return trigonometric_gradu(x).trace();};
555 
556 
557 
558  //---------------- Various regimes parametrised by global friction coefficient --------------------------------------------
559 
562 
565 
568 
571 
572  // u_D = - kappainv^{-1} grad p
574  regimes_uD = [](const VectorRd & x) -> VectorRd {
575  VectorRd uDvalue = VectorRd::Zero();
576  if (scaling_kappainv > eps) {
577  uDvalue = - std::pow(scaling_kappainv, -1) * regimes_gradp(x) ;
578  }
579  return uDvalue;
580  };
581 
583  regimes_graduD = [](const VectorRd & x) -> MatrixRd {
584  MatrixRd M = MatrixRd::Zero();
585 
586  if (scaling_kappainv>eps){
587  M.row(0) << sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
588  -cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
589  -cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2));
590  M.row(1) << -cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2)),
591  sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2)),
592  -sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2));
593  M.row(2) << -cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
594  -sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2)),
595  sin(2. * PI * x(0)) * sin(2. * PI * x(1)) * sin(2. * PI * x(2));
596  M *= 4. * std::pow(PI, 2) / scaling_kappainv;
597  }
598 
599  return M;
600  };
601 
602  // u_S = trigonometric u
605 
606  static std::function<double(const VectorRd &)>
607  XiS = [](const VectorRd & x) -> double {
608  return (scaling_mu > eps ? std::exp(-scaling_kappainv/scaling_mu) : 0.);
609  };
610 
611  static Brinkman::VelocityType
612  regimes_u = [](const VectorRd & x) -> VectorRd { return XiS(x) * regimes_uS(x) + (1.-XiS(x)) * regimes_uD(x);};
613 
615  regimes_gradu = [](const VectorRd & x) -> MatrixRd { return XiS(x) * regimes_graduS(x) + (1.-XiS(x)) * regimes_graduD(x);};
616 
618  regimes_f = [](const VectorRd & x) -> VectorRd {
619  VectorRd fvalue = XiS(x) * scaling_mu * 6. * std::pow(PI, 2) * VectorRd(
620  sin(2. * PI * x(0)) * cos(2. * PI * x(1)) * cos(2. * PI * x(2)),
621  cos(2. * PI * x(0)) * sin(2. * PI * x(1)) * cos(2. * PI * x(2)),
622  -2. * cos(2. * PI * x(0)) * cos(2. * PI * x(1)) * sin(2. * PI * x(2))
623  )
624  + regimes_gradp(x);
625 
626  if (scaling_kappainv > eps){
627  fvalue += scaling_kappainv * regimes_u(x)
628  - scaling_mu * (1-XiS(x))/scaling_kappainv * 12. * std::pow(PI, 2) * regimes_gradp(x);
629  }
630 
631  return fvalue;
632  };
633 
635  regimes_g = [](const VectorRd & x)->double { return regimes_gradu(x).trace() ;};
636 
637 
638  //---------------- V-cracked boundary conditions and data (no exact solution here)------------------------------------
639 
640  // Returns 1 in the crack, 0 outside
641  static std::function<double(const VectorRd &)>
642  XiV = [](const VectorRd & x)->double { return (x.z() > 2*std::abs(x.x()-1.)+0.5-eps); };
643 
644  // Viscosity only in crack
646  vcracked_mu = BrinkmanParameters::ViscosityType([](const VectorRd & x)->double { return 0.01 * XiV(x);});
647 
648  // Permeability only outside crack
650  vcracked_kappainv = BrinkmanParameters::PermeabilityInvType([](const VectorRd & x)->double { return 1e-5*(1.-XiV(x)) + 0.001*XiV(x);});
651 
652  // u: entry on the left (surface (0,0,0)-(0,0.2,0)-(0,0.2,1)-(0,0,1),
653  // and "wind" on the top above the crack (surface (0.75,0,1)-(1.25,0,1)-(1.25,0.2,1)-(0.75,0.2,1))
655  vcracked_u = [](const VectorRd & x) -> VectorRd {
656  VectorRd value = VectorRd::Zero();
657  if ( x.x() < .2) {
658  // Entry on the left
659  value = VectorRd(.3, 0., 0.);
660  }else if ( (x.z() > 0.8) && (x.x() > 0.75-eps) && (x.x() < 1.25+eps) ){
661  // Wind top
662  value = VectorRd(1, 0., 0.);
663  }
664 
665  return value;
666  };
667 
668 
670  vcracked_gradu = [](const VectorRd & x) -> MatrixRd {
671  return MatrixRd::Zero();
672  };
673 
675  vcracked_p = [](const VectorRd & x) -> double {
676  return 0.;
677  };
678 
680  vcracked_gradp = [](const VectorRd & x) -> VectorRd {
681  return VectorRd::Zero();
682  };
683 
684  // Source is gravity
686  vcracked_f = [](const VectorRd & x) -> VectorRd {
687  return VectorRd(0., 0., -.98);
688  };
689 
691  vcracked_g = [](const VectorRd & x) -> double {
692  return 0.;
693  };
694 
695  //---------------- BCs for the cavity in porous medium ------------------------------------
696 
697  // Characteristic functions of cavity, wedge, and surrounding box
698  static std::function<double(const VectorRd &)>
699  XiCavity = [](const VectorRd & x)->double {
700  return ( (x.x()>0) && (x.x()<1) && (x.y()>0) && (x.y()<1) && (x.z()>-1) ? 1. : 0.);
701  };
702  static std::function<double(const VectorRd &)>
703  XiWedge = [](const VectorRd & x)->double {
704  return ( (x.x()>=1) && (x.x()<2) && (x.y()>0) && (x.y()<1) && (x.z()> -.75 + .25*(x.x()-1)) ? 1. : 0.);
705  };
706  static std::function<double(const VectorRd &)>
707  XiBox = [](const VectorRd & x)->double { return 1. - XiCavity(x) - XiWedge(x); };
708 
709  // Viscosity in cavity
710  constexpr double viscosity_in_cavity = 0.01;
712 
713  // Permeability in wedge and box
714  constexpr double permeability_in_wedge = 1e-2;
715  constexpr double permeability_in_box = 1e-7;
719 
720  // u: wind above cavity (only for BCs anyway).
722  cavity_u = [](const VectorRd & x) -> VectorRd {
723  return XiCavity(x) * VectorRd(x.x()*(1.-x.x()), 0., 0.);
724  };
725 
727  cavity_gradu = [](const VectorRd & x) -> MatrixRd {
728  return MatrixRd::Zero();
729  };
730 
732  cavity_p = [](const VectorRd & x) -> double {
733  return 0.;
734  };
735 
737  cavity_gradp = [](const VectorRd & x) -> VectorRd {
738  return VectorRd::Zero();
739  };
740 
741  // Source is gravity
743  cavity_f = [](const VectorRd & x) -> VectorRd {
744  return VectorRd(0., 0., -.98);
745  };
746 
748  cavity_g = [](const VectorRd & x) -> double {
749  return 0.;
750  };
751 
752 
753  //---------------- Two regions ------------------------------------------------
754  // Two regions with different parameters, separated by x=1/2 (region S is x<1/2, D is x>1/2)
755  // u = alpha(y,z) + cos(pi x) beta_i(y,z) where i=S or D
756  // This ensures that u is continuous, and that nu \nabla u n is continuous across x=1/2
757 
758  // Characteristic functions
759  static std::function<double(const VectorRd &)>
760  Xi_S = [](const VectorRd & x)->double {
761  return ( (x.x()<0.5) ? 1. : 0.);
762  };
763  static std::function<double(const VectorRd &)>
764  Xi_D = [](const VectorRd & x)->double {
765  return 1.-Xi_S(x);
766  };
767 
768  // Viscosity and permeability
769  constexpr double viscosity_S = 1.;
770  constexpr double viscosity_D = 0.;
774 
775  constexpr double permeability_S = 1e-7;
776  constexpr double permeability_D = 1e-2;
780 
781  // Velocity (alpha and beta only depend on y,z)
782  static std::function<VectorRd(const VectorRd &)>
783  alpha = [](const VectorRd & x)->VectorRd {
784  return VectorRd(
785  exp(-x.y()-x.z()),
786  sin(PI * x.y()) * sin(PI * x.z()),
787  x.y() * x.z()
788  );
789  };
790 
791  static std::function<VectorRd(const VectorRd &)>
792  beta_S = [](const VectorRd & x)->VectorRd {
793  return cos(PI * x.x()) * (x.x() - 0.5)
794  * VectorRd(
795  x.y() + x.z(),
796  x.y() + cos(PI * x.z()),
797  sin(PI * x.y())
798  );
799  };
800 
801  static std::function<VectorRd(const VectorRd &)>
802  beta_D = [](const VectorRd & x)->VectorRd {
803  return cos(PI * x.x()) * (x.x() - 0.5)
804  * VectorRd(
805  sin(PI * x.y()) + cos(PI * x.z()),
806  std::pow(x.z(), 3),
807  std::pow(x.y(), 2) * std::pow(x.z(), 2)
808  );
809  };
810 
812  tworegions_u = [](const VectorRd & x) -> VectorRd {
813  return alpha(x) + Xi_S(x) * beta_S(x) + Xi_D(x) * beta_D(x);
814  };
815 
817  grad_alpha = [](const VectorRd & x) -> MatrixRd {
818  MatrixRd M = MatrixRd::Zero();
819  M.row(0) << 0, -exp(-x.y() - x.z()), -exp(-x.y() - x.z());
820  M.row(1) << 0, PI*sin(PI*x.z())*cos(PI*x.y()), PI*sin(PI*x.y())*cos(PI*x.z());
821  M.row(2) << 0, x.z(), x.y();
822 
823  return M;
824  };
825 
827  grad_beta_S = [](const VectorRd & x) -> MatrixRd {
828  MatrixRd M = MatrixRd::Zero();
829  M.row(0) << -PI*(x.x() - 0.5)*(x.y() + x.z())*sin(PI*x.x()) + (x.y() + x.z())*cos(PI*x.x()),
830  (x.x() - 0.5)*cos(PI*x.x()),
831  (x.x() - 0.5)*cos(PI*x.x());
832  M.row(1) << -PI*(x.x() - 0.5)*(x.y() + cos(PI*x.z()))*sin(PI*x.x()) + (x.y() + cos(PI*x.z()))*cos(PI*x.x()),
833  (x.x() - 0.5)*cos(PI*x.x()),
834  -PI*(x.x() - 0.5)*sin(PI*x.z())*cos(PI*x.x());
835  M.row(2) << -PI*(x.x() - 0.5)*sin(PI*x.x())*sin(PI*x.y()) + sin(PI*x.y())*cos(PI*x.x()),
836  PI*(x.x() - 0.5)*cos(PI*x.x())*cos(PI*x.y()),
837  0;
838 
839  return M;
840  };
841 
843  grad_beta_D = [](const VectorRd & x) -> MatrixRd {
844  MatrixRd M = MatrixRd::Zero();
845  M.row(0) << -PI*(x.x() - 0.5)*(sin(PI*x.y()) + cos(PI*x.z()))*sin(PI*x.x()) + (sin(PI*x.y()) + cos(PI*x.z()))*cos(PI*x.x()),
846  PI*(x.x() - 0.5)*cos(PI*x.x())*cos(PI*x.y()),
847  -PI*(x.x() - 0.5)*sin(PI*x.z())*cos(PI*x.x());
848  M.row(1) << -PI*pow(x.z(), 3)*(x.x() - 0.5)*sin(PI*x.x()) + pow(x.z(), 3)*cos(PI*x.x()),
849  0,
850  3*pow(x.z(), 2)*(x.x() - 0.5)*cos(PI*x.x());
851  M.row(2) << -PI*pow(x.y(), 2)*pow(x.z(), 2)*(x.x() - 0.5)*sin(PI*x.x()) + pow(x.y(), 2)*pow(x.z(), 2)*cos(PI*x.x()),
852  2*x.y()*pow(x.z(), 2)*(x.x() - 0.5)*cos(PI*x.x()),
853  2*pow(x.y(), 2)*x.z()*(x.x() - 0.5)*cos(PI*x.x());
854 
855  return M;
856  };
857 
859  tworegions_gradu = [](const VectorRd & x) -> MatrixRd {
860  return grad_alpha(x) + Xi_S(x) * grad_beta_S(x) + Xi_D(x) * grad_beta_D(x);
861  };
862 
865 
868 
869 // static Brinkman::PressureType
870 // tworegions_p = [](const VectorRd & x) -> double {
871 // return (x.x()-0.5) * x.y() * x.z();
872 // };
873 
874 // static Brinkman::PressureGradientType
875 // tworegions_gradp = [](const VectorRd & x) -> VectorRd {
876 // return VectorRd( x.y() * x.z(), (x.x()-0.5) * x.z(), (x.x()-0.5) * x.y());
877 // };
878 
880  DivGrad_alpha = [](const VectorRd & x) -> VectorRd {
881  return VectorRd(
882  2*exp(-x.y() - x.z()),
883  -2*pow(PI, 2)*sin(PI*x.y())*sin(PI*x.z()),
884  0
885  );
886  };
887 
889  DivGrad_beta_S = [](const VectorRd & x) -> VectorRd {
890  return VectorRd(
891  -pow(PI, 2)*(x.x() - 0.5)*(x.y() + x.z())*cos(PI*x.x()) - 2*PI*(x.y() + x.z())*sin(PI*x.x()),
892  -pow(PI, 2)*(x.x() - 0.5)*(x.y() + cos(PI*x.z()))*cos(PI*x.x()) - pow(PI, 2)*(x.x() - 0.5)*cos(PI*x.x())*cos(PI*x.z()) - 2*PI*(x.y() + cos(PI*x.z()))*sin(PI*x.x()),
893  -2*pow(PI, 2)*(x.x() - 0.5)*sin(PI*x.y())*cos(PI*x.x()) - 2*PI*sin(PI*x.x())*sin(PI*x.y())
894  );
895 
896  };
897 
899  DivGrad_beta_D = [](const VectorRd & x) -> VectorRd {
900  return VectorRd(
901  -pow(PI, 2)*(x.x() - 0.5)*(sin(PI*x.y()) + cos(PI*x.z()))*cos(PI*x.x()) - pow(PI, 2)*(x.x() - 0.5)*sin(PI*x.y())*cos(PI*x.x()) - pow(PI, 2)*(x.x() - 0.5)*cos(PI*x.x())*cos(PI*x.z()) - 2*PI*(sin(PI*x.y()) + cos(PI*x.z()))*sin(PI*x.x()),
902  -pow(PI, 2)*pow(x.z(), 3)*(x.x() - 0.5)*cos(PI*x.x()) - 2*PI*pow(x.z(), 3)*sin(PI*x.x()) + 6*x.z()*(x.x() - 0.5)*cos(PI*x.x()),
903  -pow(PI, 2)*pow(x.y(), 2)*pow(x.z(), 2)*(x.x() - 0.5)*cos(PI*x.x()) - 2*PI*pow(x.y(), 2)*pow(x.z(), 2)*sin(PI*x.x()) + 2*pow(x.y(), 2)*(x.x() - 0.5)*cos(PI*x.x()) + 2*pow(x.z(), 2)*(x.x() - 0.5)*cos(PI*x.x())
904  );
905  };
906 
908  tworegions_f = [](const VectorRd & x) -> VectorRd {
909  return - scaling_mu * ( viscosity_S * Xi_S(x) * (DivGrad_alpha(x) + DivGrad_beta_S(x)) +
910  viscosity_D * Xi_D(x) * (DivGrad_alpha(x) + DivGrad_beta_D(x)) )
911  + scaling_kappainv *
912  ( std::pow(permeability_S, -1) * Xi_S(x) + std::pow(permeability_D, -1) * Xi_D(x) ) * tworegions_u(x)
913  + tworegions_gradp(x);
914  };
915 
917  tworegions_g = [](const VectorRd & x)->double { return tworegions_gradu(x).trace();};
918 
919 
920 } // end of namespace HArDCore3D
921 
922 #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
Base class for global DOF spaces. Provides functions to manipulate global DOFs (the local version bei...
Definition: globaldofspace.hpp:16
Class definition: polynomial bases and operators.
Definition: vhhospace.hpp:52
Class to describe a mesh.
Definition: MeshND.hpp:17
Eigen::Matrix3d MatrixRd
Definition: basis.hpp:51
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
IntegralWeightValueType value
Definition: integralweight.hpp:69
static const double PI
Definition: ddr-magnetostatics.hpp:186
static Magnetostatics::SolutionPotentialType linear_u
Definition: ddr-magnetostatics.hpp:213
static Magnetostatics::PermeabilityType linear_mu
Definition: ddr-magnetostatics.hpp:232
static Magnetostatics::SolutionPotentialType trigonometric_u
Definition: ddr-magnetostatics.hpp:237
static Magnetostatics::PermeabilityType trigonometric_mu
Definition: ddr-magnetostatics.hpp:264
static Magnetostatics::ForcingTermType trigonometric_f
Definition: ddr-magnetostatics.hpp:255
static Magnetostatics::ForcingTermType linear_f
Definition: ddr-magnetostatics.hpp:227
static Stokes::PressureType trigonometric_p
Definition: ddr-stokes.hpp:269
static Stokes::PressureType linear_p
Definition: ddr-stokes.hpp:310
double pressure_scaling
Definition: ddr-stokes.hpp:249
const Mesh & mesh() const
Return a const reference to the mesh.
Definition: vhhospace.hpp:128
bool use_threads
Definition: HHO_DiffAdvecReac.hpp:47
static std::function< VectorRd(const VectorRd &)> beta_D
Definition: hho-brinkman.hpp:802
static Brinkman::CompressibilityForcingTermType vcracked_g
Definition: hho-brinkman.hpp:691
static Brinkman::PressureType regimes_p
Definition: hho-brinkman.hpp:567
constexpr double viscosity_in_cavity
Definition: hho-brinkman.hpp:710
size_t numDirDOFs() const
Returns the number of Dirichlet DOFs.
Definition: hho-brinkman.hpp:177
static BrinkmanParameters::PermeabilityInvType regimes_kappainv
Definition: hho-brinkman.hpp:564
static Brinkman::CompressibilityForcingTermType regimes_g
Definition: hho-brinkman.hpp:635
static Brinkman::MomentumForcingTermType cavity_f
Definition: hho-brinkman.hpp:743
Eigen::VectorXd interpolate(const VelocityType &u, const PressureType &p, const int doe_cell=-1, const int doe_face=-1) const
Interpolates velocity and pressure.
Definition: hho-brinkman.cpp:717
double energy
Global energy norm
Definition: hho-brinkman.hpp:104
static Brinkman::PressureType vcracked_p
Definition: hho-brinkman.hpp:675
static Brinkman::VelocityType cavity_u
Definition: hho-brinkman.hpp:722
std::function< VectorRd(const VectorRd &)> MomentumForcingTermType
Definition: hho-brinkman.hpp:121
static Brinkman::MomentumForcingTermType DivGrad_alpha
Definition: hho-brinkman.hpp:880
static Brinkman::PressureGradientType cubic_gradp
Definition: hho-brinkman.hpp:479
static Brinkman::VelocityType cubic_u
Definition: hho-brinkman.hpp:449
static Brinkman::VelocityGradientType regimes_graduD
Definition: hho-brinkman.hpp:583
double u
Energy norm of velocity.
Definition: hho-brinkman.hpp:100
const std::pair< double, double > & stabilizationParameter() const
Returns the stabilization parameter (scaling)
Definition: hho-brinkman.hpp:245
static const VectorRd vec_p
Definition: hho-brinkman.hpp:356
static Brinkman::VelocityType regimes_u
Definition: hho-brinkman.hpp:612
static std::function< double(const VectorRd &)> XiBox
Definition: hho-brinkman.hpp:707
IntegralWeight PermeabilityInvType
Definition: hho-brinkman.hpp:61
std::function< double(const VectorRd &)> PressureType
Definition: hho-brinkman.hpp:125
static Brinkman::PressureGradientType cavity_gradp
Definition: hho-brinkman.hpp:737
size_t sizeSystem() const
Returns the size of the final system with Lagrange multiplier, after application of SC and removal of...
Definition: hho-brinkman.hpp:201
static Brinkman::VelocityGradientType tworegions_gradu
Definition: hho-brinkman.hpp:859
constexpr double permeability_in_wedge
Definition: hho-brinkman.hpp:714
static std::function< double(const VectorRd &)> XiWedge
Definition: hho-brinkman.hpp:703
static BrinkmanParameters::ViscosityType cavity_mu
Definition: hho-brinkman.hpp:711
static Brinkman::VelocityGradientType quadratic_gradu
Definition: hho-brinkman.hpp:404
IntegralWeight ViscosityType
Definition: hho-brinkman.hpp:60
static BrinkmanParameters::PermeabilityInvType tworegions_kappainv
Definition: hho-brinkman.hpp:778
static BrinkmanParameters::PermeabilityInvType vcracked_kappainv
Definition: hho-brinkman.hpp:650
static Brinkman::VelocityGradientType regimes_gradu
Definition: hho-brinkman.hpp:615
static std::function< double(const VectorRd &)> Xi_D
Definition: hho-brinkman.hpp:764
constexpr double permeability_D
Definition: hho-brinkman.hpp:776
std::pair< double, double > computeFlux(const Eigen::VectorXd &u, const std::pair< std::vector< VectorRd >, VectorRd > &surf) const
Compute the velocity flux (integral of u.n) across a given surface, and the area of the surface,...
Definition: hho-brinkman.cpp:868
static BrinkmanParameters::PermeabilityInvType trigonometric_kappainv
Definition: hho-brinkman.hpp:500
size_t numNonSCDOFs() const
Returns the number of DOFs after SC and with Lagrange multiplier, but before eliminating Dirichlet DO...
Definition: hho-brinkman.hpp:195
static std::function< VectorRd(const VectorRd &)> beta_S
Definition: hho-brinkman.hpp:792
static std::function< double(const VectorRd &)> XiV
Definition: hho-brinkman.hpp:642
static std::function< VectorRd(const VectorRd &)> alpha
Definition: hho-brinkman.hpp:783
size_t dimVelocity() const
Returns the dimension of velocity space.
Definition: hho-brinkman.hpp:183
ViscosityType mu
Definition: hho-brinkman.hpp:86
PermeabilityInvType kappainv
Definition: hho-brinkman.hpp:87
static Brinkman::VelocityGradientType cavity_gradu
Definition: hho-brinkman.hpp:727
static Brinkman::VelocityType regimes_uD
Definition: hho-brinkman.hpp:574
constexpr double eps
Definition: hho-brinkman.hpp:344
static Brinkman::MomentumForcingTermType regimes_f
Definition: hho-brinkman.hpp:618
static BrinkmanParameters::ViscosityType tworegions_mu
Definition: hho-brinkman.hpp:772
static Brinkman::CompressibilityForcingTermType cubic_g
Definition: hho-brinkman.hpp:492
static BrinkmanParameters::PermeabilityInvType cubic_kappainv
Definition: hho-brinkman.hpp:446
constexpr double permeability_in_box
Definition: hho-brinkman.hpp:715
static Brinkman::PressureType cavity_p
Definition: hho-brinkman.hpp:732
static Brinkman::MomentumForcingTermType DivGrad_beta_S
Definition: hho-brinkman.hpp:889
static Brinkman::PressureType cubic_p
Definition: hho-brinkman.hpp:474
double Cf(const Cell &T) const
Returns the friction coefficient.
Definition: hho-brinkman.hpp:75
static BrinkmanParameters::ViscosityType cubic_mu
Definition: hho-brinkman.hpp:443
const Eigen::VectorXd & systemVector() const
Returns the linear system right-hand side vector.
Definition: hho-brinkman.hpp:235
std::function< MatrixRd(const VectorRd &)> VelocityGradientType
Definition: hho-brinkman.hpp:124
static Brinkman::PressureGradientType linear_gradp
Definition: hho-brinkman.hpp:374
static Brinkman::PressureType quadratic_p
Definition: hho-brinkman.hpp:420
std::pair< double, double > & stabilizationParameter()
Returns the stabilization parameter.
Definition: hho-brinkman.hpp:250
double scaling_mu
Definition: hho-brinkman.hpp:340
Brinkman(const VHHOSpace &vhho_space, const GlobalDOFSpace &p_space, const BoundaryConditions &BC, bool use_threads, std::ostream &output=std::cout)
Constructor.
Definition: hho-brinkman.cpp:459
std::vector< double > pressureVertexValues(const Eigen::VectorXd &p) const
Create vertex values for the pressure (from the element values), for plotting.
Definition: hho-brinkman.cpp:817
std::function< VectorRd(const VectorRd &)> VelocityType
Definition: hho-brinkman.hpp:123
static Brinkman::PressureType tworegions_p
Definition: hho-brinkman.hpp:864
static Brinkman::PressureGradientType vcracked_gradp
Definition: hho-brinkman.hpp:680
static std::function< double(const VectorRd &)> XiCavity
Definition: hho-brinkman.hpp:699
std::pair< double, double > computeConditionNum() const
Computes the condition number of the matrix.
Definition: hho-brinkman.cpp:843
static Brinkman::PressureGradientType regimes_gradp
Definition: hho-brinkman.hpp:570
static Brinkman::MomentumForcingTermType vcracked_f
Definition: hho-brinkman.hpp:686
static Brinkman::MomentumForcingTermType tworegions_f
Definition: hho-brinkman.hpp:908
static Brinkman::VelocityType regimes_uS
Definition: hho-brinkman.hpp:603
static Brinkman::CompressibilityForcingTermType cavity_g
Definition: hho-brinkman.hpp:748
size_t nloc_sc_u() const
Returns the local number of velocity statically condensed DOFs.
Definition: hho-brinkman.hpp:147
static Brinkman::MomentumForcingTermType quadratic_f
Definition: hho-brinkman.hpp:433
static BrinkmanParameters::ViscosityType quadratic_mu
Definition: hho-brinkman.hpp:389
static Brinkman::VelocityGradientType grad_beta_S
Definition: hho-brinkman.hpp:827
static BrinkmanParameters::PermeabilityInvType quadratic_kappainv
Definition: hho-brinkman.hpp:392
static Brinkman::PressureGradientType trigonometric_gradp
Definition: hho-brinkman.hpp:535
static BrinkmanParameters::ViscosityType regimes_mu
Definition: hho-brinkman.hpp:561
double scaling_kappainv
Definition: hho-brinkman.hpp:341
size_t numSCDOFs_p() const
Returns the number of pressure statically condensed DOFs.
Definition: hho-brinkman.hpp:165
std::vector< BrinkmanNorms > computeBrinkmanNorms(const std::vector< Eigen::VectorXd > &list_dofs, const BrinkmanParameters &para) const
Compute the discrete energy norm of a family of vectors representing the dofs.
Definition: hho-brinkman.cpp:753
static Brinkman::VelocityGradientType linear_gradu
Definition: hho-brinkman.hpp:369
SystemMatrixType & systemMatrix()
Returns the linear system matrix.
Definition: hho-brinkman.hpp:230
static BrinkmanParameters::ViscosityType vcracked_mu
Definition: hho-brinkman.hpp:646
void assembleLinearSystem(const MomentumForcingTermType &f, const CompressibilityForcingTermType &g, const BrinkmanParameters &para, const VelocityType &u, Eigen::VectorXd &UDir)
Assemble the global system
Definition: hho-brinkman.cpp:490
static Brinkman::MomentumForcingTermType cubic_f
Definition: hho-brinkman.hpp:487
static BrinkmanParameters::PermeabilityInvType cavity_kappainv
Definition: hho-brinkman.hpp:717
static Brinkman::VelocityType vcracked_u
Definition: hho-brinkman.hpp:655
static Brinkman::VelocityGradientType trigonometric_gradu
Definition: hho-brinkman.hpp:514
BrinkmanParameters(const ViscosityType &_mu, const PermeabilityInvType &_kappainv)
Constructor.
Definition: hho-brinkman.hpp:64
static BrinkmanParameters::PermeabilityInvType linear_kappainv
Definition: hho-brinkman.hpp:353
size_t dimPressure() const
Returns the dimension of pressure space.
Definition: hho-brinkman.hpp:189
static std::function< double(const VectorRd &)> Xi_S
Definition: hho-brinkman.hpp:760
static Brinkman::CompressibilityForcingTermType trigonometric_g
Definition: hho-brinkman.hpp:554
static const MatrixRd mat_u
Definition: hho-brinkman.hpp:355
std::function< VectorRd(const VectorRd &)> PressureGradientType
Definition: hho-brinkman.hpp:126
Eigen::VectorXd & scVector()
Returns the static condensation rhs.
Definition: hho-brinkman.hpp:260
static Brinkman::PressureGradientType tworegions_gradp
Definition: hho-brinkman.hpp:867
Eigen::SparseMatrix< double > SystemMatrixType
Definition: hho-brinkman.hpp:119
static Brinkman::MomentumForcingTermType DivGrad_beta_D
Definition: hho-brinkman.hpp:899
static Brinkman::PressureGradientType quadratic_gradp
Definition: hho-brinkman.hpp:425
static Brinkman::VelocityGradientType regimes_graduS
Definition: hho-brinkman.hpp:604
static Brinkman::VelocityType quadratic_u
Definition: hho-brinkman.hpp:395
std::function< double(const VectorRd &)> CompressibilityForcingTermType
Definition: hho-brinkman.hpp:122
constexpr double viscosity_S
Definition: hho-brinkman.hpp:769
static Brinkman::CompressibilityForcingTermType tworegions_g
Definition: hho-brinkman.hpp:917
static Brinkman::VelocityType tworegions_u
Definition: hho-brinkman.hpp:812
static Brinkman::VelocityGradientType grad_alpha
Definition: hho-brinkman.hpp:817
double p
L2 norm of p.
Definition: hho-brinkman.hpp:103
static Brinkman::VelocityGradientType cubic_gradu
Definition: hho-brinkman.hpp:458
Eigen::VectorXd & systemVector()
Returns the linear system right-hand side vector.
Definition: hho-brinkman.hpp:240
constexpr double viscosity_D
Definition: hho-brinkman.hpp:770
static Brinkman::VelocityGradientType grad_beta_D
Definition: hho-brinkman.hpp:843
const SystemMatrixType & scMatrix() const
Returns the static condensation recovery operator.
Definition: hho-brinkman.hpp:255
size_t numSCDOFs() const
Returns the number of statically condensed DOFs.
Definition: hho-brinkman.hpp:171
size_t numSCDOFs_u() const
Returns the number of velocity statically condensed DOFs.
Definition: hho-brinkman.hpp:159
constexpr double permeability_S
Definition: hho-brinkman.hpp:775
const VHHOSpace & vhhospace() const
Returns the velocity space.
Definition: hho-brinkman.hpp:207
BrinkmanNorms(double norm_u, double norm_p)
Constructor.
Definition: hho-brinkman.hpp:94
static Brinkman::CompressibilityForcingTermType linear_g
Definition: hho-brinkman.hpp:384
const GlobalDOFSpace & pspace() const
Returns the pressure space.
Definition: hho-brinkman.hpp:213
const SystemMatrixType & systemMatrix() const
Returns the linear system matrix.
Definition: hho-brinkman.hpp:225
size_t nloc_sc_p() const
Returns the local number of pressure statically condensed DOFs.
Definition: hho-brinkman.hpp:153
static Brinkman::VelocityGradientType vcracked_gradu
Definition: hho-brinkman.hpp:670
static Brinkman::CompressibilityForcingTermType quadratic_g
Definition: hho-brinkman.hpp:438
const Mesh & mesh() const
Returns the mesh.
Definition: hho-brinkman.hpp:219
static std::function< double(const VectorRd &)> XiS
Definition: hho-brinkman.hpp:607
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
MeshND::Cell< 2 > Cell
Definition: Mesh2D.hpp:13
Structure to store norm components (energy for velocity, L2 pressure, and global)
Definition: hho-brinkman.hpp:92
Structure to store physical parameter (and compute derived parameters)
Definition: hho-brinkman.hpp:59
Structure for Brinkman model.
Definition: hho-brinkman.hpp:118
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