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