MFEM  v3.0
 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(Vector &x);
43 
44 // Inflow boundary condition
45 double inflow_function(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  args.PrintUsage(cout);
127  MPI_Finalize();
128  return 1;
129  }
130  if (myid == 0)
131  args.PrintOptions(cout);
132 
133  // 3. Read the serial mesh from the given mesh file on all processors. We can
134  // handle geometrically periodic meshes in this code.
135  Mesh *mesh;
136  ifstream imesh(mesh_file);
137  if (!imesh)
138  {
139  if (myid == 0)
140  cerr << "\nCan not open mesh file: " << mesh_file << '\n' << endl;
141  MPI_Finalize();
142  return 2;
143  }
144  mesh = new Mesh(imesh, 1, 1);
145  imesh.close();
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  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
161  MPI_Finalize();
162  return 3;
163  }
164 
165  // 5. Refine the mesh in serial to increase the resolution. In this example
166  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
167  // a command-line parameter. If the mesh is of NURBS type, we convert it
168  // to a (piecewise-polynomial) high-order mesh.
169  for (int lev = 0; lev < ser_ref_levels; lev++)
170  mesh->UniformRefinement();
171 
172  if (mesh->NURBSext)
173  {
174  int mesh_order = std::max(order, 1);
175  FiniteElementCollection *mfec = new H1_FECollection(mesh_order, dim);
176  FiniteElementSpace *mfes = new FiniteElementSpace(mesh, mfec, dim);
177  mesh->SetNodalFESpace(mfes);
178  mesh->GetNodes()->MakeOwner(mfec);
179  }
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  pmesh->UniformRefinement();
188 
189  // 7. Define the parallel discontinuous DG finite element space on the
190  // parallel refined mesh of the given polynomial order.
191  DG_FECollection fec(order, dim);
192  ParFiniteElementSpace *fes = new ParFiniteElementSpace(pmesh, &fec);
193 
194  int global_vSize = fes->GlobalTrueVSize();
195  if (myid == 0)
196  cout << "Number of unknowns: " << global_vSize << endl;
197 
198  // 8. Set up and assemble the parallel bilinear and linear forms (and the
199  // parallel hypre matrices) corresponding to the DG discretization. The
200  // DGTraceIntegrator involves integrals over mesh interior faces.
204 
205  ParBilinearForm *m = new ParBilinearForm(fes);
207  ParBilinearForm *k = new ParBilinearForm(fes);
208  k->AddDomainIntegrator(new ConvectionIntegrator(velocity, -1.0));
210  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
212  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
213 
214  ParLinearForm *b = new ParLinearForm(fes);
216  new BoundaryFlowIntegrator(inflow, velocity, -1.0, -0.5));
217 
218  m->Assemble();
219  m->Finalize();
220  int skip_zeros = 0;
221  k->Assemble(skip_zeros);
222  k->Finalize(skip_zeros);
223  b->Assemble();
224 
228 
229  // 9. Define the initial conditions, save the corresponding grid function to
230  // a file and (optionally) save data in the VisIt format and initialize
231  // GLVis visualization.
232  ParGridFunction *u = new ParGridFunction(fes);
233  u->ProjectCoefficient(u0);
234  HypreParVector *U = u->GetTrueDofs();
235 
236  {
237  ostringstream mesh_name, sol_name;
238  mesh_name << "ex9-mesh." << setfill('0') << setw(6) << myid;
239  sol_name << "ex9-init." << setfill('0') << setw(6) << myid;
240  ofstream omesh(mesh_name.str().c_str());
241  omesh.precision(precision);
242  pmesh->Print(omesh);
243  ofstream osol(sol_name.str().c_str());
244  osol.precision(precision);
245  u->Save(osol);
246  }
247 
248  VisItDataCollection visit_dc("Example9-Parallel", pmesh);
249  visit_dc.RegisterField("solution", u);
250  if (visit)
251  {
252  visit_dc.SetCycle(0);
253  visit_dc.SetTime(0.0);
254  visit_dc.Save();
255  }
256 
257  socketstream sout;
258  if (visualization)
259  {
260  char vishost[] = "localhost";
261  int visport = 19916;
262  sout.open(vishost, visport);
263  if (!sout)
264  {
265  if (myid == 0)
266  cout << "Unable to connect to GLVis server at "
267  << vishost << ':' << visport << endl;
268  visualization = false;
269  if (myid == 0)
270  cout << "GLVis visualization disabled.\n";
271  }
272  else
273  {
274  sout << "parallel " << num_procs << " " << myid << "\n";
275  sout.precision(precision);
276  sout << "solution\n" << *pmesh << *u;
277  sout << "pause\n";
278  sout << flush;
279  if (myid == 0)
280  cout << "GLVis visualization paused."
281  << " Press space (in the GLVis window) to resume it.\n";
282  }
283  }
284 
285  // 10. Define the time-dependent evolution operator describing the ODE
286  // right-hand side, and perform time-integration (looping over the time
287  // iterations, ti, with a time-step dt).
288  FE_Evolution adv(*M, *K, *B);
289  ode_solver->Init(adv);
290 
291  double t = 0.0;
292  for (int ti = 0; true; )
293  {
294  if (t >= t_final - dt/2)
295  break;
296 
297  ode_solver->Step(*U, t, dt);
298  ti++;
299 
300  if (ti % vis_steps == 0)
301  {
302  if (myid == 0)
303  cout << "time step: " << ti << ", time: " << t << endl;
304 
305  // 11. Extract the parallel grid function corresponding to the finite
306  // element approximation U (the local solution on each processor).
307  *u = *U;
308 
309  if (visualization)
310  {
311  sout << "parallel " << num_procs << " " << myid << "\n";
312  sout << "solution\n" << *pmesh << *u << flush;
313  }
314 
315  if (visit)
316  {
317  visit_dc.SetCycle(ti);
318  visit_dc.SetTime(t);
319  visit_dc.Save();
320  }
321  }
322  }
323 
324  // 12. Save the final solution in parallel. This output can be viewed later
325  // using GLVis: "glvis -np <np> -m ex9-mesh -g ex9-final".
326  {
327  *u = *U;
328  ostringstream sol_name;
329  sol_name << "ex9-final." << setfill('0') << setw(6) << myid;
330  ofstream osol(sol_name.str().c_str());
331  osol.precision(precision);
332  u->Save(osol);
333  }
334 
335  // 13. Free the used memory.
336  delete U;
337  delete u;
338  delete B;
339  delete b;
340  delete K;
341  delete k;
342  delete M;
343  delete m;
344  delete fes;
345  delete pmesh;
346  delete ode_solver;
347 
348  MPI_Finalize();
349  return 0;
350 }
351 
352 
353 // Implementation of class FE_Evolution
354 FE_Evolution::FE_Evolution(HypreParMatrix &_M, HypreParMatrix &_K, const Vector &_b)
355  : TimeDependentOperator(_M.Height()),
356  M(_M), K(_K), b(_b), M_solver(M.GetComm()), z(_M.Height())
357 {
358  M_prec.SetType(HypreSmoother::Jacobi);
359  M_solver.SetPreconditioner(M_prec);
360  M_solver.SetOperator(M);
361 
362  M_solver.iterative_mode = false;
363  M_solver.SetRelTol(1e-9);
364  M_solver.SetAbsTol(0.0);
365  M_solver.SetMaxIter(100);
366  M_solver.SetPrintLevel(0);
367 }
368 
369 void FE_Evolution::Mult(const Vector &x, Vector &y) const
370 {
371  // y = M^{-1} (K x + b)
372  K.Mult(x, z);
373  z += b;
374  M_solver.Mult(z, y);
375 }
376 
377 
378 // Velocity coefficient
379 void velocity_function(const Vector &x, Vector &v)
380 {
381  int dim = x.Size();
382 
383  switch (problem)
384  {
385  case 0:
386  {
387  // Translations in 1D, 2D, and 3D
388  switch (dim)
389  {
390  case 1: v(0) = 1.0; break;
391  case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.); break;
392  case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.); break;
393  }
394  break;
395  }
396  case 1:
397  case 2:
398  {
399  // Clockwise rotation in 2D around the origin
400  const double w = M_PI/2;
401  switch (dim)
402  {
403  case 1: v(0) = 1.0; break;
404  case 2: v(0) = w*x(1); v(1) = -w*x(0); break;
405  case 3: v(0) = w*x(1); v(1) = -w*x(0); v(2) = 0.0; break;
406  }
407  break;
408  }
409  case 3:
410  {
411  // Clockwise twisting rotation in 2D around the origin
412  const double w = M_PI/2;
413  double d = max((x(0)+1.)*(1.-x(0)),0.) * max((x(1)+1.)*(1.-x(1)),0.);
414  d = d*d;
415  switch (dim)
416  {
417  case 1: v(0) = 1.0; break;
418  case 2: v(0) = d*w*x(1); v(1) = -d*w*x(0); break;
419  case 3: v(0) = d*w*x(1); v(1) = -d*w*x(0); v(2) = 0.0; break;
420  }
421  break;
422  }
423  }
424 }
425 
426 // Initial condition
427 double u0_function(Vector &x)
428 {
429  int dim = x.Size();
430 
431  switch (problem)
432  {
433  case 0:
434  case 1:
435  {
436  switch (dim)
437  {
438  case 1:
439  return exp(-40.*pow(x(0)-0.5,2));
440  case 2:
441  case 3:
442  {
443  double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
444  if (dim == 3)
445  {
446  const double s = (1. + 0.25*cos(2*M_PI*x(2)));
447  rx *= s;
448  ry *= s;
449  }
450  return ( erfc(w*(x(0)-cx-rx))*erfc(-w*(x(0)-cx+rx)) *
451  erfc(w*(x(1)-cy-ry))*erfc(-w*(x(1)-cy+ry)) )/16;
452  }
453  }
454  }
455  case 2:
456  {
457  const double r = sqrt(8.);
458  double x_ = x(0), y_ = x(1), rho, phi;
459  rho = hypot(x_, y_) / r;
460  phi = atan2(y_, x_);
461  return pow(sin(M_PI*rho),2)*sin(3*phi);
462  }
463  case 3:
464  {
465  const double f = M_PI;
466  return sin(f*x(0))*sin(f*x(1));
467  }
468  }
469  return 0.0;
470 }
471 
472 // Inflow boundary condition (zero for the problems considered in this example)
473 double inflow_function(Vector &x)
474 {
475  switch (problem)
476  {
477  case 0:
478  case 1:
479  case 2:
480  case 3: return 0.0;
481  }
482  return 0.0;
483 }
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:351
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:76
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
double inflow_function(Vector &x)
Definition: ex9.cpp:413
virtual void Save()
Save the collection and a VisIt root file.
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:250
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:229
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
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6225
Data collection with VisIt I/O routines.
void Assemble(int skip_zeros=1)
Assemble the local matrix.
Parallel smoothers in hypre.
Definition: hypre.hpp:250
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator.
Definition: linearform.cpp:29
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:3066
int problem
Definition: ex9.cpp:37
int Dimension() const
Definition: mesh.hpp:417
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:385
void SetTime(double t)
Set physical time (for time-dependent simulations)
void GetTrueDofs(Vector &tv) const
Returns the true dofs in a HypreParVector.
Definition: pgridfunc.cpp:83
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:38
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:85
int main(int argc, char *argv[])
Definition: ex1.cpp:39
void velocity_function(const Vector &x, Vector &v)
Definition: ex9.cpp:319
Abstract finite element space.
Definition: fespace.hpp:61
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
double u0_function(Vector &x)
Definition: ex9.cpp:367
NURBSExtension * NURBSext
Definition: mesh.hpp:307
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:266
Class for parallel bilinear form.
int open(const char hostname[], int port)
class for C-function coefficient
Vector data type.
Definition: vector.hpp:29
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:4948
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:52
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:103
Class for parallel meshes.
Definition: pmesh.hpp:27
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:34
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator.
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:83
virtual void Print(std::ostream &out=std::cout) const
Definition: pmesh.cpp:2444
bool Good() const
Definition: optparser.hpp:120
alpha (q . grad u, v)
Definition: bilininteg.hpp:239