MFEM  v3.3
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex9p.cpp
Go to the documentation of this file.
1 // MFEM Example 9 - Parallel Version
2 //
3 // Compile with: make ex9p
4 //
5 // Sample runs:
6 // mpirun -np 4 ex9p -m ../data/periodic-segment.mesh -p 0 -dt 0.005
7 // mpirun -np 4 ex9p -m ../data/periodic-square.mesh -p 0 -dt 0.01
8 // mpirun -np 4 ex9p -m ../data/periodic-hexagon.mesh -p 0 -dt 0.01
9 // mpirun -np 4 ex9p -m ../data/periodic-square.mesh -p 1 -dt 0.005 -tf 9
10 // mpirun -np 4 ex9p -m ../data/periodic-hexagon.mesh -p 1 -dt 0.005 -tf 9
11 // mpirun -np 4 ex9p -m ../data/amr-quad.mesh -p 1 -rp 1 -dt 0.002 -tf 9
12 // mpirun -np 4 ex9p -m ../data/star-q3.mesh -p 1 -rp 1 -dt 0.004 -tf 9
13 // mpirun -np 4 ex9p -m ../data/disc-nurbs.mesh -p 1 -rp 1 -dt 0.005 -tf 9
14 // mpirun -np 4 ex9p -m ../data/disc-nurbs.mesh -p 2 -rp 1 -dt 0.005 -tf 9
15 // mpirun -np 4 ex9p -m ../data/periodic-square.mesh -p 3 -rp 2 -dt 0.0025 -tf 9 -vs 20
16 // mpirun -np 4 ex9p -m ../data/periodic-cube.mesh -p 0 -o 2 -rp 1 -dt 0.01 -tf 8
17 //
18 // Description: This example code solves the time-dependent advection equation
19 // du/dt + v.grad(u) = 0, where v is a given fluid velocity, and
20 // u0(x)=u(0,x) is a given initial condition.
21 //
22 // The example demonstrates the use of Discontinuous Galerkin (DG)
23 // bilinear forms in MFEM (face integrators), the use of explicit
24 // ODE time integrators, the definition of periodic boundary
25 // conditions through periodic meshes, as well as the use of GLVis
26 // for persistent visualization of a time-evolving solution. The
27 // saving of time-dependent data files for external visualization
28 // with VisIt (visit.llnl.gov) is also illustrated.
29 
30 #include "mfem.hpp"
31 #include <fstream>
32 #include <iostream>
33 
34 using namespace std;
35 using namespace mfem;
36 
37 // Choice for the problem setup. The fluid velocity, initial condition and
38 // inflow boundary condition are chosen based on this parameter.
39 int problem;
40 
41 // Velocity coefficient
42 void velocity_function(const Vector &x, Vector &v);
43 
44 // Initial condition
45 double u0_function(const Vector &x);
46 
47 // Inflow boundary condition
48 double inflow_function(const Vector &x);
49 
50 // Mesh bounding box
52 
53 
54 /** A time-dependent operator for the right-hand side of the ODE. The DG weak
55  form of du/dt = -v.grad(u) is M du/dt = K u + b, where M and K are the mass
56  and advection matrices, and b describes the flow on the boundary. This can
57  be written as a general ODE, du/dt = M^{-1} (K u + b), and this class is
58  used to evaluate the right-hand side. */
59 class FE_Evolution : public TimeDependentOperator
60 {
61 private:
62  HypreParMatrix &M, &K;
63  const Vector &b;
64  HypreSmoother M_prec;
65  CGSolver M_solver;
66 
67  mutable Vector z;
68 
69 public:
70  FE_Evolution(HypreParMatrix &_M, HypreParMatrix &_K, const Vector &_b);
71 
72  virtual void Mult(const Vector &x, Vector &y) const;
73 
74  virtual ~FE_Evolution() { }
75 };
76 
77 
78 int main(int argc, char *argv[])
79 {
80  // 1. Initialize MPI.
81  int num_procs, myid;
82  MPI_Init(&argc, &argv);
83  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
84  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
85 
86  // 2. Parse command-line options.
87  problem = 0;
88  const char *mesh_file = "../data/periodic-hexagon.mesh";
89  int ser_ref_levels = 2;
90  int par_ref_levels = 0;
91  int order = 3;
92  int ode_solver_type = 4;
93  double t_final = 10.0;
94  double dt = 0.01;
95  bool visualization = true;
96  bool visit = false;
97  bool binary = false;
98  int vis_steps = 5;
99 
100  int precision = 8;
101  cout.precision(precision);
102 
103  OptionsParser args(argc, argv);
104  args.AddOption(&mesh_file, "-m", "--mesh",
105  "Mesh file to use.");
106  args.AddOption(&problem, "-p", "--problem",
107  "Problem setup to use. See options in velocity_function().");
108  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
109  "Number of times to refine the mesh uniformly in serial.");
110  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
111  "Number of times to refine the mesh uniformly in parallel.");
112  args.AddOption(&order, "-o", "--order",
113  "Order (degree) of the finite elements.");
114  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
115  "ODE solver: 1 - Forward Euler,\n\t"
116  " 2 - RK2 SSP, 3 - RK3 SSP, 4 - RK4, 6 - RK6.");
117  args.AddOption(&t_final, "-tf", "--t-final",
118  "Final time; start time is 0.");
119  args.AddOption(&dt, "-dt", "--time-step",
120  "Time step.");
121  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
122  "--no-visualization",
123  "Enable or disable GLVis visualization.");
124  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
125  "--no-visit-datafiles",
126  "Save data files for VisIt (visit.llnl.gov) visualization.");
127  args.AddOption(&binary, "-binary", "--binary-datafiles", "-ascii",
128  "--ascii-datafiles",
129  "Use binary (Sidre) or ascii format for VisIt data files.");
130  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
131  "Visualize every n-th timestep.");
132  args.Parse();
133  if (!args.Good())
134  {
135  if (myid == 0)
136  {
137  args.PrintUsage(cout);
138  }
139  MPI_Finalize();
140  return 1;
141  }
142  if (myid == 0)
143  {
144  args.PrintOptions(cout);
145  }
146 
147  // 3. Read the serial mesh from the given mesh file on all processors. We can
148  // handle geometrically periodic meshes in this code.
149  Mesh *mesh = new Mesh(mesh_file, 1, 1);
150  int dim = mesh->Dimension();
151 
152  // 4. Define the ODE solver used for time integration. Several explicit
153  // Runge-Kutta methods are available.
154  ODESolver *ode_solver = NULL;
155  switch (ode_solver_type)
156  {
157  case 1: ode_solver = new ForwardEulerSolver; break;
158  case 2: ode_solver = new RK2Solver(1.0); break;
159  case 3: ode_solver = new RK3SSPSolver; break;
160  case 4: ode_solver = new RK4Solver; break;
161  case 6: ode_solver = new RK6Solver; break;
162  default:
163  if (myid == 0)
164  {
165  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
166  }
167  MPI_Finalize();
168  return 3;
169  }
170 
171  // 5. Refine the mesh in serial to increase the resolution. In this example
172  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
173  // a command-line parameter. If the mesh is of NURBS type, we convert it
174  // to a (piecewise-polynomial) high-order mesh.
175  for (int lev = 0; lev < ser_ref_levels; lev++)
176  {
177  mesh->UniformRefinement();
178  }
179  if (mesh->NURBSext)
180  {
181  mesh->SetCurvature(max(order, 1));
182  }
183  mesh->GetBoundingBox(bb_min, bb_max, max(order, 1));
184 
185  // 6. Define the parallel mesh by a partitioning of the serial mesh. Refine
186  // this mesh further in parallel to increase the resolution. Once the
187  // parallel mesh is defined, the serial mesh can be deleted.
188  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
189  delete mesh;
190  for (int lev = 0; lev < par_ref_levels; lev++)
191  {
192  pmesh->UniformRefinement();
193  }
194 
195  // 7. Define the parallel discontinuous DG finite element space on the
196  // parallel refined mesh of the given polynomial order.
197  DG_FECollection fec(order, dim);
198  ParFiniteElementSpace *fes = new ParFiniteElementSpace(pmesh, &fec);
199 
200  HYPRE_Int global_vSize = fes->GlobalTrueVSize();
201  if (myid == 0)
202  {
203  cout << "Number of unknowns: " << global_vSize << endl;
204  }
205 
206  // 8. Set up and assemble the parallel bilinear and linear forms (and the
207  // parallel hypre matrices) corresponding to the DG discretization. The
208  // DGTraceIntegrator involves integrals over mesh interior faces.
212 
213  ParBilinearForm *m = new ParBilinearForm(fes);
215  ParBilinearForm *k = new ParBilinearForm(fes);
216  k->AddDomainIntegrator(new ConvectionIntegrator(velocity, -1.0));
218  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
220  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
221 
222  ParLinearForm *b = new ParLinearForm(fes);
224  new BoundaryFlowIntegrator(inflow, velocity, -1.0, -0.5));
225 
226  m->Assemble();
227  m->Finalize();
228  int skip_zeros = 0;
229  k->Assemble(skip_zeros);
230  k->Finalize(skip_zeros);
231  b->Assemble();
232 
236 
237  // 9. Define the initial conditions, save the corresponding grid function to
238  // a file and (optionally) save data in the VisIt format and initialize
239  // GLVis visualization.
240  ParGridFunction *u = new ParGridFunction(fes);
241  u->ProjectCoefficient(u0);
242  HypreParVector *U = u->GetTrueDofs();
243 
244  {
245  ostringstream mesh_name, sol_name;
246  mesh_name << "ex9-mesh." << setfill('0') << setw(6) << myid;
247  sol_name << "ex9-init." << setfill('0') << setw(6) << myid;
248  ofstream omesh(mesh_name.str().c_str());
249  omesh.precision(precision);
250  pmesh->Print(omesh);
251  ofstream osol(sol_name.str().c_str());
252  osol.precision(precision);
253  u->Save(osol);
254  }
255 
256  // Create data collection for solution output: either VisItDataCollection for
257  // ascii data files, or SidreDataCollection for binary data files.
258  DataCollection *dc = NULL;
259  if (visit)
260  {
261  if (binary)
262  {
263 #ifdef MFEM_USE_SIDRE
264  dc = new SidreDataCollection("Example9-Parallel", pmesh);
265 #else
266  MFEM_ABORT("Must build with MFEM_USE_SIDRE=YES for binary output.");
267 #endif
268  }
269  else
270  {
271  dc = new VisItDataCollection("Example9-Parallel", pmesh);
272  dc->SetPrecision(precision);
273  }
274  dc->RegisterField("solution", u);
275  dc->SetCycle(0);
276  dc->SetTime(0.0);
277  dc->Save();
278  }
279 
280  socketstream sout;
281  if (visualization)
282  {
283  char vishost[] = "localhost";
284  int visport = 19916;
285  sout.open(vishost, visport);
286  if (!sout)
287  {
288  if (myid == 0)
289  cout << "Unable to connect to GLVis server at "
290  << vishost << ':' << visport << endl;
291  visualization = false;
292  if (myid == 0)
293  {
294  cout << "GLVis visualization disabled.\n";
295  }
296  }
297  else
298  {
299  sout << "parallel " << num_procs << " " << myid << "\n";
300  sout.precision(precision);
301  sout << "solution\n" << *pmesh << *u;
302  sout << "pause\n";
303  sout << flush;
304  if (myid == 0)
305  cout << "GLVis visualization paused."
306  << " Press space (in the GLVis window) to resume it.\n";
307  }
308  }
309 
310  // 10. Define the time-dependent evolution operator describing the ODE
311  // right-hand side, and perform time-integration (looping over the time
312  // iterations, ti, with a time-step dt).
313  FE_Evolution adv(*M, *K, *B);
314 
315  double t = 0.0;
316  adv.SetTime(t);
317  ode_solver->Init(adv);
318 
319  bool done = false;
320  for (int ti = 0; !done; )
321  {
322  double dt_real = min(dt, t_final - t);
323  ode_solver->Step(*U, t, dt_real);
324  ti++;
325 
326  done = (t >= t_final - 1e-8*dt);
327 
328  if (done || ti % vis_steps == 0)
329  {
330  if (myid == 0)
331  {
332  cout << "time step: " << ti << ", time: " << t << endl;
333  }
334 
335  // 11. Extract the parallel grid function corresponding to the finite
336  // element approximation U (the local solution on each processor).
337  *u = *U;
338 
339  if (visualization)
340  {
341  sout << "parallel " << num_procs << " " << myid << "\n";
342  sout << "solution\n" << *pmesh << *u << flush;
343  }
344 
345  if (visit)
346  {
347  dc->SetCycle(ti);
348  dc->SetTime(t);
349  dc->Save();
350  }
351  }
352  }
353 
354  // 12. Save the final solution in parallel. This output can be viewed later
355  // using GLVis: "glvis -np <np> -m ex9-mesh -g ex9-final".
356  {
357  *u = *U;
358  ostringstream sol_name;
359  sol_name << "ex9-final." << setfill('0') << setw(6) << myid;
360  ofstream osol(sol_name.str().c_str());
361  osol.precision(precision);
362  u->Save(osol);
363  }
364 
365  // 13. Free the used memory.
366  delete U;
367  delete u;
368  delete B;
369  delete b;
370  delete K;
371  delete k;
372  delete M;
373  delete m;
374  delete fes;
375  delete pmesh;
376  delete ode_solver;
377  delete dc;
378 
379  MPI_Finalize();
380  return 0;
381 }
382 
383 
384 // Implementation of class FE_Evolution
385 FE_Evolution::FE_Evolution(HypreParMatrix &_M, HypreParMatrix &_K,
386  const Vector &_b)
387  : TimeDependentOperator(_M.Height()),
388  M(_M), K(_K), b(_b), M_solver(M.GetComm()), z(_M.Height())
389 {
390  M_prec.SetType(HypreSmoother::Jacobi);
391  M_solver.SetPreconditioner(M_prec);
392  M_solver.SetOperator(M);
393 
394  M_solver.iterative_mode = false;
395  M_solver.SetRelTol(1e-9);
396  M_solver.SetAbsTol(0.0);
397  M_solver.SetMaxIter(100);
398  M_solver.SetPrintLevel(0);
399 }
400 
401 void FE_Evolution::Mult(const Vector &x, Vector &y) const
402 {
403  // y = M^{-1} (K x + b)
404  K.Mult(x, z);
405  z += b;
406  M_solver.Mult(z, y);
407 }
408 
409 
410 // Velocity coefficient
411 void velocity_function(const Vector &x, Vector &v)
412 {
413  int dim = x.Size();
414 
415  // map to the reference [-1,1] domain
416  Vector X(dim);
417  for (int i = 0; i < dim; i++)
418  {
419  double center = (bb_min[i] + bb_max[i]) * 0.5;
420  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
421  }
422 
423  switch (problem)
424  {
425  case 0:
426  {
427  // Translations in 1D, 2D, and 3D
428  switch (dim)
429  {
430  case 1: v(0) = 1.0; break;
431  case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.); break;
432  case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
433  break;
434  }
435  break;
436  }
437  case 1:
438  case 2:
439  {
440  // Clockwise rotation in 2D around the origin
441  const double w = M_PI/2;
442  switch (dim)
443  {
444  case 1: v(0) = 1.0; break;
445  case 2: v(0) = w*X(1); v(1) = -w*X(0); break;
446  case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0; break;
447  }
448  break;
449  }
450  case 3:
451  {
452  // Clockwise twisting rotation in 2D around the origin
453  const double w = M_PI/2;
454  double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
455  d = d*d;
456  switch (dim)
457  {
458  case 1: v(0) = 1.0; break;
459  case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0); break;
460  case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0; break;
461  }
462  break;
463  }
464  }
465 }
466 
467 // Initial condition
468 double u0_function(const Vector &x)
469 {
470  int dim = x.Size();
471 
472  // map to the reference [-1,1] domain
473  Vector X(dim);
474  for (int i = 0; i < dim; i++)
475  {
476  double center = (bb_min[i] + bb_max[i]) * 0.5;
477  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
478  }
479 
480  switch (problem)
481  {
482  case 0:
483  case 1:
484  {
485  switch (dim)
486  {
487  case 1:
488  return exp(-40.*pow(X(0)-0.5,2));
489  case 2:
490  case 3:
491  {
492  double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
493  if (dim == 3)
494  {
495  const double s = (1. + 0.25*cos(2*M_PI*X(2)));
496  rx *= s;
497  ry *= s;
498  }
499  return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
500  erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
501  }
502  }
503  }
504  case 2:
505  {
506  double x_ = X(0), y_ = X(1), rho, phi;
507  rho = hypot(x_, y_);
508  phi = atan2(y_, x_);
509  return pow(sin(M_PI*rho),2)*sin(3*phi);
510  }
511  case 3:
512  {
513  const double f = M_PI;
514  return sin(f*X(0))*sin(f*X(1));
515  }
516  }
517  return 0.0;
518 }
519 
520 // Inflow boundary condition (zero for the problems considered in this example)
521 double inflow_function(const Vector &x)
522 {
523  switch (problem)
524  {
525  case 0:
526  case 1:
527  case 2:
528  case 3: return 0.0;
529  }
530  return 0.0;
531 }
Vector bb_max
Definition: ex9.cpp:52
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
Conjugate gradient method.
Definition: solvers.hpp:111
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
virtual void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.hpp:37
Base abstract class for time dependent operators.
Definition: operator.hpp:140
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:468
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Definition: mesh.cpp:88
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:42
virtual void Step(Vector &x, double &t, double &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
int Size() const
Returns the size of the vector.
Definition: vector.hpp:106
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:314
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:249
int main(int argc, char *argv[])
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
Data collection with Sidre routines following the Conduit mesh blueprint specification.
Class for parallel linear form.
Definition: plinearform.hpp:26
int dim
Definition: ex3.cpp:47
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6684
Data collection with VisIt I/O routines.
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3278
void Assemble(int skip_zeros=1)
Assemble the local matrix.
Parallel smoothers in hypre.
Definition: hypre.hpp:491
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator.
Definition: linearform.cpp:29
int Dimension() const
Definition: mesh.hpp:611
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetTime(double t)
Set physical time (for time-dependent simulations)
Vector bb_min
Definition: ex9.cpp:52
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:72
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:147
void velocity_function(const Vector &x, Vector &v)
Definition: ex9.cpp:339
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)
Definition: optparser.hpp:74
int problem
Definition: ex15.cpp:57
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:134
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:144
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:102
Class for parallel bilinear form.
int open(const char hostname[], int port)
double u0_function(const Vector &x)
Definition: ex9.cpp:396
class for C-function coefficient
Vector data type.
Definition: vector.hpp:36
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
Class for parallel grid function.
Definition: pgridfunc.hpp:31
The classical forward Euler method.
Definition: ode.hpp:101
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:174
Class for parallel meshes.
Definition: pmesh.hpp:28
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:34
double inflow_function(const Vector &x)
Definition: ex9.cpp:449
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator.
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:195
virtual void Print(std::ostream &out=std::cout) const
Definition: pmesh.cpp:3395
bool Good() const
Definition: optparser.hpp:120
alpha (q . grad u, v)