MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
navier_solver.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_NAVIER_SOLVER_HPP
13 #define MFEM_NAVIER_SOLVER_HPP
14 
15 #define NAVIER_VERSION 0.1
16 
17 #include "mfem.hpp"
18 #include "ortho_solver.hpp"
19 
20 namespace mfem
21 {
22 namespace navier
23 {
24 using VecFuncT = void(const Vector &x, double t, Vector &u);
25 using ScalarFuncT = double(const Vector &x, double t);
26 
27 /// Container for a Dirichlet boundary condition of the velocity field.
29 {
30 public:
32  : attr(attr), coeff(coeff)
33  {}
34 
36  {
37  // Deep copy the attribute array
38  this->attr = obj.attr;
39 
40  // Move the coefficient pointer
41  this->coeff = obj.coeff;
42  obj.coeff = nullptr;
43  }
44 
45  ~VelDirichletBC_T() { delete coeff; }
46 
49 };
50 
51 /// Container for a Dirichlet boundary condition of the pressure field.
53 {
54 public:
56  : attr(attr), coeff(coeff)
57  {}
58 
60  {
61  // Deep copy the attribute array
62  this->attr = obj.attr;
63 
64  // Move the coefficient pointer
65  this->coeff = obj.coeff;
66  obj.coeff = nullptr;
67  }
68 
69  ~PresDirichletBC_T() { delete coeff; }
70 
73 };
74 
75 /// Container for an acceleration term.
77 {
78 public:
80  : attr(attr), coeff(coeff)
81  {}
82 
84  {
85  // Deep copy the attribute array
86  this->attr = obj.attr;
87 
88  // Move the coefficient pointer
89  this->coeff = obj.coeff;
90  obj.coeff = nullptr;
91  }
92 
93  ~AccelTerm_T() { delete coeff; }
94 
97 };
98 
99 /// Transient incompressible Navier Stokes solver in a split scheme formulation.
100 /**
101  * This implementation of a transient incompressible Navier Stokes solver uses
102  * the non-dimensionalized formulation. The coupled momentum and
103  * incompressibilty equations are decoupled using the split scheme described in
104  * [1]. This leads to three solving steps.
105  *
106  * 1. An extrapolation step for all nonlinear terms which are treated
107  * explicitly. This step avoids a fully coupled nonlinear solve and only
108  * requires a solve of the mass matrix in velocity space \f$M_v^{-1}\f$. On
109  * the other hand this introduces a CFL stability condition on the maximum
110  * timestep.
111  *
112  * 2. A Poisson solve \f$S_p^{-1}\f$.
113  *
114  * 3. A Helmholtz like solve \f$(M_v - \partial t K_v)^{-1}\f$.
115  *
116  * The numerical solver setup for each step are as follows.
117  *
118  * \f$M_v^{-1}\f$ is solved using CG with Jacobi as preconditioner.
119  *
120  * \f$S_p^{-1}\f$ is solved using CG with AMG applied to the low order refined
121  * (LOR) assembled pressure Poisson matrix. To avoid assembling a matrix for
122  * preconditioning, one can use p-MG as an alternative (NYI).
123  *
124  * \f$(M_v - \partial t K_v)^{-1}\f$ due to the CFL condition we expect the time
125  * step to be small. Therefore this is solved using CG with Jacobi as
126  * preconditioner. For large time steps a preconditioner like AMG or p-MG should
127  * be used (NYI).
128  *
129  * Statements marked with NYI mean this feature is planned but Not Yet
130  * Implemented.
131  *
132  * A detailed description is available in [1] section 4.2. The algorithm is
133  * originated from [2].
134  *
135  * [1] Michael Franco, Jean-Sylvain Camier, Julian Andrej, Will Pazner (2020)
136  * High-order matrix-free incompressible flow solvers with GPU acceleration and
137  * low-order refined preconditioners (https://arxiv.org/abs/1910.03032)
138  *
139  * [2] A. G. Tomboulides, J. C. Y. Lee & S. A. Orszag (1997) Numerical
140  * Simulation of Low Mach Number Reactive Flows
141  */
143 {
144 public:
145  /// Initialize data structures, set FE space order and kinematic viscosity.
146  /**
147  * The ParMesh @a mesh can be a linear or curved parallel mesh. The @a order
148  * of the finite element spaces is this algorithm is of equal order
149  * \f$(P_N)^d P_N\f$ for velocity and pressure respectively. This means the
150  * pressure is in discretized in the same space (just scalar instead of a
151  * vector space) as the velocity.
152  *
153  * Kinematic viscosity (dimensionless) is set using @a kin_vis and
154  * automatically converted to the Reynolds number. If you want to set the
155  * Reynolds number directly, you can provide the inverse.
156  */
157  NavierSolver(ParMesh *mesh, int order, double kin_vis);
158 
159  /// Initialize forms, solvers and preconditioners.
160  void Setup(double dt);
161 
162  /// Compute solution at the next time step t+dt.
163  void Step(double &time, double dt, int cur_step);
164 
165  /// Return a pointer to the current velocity ParGridFunction.
167 
168  /// Return a pointer to the current pressure ParGridFunction.
170 
171  /// Add a Dirichlet boundary condition to the velocity field.
172  void AddVelDirichletBC(VectorCoefficient *coeff, Array<int> &attr);
173 
174  void AddVelDirichletBC(VecFuncT *f, Array<int> &attr);
175 
176  /// Add a Dirichlet boundary condition to the pressure field.
177  void AddPresDirichletBC(Coefficient *coeff, Array<int> &attr);
178 
180 
181  /// Add an acceleration term to the RHS of the equation.
182  /**
183  * The VecFuncT @a f is evaluated at the current time t and extrapolated
184  * together with the nonlinear parts of the Navier Stokes equation.
185  */
186  void AddAccelTerm(VectorCoefficient *coeff, Array<int> &attr);
187 
188  void AddAccelTerm(VecFuncT *f, Array<int> &attr);
189 
190  /// Enable partial assembly for every operator.
191  void EnablePA(bool pa) { partial_assembly = pa; }
192 
193  /// Enable numerical integration rules. This means collocated quadrature at
194  /// the nodal points.
195  void EnableNI(bool ni) { numerical_integ = ni; }
196 
197  /// Print timing summary of the solving routine.
198  /**
199  * The summary shows the timing in seconds in the first row of
200  *
201  * 1. SETUP: Time spent for the setup of all forms, solvers and
202  * preconditioners.
203  * 2. STEP: Time spent computing a full time step. It includes all three
204  * solves.
205  * 3. EXTRAP: Time spent for extrapolation of all forcing and nonlinear
206  * terms.
207  * 4. CURLCURL: Time spent for computing the curl curl term in the pressure
208  * Poisson equation (see references for detailed explanation).
209  * 5. PSOLVE: Time spent in the pressure Poisson solve.
210  * 6. HSOLVE: Time spent in the Helmholtz solve.
211  *
212  * The second row shows a proportion of a column relative to the whole
213  * time step.
214  */
215  void PrintTimingData();
216 
217  ~NavierSolver();
218 
219  /// Compute \f$\nabla \times \nabla \times u\f$ for \f$u \in (H^1)^2\f$.
221  ParGridFunction &cu,
222  bool assume_scalar = false);
223 
224  /// Compute \f$\nabla \times \nabla \times u\f$ for \f$u \in (H^1)^3\f$.
226 
227  /// Remove mean from a Vector.
228  /**
229  * Modify the Vector @a v by subtracting its mean using
230  * \f$v = v - \frac{\sum_i^N v_i}{N} \f$
231  */
232  void Orthogonalize(Vector &v);
233 
234  /// Remove the mean from a ParGridFunction.
235  /**
236  * Modify the ParGridFunction @a v by subtracting its mean using
237  * \f$ v = v - \int_\Omega \frac{v}{vol(\Omega)} dx \f$.
238  */
239  void MeanZero(ParGridFunction &v);
240 
241  /// Compute CFL
242  double ComputeCFL(ParGridFunction &u, double dt);
243 
244 protected:
245  /// Print informations about the Navier version.
246  void PrintInfo();
247 
248  /// Update the EXTk/BDF time integration coefficient.
249  /**
250  * Depending on which time step the computation is in, the EXTk/BDF time
251  * integration coefficients have to be set accordingly. This allows
252  * bootstrapping with a BDF scheme of order 1 and increasing the order each
253  * following time step, up to order 3.
254  */
255  void SetTimeIntegrationCoefficients(int step);
256 
257  /// Eliminate essential BCs in an Operator and apply to RHS.
258  void EliminateRHS(Operator &A,
259  ConstrainedOperator &constrainedA,
260  const Array<int> &ess_tdof_list,
261  Vector &x,
262  Vector &b,
263  Vector &X,
264  Vector &B,
265  int copy_interior = 0);
266 
267  /// Enable/disable debug output.
268  bool debug = false;
269 
270  /// Enable/disable verbose output.
271  bool verbose = true;
272 
273  /// Enable/disable partial assembly of forms.
274  bool partial_assembly = false;
275 
276  /// Enable/disable numerical integration rules of forms.
277  bool numerical_integ = false;
278 
279  /// The parallel mesh.
280  ParMesh *pmesh = nullptr;
281 
282  /// The order of the velocity and pressure space.
283  int order;
284 
285  /// Kinematic viscosity (dimensionless).
286  double kin_vis;
287 
288  /// Velocity \f$H^1\f$ finite element collection.
290 
291  /// Pressure \f$H^1\f$ finite element collection.
293 
294  /// Velocity \f$(H^1)^d\f$ finite element space.
296 
297  /// Pressure \f$H^1\f$ finite element space.
299 
300  ParNonlinearForm *N = nullptr;
301 
303 
305 
307 
309 
311 
313 
315 
316  ParLinearForm *f_form = nullptr;
317 
319 
320  /// Linear form to compute the mass matrix in various subroutines.
321  ParLinearForm *mass_lf = nullptr;
323  double volume = 0.0;
324 
329 
335 
336  Solver *MvInvPC = nullptr;
337  CGSolver *MvInv = nullptr;
338 
341  CGSolver *SpInv = nullptr;
342 
343  Solver *HInvPC = nullptr;
344  CGSolver *HInv = nullptr;
345 
348 
350 
352 
354 
355  // All essential attributes.
358 
359  // All essential true dofs.
362 
363  // Bookkeeping for velocity dirichlet bcs.
364  std::vector<VelDirichletBC_T> vel_dbcs;
365 
366  // Bookkeeping for pressure dirichlet bcs.
367  std::vector<PresDirichletBC_T> pres_dbcs;
368 
369  // Bookkeeping for acceleration (forcing) terms.
370  std::vector<AccelTerm_T> accel_terms;
371 
372  int cur_step = 0;
373 
374  // BDFk/EXTk coefficients.
375  double bd0 = 0.0;
376  double bd1 = 0.0;
377  double bd2 = 0.0;
378  double bd3 = 0.0;
379  double ab1 = 0.0;
380  double ab2 = 0.0;
381  double ab3 = 0.0;
382 
383  // Timers.
385 
386  // Print levels.
387  int pl_mvsolve = 0;
388  int pl_spsolve = 0;
389  int pl_hsolve = 0;
390  int pl_amg = 0;
391 
392  // Relative tolerances.
393  double rtol_spsolve = 1e-6;
394  double rtol_hsolve = 1e-8;
395 
396  // Iteration counts.
398 
399  // Residuals.
400  double res_mvsolve = 0.0, res_spsolve = 0.0, res_hsolve = 0.0;
401 
402  // LOR related.
403  ParMesh *pmesh_lor = nullptr;
406  InterpolationGridTransfer *vgt = nullptr, *pgt = nullptr;
407 
411 
415 };
416 
417 } // namespace navier
418 
419 } // namespace mfem
420 
421 #endif
VectorCoefficient * coeff
Container for a Dirichlet boundary condition of the velocity field.
void PrintTimingData()
Print timing summary of the solving routine.
void Orthogonalize(Vector &v)
Remove mean from a Vector.
Conjugate gradient method.
Definition: solvers.hpp:258
FiniteElementCollection * pfec
Pressure finite element collection.
void EliminateRHS(Operator &A, ConstrainedOperator &constrainedA, const Array< int > &ess_tdof_list, Vector &x, Vector &b, Vector &X, Vector &B, int copy_interior=0)
Eliminate essential BCs in an Operator and apply to RHS.
ParBilinearForm * Mv_form_lor
Container for a Dirichlet boundary condition of the pressure field.
Base class for vector Coefficients that optionally depend on time and space.
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
void PrintInfo()
Print informations about the Navier version.
FiniteElementCollection * pfec_lor
ConstantCoefficient nlcoeff
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
void(const Vector &x, double t, Vector &u) VecFuncT
Parallel non-linear operator on the true dofs.
double(const Vector &x, double t) ScalarFuncT
ParLinearForm * FText_bdr_form
std::vector< PresDirichletBC_T > pres_dbcs
bool numerical_integ
Enable/disable numerical integration rules of forms.
Abstract parallel finite element space.
Definition: pfespace.hpp:28
ConstantCoefficient H_lincoeff
Container for an acceleration term.
NavierSolver(ParMesh *mesh, int order, double kin_vis)
Initialize data structures, set FE space order and kinematic viscosity.
AccelTerm_T(Array< int > attr, VectorCoefficient *coeff)
PresDirichletBC_T(PresDirichletBC_T &&obj)
ParMixedBilinearForm * G_form
void Setup(double dt)
Initialize forms, solvers and preconditioners.
void MeanZero(ParGridFunction &v)
Remove the mean from a ParGridFunction.
VelDirichletBC_T(VelDirichletBC_T &&obj)
bool partial_assembly
Enable/disable partial assembly of forms.
ConstantCoefficient Sp_coeff
void Step(double &time, double dt, int cur_step)
Compute solution at the next time step t+dt.
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1079
void ComputeCurl3D(ParGridFunction &u, ParGridFunction &cu)
Compute for .
Class for parallel linear form.
Definition: plinearform.hpp:26
FiniteElementCollection * vfec
Velocity finite element collection.
ParFiniteElementSpace * pfes_lor
InterpolationGridTransfer * vgt
double ComputeCFL(ParGridFunction &u, double dt)
Compute CFL.
double b
Definition: lissajous.cpp:42
ParGridFunction * GetCurrentPressure()
Return a pointer to the current pressure ParGridFunction.
Timing object.
Definition: tic_toc.hpp:34
ConstantCoefficient onecoeff
void SetTimeIntegrationCoefficients(int step)
Update the EXTk/BDF time integration coefficient.
void EnablePA(bool pa)
Enable partial assembly for every operator.
ParBilinearForm * Sp_form_lor
void AddVelDirichletBC(VectorCoefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the velocity field.
void ComputeCurl2D(ParGridFunction &u, ParGridFunction &cu, bool assume_scalar=false)
Compute for .
ParBilinearForm * H_form_lor
ParFiniteElementSpace * vfes
Velocity finite element space.
Transfer data between a coarse mesh and an embedded refined mesh using interpolation.
Definition: fespace.hpp:870
ParLinearForm * mass_lf
Linear form to compute the mass matrix in various subroutines.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
ParGridFunction * GetCurrentVelocity()
Return a pointer to the current velocity ParGridFunction.
std::vector< AccelTerm_T > accel_terms
bool verbose
Enable/disable verbose output.
InterpolationGridTransfer * pgt
double kin_vis
Kinematic viscosity (dimensionless).
Class for parallel bilinear form using different test and trial FE spaces.
std::vector< VelDirichletBC_T > vel_dbcs
VectorGridFunctionCoefficient * FText_gfcoeff
ParFiniteElementSpace * pfes
Pressure finite element space.
PresDirichletBC_T(Array< int > attr, Coefficient *coeff)
Class for parallel bilinear form.
AccelTerm_T(AccelTerm_T &&obj)
ParMixedBilinearForm * D_form
Vector data type.
Definition: vector.hpp:51
Vector coefficient defined by a vector GridFunction.
void AddPresDirichletBC(Coefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the pressure field.
Base class for solvers.
Definition: operator.hpp:634
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Square Operator for imposing essential boundary conditions using only the action, Mult()...
Definition: operator.hpp:800
Abstract operator.
Definition: operator.hpp:24
bool debug
Enable/disable debug output.
Class for parallel meshes.
Definition: pmesh.hpp:32
Transient incompressible Navier Stokes solver in a split scheme formulation.
Solver wrapper which orthogonalizes the input and output vector.
ParMesh * pmesh
The parallel mesh.
void AddAccelTerm(VectorCoefficient *coeff, Array< int > &attr)
Add an acceleration term to the RHS of the equation.
int order
The order of the velocity and pressure space.
ConstantCoefficient H_bdfcoeff
double f(const Vector &p)
VelDirichletBC_T(Array< int > attr, VectorCoefficient *coeff)