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