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