MFEM  v3.2
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 
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  int vis_steps = 5;
98 
99  int precision = 8;
100  cout.precision(precision);
101 
102  OptionsParser args(argc, argv);
103  args.AddOption(&mesh_file, "-m", "--mesh",
104  "Mesh file to use.");
105  args.AddOption(&problem, "-p", "--problem",
106  "Problem setup to use. See options in velocity_function().");
107  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
108  "Number of times to refine the mesh uniformly in serial.");
109  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
110  "Number of times to refine the mesh uniformly in parallel.");
111  args.AddOption(&order, "-o", "--order",
112  "Order (degree) of the finite elements.");
113  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
114  "ODE solver: 1 - Forward Euler, 2 - RK2 SSP, 3 - RK3 SSP,"
115  " 4 - RK4, 6 - RK6.");
116  args.AddOption(&t_final, "-tf", "--t-final",
117  "Final time; start time is 0.");
118  args.AddOption(&dt, "-dt", "--time-step",
119  "Time step.");
120  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
121  "--no-visualization",
122  "Enable or disable GLVis visualization.");
123  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
124  "--no-visit-datafiles",
125  "Save data files for VisIt (visit.llnl.gov) visualization.");
126  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
127  "Visualize every n-th timestep.");
128  args.Parse();
129  if (!args.Good())
130  {
131  if (myid == 0)
132  {
133  args.PrintUsage(cout);
134  }
135  MPI_Finalize();
136  return 1;
137  }
138  if (myid == 0)
139  {
140  args.PrintOptions(cout);
141  }
142 
143  // 3. Read the serial mesh from the given mesh file on all processors. We can
144  // handle geometrically periodic meshes in this code.
145  Mesh *mesh = new Mesh(mesh_file, 1, 1);
146  int dim = mesh->Dimension();
147 
148  // 4. Define the ODE solver used for time integration. Several explicit
149  // Runge-Kutta methods are available.
150  ODESolver *ode_solver = NULL;
151  switch (ode_solver_type)
152  {
153  case 1: ode_solver = new ForwardEulerSolver; break;
154  case 2: ode_solver = new RK2Solver(1.0); break;
155  case 3: ode_solver = new RK3SSPSolver; break;
156  case 4: ode_solver = new RK4Solver; break;
157  case 6: ode_solver = new RK6Solver; break;
158  default:
159  if (myid == 0)
160  {
161  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
162  }
163  MPI_Finalize();
164  return 3;
165  }
166 
167  // 5. Refine the mesh in serial to increase the resolution. In this example
168  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
169  // a command-line parameter. If the mesh is of NURBS type, we convert it
170  // to a (piecewise-polynomial) high-order mesh.
171  for (int lev = 0; lev < ser_ref_levels; lev++)
172  {
173  mesh->UniformRefinement();
174  }
175  if (mesh->NURBSext)
176  {
177  mesh->SetCurvature(max(order, 1));
178  }
179  mesh->GetBoundingBox(bb_min, bb_max, max(order, 1));
180 
181  // 6. Define the parallel mesh by a partitioning of the serial mesh. Refine
182  // this mesh further in parallel to increase the resolution. Once the
183  // parallel mesh is defined, the serial mesh can be deleted.
184  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
185  delete mesh;
186  for (int lev = 0; lev < par_ref_levels; lev++)
187  {
188  pmesh->UniformRefinement();
189  }
190 
191  // 7. Define the parallel discontinuous DG finite element space on the
192  // parallel refined mesh of the given polynomial order.
193  DG_FECollection fec(order, dim);
194  ParFiniteElementSpace *fes = new ParFiniteElementSpace(pmesh, &fec);
195 
196  HYPRE_Int global_vSize = fes->GlobalTrueVSize();
197  if (myid == 0)
198  {
199  cout << "Number of unknowns: " << global_vSize << endl;
200  }
201 
202  // 8. Set up and assemble the parallel bilinear and linear forms (and the
203  // parallel hypre matrices) corresponding to the DG discretization. The
204  // DGTraceIntegrator involves integrals over mesh interior faces.
208 
209  ParBilinearForm *m = new ParBilinearForm(fes);
211  ParBilinearForm *k = new ParBilinearForm(fes);
212  k->AddDomainIntegrator(new ConvectionIntegrator(velocity, -1.0));
214  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
216  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
217 
218  ParLinearForm *b = new ParLinearForm(fes);
220  new BoundaryFlowIntegrator(inflow, velocity, -1.0, -0.5));
221 
222  m->Assemble();
223  m->Finalize();
224  int skip_zeros = 0;
225  k->Assemble(skip_zeros);
226  k->Finalize(skip_zeros);
227  b->Assemble();
228 
232 
233  // 9. Define the initial conditions, save the corresponding grid function to
234  // a file and (optionally) save data in the VisIt format and initialize
235  // GLVis visualization.
236  ParGridFunction *u = new ParGridFunction(fes);
237  u->ProjectCoefficient(u0);
238  HypreParVector *U = u->GetTrueDofs();
239 
240  {
241  ostringstream mesh_name, sol_name;
242  mesh_name << "ex9-mesh." << setfill('0') << setw(6) << myid;
243  sol_name << "ex9-init." << setfill('0') << setw(6) << myid;
244  ofstream omesh(mesh_name.str().c_str());
245  omesh.precision(precision);
246  pmesh->Print(omesh);
247  ofstream osol(sol_name.str().c_str());
248  osol.precision(precision);
249  u->Save(osol);
250  }
251 
252  VisItDataCollection visit_dc("Example9-Parallel", pmesh);
253  visit_dc.RegisterField("solution", u);
254  if (visit)
255  {
256  visit_dc.SetCycle(0);
257  visit_dc.SetTime(0.0);
258  visit_dc.Save();
259  }
260 
261  socketstream sout;
262  if (visualization)
263  {
264  char vishost[] = "localhost";
265  int visport = 19916;
266  sout.open(vishost, visport);
267  if (!sout)
268  {
269  if (myid == 0)
270  cout << "Unable to connect to GLVis server at "
271  << vishost << ':' << visport << endl;
272  visualization = false;
273  if (myid == 0)
274  {
275  cout << "GLVis visualization disabled.\n";
276  }
277  }
278  else
279  {
280  sout << "parallel " << num_procs << " " << myid << "\n";
281  sout.precision(precision);
282  sout << "solution\n" << *pmesh << *u;
283  sout << "pause\n";
284  sout << flush;
285  if (myid == 0)
286  cout << "GLVis visualization paused."
287  << " Press space (in the GLVis window) to resume it.\n";
288  }
289  }
290 
291  // 10. Define the time-dependent evolution operator describing the ODE
292  // right-hand side, and perform time-integration (looping over the time
293  // iterations, ti, with a time-step dt).
294  FE_Evolution adv(*M, *K, *B);
295  ode_solver->Init(adv);
296 
297  double t = 0.0;
298  for (int ti = 0; true; )
299  {
300  if (t >= t_final - dt/2)
301  {
302  break;
303  }
304 
305  ode_solver->Step(*U, t, dt);
306  ti++;
307 
308  if (ti % vis_steps == 0)
309  {
310  if (myid == 0)
311  {
312  cout << "time step: " << ti << ", time: " << t << endl;
313  }
314 
315  // 11. Extract the parallel grid function corresponding to the finite
316  // element approximation U (the local solution on each processor).
317  *u = *U;
318 
319  if (visualization)
320  {
321  sout << "parallel " << num_procs << " " << myid << "\n";
322  sout << "solution\n" << *pmesh << *u << flush;
323  }
324 
325  if (visit)
326  {
327  visit_dc.SetCycle(ti);
328  visit_dc.SetTime(t);
329  visit_dc.Save();
330  }
331  }
332  }
333 
334  // 12. Save the final solution in parallel. This output can be viewed later
335  // using GLVis: "glvis -np <np> -m ex9-mesh -g ex9-final".
336  {
337  *u = *U;
338  ostringstream sol_name;
339  sol_name << "ex9-final." << setfill('0') << setw(6) << myid;
340  ofstream osol(sol_name.str().c_str());
341  osol.precision(precision);
342  u->Save(osol);
343  }
344 
345  // 13. Free the used memory.
346  delete U;
347  delete u;
348  delete B;
349  delete b;
350  delete K;
351  delete k;
352  delete M;
353  delete m;
354  delete fes;
355  delete pmesh;
356  delete ode_solver;
357 
358  MPI_Finalize();
359  return 0;
360 }
361 
362 
363 // Implementation of class FE_Evolution
364 FE_Evolution::FE_Evolution(HypreParMatrix &_M, HypreParMatrix &_K,
365  const Vector &_b)
366  : TimeDependentOperator(_M.Height()),
367  M(_M), K(_K), b(_b), M_solver(M.GetComm()), z(_M.Height())
368 {
369  M_prec.SetType(HypreSmoother::Jacobi);
370  M_solver.SetPreconditioner(M_prec);
371  M_solver.SetOperator(M);
372 
373  M_solver.iterative_mode = false;
374  M_solver.SetRelTol(1e-9);
375  M_solver.SetAbsTol(0.0);
376  M_solver.SetMaxIter(100);
377  M_solver.SetPrintLevel(0);
378 }
379 
380 void FE_Evolution::Mult(const Vector &x, Vector &y) const
381 {
382  // y = M^{-1} (K x + b)
383  K.Mult(x, z);
384  z += b;
385  M_solver.Mult(z, y);
386 }
387 
388 
389 // Velocity coefficient
390 void velocity_function(const Vector &x, Vector &v)
391 {
392  int dim = x.Size();
393 
394  // map to the reference [-1,1] domain
395  Vector X(dim);
396  for (int i = 0; i < dim; i++)
397  {
398  double center = (bb_min[i] + bb_max[i]) * 0.5;
399  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
400  }
401 
402  switch (problem)
403  {
404  case 0:
405  {
406  // Translations in 1D, 2D, and 3D
407  switch (dim)
408  {
409  case 1: v(0) = 1.0; break;
410  case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.); break;
411  case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
412  break;
413  }
414  break;
415  }
416  case 1:
417  case 2:
418  {
419  // Clockwise rotation in 2D around the origin
420  const double w = M_PI/2;
421  switch (dim)
422  {
423  case 1: v(0) = 1.0; break;
424  case 2: v(0) = w*X(1); v(1) = -w*X(0); break;
425  case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0; break;
426  }
427  break;
428  }
429  case 3:
430  {
431  // Clockwise twisting rotation in 2D around the origin
432  const double w = M_PI/2;
433  double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
434  d = d*d;
435  switch (dim)
436  {
437  case 1: v(0) = 1.0; break;
438  case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0); break;
439  case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0; break;
440  }
441  break;
442  }
443  }
444 }
445 
446 // Initial condition
447 double u0_function(const Vector &x)
448 {
449  int dim = x.Size();
450 
451  // map to the reference [-1,1] domain
452  Vector X(dim);
453  for (int i = 0; i < dim; i++)
454  {
455  double center = (bb_min[i] + bb_max[i]) * 0.5;
456  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
457  }
458 
459  switch (problem)
460  {
461  case 0:
462  case 1:
463  {
464  switch (dim)
465  {
466  case 1:
467  return exp(-40.*pow(X(0)-0.5,2));
468  case 2:
469  case 3:
470  {
471  double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
472  if (dim == 3)
473  {
474  const double s = (1. + 0.25*cos(2*M_PI*X(2)));
475  rx *= s;
476  ry *= s;
477  }
478  return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
479  erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
480  }
481  }
482  }
483  case 2:
484  {
485  double x_ = X(0), y_ = X(1), rho, phi;
486  rho = hypot(x_, y_);
487  phi = atan2(y_, x_);
488  return pow(sin(M_PI*rho),2)*sin(3*phi);
489  }
490  case 3:
491  {
492  const double f = M_PI;
493  return sin(f*X(0))*sin(f*X(1));
494  }
495  }
496  return 0.0;
497 }
498 
499 // Inflow boundary condition (zero for the problems considered in this example)
500 double inflow_function(const Vector &x)
501 {
502  switch (problem)
503  {
504  case 0:
505  case 1:
506  case 2:
507  case 3: return 0.0;
508  }
509  return 0.0;
510 }
Vector bb_max
Definition: ex9.cpp:52
Conjugate gradient method.
Definition: solvers.hpp:110
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
virtual void RegisterField(const char *field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
Base abstract class for time dependent operators: (x,t) -&gt; f(x,t)
Definition: operator.hpp:68
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:451
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Definition: mesh.cpp:87
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:34
virtual void Step(Vector &x, double &t, double &dt)=0
int Size() const
Returns the size of the vector.
Definition: vector.hpp:86
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
virtual void Save()
Save the collection and a VisIt root file.
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:312
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:247
virtual void Init(TimeDependentOperator &_f)
Definition: ode.hpp:30
int main(int argc, char *argv[])
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
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:6402
Data collection with VisIt I/O routines.
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3020
void Assemble(int skip_zeros=1)
Assemble the local matrix.
Parallel smoothers in hypre.
Definition: hypre.hpp:433
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator.
Definition: linearform.cpp:29
int Dimension() const
Definition: mesh.hpp:523
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
void GetTrueDofs(Vector &tv) const
Returns the true dofs in a Vector.
Definition: pgridfunc.cpp:95
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:58
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:85
void velocity_function(const Vector &x, Vector &v)
Definition: ex9.cpp:319
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:72
NURBSExtension * NURBSext
Definition: mesh.hpp:143
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
Class for parallel bilinear form.
int open(const char hostname[], int port)
double u0_function(const Vector &x)
Definition: ex9.cpp:375
class for C-function coefficient
Vector data type.
Definition: vector.hpp:33
Class for parallel grid function.
Definition: pgridfunc.hpp:31
The classical forward Euler method.
Definition: ode.hpp:39
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:150
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:426
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator.
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:130
virtual void Print(std::ostream &out=std::cout) const
Definition: pmesh.cpp:2927
bool Good() const
Definition: optparser.hpp:120
alpha (q . grad u, v)
Definition: bilininteg.hpp:259