MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
navier_kovasznay.cpp
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 // Navier Kovasznay example
13 //
14 // Solve for the steady Kovasznay flow at Re = 40 defined by
15 //
16 // u = [1 - exp(L * x) * cos(2 * pi * y),
17 // L / (2 * pi) * exp(L * x) * sin(2 * pi * y)],
18 //
19 // p = 1/2 * (1 - exp(2 * L * x)),
20 //
21 // with L = Re/2 - sqrt(Re^2/4 + 4 * pi^2).
22 //
23 // The problem domain is set up like this
24 //
25 // +-------------+
26 // | |
27 // | |
28 // | |
29 // | |
30 // Inflow -> | | -> Outflow
31 // | |
32 // | |
33 // | |
34 // | |
35 // | |
36 // +-------------+
37 //
38 // and Dirichlet boundary conditions are applied for the velocity on every
39 // boundary. The problem, although steady state, is time integrated up to the
40 // final time and the solution is compared with the known exact solution.
41 
42 #include "navier_solver.hpp"
43 #include <fstream>
44 
45 using namespace mfem;
46 using namespace navier;
47 
48 struct s_NavierContext
49 {
50  int ser_ref_levels = 1;
51  int order = 6;
52  double kinvis = 1.0 / 40.0;
53  double t_final = 10 * 0.001;
54  double dt = 0.001;
55  double reference_pressure = 0.0;
56  double reynolds = 1.0 / kinvis;
57  double lam = 0.5 * reynolds
58  - sqrt(0.25 * reynolds * reynolds + 4.0 * M_PI * M_PI);
59  bool pa = true;
60  bool ni = false;
61  bool visualization = false;
62  bool checkres = false;
63 } ctx;
64 
65 void vel_kovasznay(const Vector &x, double t, Vector &u)
66 {
67  double xi = x(0);
68  double yi = x(1);
69 
70  u(0) = 1.0 - exp(ctx.lam * xi) * cos(2.0 * M_PI * yi);
71  u(1) = ctx.lam / (2.0 * M_PI) * exp(ctx.lam * xi) * sin(2.0 * M_PI * yi);
72 }
73 
74 double pres_kovasznay(const Vector &x, double t)
75 {
76  double xi = x(0);
77 
78  return 0.5 * (1.0 - exp(2.0 * ctx.lam * xi)) + ctx.reference_pressure;
79 }
80 
81 int main(int argc, char *argv[])
82 {
83  MPI_Session mpi(argc, argv);
84 
85  OptionsParser args(argc, argv);
86  args.AddOption(&ctx.ser_ref_levels,
87  "-rs",
88  "--refine-serial",
89  "Number of times to refine the mesh uniformly in serial.");
90  args.AddOption(&ctx.order,
91  "-o",
92  "--order",
93  "Order (degree) of the finite elements.");
94  args.AddOption(&ctx.dt, "-dt", "--time-step", "Time step.");
95  args.AddOption(&ctx.t_final, "-tf", "--final-time", "Final time.");
96  args.AddOption(&ctx.pa,
97  "-pa",
98  "--enable-pa",
99  "-no-pa",
100  "--disable-pa",
101  "Enable partial assembly.");
102  args.AddOption(&ctx.ni,
103  "-ni",
104  "--enable-ni",
105  "-no-ni",
106  "--disable-ni",
107  "Enable numerical integration rules.");
108  args.AddOption(&ctx.visualization,
109  "-vis",
110  "--visualization",
111  "-no-vis",
112  "--no-visualization",
113  "Enable or disable GLVis visualization.");
114  args.AddOption(
115  &ctx.checkres,
116  "-cr",
117  "--checkresult",
118  "-no-cr",
119  "--no-checkresult",
120  "Enable or disable checking of the result. Returns -1 on failure.");
121  args.Parse();
122  if (!args.Good())
123  {
124  if (mpi.Root())
125  {
126  args.PrintUsage(mfem::out);
127  }
128  return 1;
129  }
130  if (mpi.Root())
131  {
132  args.PrintOptions(mfem::out);
133  }
134 
135  Mesh *mesh = new Mesh(2, 4, Element::QUADRILATERAL, false, 1.5, 2.0);
136 
137  mesh->EnsureNodes();
138  GridFunction *nodes = mesh->GetNodes();
139  *nodes -= 0.5;
140 
141  for (int i = 0; i < ctx.ser_ref_levels; ++i)
142  {
143  mesh->UniformRefinement();
144  }
145 
146  if (mpi.Root())
147  {
148  std::cout << "Number of elements: " << mesh->GetNE() << std::endl;
149  }
150 
151  auto *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
152  delete mesh;
153 
154  // Create the flow solver.
155  NavierSolver flowsolver(pmesh, ctx.order, ctx.kinvis);
156  flowsolver.EnablePA(ctx.pa);
157  flowsolver.EnableNI(ctx.ni);
158 
159  // Set the initial condition.
160  ParGridFunction *u_ic = flowsolver.GetCurrentVelocity();
161  VectorFunctionCoefficient u_excoeff(pmesh->Dimension(), vel_kovasznay);
162  u_ic->ProjectCoefficient(u_excoeff);
163 
165 
166  // Add Dirichlet boundary conditions to velocity space restricted to
167  // selected attributes on the mesh.
168  Array<int> attr(pmesh->bdr_attributes.Max());
169  attr = 1;
170  flowsolver.AddVelDirichletBC(vel_kovasznay, attr);
171 
172  double t = 0.0;
173  double dt = ctx.dt;
174  double t_final = ctx.t_final;
175  bool last_step = false;
176 
177  flowsolver.Setup(dt);
178 
179  double err_u = 0.0;
180  double err_p = 0.0;
181  ParGridFunction *u_gf = nullptr;
182  ParGridFunction *p_gf = nullptr;
183 
184  ParGridFunction p_ex_gf(flowsolver.GetCurrentPressure()->ParFESpace());
185  GridFunctionCoefficient p_ex_gf_coeff(&p_ex_gf);
186 
187  for (int step = 0; !last_step; ++step)
188  {
189  if (t + dt >= t_final - dt / 2)
190  {
191  last_step = true;
192  }
193 
194  flowsolver.Step(t, dt, step);
195 
196  // Compare against exact solution of velocity and pressure.
197  u_gf = flowsolver.GetCurrentVelocity();
198  p_gf = flowsolver.GetCurrentPressure();
199 
200  u_excoeff.SetTime(t);
201  p_excoeff.SetTime(t);
202 
203  // Remove mean value from exact pressure solution.
204  p_ex_gf.ProjectCoefficient(p_excoeff);
205  flowsolver.MeanZero(p_ex_gf);
206 
207  err_u = u_gf->ComputeL2Error(u_excoeff);
208  err_p = p_gf->ComputeL2Error(p_ex_gf_coeff);
209 
210  double cfl = flowsolver.ComputeCFL(*u_gf, dt);
211 
212  if (mpi.Root())
213  {
214  printf("%5s %8s %8s %8s %11s %11s\n",
215  "Order",
216  "CFL",
217  "Time",
218  "dt",
219  "err_u",
220  "err_p");
221  printf("%5.2d %8.2E %.2E %.2E %.5E %.5E err\n",
222  ctx.order,
223  cfl,
224  t,
225  dt,
226  err_u,
227  err_p);
228  fflush(stdout);
229  }
230  }
231 
232  if (ctx.visualization)
233  {
234  char vishost[] = "localhost";
235  int visport = 19916;
236  socketstream sol_sock(vishost, visport);
237  sol_sock.precision(8);
238  sol_sock << "parallel " << mpi.WorldSize() << " " << mpi.WorldRank()
239  << "\n";
240  sol_sock << "solution\n" << *pmesh << *u_ic << std::flush;
241  }
242 
243  flowsolver.PrintTimingData();
244 
245  // Test if the result for the test run is as expected.
246  if (ctx.checkres)
247  {
248  double tol_u = 1e-6;
249  double tol_p = 1e-5;
250  if (err_u > tol_u || err_p > tol_p)
251  {
252  if (mpi.Root())
253  {
254  mfem::out << "Result has a larger error than expected."
255  << std::endl;
256  }
257  return -1;
258  }
259  }
260 
261  delete pmesh;
262 
263  return 0;
264 }
int WorldSize() const
Return MPI_COMM_WORLD&#39;s size.
void PrintTimingData()
Print timing summary of the solving routine.
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:474
int main(int argc, char *argv[])
Definition: ex1.cpp:66
void Setup(double dt)
Initialize forms, solvers and preconditioners.
void MeanZero(ParGridFunction &v)
Remove the mean from a ParGridFunction.
void Step(double &time, double dt, int cur_step)
Compute solution at the next time step t+dt.
A simple convenience class that calls MPI_Init() at construction and MPI_Finalize() at destruction...
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
double ComputeCFL(ParGridFunction &u, double dt)
Compute CFL.
constexpr char vishost[]
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:8382
constexpr int visport
ParGridFunction * GetCurrentPressure()
Return a pointer to the current pressure ParGridFunction.
bool Root() const
Return true if WorldRank() == 0.
void EnablePA(bool pa)
Enable partial assembly for every operator.
void AddVelDirichletBC(VectorCoefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the velocity field.
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:434
void SetTime(double t)
Set the time for time dependent coefficients.
Definition: coefficient.hpp:48
A general vector function coefficient.
int WorldRank() const
Return MPI_COMM_WORLD&#39;s rank.
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
ParGridFunction * GetCurrentVelocity()
Return a pointer to the current velocity ParGridFunction.
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:267
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:304
A general function coefficient.
Vector data type.
Definition: vector.hpp:51
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:6603
Class for parallel grid function.
Definition: pgridfunc.hpp:32
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
Class for parallel meshes.
Definition: pmesh.hpp:32
Transient incompressible Navier Stokes solver in a split scheme formulation.
void EnsureNodes()
Definition: mesh.cpp:4132
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:108
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:145