MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex9.cpp
Go to the documentation of this file.
1 // MFEM Example 9
2 //
3 // Compile with: make ex9
4 //
5 // Sample runs:
6 // ex9 -m ../data/periodic-segment.mesh -p 0 -r 2 -dt 0.005
7 // ex9 -m ../data/periodic-square.mesh -p 0 -r 2 -dt 0.01 -tf 10
8 // ex9 -m ../data/periodic-hexagon.mesh -p 0 -r 2 -dt 0.01 -tf 10
9 // ex9 -m ../data/periodic-square.mesh -p 1 -r 2 -dt 0.005 -tf 9
10 // ex9 -m ../data/periodic-hexagon.mesh -p 1 -r 2 -dt 0.005 -tf 9
11 // ex9 -m ../data/amr-quad.mesh -p 1 -r 2 -dt 0.002 -tf 9
12 // ex9 -m ../data/amr-quad.mesh -p 1 -r 2 -dt 0.02 -s 13 -tf 9
13 // ex9 -m ../data/star-q3.mesh -p 1 -r 2 -dt 0.005 -tf 9
14 // ex9 -m ../data/star-mixed.mesh -p 1 -r 2 -dt 0.005 -tf 9
15 // ex9 -m ../data/disc-nurbs.mesh -p 1 -r 3 -dt 0.005 -tf 9
16 // ex9 -m ../data/disc-nurbs.mesh -p 2 -r 3 -dt 0.005 -tf 9
17 // ex9 -m ../data/periodic-square.mesh -p 3 -r 4 -dt 0.0025 -tf 9 -vs 20
18 // ex9 -m ../data/periodic-cube.mesh -p 0 -r 2 -o 2 -dt 0.02 -tf 8
19 // ex9 -m ../data/periodic-square.msh -p 0 -r 2 -dt 0.005 -tf 2
20 // ex9 -m ../data/periodic-cube.msh -p 0 -r 1 -o 2 -tf 2
21 //
22 // Device sample runs:
23 // ex9 -pa
24 // ex9 -ea
25 // ex9 -fa
26 // ex9 -pa -m ../data/periodic-cube.mesh
27 // ex9 -pa -m ../data/periodic-cube.mesh -d cuda
28 // ex9 -ea -m ../data/periodic-cube.mesh -d cuda
29 // ex9 -fa -m ../data/periodic-cube.mesh -d cuda
30 //
31 // Description: This example code solves the time-dependent advection equation
32 // du/dt + v.grad(u) = 0, where v is a given fluid velocity, and
33 // u0(x)=u(0,x) is a given initial condition.
34 //
35 // The example demonstrates the use of Discontinuous Galerkin (DG)
36 // bilinear forms in MFEM (face integrators), the use of implicit
37 // and explicit ODE time integrators, the definition of periodic
38 // boundary conditions through periodic meshes, as well as the use
39 // of GLVis for persistent visualization of a time-evolving
40 // solution. The saving of time-dependent data files for external
41 // visualization with VisIt (visit.llnl.gov) and ParaView
42 // (paraview.org) is also illustrated.
43 
44 #include "mfem.hpp"
45 #include <fstream>
46 #include <iostream>
47 #include <algorithm>
48 
49 using namespace std;
50 using namespace mfem;
51 
52 // Choice for the problem setup. The fluid velocity, initial condition and
53 // inflow boundary condition are chosen based on this parameter.
54 int problem;
55 
56 // Velocity coefficient
57 void velocity_function(const Vector &x, Vector &v);
58 
59 // Initial condition
60 double u0_function(const Vector &x);
61 
62 // Inflow boundary condition
63 double inflow_function(const Vector &x);
64 
65 // Mesh bounding box
67 
68 class DG_Solver : public Solver
69 {
70 private:
71  SparseMatrix &M, &K, A;
72  GMRESSolver linear_solver;
73  BlockILU prec;
74  double dt;
75 public:
76  DG_Solver(SparseMatrix &M_, SparseMatrix &K_, const FiniteElementSpace &fes)
77  : M(M_),
78  K(K_),
79  prec(fes.GetFE(0)->GetDof(),
80  BlockILU::Reordering::MINIMUM_DISCARDED_FILL),
81  dt(-1.0)
82  {
83  linear_solver.iterative_mode = false;
84  linear_solver.SetRelTol(1e-9);
85  linear_solver.SetAbsTol(0.0);
86  linear_solver.SetMaxIter(100);
87  linear_solver.SetPrintLevel(0);
88  linear_solver.SetPreconditioner(prec);
89  }
90 
91  void SetTimeStep(double dt_)
92  {
93  if (dt_ != dt)
94  {
95  dt = dt_;
96  // Form operator A = M - dt*K
97  A = K;
98  A *= -dt;
99  A += M;
100 
101  // this will also call SetOperator on the preconditioner
102  linear_solver.SetOperator(A);
103  }
104  }
105 
106  void SetOperator(const Operator &op)
107  {
108  linear_solver.SetOperator(op);
109  }
110 
111  virtual void Mult(const Vector &x, Vector &y) const
112  {
113  linear_solver.Mult(x, y);
114  }
115 };
116 
117 /** A time-dependent operator for the right-hand side of the ODE. The DG weak
118  form of du/dt = -v.grad(u) is M du/dt = K u + b, where M and K are the mass
119  and advection matrices, and b describes the flow on the boundary. This can
120  be written as a general ODE, du/dt = M^{-1} (K u + b), and this class is
121  used to evaluate the right-hand side. */
122 class FE_Evolution : public TimeDependentOperator
123 {
124 private:
125  BilinearForm &M, &K;
126  const Vector &b;
127  Solver *M_prec;
128  CGSolver M_solver;
129  DG_Solver *dg_solver;
130 
131  mutable Vector z;
132 
133 public:
134  FE_Evolution(BilinearForm &M_, BilinearForm &K_, const Vector &b_);
135 
136  virtual void Mult(const Vector &x, Vector &y) const;
137  virtual void ImplicitSolve(const double dt, const Vector &x, Vector &k);
138 
139  virtual ~FE_Evolution();
140 };
141 
142 
143 int main(int argc, char *argv[])
144 {
145  // 1. Parse command-line options.
146  problem = 0;
147  const char *mesh_file = "../data/periodic-hexagon.mesh";
148  int ref_levels = 2;
149  int order = 3;
150  bool pa = false;
151  bool ea = false;
152  bool fa = false;
153  const char *device_config = "cpu";
154  int ode_solver_type = 4;
155  double t_final = 10.0;
156  double dt = 0.01;
157  bool visualization = true;
158  bool visit = false;
159  bool paraview = false;
160  bool binary = false;
161  int vis_steps = 5;
162 
163  int precision = 8;
164  cout.precision(precision);
165 
166  OptionsParser args(argc, argv);
167  args.AddOption(&mesh_file, "-m", "--mesh",
168  "Mesh file to use.");
169  args.AddOption(&problem, "-p", "--problem",
170  "Problem setup to use. See options in velocity_function().");
171  args.AddOption(&ref_levels, "-r", "--refine",
172  "Number of times to refine the mesh uniformly.");
173  args.AddOption(&order, "-o", "--order",
174  "Order (degree) of the finite elements.");
175  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
176  "--no-partial-assembly", "Enable Partial Assembly.");
177  args.AddOption(&ea, "-ea", "--element-assembly", "-no-ea",
178  "--no-element-assembly", "Enable Element Assembly.");
179  args.AddOption(&fa, "-fa", "--full-assembly", "-no-fa",
180  "--no-full-assembly", "Enable Full Assembly.");
181  args.AddOption(&device_config, "-d", "--device",
182  "Device configuration string, see Device::Configure().");
183  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
184  "ODE solver: 1 - Forward Euler,\n\t"
185  " 2 - RK2 SSP, 3 - RK3 SSP, 4 - RK4, 6 - RK6,\n\t"
186  " 11 - Backward Euler,\n\t"
187  " 12 - SDIRK23 (L-stable), 13 - SDIRK33,\n\t"
188  " 22 - Implicit Midpoint Method,\n\t"
189  " 23 - SDIRK23 (A-stable), 24 - SDIRK34");
190  args.AddOption(&t_final, "-tf", "--t-final",
191  "Final time; start time is 0.");
192  args.AddOption(&dt, "-dt", "--time-step",
193  "Time step.");
194  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
195  "--no-visualization",
196  "Enable or disable GLVis visualization.");
197  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
198  "--no-visit-datafiles",
199  "Save data files for VisIt (visit.llnl.gov) visualization.");
200  args.AddOption(&paraview, "-paraview", "--paraview-datafiles", "-no-paraview",
201  "--no-paraview-datafiles",
202  "Save data files for ParaView (paraview.org) visualization.");
203  args.AddOption(&binary, "-binary", "--binary-datafiles", "-ascii",
204  "--ascii-datafiles",
205  "Use binary (Sidre) or ascii format for VisIt data files.");
206  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
207  "Visualize every n-th timestep.");
208  args.Parse();
209  if (!args.Good())
210  {
211  args.PrintUsage(cout);
212  return 1;
213  }
214  args.PrintOptions(cout);
215 
216  Device device(device_config);
217  device.Print();
218 
219  // 2. Read the mesh from the given mesh file. We can handle geometrically
220  // periodic meshes in this code.
221  Mesh mesh(mesh_file, 1, 1);
222  int dim = mesh.Dimension();
223 
224  // 3. Define the ODE solver used for time integration. Several explicit
225  // Runge-Kutta methods are available.
226  ODESolver *ode_solver = NULL;
227  switch (ode_solver_type)
228  {
229  // Explicit methods
230  case 1: ode_solver = new ForwardEulerSolver; break;
231  case 2: ode_solver = new RK2Solver(1.0); break;
232  case 3: ode_solver = new RK3SSPSolver; break;
233  case 4: ode_solver = new RK4Solver; break;
234  case 6: ode_solver = new RK6Solver; break;
235  // Implicit (L-stable) methods
236  case 11: ode_solver = new BackwardEulerSolver; break;
237  case 12: ode_solver = new SDIRK23Solver(2); break;
238  case 13: ode_solver = new SDIRK33Solver; break;
239  // Implicit A-stable methods (not L-stable)
240  case 22: ode_solver = new ImplicitMidpointSolver; break;
241  case 23: ode_solver = new SDIRK23Solver; break;
242  case 24: ode_solver = new SDIRK34Solver; break;
243 
244  default:
245  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
246  return 3;
247  }
248 
249  // 4. Refine the mesh to increase the resolution. In this example we do
250  // 'ref_levels' of uniform refinement, where 'ref_levels' is a
251  // command-line parameter. If the mesh is of NURBS type, we convert it to
252  // a (piecewise-polynomial) high-order mesh.
253  for (int lev = 0; lev < ref_levels; lev++)
254  {
255  mesh.UniformRefinement();
256  }
257  if (mesh.NURBSext)
258  {
259  mesh.SetCurvature(max(order, 1));
260  }
261  mesh.GetBoundingBox(bb_min, bb_max, max(order, 1));
262 
263  // 5. Define the discontinuous DG finite element space of the given
264  // polynomial order on the refined mesh.
265  DG_FECollection fec(order, dim, BasisType::GaussLobatto);
266  FiniteElementSpace fes(&mesh, &fec);
267 
268  cout << "Number of unknowns: " << fes.GetVSize() << endl;
269 
270  // 6. Set up and assemble the bilinear and linear forms corresponding to the
271  // DG discretization. The DGTraceIntegrator involves integrals over mesh
272  // interior faces.
276 
277  BilinearForm m(&fes);
278  BilinearForm k(&fes);
279  if (pa)
280  {
281  m.SetAssemblyLevel(AssemblyLevel::PARTIAL);
282  k.SetAssemblyLevel(AssemblyLevel::PARTIAL);
283  }
284  else if (ea)
285  {
286  m.SetAssemblyLevel(AssemblyLevel::ELEMENT);
287  k.SetAssemblyLevel(AssemblyLevel::ELEMENT);
288  }
289  else if (fa)
290  {
291  m.SetAssemblyLevel(AssemblyLevel::FULL);
292  k.SetAssemblyLevel(AssemblyLevel::FULL);
293  }
295  constexpr double alpha = -1.0;
296  k.AddDomainIntegrator(new ConvectionIntegrator(velocity, alpha));
298  new NonconservativeDGTraceIntegrator(velocity, alpha));
300  new NonconservativeDGTraceIntegrator(velocity, alpha));
301 
302  LinearForm b(&fes);
304  new BoundaryFlowIntegrator(inflow, velocity, alpha));
305 
306  m.Assemble();
307  int skip_zeros = 0;
308  k.Assemble(skip_zeros);
309  b.Assemble();
310  m.Finalize();
311  k.Finalize(skip_zeros);
312 
313  // 7. Define the initial conditions, save the corresponding grid function to
314  // a file and (optionally) save data in the VisIt format and initialize
315  // GLVis visualization.
316  GridFunction u(&fes);
317  u.ProjectCoefficient(u0);
318 
319  {
320  ofstream omesh("ex9.mesh");
321  omesh.precision(precision);
322  mesh.Print(omesh);
323  ofstream osol("ex9-init.gf");
324  osol.precision(precision);
325  u.Save(osol);
326  }
327 
328  // Create data collection for solution output: either VisItDataCollection for
329  // ascii data files, or SidreDataCollection for binary data files.
330  DataCollection *dc = NULL;
331  if (visit)
332  {
333  if (binary)
334  {
335 #ifdef MFEM_USE_SIDRE
336  dc = new SidreDataCollection("Example9", &mesh);
337 #else
338  MFEM_ABORT("Must build with MFEM_USE_SIDRE=YES for binary output.");
339 #endif
340  }
341  else
342  {
343  dc = new VisItDataCollection("Example9", &mesh);
344  dc->SetPrecision(precision);
345  }
346  dc->RegisterField("solution", &u);
347  dc->SetCycle(0);
348  dc->SetTime(0.0);
349  dc->Save();
350  }
351 
352  ParaViewDataCollection *pd = NULL;
353  if (paraview)
354  {
355  pd = new ParaViewDataCollection("Example9", &mesh);
356  pd->SetPrefixPath("ParaView");
357  pd->RegisterField("solution", &u);
358  pd->SetLevelsOfDetail(order);
359  pd->SetDataFormat(VTKFormat::BINARY);
360  pd->SetHighOrderOutput(true);
361  pd->SetCycle(0);
362  pd->SetTime(0.0);
363  pd->Save();
364  }
365 
366  socketstream sout;
367  if (visualization)
368  {
369  char vishost[] = "localhost";
370  int visport = 19916;
371  sout.open(vishost, visport);
372  if (!sout)
373  {
374  cout << "Unable to connect to GLVis server at "
375  << vishost << ':' << visport << endl;
376  visualization = false;
377  cout << "GLVis visualization disabled.\n";
378  }
379  else
380  {
381  sout.precision(precision);
382  sout << "solution\n" << mesh << u;
383  sout << "pause\n";
384  sout << flush;
385  cout << "GLVis visualization paused."
386  << " Press space (in the GLVis window) to resume it.\n";
387  }
388  }
389 
390  // 8. Define the time-dependent evolution operator describing the ODE
391  // right-hand side, and perform time-integration (looping over the time
392  // iterations, ti, with a time-step dt).
393  FE_Evolution adv(m, k, b);
394 
395  double t = 0.0;
396  adv.SetTime(t);
397  ode_solver->Init(adv);
398 
399  bool done = false;
400  for (int ti = 0; !done; )
401  {
402  double dt_real = min(dt, t_final - t);
403  ode_solver->Step(u, t, dt_real);
404  ti++;
405 
406  done = (t >= t_final - 1e-8*dt);
407 
408  if (done || ti % vis_steps == 0)
409  {
410  cout << "time step: " << ti << ", time: " << t << endl;
411 
412  if (visualization)
413  {
414  sout << "solution\n" << mesh << u << flush;
415  }
416 
417  if (visit)
418  {
419  dc->SetCycle(ti);
420  dc->SetTime(t);
421  dc->Save();
422  }
423 
424  if (paraview)
425  {
426  pd->SetCycle(ti);
427  pd->SetTime(t);
428  pd->Save();
429  }
430  }
431  }
432 
433  // 9. Save the final solution. This output can be viewed later using GLVis:
434  // "glvis -m ex9.mesh -g ex9-final.gf".
435  {
436  ofstream osol("ex9-final.gf");
437  osol.precision(precision);
438  u.Save(osol);
439  }
440 
441  // 10. Free the used memory.
442  delete ode_solver;
443  delete pd;
444  delete dc;
445 
446  return 0;
447 }
448 
449 
450 // Implementation of class FE_Evolution
452  : TimeDependentOperator(M_.Height()), M(M_), K(K_), b(b_), z(M_.Height())
453 {
454  Array<int> ess_tdof_list;
455  if (M.GetAssemblyLevel() == AssemblyLevel::LEGACY)
456  {
457  M_prec = new DSmoother(M.SpMat());
458  M_solver.SetOperator(M.SpMat());
459  dg_solver = new DG_Solver(M.SpMat(), K.SpMat(), *M.FESpace());
460  }
461  else
462  {
463  M_prec = new OperatorJacobiSmoother(M, ess_tdof_list);
464  M_solver.SetOperator(M);
465  dg_solver = NULL;
466  }
467  M_solver.SetPreconditioner(*M_prec);
468  M_solver.iterative_mode = false;
469  M_solver.SetRelTol(1e-9);
470  M_solver.SetAbsTol(0.0);
471  M_solver.SetMaxIter(100);
472  M_solver.SetPrintLevel(0);
473 }
474 
475 void FE_Evolution::Mult(const Vector &x, Vector &y) const
476 {
477  // y = M^{-1} (K x + b)
478  K.Mult(x, z);
479  z += b;
480  M_solver.Mult(z, y);
481 }
482 
483 void FE_Evolution::ImplicitSolve(const double dt, const Vector &x, Vector &k)
484 {
485  MFEM_VERIFY(dg_solver != NULL,
486  "Implicit time integration is not supported with partial assembly");
487  K.Mult(x, z);
488  z += b;
489  dg_solver->SetTimeStep(dt);
490  dg_solver->Mult(z, k);
491 }
492 
494 {
495  delete M_prec;
496  delete dg_solver;
497 }
498 
499 // Velocity coefficient
500 void velocity_function(const Vector &x, Vector &v)
501 {
502  int dim = x.Size();
503 
504  // map to the reference [-1,1] domain
505  Vector X(dim);
506  for (int i = 0; i < dim; i++)
507  {
508  double center = (bb_min[i] + bb_max[i]) * 0.5;
509  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
510  }
511 
512  switch (problem)
513  {
514  case 0:
515  {
516  // Translations in 1D, 2D, and 3D
517  switch (dim)
518  {
519  case 1: v(0) = 1.0; break;
520  case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.); break;
521  case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
522  break;
523  }
524  break;
525  }
526  case 1:
527  case 2:
528  {
529  // Clockwise rotation in 2D around the origin
530  const double w = M_PI/2;
531  switch (dim)
532  {
533  case 1: v(0) = 1.0; break;
534  case 2: v(0) = w*X(1); v(1) = -w*X(0); break;
535  case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0; break;
536  }
537  break;
538  }
539  case 3:
540  {
541  // Clockwise twisting rotation in 2D around the origin
542  const double w = M_PI/2;
543  double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
544  d = d*d;
545  switch (dim)
546  {
547  case 1: v(0) = 1.0; break;
548  case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0); break;
549  case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0; break;
550  }
551  break;
552  }
553  }
554 }
555 
556 // Initial condition
557 double u0_function(const Vector &x)
558 {
559  int dim = x.Size();
560 
561  // map to the reference [-1,1] domain
562  Vector X(dim);
563  for (int i = 0; i < dim; i++)
564  {
565  double center = (bb_min[i] + bb_max[i]) * 0.5;
566  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
567  }
568 
569  switch (problem)
570  {
571  case 0:
572  case 1:
573  {
574  switch (dim)
575  {
576  case 1:
577  return exp(-40.*pow(X(0)-0.5,2));
578  case 2:
579  case 3:
580  {
581  double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
582  if (dim == 3)
583  {
584  const double s = (1. + 0.25*cos(2*M_PI*X(2)));
585  rx *= s;
586  ry *= s;
587  }
588  return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
589  erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
590  }
591  }
592  }
593  case 2:
594  {
595  double x_ = X(0), y_ = X(1), rho, phi;
596  rho = hypot(x_, y_);
597  phi = atan2(y_, x_);
598  return pow(sin(M_PI*rho),2)*sin(3*phi);
599  }
600  case 3:
601  {
602  const double f = M_PI;
603  return sin(f*X(0))*sin(f*X(1));
604  }
605  }
606  return 0.0;
607 }
608 
609 // Inflow boundary condition (zero for the problems considered in this example)
610 double inflow_function(const Vector &x)
611 {
612  switch (problem)
613  {
614  case 0:
615  case 1:
616  case 2:
617  case 3: return 0.0;
618  }
619  return 0.0;
620 }
Vector bb_max
Definition: ex9.cpp:66
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:540
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1387
Conjugate gradient method.
Definition: solvers.hpp:316
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication: .
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
Data type for scaled Jacobi-type smoother of sparse matrix.
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
void SetDataFormat(VTKFormat fmt)
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:616
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: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
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:102
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]...
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
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
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:652
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:389
Data collection with Sidre routines following the Conduit mesh blueprint specification.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3484
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
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:41
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
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 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: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.
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:162
void SetAbsTol(double atol)
Definition: solvers.hpp:99
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:483
void SetRelTol(double rtol)
Definition: solvers.hpp:98
void velocity_function(const Vector &x, Vector &v)
Definition: ex9.cpp:500
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: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
GMRES method.
Definition: solvers.hpp:348
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
FE_Evolution(FiniteElementSpace &vfes_, Operator &A_, SparseMatrix &Aflux_)
Definition: ex18.hpp:102
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:206
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
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.
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:2388
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:327
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: gridfunc.cpp:2278
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:557
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:330
const double alpha
Definition: ex15.cpp:369
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
Vector with associated FE space and LinearFormIntegrators.
Definition: linearform.hpp:23
virtual void Save() override
RefCoord s[3]
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix.
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
Base class for solvers.
Definition: operator.hpp:648
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
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: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)