MFEM  v4.3.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 // PETSc Modification
3 //
4 // Compile with: make ex9p
5 //
6 // Sample runs:
7 // mpirun -np 4 ex9p -m ../../data/periodic-hexagon.mesh --petscopts rc_ex9p_expl
8 // mpirun -np 4 ex9p -m ../../data/periodic-hexagon.mesh --petscopts rc_ex9p_impl -implicit
9 //
10 // Description: This example code solves the time-dependent advection equation
11 // du/dt + v.grad(u) = 0, where v is a given fluid velocity, and
12 // u0(x)=u(0,x) is a given initial condition.
13 //
14 // The example demonstrates the use of Discontinuous Galerkin (DG)
15 // bilinear forms in MFEM (face integrators), the use of explicit
16 // ODE time integrators, the definition of periodic boundary
17 // conditions through periodic meshes, as well as the use of GLVis
18 // for persistent visualization of a time-evolving solution. The
19 // saving of time-dependent data files for external visualization
20 // with VisIt (visit.llnl.gov) is also illustrated.
21 //
22 // The example also demonstrates how to use PETSc ODE solvers and
23 // customize them by command line (see rc_ex9p_expl and
24 // rc_ex9p_impl). The split in left-hand side and right-hand side
25 // of the TimeDependentOperator is amenable for IMEX methods.
26 // When using fully implicit methods, just the left-hand side of
27 // the operator should be provided for efficiency reasons when
28 // assembling the Jacobians. Here, we provide two Jacobian
29 // routines just to illustrate the capabilities of the
30 // PetscODESolver class. We also show how to monitor the time
31 // dependent solution inside a call to PetscODESolver:Mult.
32 
33 #include "mfem.hpp"
34 #include <fstream>
35 #include <iostream>
36 
37 #ifndef MFEM_USE_PETSC
38 #error This example requires that MFEM is built with MFEM_USE_PETSC=YES
39 #endif
40 
41 using namespace std;
42 using namespace mfem;
43 
44 // Choice for the problem setup. The fluid velocity, initial condition and
45 // inflow boundary condition are chosen based on this parameter.
46 int problem;
47 
48 // Velocity coefficient
49 void velocity_function(const Vector &x, Vector &v);
50 
51 // Initial condition
52 double u0_function(const Vector &x);
53 
54 // Inflow boundary condition
55 double inflow_function(const Vector &x);
56 
57 // Mesh bounding box
59 
60 
61 /** A time-dependent operator for the ODE as F(u,du/dt,t) = G(u,t)
62  The DG weak form of du/dt = -v.grad(u) is M du/dt = K u + b, where M and K are the mass
63  and advection matrices, and b describes the flow on the boundary. This can
64  be also written as a general ODE with the right-hand side only as
65  du/dt = M^{-1} (K u + b).
66  This class is used to evaluate the right-hand side and the left-hand side. */
68 {
69 private:
70  OperatorHandle M, K;
71  const Vector &b;
72  MPI_Comm comm;
73  Solver *M_prec;
74  CGSolver M_solver;
75  AssemblyLevel MAlev,KAlev;
76 
77  mutable Vector z;
78  mutable PetscParMatrix* iJacobian;
79  mutable PetscParMatrix* rJacobian;
80 
81 public:
83  bool implicit);
84 
85  virtual void ExplicitMult(const Vector &x, Vector &y) const;
86  virtual void ImplicitMult(const Vector &x, const Vector &xp, Vector &y) const;
87  virtual void Mult(const Vector &x, Vector &y) const;
88  virtual Operator& GetExplicitGradient(const Vector &x) const;
89  virtual Operator& GetImplicitGradient(const Vector &x, const Vector &xp,
90  double shift) const;
91  virtual ~FE_Evolution() { delete iJacobian; delete rJacobian; }
92 };
93 
94 
95 // Monitor the solution at time step "step", explicitly in the time loop
96 class UserMonitor : public PetscSolverMonitor
97 {
98 private:
99  socketstream& sout;
100  ParMesh* pmesh;
102  int vt;
103  bool pause;
104 
105 public:
106  UserMonitor(socketstream& s_, ParMesh* m_, ParGridFunction* u_, int vt_) :
107  PetscSolverMonitor(true,false), sout(s_), pmesh(m_), u(u_), vt(vt_),
108  pause(true) {}
109 
110  void MonitorSolution(PetscInt step, PetscReal norm, const Vector &X)
111  {
112  if (step % vt == 0)
113  {
114  int num_procs, myid;
115 
116  *u = X;
117  MPI_Comm_size(pmesh->GetComm(),&num_procs);
118  MPI_Comm_rank(pmesh->GetComm(),&myid);
119  sout << "parallel " << num_procs << " " << myid << "\n";
120  sout << "solution\n" << *pmesh << *u;
121  if (pause) { sout << "pause\n"; }
122  sout << flush;
123  if (pause)
124  {
125  pause = false;
126  if (myid == 0)
127  {
128  cout << "GLVis visualization paused."
129  << " Press space (in the GLVis window) to resume it.\n";
130  }
131  }
132  }
133  }
134 };
135 
136 
137 int main(int argc, char *argv[])
138 {
139  // 1. Initialize MPI.
140  int num_procs, myid;
141  MPI_Init(&argc, &argv);
142  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
143  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
144 
145  // 2. Parse command-line options.
146  problem = 0;
147  const char *mesh_file = "../../data/periodic-hexagon.mesh";
148  int ser_ref_levels = 2;
149  int par_ref_levels = 0;
150  int order = 3;
151  bool pa = false;
152  bool ea = false;
153  bool fa = false;
154  const char *device_config = "cpu";
155  int ode_solver_type = 4;
156  double t_final = 10.0;
157  double dt = 0.01;
158  bool visualization = true;
159  bool visit = false;
160  bool binary = false;
161  int vis_steps = 5;
162  bool use_petsc = true;
163  bool implicit = false;
164  bool use_step = true;
165  const char *petscrc_file = "";
166 
167  int precision = 8;
168  cout.precision(precision);
169 
170  OptionsParser args(argc, argv);
171  args.AddOption(&mesh_file, "-m", "--mesh",
172  "Mesh file to use.");
173  args.AddOption(&problem, "-p", "--problem",
174  "Problem setup to use. See options in velocity_function().");
175  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
176  "Number of times to refine the mesh uniformly in serial.");
177  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
178  "Number of times to refine the mesh uniformly in parallel.");
179  args.AddOption(&order, "-o", "--order",
180  "Order (degree) of the finite elements.");
181  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
182  "--no-partial-assembly", "Enable Partial Assembly.");
183  args.AddOption(&ea, "-ea", "--element-assembly", "-no-ea",
184  "--no-element-assembly", "Enable Element Assembly.");
185  args.AddOption(&fa, "-fa", "--full-assembly", "-no-fa",
186  "--no-full-assembly", "Enable Full Assembly.");
187  args.AddOption(&device_config, "-d", "--device",
188  "Device configuration string, see Device::Configure().");
189  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
190  "ODE solver: 1 - Forward Euler,\n\t"
191  " 2 - RK2 SSP, 3 - RK3 SSP, 4 - RK4, 6 - RK6.");
192  args.AddOption(&t_final, "-tf", "--t-final",
193  "Final time; start time is 0.");
194  args.AddOption(&dt, "-dt", "--time-step",
195  "Time step.");
196  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
197  "--no-visualization",
198  "Enable or disable GLVis visualization.");
199  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
200  "--no-visit-datafiles",
201  "Save data files for VisIt (visit.llnl.gov) visualization.");
202  args.AddOption(&binary, "-binary", "--binary-datafiles", "-ascii",
203  "--ascii-datafiles",
204  "Use binary (Sidre) or ascii format for VisIt data files.");
205  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
206  "Visualize every n-th timestep.");
207  args.AddOption(&use_petsc, "-usepetsc", "--usepetsc", "-no-petsc",
208  "--no-petsc",
209  "Use or not PETSc to solve the ODE system.");
210  args.AddOption(&petscrc_file, "-petscopts", "--petscopts",
211  "PetscOptions file to use.");
212  args.AddOption(&use_step, "-usestep", "--usestep", "-no-step",
213  "--no-step",
214  "Use the Step() or Run() method to solve the ODE system.");
215  args.AddOption(&implicit, "-implicit", "--implicit", "-no-implicit",
216  "--no-implicit",
217  "Use or not an implicit method in PETSc to solve the ODE system.");
218  args.Parse();
219  if (!args.Good())
220  {
221  if (myid == 0)
222  {
223  args.PrintUsage(cout);
224  }
225  MPI_Finalize();
226  return 1;
227  }
228  if (myid == 0)
229  {
230  args.PrintOptions(cout);
231  }
232 
233  Device device(device_config);
234  if (myid == 0) { device.Print(); }
235 
236  // 3. Read the serial mesh from the given mesh file on all processors. We can
237  // handle geometrically periodic meshes in this code.
238  Mesh *mesh = new Mesh(mesh_file, 1, 1);
239  int dim = mesh->Dimension();
240 
241  // 4. Define the ODE solver used for time integration. Several explicit
242  // Runge-Kutta methods are available.
243  ODESolver *ode_solver = NULL;
244  PetscODESolver *pode_solver = NULL;
245  UserMonitor *pmon = NULL;
246  if (!use_petsc)
247  {
248  switch (ode_solver_type)
249  {
250  case 1: ode_solver = new ForwardEulerSolver; break;
251  case 2: ode_solver = new RK2Solver(1.0); break;
252  case 3: ode_solver = new RK3SSPSolver; break;
253  case 4: ode_solver = new RK4Solver; break;
254  case 6: ode_solver = new RK6Solver; break;
255  default:
256  if (myid == 0)
257  {
258  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
259  }
260  MPI_Finalize();
261  return 3;
262  }
263  }
264  else
265  {
266  // When using PETSc, we just create the ODE solver. We use command line
267  // customization to select a specific solver.
268  MFEMInitializePetsc(NULL, NULL, petscrc_file, NULL);
269  ode_solver = pode_solver = new PetscODESolver(MPI_COMM_WORLD);
270  }
271 
272  // 5. Refine the mesh in serial to increase the resolution. In this example
273  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
274  // a command-line parameter. If the mesh is of NURBS type, we convert it
275  // to a (piecewise-polynomial) high-order mesh.
276  for (int lev = 0; lev < ser_ref_levels; lev++)
277  {
278  mesh->UniformRefinement();
279  }
280  if (mesh->NURBSext)
281  {
282  mesh->SetCurvature(max(order, 1));
283  }
284  mesh->GetBoundingBox(bb_min, bb_max, max(order, 1));
285 
286  // 6. Define the parallel mesh by a partitioning of the serial mesh. Refine
287  // this mesh further in parallel to increase the resolution. Once the
288  // parallel mesh is defined, the serial mesh can be deleted.
289  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
290  delete mesh;
291  for (int lev = 0; lev < par_ref_levels; lev++)
292  {
293  pmesh->UniformRefinement();
294  }
295 
296  // 7. Define the parallel discontinuous DG finite element space on the
297  // parallel refined mesh of the given polynomial order.
298  DG_FECollection fec(order, dim, BasisType::GaussLobatto);
299  ParFiniteElementSpace *fes = new ParFiniteElementSpace(pmesh, &fec);
300 
301  HYPRE_BigInt global_vSize = fes->GlobalTrueVSize();
302  if (myid == 0)
303  {
304  cout << "Number of unknowns: " << global_vSize << endl;
305  }
306 
307  // 8. Set up and assemble the parallel bilinear and linear forms (and the
308  // parallel hypre matrices) corresponding to the DG discretization. The
309  // DGTraceIntegrator involves integrals over mesh interior faces.
313 
314  ParBilinearForm *m = new ParBilinearForm(fes);
315  ParBilinearForm *k = new ParBilinearForm(fes);
316  if (pa)
317  {
318  m->SetAssemblyLevel(AssemblyLevel::PARTIAL);
319  k->SetAssemblyLevel(AssemblyLevel::PARTIAL);
320  }
321  else if (ea)
322  {
323  m->SetAssemblyLevel(AssemblyLevel::ELEMENT);
324  k->SetAssemblyLevel(AssemblyLevel::ELEMENT);
325  }
326  else if (fa)
327  {
328  m->SetAssemblyLevel(AssemblyLevel::FULL);
329  k->SetAssemblyLevel(AssemblyLevel::FULL);
330  }
331 
333  k->AddDomainIntegrator(new ConvectionIntegrator(velocity, -1.0));
335  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
337  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
338 
339  ParLinearForm *b = new ParLinearForm(fes);
341  new BoundaryFlowIntegrator(inflow, velocity, -1.0, -0.5));
342 
343  int skip_zeros = 0;
344  m->Assemble();
345  k->Assemble(skip_zeros);
346  b->Assemble();
347  m->Finalize();
348  k->Finalize(skip_zeros);
349 
350 
352 
353  // 9. Define the initial conditions, save the corresponding grid function to
354  // a file and (optionally) save data in the VisIt format and initialize
355  // GLVis visualization.
356  ParGridFunction *u = new ParGridFunction(fes);
357  u->ProjectCoefficient(u0);
358  HypreParVector *U = u->GetTrueDofs();
359 
360  {
361  ostringstream mesh_name, sol_name;
362  mesh_name << "ex9-mesh." << setfill('0') << setw(6) << myid;
363  sol_name << "ex9-init." << setfill('0') << setw(6) << myid;
364  ofstream omesh(mesh_name.str().c_str());
365  omesh.precision(precision);
366  pmesh->Print(omesh);
367  ofstream osol(sol_name.str().c_str());
368  osol.precision(precision);
369  u->Save(osol);
370  }
371 
372  // Create data collection for solution output: either VisItDataCollection for
373  // ascii data files, or SidreDataCollection for binary data files.
374  DataCollection *dc = NULL;
375  if (visit)
376  {
377  if (binary)
378  {
379 #ifdef MFEM_USE_SIDRE
380  dc = new SidreDataCollection("Example9-Parallel", pmesh);
381 #else
382  MFEM_ABORT("Must build with MFEM_USE_SIDRE=YES for binary output.");
383 #endif
384  }
385  else
386  {
387  dc = new VisItDataCollection("Example9-Parallel", pmesh);
388  dc->SetPrecision(precision);
389  }
390  dc->RegisterField("solution", u);
391  dc->SetCycle(0);
392  dc->SetTime(0.0);
393  dc->Save();
394  }
395 
396  socketstream sout;
397  if (visualization)
398  {
399  char vishost[] = "localhost";
400  int visport = 19916;
401  sout.open(vishost, visport);
402  if (!sout)
403  {
404  if (myid == 0)
405  cout << "Unable to connect to GLVis server at "
406  << vishost << ':' << visport << endl;
407  visualization = false;
408  if (myid == 0)
409  {
410  cout << "GLVis visualization disabled.\n";
411  }
412  }
413  else if (use_step)
414  {
415  sout << "parallel " << num_procs << " " << myid << "\n";
416  sout.precision(precision);
417  sout << "solution\n" << *pmesh << *u;
418  sout << "pause\n";
419  sout << flush;
420  if (myid == 0)
421  cout << "GLVis visualization paused."
422  << " Press space (in the GLVis window) to resume it.\n";
423  }
424  else if (use_petsc)
425  {
426  // Set the monitoring routine for the PetscODESolver.
427  sout.precision(precision);
428  pmon = new UserMonitor(sout,pmesh,u,vis_steps);
429  pode_solver->SetMonitor(pmon);
430  }
431  }
432 
433  // 10. Define the time-dependent evolution operator describing the ODE
434  FE_Evolution *adv = new FE_Evolution(*m, *k, *B, implicit);
435 
436  double t = 0.0;
437  adv->SetTime(t);
438  if (use_petsc)
439  {
440  pode_solver->Init(*adv,PetscODESolver::ODE_SOLVER_LINEAR);
441  }
442  else
443  {
444  ode_solver->Init(*adv);
445  }
446 
447  // Explicitly perform time-integration (looping over the time iterations, ti,
448  // with a time-step dt), or use the Run method of the ODE solver class.
449  if (use_step)
450  {
451  bool done = false;
452  for (int ti = 0; !done; )
453  {
454  // We cannot match exactly the time history of the Run method
455  // since we are explictly telling PETSc to use a time step
456  double dt_real = min(dt, t_final - t);
457  ode_solver->Step(*U, t, dt_real);
458  ti++;
459 
460  done = (t >= t_final - 1e-8*dt);
461 
462  if (done || ti % vis_steps == 0)
463  {
464  if (myid == 0)
465  {
466  cout << "time step: " << ti << ", time: " << t << endl;
467  }
468  // 11. Extract the parallel grid function corresponding to the finite
469  // element approximation U (the local solution on each processor).
470  *u = *U;
471 
472  if (visualization)
473  {
474  sout << "parallel " << num_procs << " " << myid << "\n";
475  sout << "solution\n" << *pmesh << *u << flush;
476  }
477 
478  if (visit)
479  {
480  dc->SetCycle(ti);
481  dc->SetTime(t);
482  dc->Save();
483  }
484  }
485  }
486  }
487  else { ode_solver->Run(*U, t, dt, t_final); }
488 
489  // 12. Save the final solution in parallel. This output can be viewed later
490  // using GLVis: "glvis -np <np> -m ex9-mesh -g ex9-final".
491  {
492  *u = *U;
493  ostringstream sol_name;
494  sol_name << "ex9-final." << setfill('0') << setw(6) << myid;
495  ofstream osol(sol_name.str().c_str());
496  osol.precision(precision);
497  u->Save(osol);
498  }
499 
500  // 13. Free the used memory.
501  delete U;
502  delete u;
503  delete B;
504  delete b;
505  delete k;
506  delete m;
507  delete fes;
508  delete pmesh;
509  delete ode_solver;
510  delete dc;
511  delete adv;
512 
513  delete pmon;
514 
515  // We finalize PETSc
516  if (use_petsc) { MFEMFinalizePetsc(); }
517 
518  MPI_Finalize();
519  return 0;
520 }
521 
522 
523 // Implementation of class FE_Evolution
525  const Vector &b_,bool M_in_lhs)
526  : TimeDependentOperator(M_.Height(), 0.0,
527  M_in_lhs ? TimeDependentOperator::IMPLICIT
528  : TimeDependentOperator::EXPLICIT),
529  b(b_), comm(M_.ParFESpace()->GetComm()), M_solver(comm), z(M_.Height()),
530  iJacobian(NULL), rJacobian(NULL)
531 {
532  MAlev = M_.GetAssemblyLevel();
533  KAlev = K_.GetAssemblyLevel();
534  if (M_.GetAssemblyLevel()==AssemblyLevel::LEGACY)
535  {
536  M.Reset(M_.ParallelAssemble(), true);
537  K.Reset(K_.ParallelAssemble(), true);
538  }
539  else
540  {
541  M.Reset(&M_, false);
542  K.Reset(&K_, false);
543  }
544 
545  M_solver.SetOperator(*M);
546 
547  Array<int> ess_tdof_list;
548  if (M_.GetAssemblyLevel()==AssemblyLevel::LEGACY)
549  {
550  HypreParMatrix &M_mat = *M.As<HypreParMatrix>();
551 
552  HypreSmoother *hypre_prec = new HypreSmoother(M_mat, HypreSmoother::Jacobi);
553  M_prec = hypre_prec;
554  }
555  else
556  {
557  M_prec = new OperatorJacobiSmoother(M_, ess_tdof_list);
558  }
559 
560  M_solver.SetPreconditioner(*M_prec);
561  M_solver.iterative_mode = false;
562  M_solver.SetRelTol(1e-9);
563  M_solver.SetAbsTol(0.0);
564  M_solver.SetMaxIter(100);
565  M_solver.SetPrintLevel(0);
566 }
567 
568 // RHS evaluation
569 void FE_Evolution::ExplicitMult(const Vector &x, Vector &y) const
570 {
571  if (isExplicit())
572  {
573  // y = M^{-1} (K x + b)
574  K->Mult(x, z);
575  z += b;
576  M_solver.Mult(z, y);
577  }
578  else
579  {
580  // y = K x + b
581  K->Mult(x, y);
582  y += b;
583  }
584 }
585 
586 // LHS evaluation
587 void FE_Evolution::ImplicitMult(const Vector &x, const Vector &xp,
588  Vector &y) const
589 {
590  if (isImplicit())
591  {
592  M->Mult(xp, y);
593  }
594  else
595  {
596  y = xp;
597  }
598 }
599 
600 void FE_Evolution::Mult(const Vector &x, Vector &y) const
601 {
602  // y = M^{-1} (K x + b)
603  K->Mult(x, z);
604  z += b;
605  M_solver.Mult(z, y);
606 }
607 
608 // RHS Jacobian
610 {
611  delete rJacobian;
612  Operator::Type otype = (KAlev == AssemblyLevel::LEGACY ?
613  Operator::PETSC_MATAIJ : Operator::ANY_TYPE);
614  if (isImplicit())
615  {
616  rJacobian = new PetscParMatrix(comm, K.Ptr(), otype);
617  }
618  else
619  {
620  mfem_error("FE_Evolution::GetExplicitGradient(x): Capability not coded!");
621  }
622  return *rJacobian;
623 }
624 
625 // LHS Jacobian, evaluated as shift*F_du/dt + F_u
627  double shift) const
628 {
629  Operator::Type otype = (MAlev == AssemblyLevel::LEGACY ?
630  Operator::PETSC_MATAIJ : Operator::ANY_TYPE);
631  delete iJacobian;
632  if (isImplicit())
633  {
634  iJacobian = new PetscParMatrix(comm, M.Ptr(), otype);
635  *iJacobian *= shift;
636  }
637  else
638  {
639  mfem_error("FE_Evolution::GetImplicitGradient(x,xp,shift):"
640  " Capability not coded!");
641  }
642  return *iJacobian;
643 }
644 
645 // Velocity coefficient
646 void velocity_function(const Vector &x, Vector &v)
647 {
648  int dim = x.Size();
649 
650  // map to the reference [-1,1] domain
651  Vector X(dim);
652  for (int i = 0; i < dim; i++)
653  {
654  double center = (bb_min[i] + bb_max[i]) * 0.5;
655  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
656  }
657 
658  switch (problem)
659  {
660  case 0:
661  {
662  // Translations in 1D, 2D, and 3D
663  switch (dim)
664  {
665  case 1: v(0) = 1.0; break;
666  case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.); break;
667  case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
668  break;
669  }
670  break;
671  }
672  case 1:
673  case 2:
674  {
675  // Clockwise rotation in 2D around the origin
676  const double w = M_PI/2;
677  switch (dim)
678  {
679  case 1: v(0) = 1.0; break;
680  case 2: v(0) = w*X(1); v(1) = -w*X(0); break;
681  case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0; break;
682  }
683  break;
684  }
685  case 3:
686  {
687  // Clockwise twisting rotation in 2D around the origin
688  const double w = M_PI/2;
689  double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
690  d = d*d;
691  switch (dim)
692  {
693  case 1: v(0) = 1.0; break;
694  case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0); break;
695  case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0; break;
696  }
697  break;
698  }
699  }
700 }
701 
702 // Initial condition
703 double u0_function(const Vector &x)
704 {
705  int dim = x.Size();
706 
707  // map to the reference [-1,1] domain
708  Vector X(dim);
709  for (int i = 0; i < dim; i++)
710  {
711  double center = (bb_min[i] + bb_max[i]) * 0.5;
712  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
713  }
714 
715  switch (problem)
716  {
717  case 0:
718  case 1:
719  {
720  switch (dim)
721  {
722  case 1:
723  return exp(-40.*pow(X(0)-0.5,2));
724  case 2:
725  case 3:
726  {
727  double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
728  if (dim == 3)
729  {
730  const double s = (1. + 0.25*cos(2*M_PI*X(2)));
731  rx *= s;
732  ry *= s;
733  }
734  return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
735  erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
736  }
737  }
738  }
739  case 2:
740  {
741  double x_ = X(0), y_ = X(1), rho, phi;
742  rho = hypot(x_, y_);
743  phi = atan2(y_, x_);
744  return pow(sin(M_PI*rho),2)*sin(3*phi);
745  }
746  case 3:
747  {
748  const double f = M_PI;
749  return sin(f*X(0))*sin(f*X(1));
750  }
751  }
752  return 0.0;
753 }
754 
755 // Inflow boundary condition (zero for the problems considered in this example)
756 double inflow_function(const Vector &x)
757 {
758  switch (problem)
759  {
760  case 0:
761  case 1:
762  case 2:
763  case 3: return 0.0;
764  }
765  return 0.0;
766 }
virtual void ImplicitMult(const Vector &x, const Vector &xp, Vector &y) const
Perform the action of the implicit part of the operator, F: y = F(x, k, t) where t is the current tim...
Definition: ex9p.cpp:587
Vector bb_max
Definition: ex9.cpp:66
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
Conjugate gradient method.
Definition: solvers.hpp:316
AssemblyLevel
Enumeration defining the assembly level for bilinear and nonlinear form classes derived from Operator...
bool isImplicit() const
True if type is IMPLICIT or HOMOGENEOUS.
Definition: operator.hpp:331
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication: .
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
HYPRE_Int PetscInt
Definition: petsc.hpp:47
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:616
Wrapper for PETSc&#39;s matrix class.
Definition: petsc.hpp:307
Base abstract class for first order time dependent operators.
Definition: operator.hpp:282
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: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:275
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
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:841
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:326
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:493
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:652
void SetMonitor(PetscSolverMonitor *ctx)
Sets user-defined monitoring routine.
Definition: petsc.cpp:2459
Abstract class for PETSc&#39;s ODE solvers.
Definition: petsc.hpp:920
double PetscReal
Definition: petsc.hpp:49
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
Data collection with Sidre routines following the Conduit mesh blueprint specification.
Class for parallel linear form.
Definition: plinearform.hpp:26
virtual void Init(TimeDependentOperator &f_, enum PetscODESolver::Type type)
Initialize the ODE solver.
Definition: petsc.cpp:4074
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:71
Abstract class for monitoring PETSc&#39;s solvers.
Definition: petsc.hpp:955
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
constexpr char vishost[]
virtual void Save()
Save the collection to disk.
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:130
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:153
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9143
constexpr int visport
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
Definition: solvers.hpp:100
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:4882
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4382
Parallel smoothers in hypre.
Definition: hypre.hpp:840
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:83
bool isExplicit() const
True if type is EXPLICIT.
Definition: operator.hpp:329
int Dimension() const
Definition: mesh.hpp:911
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
void SetTime(double t)
Set physical time (for time-dependent simulations)
Vector bb_min
Definition: ex9.cpp:66
A general vector function coefficient.
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:99
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:252
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:162
void SetAbsTol(double atol)
Definition: solvers.hpp:99
void SetRelTol(double rtol)
Definition: solvers.hpp:98
void velocity_function(const Vector &x, Vector &v)
Definition: ex9.cpp:500
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
virtual Operator & GetExplicitGradient(const Vector &x) const
Return an Operator representing dG/dx at the given point x and the currently set time.
Definition: ex9p.cpp:609
virtual ~FE_Evolution()
Definition: ex9p.cpp:91
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:206
void MFEMFinalizePetsc()
Definition: petsc.cpp:192
virtual Operator & GetImplicitGradient(const Vector &x, const Vector &xp, double shift) const
Return an Operator representing (dF/dk shift + dF/dx) at the given x, k, and the currently set time...
Definition: ex9p.cpp:626
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
Definition: petsc.cpp:168
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.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:327
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:143
Class for parallel bilinear form.
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:557
RefCoord t[3]
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:330
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:92
virtual void Run(Vector &x, double &t, double &dt, double tf)
Perform time integration from time t [in] to time tf [in].
Definition: ode.hpp:89
RefCoord s[3]
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
Base class for solvers.
Definition: operator.hpp:648
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:116
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:277
Class for parallel meshes.
Definition: pmesh.hpp:32
virtual void ExplicitMult(const Vector &x, Vector &y) const
Perform the action of the explicit part of the operator, G: y = G(x, t) where t is the current time...
Definition: ex9p.cpp:569
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()
double inflow_function(const Vector &x)
Definition: ex9.cpp:610
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:285
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)