MFEM v4.8.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-2025, 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
45using namespace mfem;
46using namespace navier;
47
48struct s_NavierContext
49{
50 int ser_ref_levels = 1;
51 int order = 6;
52 real_t kinvis = 1.0 / 40.0;
53 real_t t_final = 10 * 0.001;
54 real_t dt = 0.001;
55 real_t reference_pressure = 0.0;
56 real_t reynolds = 1.0 / kinvis;
57 real_t 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;
64
65void vel_kovasznay(const Vector &x, real_t t, Vector &u)
66{
67 real_t xi = x(0);
68 real_t 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
75{
76 real_t xi = x(0);
77
78 return 0.5 * (1.0 - exp(2.0 * ctx.lam * xi)) + ctx.reference_pressure;
79}
80
81int main(int argc, char *argv[])
82{
83 Mpi::Init(argc, argv);
85 int visport = 19916;
86
87 OptionsParser args(argc, argv);
88 args.AddOption(&ctx.ser_ref_levels,
89 "-rs",
90 "--refine-serial",
91 "Number of times to refine the mesh uniformly in serial.");
92 args.AddOption(&ctx.order,
93 "-o",
94 "--order",
95 "Order (degree) of the finite elements.");
96 args.AddOption(&ctx.dt, "-dt", "--time-step", "Time step.");
97 args.AddOption(&ctx.t_final, "-tf", "--final-time", "Final time.");
98 args.AddOption(&ctx.pa,
99 "-pa",
100 "--enable-pa",
101 "-no-pa",
102 "--disable-pa",
103 "Enable partial assembly.");
104 args.AddOption(&ctx.ni,
105 "-ni",
106 "--enable-ni",
107 "-no-ni",
108 "--disable-ni",
109 "Enable numerical integration rules.");
110 args.AddOption(&ctx.visualization,
111 "-vis",
112 "--visualization",
113 "-no-vis",
114 "--no-visualization",
115 "Enable or disable GLVis visualization.");
116 args.AddOption(
117 &ctx.checkres,
118 "-cr",
119 "--checkresult",
120 "-no-cr",
121 "--no-checkresult",
122 "Enable or disable checking of the result. Returns -1 on failure.");
123 args.AddOption(&visport, "-p", "--send-port", "Socket for GLVis.");
124 args.Parse();
125 if (!args.Good())
126 {
127 if (Mpi::Root())
128 {
129 args.PrintUsage(mfem::out);
130 }
131 return 1;
132 }
133 if (Mpi::Root())
134 {
136 }
137
138 Mesh mesh = Mesh::MakeCartesian2D(2, 4, Element::QUADRILATERAL, false, 1.5,
139 2.0);
140
141 mesh.EnsureNodes();
142 GridFunction *nodes = mesh.GetNodes();
143 *nodes -= 0.5;
144
145 for (int i = 0; i < ctx.ser_ref_levels; ++i)
146 {
147 mesh.UniformRefinement();
148 }
149
150 if (Mpi::Root())
151 {
152 std::cout << "Number of elements: " << mesh.GetNE() << std::endl;
153 }
154
155 auto *pmesh = new ParMesh(MPI_COMM_WORLD, mesh);
156 mesh.Clear();
157
158 // Create the flow solver.
159 NavierSolver flowsolver(pmesh, ctx.order, ctx.kinvis);
160 flowsolver.EnablePA(ctx.pa);
161 flowsolver.EnableNI(ctx.ni);
162
163 // Set the initial condition.
164 ParGridFunction *u_ic = flowsolver.GetCurrentVelocity();
165 VectorFunctionCoefficient u_excoeff(pmesh->Dimension(), vel_kovasznay);
166 u_ic->ProjectCoefficient(u_excoeff);
167
169
170 // Add Dirichlet boundary conditions to velocity space restricted to
171 // selected attributes on the mesh.
172 Array<int> attr(pmesh->bdr_attributes.Max());
173 attr = 1;
174 flowsolver.AddVelDirichletBC(vel_kovasznay, attr);
175
176 real_t t = 0.0;
177 real_t dt = ctx.dt;
178 real_t t_final = ctx.t_final;
179 bool last_step = false;
180
181 flowsolver.Setup(dt);
182
183 real_t err_u = 0.0;
184 real_t err_p = 0.0;
185 ParGridFunction *u_gf = nullptr;
186 ParGridFunction *p_gf = nullptr;
187
188 ParGridFunction p_ex_gf(flowsolver.GetCurrentPressure()->ParFESpace());
189 GridFunctionCoefficient p_ex_gf_coeff(&p_ex_gf);
190
191 for (int step = 0; !last_step; ++step)
192 {
193 if (t + dt >= t_final - dt / 2)
194 {
195 last_step = true;
196 }
197
198 flowsolver.Step(t, dt, step);
199
200 // Compare against exact solution of velocity and pressure.
201 u_gf = flowsolver.GetCurrentVelocity();
202 p_gf = flowsolver.GetCurrentPressure();
203
204 u_excoeff.SetTime(t);
205 p_excoeff.SetTime(t);
206
207 // Remove mean value from exact pressure solution.
208 p_ex_gf.ProjectCoefficient(p_excoeff);
209 flowsolver.MeanZero(p_ex_gf);
210
211 err_u = u_gf->ComputeL2Error(u_excoeff);
212 err_p = p_gf->ComputeL2Error(p_ex_gf_coeff);
213
214 real_t cfl = flowsolver.ComputeCFL(*u_gf, dt);
215
216 if (Mpi::Root())
217 {
218 printf("%5s %8s %8s %8s %11s %11s\n",
219 "Order",
220 "CFL",
221 "Time",
222 "dt",
223 "err_u",
224 "err_p");
225 printf("%5.2d %8.2E %.2E %.2E %.5E %.5E err\n",
226 ctx.order,
227 cfl,
228 t,
229 dt,
230 err_u,
231 err_p);
232 fflush(stdout);
233 }
234 }
235
236 if (ctx.visualization)
237 {
238 char vishost[] = "localhost";
239 socketstream sol_sock(vishost, visport);
240 sol_sock.precision(8);
241 sol_sock << "parallel " << Mpi::WorldSize() << " "
242 << Mpi::WorldRank() << "\n";
243 sol_sock << "solution\n" << *pmesh << *u_ic << std::flush;
244 }
245
246 flowsolver.PrintTimingData();
247
248 // Test if the result for the test run is as expected.
249 if (ctx.checkres)
250 {
251#if defined(MFEM_USE_DOUBLE)
252 real_t tol_u = 1e-6;
253 real_t tol_p = 1e-5;
254#elif defined(MFEM_USE_SINGLE)
255 real_t tol_u = 1e-5;
256 real_t tol_p = 2e-4;
257#else
258#error "Only single and double precision are supported!"
259 real_t tol_u = 0;
260 real_t tol_p = 0;
261#endif
262 if (err_u > tol_u || err_p > tol_p)
263 {
264 if (Mpi::Root())
265 {
266 mfem::out << "Result has a larger error than expected."
267 << std::endl;
268 }
269 return -1;
270 }
271 }
272
273 delete pmesh;
274
275 return 0;
276}
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
A general function coefficient.
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.cpp:33
Mesh data type.
Definition mesh.hpp:64
void EnsureNodes()
Make sure that the mesh has valid nodes, i.e. its geometry is described by a vector finite element gr...
Definition mesh.cpp:6432
void Clear()
Clear the contents of the Mesh.
Definition mesh.hpp:761
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1282
void GetNodes(Vector &node_coord) const
Definition mesh.cpp:9294
static Mesh MakeCartesian2D(int nx, int ny, Element::Type type, bool generate_edges=false, real_t sx=1.0, real_t sy=1.0, bool sfc_ordering=true)
Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny quadrilaterals if type = QUADRILATER...
Definition mesh.cpp:4471
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:11295
static bool Root()
Return true if the rank in MPI_COMM_WORLD is zero.
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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 'var' to receive the value. Enable/disable tags are used to set the bool...
Definition optparser.hpp:82
bool Good() const
Return true if the command line options were parsed successfully.
Class for parallel grid function.
Definition pgridfunc.hpp:50
real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
Returns ||u_ex - u_h||_L2 in parallel for H1 or L2 elements.
ParFiniteElementSpace * ParFESpace() const
void ProjectCoefficient(Coefficient &coeff) override
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Class for parallel meshes.
Definition pmesh.hpp:34
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
A general vector function coefficient.
Vector data type.
Definition vector.hpp:82
Transient incompressible Navier Stokes solver in a split scheme formulation.
void Setup(real_t dt)
Initialize forms, solvers and preconditioners.
real_t ComputeCFL(ParGridFunction &u, real_t dt)
Compute CFL.
void AddVelDirichletBC(VectorCoefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the velocity field.
ParGridFunction * GetCurrentPressure()
Return a pointer to the current pressure ParGridFunction.
void PrintTimingData()
Print timing summary of the solving routine.
void EnablePA(bool pa)
Enable partial assembly for every operator.
void MeanZero(ParGridFunction &v)
Remove the mean from a ParGridFunction.
ParGridFunction * GetCurrentVelocity()
Return a pointer to the current velocity ParGridFunction.
void Step(real_t &time, real_t dt, int cur_step, bool provisional=false)
Compute solution at the next time step t+dt.
int main()
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
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
float real_t
Definition config.hpp:43
const char vishost[]
real_t pres_kovasznay(const Vector &x, real_t t)
struct s_NavierContext ctx
void vel_kovasznay(const Vector &x, real_t t, Vector &u)
std::array< int, NCMesh::MaxFaceNodes > nodes