MFEM  v4.6.0
Finite element discretization library
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/amr-quad.mesh -p 1 -rp 1 -dt 0.02 -s 13 -tf 9
13 // mpirun -np 4 ex9p -m ../data/star-q3.mesh -p 1 -rp 1 -dt 0.004 -tf 9
14 // mpirun -np 4 ex9p -m ../data/star-mixed.mesh -p 1 -rp 1 -dt 0.004 -tf 9
15 // mpirun -np 4 ex9p -m ../data/disc-nurbs.mesh -p 1 -rp 1 -dt 0.005 -tf 9
16 // mpirun -np 4 ex9p -m ../data/disc-nurbs.mesh -p 2 -rp 1 -dt 0.005 -tf 9
17 // mpirun -np 4 ex9p -m ../data/periodic-square.mesh -p 3 -rp 2 -dt 0.0025 -tf 9 -vs 20
18 // mpirun -np 4 ex9p -m ../data/periodic-cube.mesh -p 0 -o 2 -rp 1 -dt 0.01 -tf 8
19 // mpirun -np 4 ex9p -m ../data/periodic-square.msh -p 0 -rs 2 -dt 0.005 -tf 2
20 // mpirun -np 4 ex9p -m ../data/periodic-cube.msh -p 0 -rs 1 -o 2 -tf 2
21 // mpirun -np 3 ex9p -m ../data/amr-hex.mesh -p 1 -rs 1 -rp 0 -dt 0.005 -tf 0.5
22 //
23 // Device sample runs:
24 // mpirun -np 4 ex9p -pa
25 // mpirun -np 4 ex9p -ea
26 // mpirun -np 4 ex9p -fa
27 // mpirun -np 4 ex9p -pa -m ../data/periodic-cube.mesh
28 // mpirun -np 4 ex9p -pa -m ../data/periodic-cube.mesh -d cuda
29 // mpirun -np 4 ex9p -ea -m ../data/periodic-cube.mesh -d cuda
30 // mpirun -np 4 ex9p -fa -m ../data/periodic-cube.mesh -d cuda
31 // mpirun -np 4 ex9p -pa -m ../data/amr-quad.mesh -p 1 -rp 1 -dt 0.002 -tf 9 -d cuda
32 //
33 // Description: This example code solves the time-dependent advection equation
34 // du/dt + v.grad(u) = 0, where v is a given fluid velocity, and
35 // u0(x)=u(0,x) is a given initial condition.
36 //
37 // The example demonstrates the use of Discontinuous Galerkin (DG)
38 // bilinear forms in MFEM (face integrators), the use of implicit
39 // and explicit ODE time integrators, the definition of periodic
40 // boundary conditions through periodic meshes, as well as the use
41 // of GLVis for persistent visualization of a time-evolving
42 // solution. Saving of time-dependent data files for visualization
43 // with VisIt (visit.llnl.gov) and ParaView (paraview.org), as
44 // well as the optional saving with ADIOS2 (adios2.readthedocs.io)
45 // are also illustrated.
46 
47 #include "mfem.hpp"
48 #include <fstream>
49 #include <iostream>
50 
51 using namespace std;
52 using namespace mfem;
53 
54 // Choice for the problem setup. The fluid velocity, initial condition and
55 // inflow boundary condition are chosen based on this parameter.
56 int problem;
57 
58 // Velocity coefficient
59 void velocity_function(const Vector &x, Vector &v);
60 
61 // Initial condition
62 double u0_function(const Vector &x);
63 
64 // Inflow boundary condition
65 double inflow_function(const Vector &x);
66 
67 // Mesh bounding box
69 
70 // Type of preconditioner for implicit time integrator
71 enum class PrecType : int
72 {
73  ILU = 0,
74  AIR = 1
75 };
76 
77 #if MFEM_HYPRE_VERSION >= 21800
78 // Algebraic multigrid preconditioner for advective problems based on
79 // approximate ideal restriction (AIR). Most effective when matrix is
80 // first scaled by DG block inverse, and AIR applied to scaled matrix.
81 // See https://doi.org/10.1137/17M1144350.
82 class AIR_prec : public Solver
83 {
84 private:
85  const HypreParMatrix *A;
86  // Copy of A scaled by block-diagonal inverse
87  HypreParMatrix A_s;
88 
89  HypreBoomerAMG *AIR_solver;
90  int blocksize;
91 
92 public:
93  AIR_prec(int blocksize_) : AIR_solver(NULL), blocksize(blocksize_) { }
94 
95  void SetOperator(const Operator &op)
96  {
97  width = op.Width();
98  height = op.Height();
99 
100  A = dynamic_cast<const HypreParMatrix *>(&op);
101  MFEM_VERIFY(A != NULL, "AIR_prec requires a HypreParMatrix.")
102 
103  // Scale A by block-diagonal inverse
104  BlockInverseScale(A, &A_s, NULL, NULL, blocksize,
105  BlockInverseScaleJob::MATRIX_ONLY);
106  delete AIR_solver;
107  AIR_solver = new HypreBoomerAMG(A_s);
108  AIR_solver->SetAdvectiveOptions(1, "", "FA");
109  AIR_solver->SetPrintLevel(0);
110  AIR_solver->SetMaxLevels(50);
111  }
112 
113  virtual void Mult(const Vector &x, Vector &y) const
114  {
115  // Scale the rhs by block inverse and solve system
116  HypreParVector z_s;
117  BlockInverseScale(A, NULL, &x, &z_s, blocksize,
118  BlockInverseScaleJob::RHS_ONLY);
119  AIR_solver->Mult(z_s, y);
120  }
121 
122  ~AIR_prec()
123  {
124  delete AIR_solver;
125  }
126 };
127 #endif
128 
129 
130 class DG_Solver : public Solver
131 {
132 private:
133  HypreParMatrix &M, &K;
134  SparseMatrix M_diag;
135  HypreParMatrix *A;
136  GMRESSolver linear_solver;
137  Solver *prec;
138  double dt;
139 public:
140  DG_Solver(HypreParMatrix &M_, HypreParMatrix &K_, const FiniteElementSpace &fes,
141  PrecType prec_type)
142  : M(M_),
143  K(K_),
144  A(NULL),
145  linear_solver(M.GetComm()),
146  dt(-1.0)
147  {
148  int block_size = fes.GetFE(0)->GetDof();
149  if (prec_type == PrecType::ILU)
150  {
151  prec = new BlockILU(block_size,
152  BlockILU::Reordering::MINIMUM_DISCARDED_FILL);
153  }
154  else if (prec_type == PrecType::AIR)
155  {
156 #if MFEM_HYPRE_VERSION >= 21800
157  prec = new AIR_prec(block_size);
158 #else
159  MFEM_ABORT("Must have MFEM_HYPRE_VERSION >= 21800 to use AIR.\n");
160 #endif
161  }
162  linear_solver.iterative_mode = false;
163  linear_solver.SetRelTol(1e-9);
164  linear_solver.SetAbsTol(0.0);
165  linear_solver.SetMaxIter(100);
166  linear_solver.SetPrintLevel(0);
167  linear_solver.SetPreconditioner(*prec);
168 
169  M.GetDiag(M_diag);
170  }
171 
172  void SetTimeStep(double dt_)
173  {
174  if (dt_ != dt)
175  {
176  dt = dt_;
177  // Form operator A = M - dt*K
178  delete A;
179  A = Add(-dt, K, 0.0, K);
180  SparseMatrix A_diag;
181  A->GetDiag(A_diag);
182  A_diag.Add(1.0, M_diag);
183  // this will also call SetOperator on the preconditioner
184  linear_solver.SetOperator(*A);
185  }
186  }
187 
188  void SetOperator(const Operator &op)
189  {
190  linear_solver.SetOperator(op);
191  }
192 
193  virtual void Mult(const Vector &x, Vector &y) const
194  {
195  linear_solver.Mult(x, y);
196  }
197 
198  ~DG_Solver()
199  {
200  delete prec;
201  delete A;
202  }
203 };
204 
205 
206 /** A time-dependent operator for the right-hand side of the ODE. The DG weak
207  form of du/dt = -v.grad(u) is M du/dt = K u + b, where M and K are the mass
208  and advection matrices, and b describes the flow on the boundary. This can
209  be written as a general ODE, du/dt = M^{-1} (K u + b), and this class is
210  used to evaluate the right-hand side. */
211 class FE_Evolution : public TimeDependentOperator
212 {
213 private:
214  OperatorHandle M, K;
215  const Vector &b;
216  Solver *M_prec;
217  CGSolver M_solver;
218  DG_Solver *dg_solver;
219 
220  mutable Vector z;
221 
222 public:
224  PrecType prec_type);
225 
226  virtual void Mult(const Vector &x, Vector &y) const;
227  virtual void ImplicitSolve(const double dt, const Vector &x, Vector &k);
228 
229  virtual ~FE_Evolution();
230 };
231 
232 
233 int main(int argc, char *argv[])
234 {
235  // 1. Initialize MPI and HYPRE.
236  Mpi::Init();
237  int num_procs = Mpi::WorldSize();
238  int myid = Mpi::WorldRank();
239  Hypre::Init();
240 
241  // 2. Parse command-line options.
242  problem = 0;
243  const char *mesh_file = "../data/periodic-hexagon.mesh";
244  int ser_ref_levels = 2;
245  int par_ref_levels = 0;
246  int order = 3;
247  bool pa = false;
248  bool ea = false;
249  bool fa = false;
250  const char *device_config = "cpu";
251  int ode_solver_type = 4;
252  double t_final = 10.0;
253  double dt = 0.01;
254  bool visualization = true;
255  bool visit = false;
256  bool paraview = false;
257  bool adios2 = false;
258  bool binary = false;
259  int vis_steps = 5;
260 #if MFEM_HYPRE_VERSION >= 21800
261  PrecType prec_type = PrecType::AIR;
262 #else
263  PrecType prec_type = PrecType::ILU;
264 #endif
265  int precision = 8;
266  cout.precision(precision);
267 
268  OptionsParser args(argc, argv);
269  args.AddOption(&mesh_file, "-m", "--mesh",
270  "Mesh file to use.");
271  args.AddOption(&problem, "-p", "--problem",
272  "Problem setup to use. See options in velocity_function().");
273  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
274  "Number of times to refine the mesh uniformly in serial.");
275  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
276  "Number of times to refine the mesh uniformly in parallel.");
277  args.AddOption(&order, "-o", "--order",
278  "Order (degree) of the finite elements.");
279  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
280  "--no-partial-assembly", "Enable Partial Assembly.");
281  args.AddOption(&ea, "-ea", "--element-assembly", "-no-ea",
282  "--no-element-assembly", "Enable Element Assembly.");
283  args.AddOption(&fa, "-fa", "--full-assembly", "-no-fa",
284  "--no-full-assembly", "Enable Full Assembly.");
285  args.AddOption(&device_config, "-d", "--device",
286  "Device configuration string, see Device::Configure().");
287  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
288  "ODE solver: 1 - Forward Euler,\n\t"
289  " 2 - RK2 SSP, 3 - RK3 SSP, 4 - RK4, 6 - RK6,\n\t"
290  " 11 - Backward Euler,\n\t"
291  " 12 - SDIRK23 (L-stable), 13 - SDIRK33,\n\t"
292  " 22 - Implicit Midpoint Method,\n\t"
293  " 23 - SDIRK23 (A-stable), 24 - SDIRK34");
294  args.AddOption(&t_final, "-tf", "--t-final",
295  "Final time; start time is 0.");
296  args.AddOption(&dt, "-dt", "--time-step",
297  "Time step.");
298  args.AddOption((int *)&prec_type, "-pt", "--prec-type", "Preconditioner for "
299  "implicit solves. 0 for ILU, 1 for pAIR-AMG.");
300  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
301  "--no-visualization",
302  "Enable or disable GLVis visualization.");
303  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
304  "--no-visit-datafiles",
305  "Save data files for VisIt (visit.llnl.gov) visualization.");
306  args.AddOption(&paraview, "-paraview", "--paraview-datafiles", "-no-paraview",
307  "--no-paraview-datafiles",
308  "Save data files for ParaView (paraview.org) visualization.");
309  args.AddOption(&adios2, "-adios2", "--adios2-streams", "-no-adios2",
310  "--no-adios2-streams",
311  "Save data using adios2 streams.");
312  args.AddOption(&binary, "-binary", "--binary-datafiles", "-ascii",
313  "--ascii-datafiles",
314  "Use binary (Sidre) or ascii format for VisIt data files.");
315  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
316  "Visualize every n-th timestep.");
317  args.Parse();
318  if (!args.Good())
319  {
320  if (Mpi::Root())
321  {
322  args.PrintUsage(cout);
323  }
324  return 1;
325  }
326  if (Mpi::Root())
327  {
328  args.PrintOptions(cout);
329  }
330 
331  Device device(device_config);
332  if (Mpi::Root()) { device.Print(); }
333 
334  // 3. Read the serial mesh from the given mesh file on all processors. We can
335  // handle geometrically periodic meshes in this code.
336  Mesh *mesh = new Mesh(mesh_file, 1, 1);
337  int dim = mesh->Dimension();
338 
339  // 4. Define the ODE solver used for time integration. Several explicit
340  // Runge-Kutta methods are available.
341  ODESolver *ode_solver = NULL;
342  switch (ode_solver_type)
343  {
344  // Explicit methods
345  case 1: ode_solver = new ForwardEulerSolver; break;
346  case 2: ode_solver = new RK2Solver(1.0); break;
347  case 3: ode_solver = new RK3SSPSolver; break;
348  case 4: ode_solver = new RK4Solver; break;
349  case 6: ode_solver = new RK6Solver; break;
350  // Implicit (L-stable) methods
351  case 11: ode_solver = new BackwardEulerSolver; break;
352  case 12: ode_solver = new SDIRK23Solver(2); break;
353  case 13: ode_solver = new SDIRK33Solver; break;
354  // Implicit A-stable methods (not L-stable)
355  case 22: ode_solver = new ImplicitMidpointSolver; break;
356  case 23: ode_solver = new SDIRK23Solver; break;
357  case 24: ode_solver = new SDIRK34Solver; break;
358  default:
359  if (Mpi::Root())
360  {
361  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
362  }
363  delete mesh;
364  return 3;
365  }
366 
367  // 5. Refine the mesh in serial to increase the resolution. In this example
368  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
369  // a command-line parameter. If the mesh is of NURBS type, we convert it
370  // to a (piecewise-polynomial) high-order mesh.
371  for (int lev = 0; lev < ser_ref_levels; lev++)
372  {
373  mesh->UniformRefinement();
374  }
375  if (mesh->NURBSext)
376  {
377  mesh->SetCurvature(max(order, 1));
378  }
379  mesh->GetBoundingBox(bb_min, bb_max, max(order, 1));
380 
381  // 6. Define the parallel mesh by a partitioning of the serial mesh. Refine
382  // this mesh further in parallel to increase the resolution. Once the
383  // parallel mesh is defined, the serial mesh can be deleted.
384  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
385  delete mesh;
386  for (int lev = 0; lev < par_ref_levels; lev++)
387  {
388  pmesh->UniformRefinement();
389  }
390 
391  // 7. Define the parallel discontinuous DG finite element space on the
392  // parallel refined mesh of the given polynomial order.
393  DG_FECollection fec(order, dim, BasisType::GaussLobatto);
394  ParFiniteElementSpace *fes = new ParFiniteElementSpace(pmesh, &fec);
395 
396  HYPRE_BigInt global_vSize = fes->GlobalTrueVSize();
397  if (Mpi::Root())
398  {
399  cout << "Number of unknowns: " << global_vSize << endl;
400  }
401 
402  // 8. Set up and assemble the parallel bilinear and linear forms (and the
403  // parallel hypre matrices) corresponding to the DG discretization. The
404  // DGTraceIntegrator involves integrals over mesh interior faces.
408 
409  ParBilinearForm *m = new ParBilinearForm(fes);
410  ParBilinearForm *k = new ParBilinearForm(fes);
411  if (pa)
412  {
413  m->SetAssemblyLevel(AssemblyLevel::PARTIAL);
414  k->SetAssemblyLevel(AssemblyLevel::PARTIAL);
415  }
416  else if (ea)
417  {
418  m->SetAssemblyLevel(AssemblyLevel::ELEMENT);
419  k->SetAssemblyLevel(AssemblyLevel::ELEMENT);
420  }
421  else if (fa)
422  {
423  m->SetAssemblyLevel(AssemblyLevel::FULL);
424  k->SetAssemblyLevel(AssemblyLevel::FULL);
425  }
426 
428  constexpr double alpha = -1.0;
431  new NonconservativeDGTraceIntegrator(velocity, alpha));
433  new NonconservativeDGTraceIntegrator(velocity, alpha));
434 
435  ParLinearForm *b = new ParLinearForm(fes);
436  b->AddBdrFaceIntegrator(
437  new BoundaryFlowIntegrator(inflow, velocity, alpha));
438 
439  int skip_zeros = 0;
440  m->Assemble();
441  k->Assemble(skip_zeros);
442  b->Assemble();
443  m->Finalize();
444  k->Finalize(skip_zeros);
445 
446 
447  HypreParVector *B = b->ParallelAssemble();
448 
449  // 9. Define the initial conditions, save the corresponding grid function to
450  // a file and (optionally) save data in the VisIt format and initialize
451  // GLVis visualization.
452  ParGridFunction *u = new ParGridFunction(fes);
453  u->ProjectCoefficient(u0);
454  HypreParVector *U = u->GetTrueDofs();
455 
456  {
457  ostringstream mesh_name, sol_name;
458  mesh_name << "ex9-mesh." << setfill('0') << setw(6) << myid;
459  sol_name << "ex9-init." << setfill('0') << setw(6) << myid;
460  ofstream omesh(mesh_name.str().c_str());
461  omesh.precision(precision);
462  pmesh->Print(omesh);
463  ofstream osol(sol_name.str().c_str());
464  osol.precision(precision);
465  u->Save(osol);
466  }
467 
468  // Create data collection for solution output: either VisItDataCollection for
469  // ascii data files, or SidreDataCollection for binary data files.
470  DataCollection *dc = NULL;
471  if (visit)
472  {
473  if (binary)
474  {
475 #ifdef MFEM_USE_SIDRE
476  dc = new SidreDataCollection("Example9-Parallel", pmesh);
477 #else
478  MFEM_ABORT("Must build with MFEM_USE_SIDRE=YES for binary output.");
479 #endif
480  }
481  else
482  {
483  dc = new VisItDataCollection("Example9-Parallel", pmesh);
484  dc->SetPrecision(precision);
485  // To save the mesh using MFEM's parallel mesh format:
486  // dc->SetFormat(DataCollection::PARALLEL_FORMAT);
487  }
488  dc->RegisterField("solution", u);
489  dc->SetCycle(0);
490  dc->SetTime(0.0);
491  dc->Save();
492  }
493 
494  ParaViewDataCollection *pd = NULL;
495  if (paraview)
496  {
497  pd = new ParaViewDataCollection("Example9P", pmesh);
498  pd->SetPrefixPath("ParaView");
499  pd->RegisterField("solution", u);
500  pd->SetLevelsOfDetail(order);
501  pd->SetDataFormat(VTKFormat::BINARY);
502  pd->SetHighOrderOutput(true);
503  pd->SetCycle(0);
504  pd->SetTime(0.0);
505  pd->Save();
506  }
507 
508  // Optionally output a BP (binary pack) file using ADIOS2. This can be
509  // visualized with the ParaView VTX reader.
510 #ifdef MFEM_USE_ADIOS2
511  ADIOS2DataCollection *adios2_dc = NULL;
512  if (adios2)
513  {
514  std::string postfix(mesh_file);
515  postfix.erase(0, std::string("../data/").size() );
516  postfix += "_o" + std::to_string(order);
517  const std::string collection_name = "ex9-p-" + postfix + ".bp";
518 
519  adios2_dc = new ADIOS2DataCollection(MPI_COMM_WORLD, collection_name, pmesh);
520  // output data substreams are half the number of mpi processes
521  adios2_dc->SetParameter("SubStreams", std::to_string(num_procs/2) );
522  // adios2_dc->SetLevelsOfDetail(2);
523  adios2_dc->RegisterField("solution", u);
524  adios2_dc->SetCycle(0);
525  adios2_dc->SetTime(0.0);
526  adios2_dc->Save();
527  }
528 #endif
529 
530  socketstream sout;
531  if (visualization)
532  {
533  char vishost[] = "localhost";
534  int visport = 19916;
535  sout.open(vishost, visport);
536  if (!sout)
537  {
538  if (Mpi::Root())
539  {
540  cout << "Unable to connect to GLVis server at "
541  << vishost << ':' << visport << endl;
542  }
543  visualization = false;
544  if (Mpi::Root())
545  {
546  cout << "GLVis visualization disabled.\n";
547  }
548  }
549  else
550  {
551  sout << "parallel " << num_procs << " " << myid << "\n";
552  sout.precision(precision);
553  sout << "solution\n" << *pmesh << *u;
554  sout << "pause\n";
555  sout << flush;
556  if (Mpi::Root())
557  {
558  cout << "GLVis visualization paused."
559  << " Press space (in the GLVis window) to resume it.\n";
560  }
561  }
562  }
563 
564  // 10. Define the time-dependent evolution operator describing the ODE
565  // right-hand side, and perform time-integration (looping over the time
566  // iterations, ti, with a time-step dt).
567  FE_Evolution adv(*m, *k, *B, prec_type);
568 
569  double t = 0.0;
570  adv.SetTime(t);
571  ode_solver->Init(adv);
572 
573  bool done = false;
574  for (int ti = 0; !done; )
575  {
576  double dt_real = min(dt, t_final - t);
577  ode_solver->Step(*U, t, dt_real);
578  ti++;
579 
580  done = (t >= t_final - 1e-8*dt);
581 
582  if (done || ti % vis_steps == 0)
583  {
584  if (Mpi::Root())
585  {
586  cout << "time step: " << ti << ", time: " << t << endl;
587  }
588 
589  // 11. Extract the parallel grid function corresponding to the finite
590  // element approximation U (the local solution on each processor).
591  *u = *U;
592 
593  if (visualization)
594  {
595  sout << "parallel " << num_procs << " " << myid << "\n";
596  sout << "solution\n" << *pmesh << *u << flush;
597  }
598 
599  if (visit)
600  {
601  dc->SetCycle(ti);
602  dc->SetTime(t);
603  dc->Save();
604  }
605 
606  if (paraview)
607  {
608  pd->SetCycle(ti);
609  pd->SetTime(t);
610  pd->Save();
611  }
612 
613 #ifdef MFEM_USE_ADIOS2
614  // transient solutions can be visualized with ParaView
615  if (adios2)
616  {
617  adios2_dc->SetCycle(ti);
618  adios2_dc->SetTime(t);
619  adios2_dc->Save();
620  }
621 #endif
622  }
623  }
624 
625  // 12. Save the final solution in parallel. This output can be viewed later
626  // using GLVis: "glvis -np <np> -m ex9-mesh -g ex9-final".
627  {
628  *u = *U;
629  ostringstream sol_name;
630  sol_name << "ex9-final." << setfill('0') << setw(6) << myid;
631  ofstream osol(sol_name.str().c_str());
632  osol.precision(precision);
633  u->Save(osol);
634  }
635 
636  // 13. Free the used memory.
637  delete U;
638  delete u;
639  delete B;
640  delete b;
641  delete k;
642  delete m;
643  delete fes;
644  delete pmesh;
645  delete ode_solver;
646  delete pd;
647 #ifdef MFEM_USE_ADIOS2
648  if (adios2)
649  {
650  delete adios2_dc;
651  }
652 #endif
653  delete dc;
654 
655  return 0;
656 }
657 
658 
659 // Implementation of class FE_Evolution
661  const Vector &b_, PrecType prec_type)
662  : TimeDependentOperator(M_.Height()), b(b_),
663  M_solver(M_.ParFESpace()->GetComm()),
664  z(M_.Height())
665 {
666  if (M_.GetAssemblyLevel()==AssemblyLevel::LEGACY)
667  {
668  M.Reset(M_.ParallelAssemble(), true);
669  K.Reset(K_.ParallelAssemble(), true);
670  }
671  else
672  {
673  M.Reset(&M_, false);
674  K.Reset(&K_, false);
675  }
676 
677  M_solver.SetOperator(*M);
678 
679  Array<int> ess_tdof_list;
680  if (M_.GetAssemblyLevel()==AssemblyLevel::LEGACY)
681  {
682  HypreParMatrix &M_mat = *M.As<HypreParMatrix>();
683  HypreParMatrix &K_mat = *K.As<HypreParMatrix>();
684  HypreSmoother *hypre_prec = new HypreSmoother(M_mat, HypreSmoother::Jacobi);
685  M_prec = hypre_prec;
686 
687  dg_solver = new DG_Solver(M_mat, K_mat, *M_.FESpace(), prec_type);
688  }
689  else
690  {
691  M_prec = new OperatorJacobiSmoother(M_, ess_tdof_list);
692  dg_solver = NULL;
693  }
694 
695  M_solver.SetPreconditioner(*M_prec);
696  M_solver.iterative_mode = false;
697  M_solver.SetRelTol(1e-9);
698  M_solver.SetAbsTol(0.0);
699  M_solver.SetMaxIter(100);
700  M_solver.SetPrintLevel(0);
701 }
702 
703 // Solve the equation:
704 // u_t = M^{-1}(Ku + b),
705 // by solving associated linear system
706 // (M - dt*K) d = K*u + b
707 void FE_Evolution::ImplicitSolve(const double dt, const Vector &x, Vector &k)
708 {
709  K->Mult(x, z);
710  z += b;
711  dg_solver->SetTimeStep(dt);
712  dg_solver->Mult(z, k);
713 }
714 
715 void FE_Evolution::Mult(const Vector &x, Vector &y) const
716 {
717  // y = M^{-1} (K x + b)
718  K->Mult(x, z);
719  z += b;
720  M_solver.Mult(z, y);
721 }
722 
724 {
725  delete M_prec;
726  delete dg_solver;
727 }
728 
729 
730 // Velocity coefficient
731 void velocity_function(const Vector &x, Vector &v)
732 {
733  int dim = x.Size();
734 
735  // map to the reference [-1,1] domain
736  Vector X(dim);
737  for (int i = 0; i < dim; i++)
738  {
739  double center = (bb_min[i] + bb_max[i]) * 0.5;
740  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
741  }
742 
743  switch (problem)
744  {
745  case 0:
746  {
747  // Translations in 1D, 2D, and 3D
748  switch (dim)
749  {
750  case 1: v(0) = 1.0; break;
751  case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.); break;
752  case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
753  break;
754  }
755  break;
756  }
757  case 1:
758  case 2:
759  {
760  // Clockwise rotation in 2D around the origin
761  const double w = M_PI/2;
762  switch (dim)
763  {
764  case 1: v(0) = 1.0; break;
765  case 2: v(0) = w*X(1); v(1) = -w*X(0); break;
766  case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0; break;
767  }
768  break;
769  }
770  case 3:
771  {
772  // Clockwise twisting rotation in 2D around the origin
773  const double w = M_PI/2;
774  double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
775  d = d*d;
776  switch (dim)
777  {
778  case 1: v(0) = 1.0; break;
779  case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0); break;
780  case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0; break;
781  }
782  break;
783  }
784  }
785 }
786 
787 // Initial condition
788 double u0_function(const Vector &x)
789 {
790  int dim = x.Size();
791 
792  // map to the reference [-1,1] domain
793  Vector X(dim);
794  for (int i = 0; i < dim; i++)
795  {
796  double center = (bb_min[i] + bb_max[i]) * 0.5;
797  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
798  }
799 
800  switch (problem)
801  {
802  case 0:
803  case 1:
804  {
805  switch (dim)
806  {
807  case 1:
808  return exp(-40.*pow(X(0)-0.5,2));
809  case 2:
810  case 3:
811  {
812  double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
813  if (dim == 3)
814  {
815  const double s = (1. + 0.25*cos(2*M_PI*X(2)));
816  rx *= s;
817  ry *= s;
818  }
819  return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
820  erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
821  }
822  }
823  }
824  case 2:
825  {
826  double x_ = X(0), y_ = X(1), rho, phi;
827  rho = hypot(x_, y_);
828  phi = atan2(y_, x_);
829  return pow(sin(M_PI*rho),2)*sin(3*phi);
830  }
831  case 3:
832  {
833  const double f = M_PI;
834  return sin(f*X(0))*sin(f*X(1));
835  }
836  }
837  return 0.0;
838 }
839 
840 // Inflow boundary condition (zero for the problems considered in this example)
841 double inflow_function(const Vector &x)
842 {
843  switch (problem)
844  {
845  case 0:
846  case 1:
847  case 2:
848  case 3: return 0.0;
849  }
850  return 0.0;
851 }
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
int visport
Conjugate gradient method.
Definition: solvers.hpp:493
double u0_function(const Vector &x)
Definition: ex9p.cpp:788
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetDataFormat(VTKFormat fmt)
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:331
FiniteElementSpace * FESpace()
Return the FE space associated with the BilinearForm.
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1020
Helper class for ParaView visualization data.
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:462
Base abstract class for first order time dependent operators.
Definition: operator.hpp:316
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:475
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
Definition: mesh.cpp:136
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
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]...
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:718
int Size() const
Returns the size of the vector.
Definition: vector.hpp:197
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:980
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:23
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:279
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
virtual void SetTime(const double t_)
Set the current time.
Definition: operator.hpp:360
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:159
Abstract parallel finite element space.
Definition: pfespace.hpp:28
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:686
STL namespace.
void BlockInverseScale(const HypreParMatrix *A, HypreParMatrix *C, const Vector *b, HypreParVector *d, int blocksize, BlockInverseScaleJob job)
Definition: hypre.cpp:2766
AssemblyLevel GetAssemblyLevel() const
Returns the assembly level.
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:411
Vector bb_max
Definition: ex9p.cpp:68
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:2841
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
virtual void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:179
Data collection with Sidre routines following the Conduit mesh blueprint specification.
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:2320
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1590
Class for parallel linear form.
Definition: plinearform.hpp:26
PrecType
Definition: ex9p.cpp:71
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
char vishost[]
Data type sparse matrix.
Definition: sparsemat.hpp:50
BlockInverseScaleJob
Definition: hypre.hpp:915
virtual void Save()
Save the collection to disk.
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:302
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:10232
Vector bb_min
Definition: ex9p.cpp:68
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
Definition: solvers.hpp:201
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
Definition: mesh.cpp:5635
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void Mult(const Vector &x, Vector &y) const
Perform the action of the operator: y = k = f(x, t), where k solves the algebraic equation F(x...
Definition: ex18.hpp:107
void SetParameter(const std::string key, const std::string value) noexcept
Parallel smoothers in hypre.
Definition: hypre.hpp:971
void SetHighOrderOutput(bool high_order_output_)
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:1481
void SetTime(double t)
Set physical time (for time-dependent simulations)
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:281
void Add(const int i, const int j, const double val)
Definition: sparsemat.cpp:2787
A general vector function coefficient.
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:161
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:163
void SetAbsTol(double atol)
Definition: solvers.hpp:200
virtual void ImplicitSolve(const double dt, const Vector &x, Vector &k)
Solve the equation: k = f(x + dt k, t), for the unknown k at the current time t.
Definition: ex9.cpp:484
void SetRelTol(double rtol)
Definition: solvers.hpp:199
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
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 &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:150
HYPRE_Int HYPRE_BigInt
GMRES method.
Definition: solvers.hpp:525
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
Implicit midpoint method. A-stable, not L-stable.
Definition: ode.hpp:424
FE_Evolution(FiniteElementSpace &vfes_, Operator &A_, SparseMatrix &Aflux_)
Definition: ex18.hpp:81
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:277
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
int dim
Definition: ex24.cpp:53
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator. Assumes ownership of bfi.
int main(int argc, char *argv[])
Definition: ex9p.cpp:233
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
Class for parallel bilinear form.
virtual ~FE_Evolution()
Definition: ex18.hpp:43
int open(const char hostname[], int port)
Open the socket stream on &#39;port&#39; at &#39;hostname&#39;.
RefCoord t[3]
double inflow_function(const Vector &x)
Definition: ex9p.cpp:841
void SetLevelsOfDetail(int levels_of_detail_)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:507
const double alpha
Definition: ex15.cpp:369
A general function coefficient.
Vector data type.
Definition: vector.hpp:58
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:173
void velocity_function(const Vector &x, Vector &v)
Definition: ex9p.cpp:731
virtual void Save() override
RefCoord s[3]
void Print(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4825
double u(const Vector &xvec)
Definition: lor_mms.hpp:22
Base class for solvers.
Definition: operator.hpp:682
Class for parallel grid function.
Definition: pgridfunc.hpp:32
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:121
The classical forward Euler method.
Definition: ode.hpp:117
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:343
Class for parallel meshes.
Definition: pmesh.hpp:32
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:18
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0) const
Computes y = alpha * A * x + beta * y.
Definition: hypre.cpp:1736
int problem
Definition: ex9p.cpp:56
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator. Assumes ownership of bfi.
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition: fe_coll.hpp:327
double f(const Vector &p)
alpha (q . grad u, v)