MFEM  v4.5.1
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 // SUNDIALS Modification
3 //
4 // Compile with: make ex9
5 //
6 // Sample runs:
7 // ex9 -m ../../data/periodic-segment.mesh -p 0 -r 2 -s 7 -dt 0.005
8 // ex9 -m ../../data/periodic-square.mesh -p 1 -r 2 -s 8 -dt 0.005 -tf 9
9 // ex9 -m ../../data/periodic-hexagon.mesh -p 0 -r 2 -s 7 -dt 0.0018 -vs 25
10 // ex9 -m ../../data/periodic-hexagon.mesh -p 0 -r 2 -s 9 -dt 0.01 -vs 15
11 // ex9 -m ../../data/amr-quad.mesh -p 1 -r 2 -s 9 -dt 0.002 -tf 9
12 // ex9 -m ../../data/star-q3.mesh -p 1 -r 2 -s 9 -dt 0.005 -tf 9
13 // ex9 -m ../../data/disc-nurbs.mesh -p 1 -r 3 -s 7 -dt 0.005 -tf 9
14 // ex9 -m ../../data/periodic-cube.mesh -p 0 -r 2 -s 8 -dt 0.02 -tf 8 -o 2
15 //
16 // Device sample runs:
17 // ex9 -pa
18 // ex9 -ea
19 // ex9 -fa
20 // ex9 -pa -m ../../data/periodic-cube.mesh
21 // ex9 -pa -m ../../data/periodic-cube.mesh -d cuda
22 // ex9 -ea -m ../../data/periodic-cube.mesh -d cuda
23 // ex9 -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. The saving of time-dependent data files for external
35 // visualization with VisIt (visit.llnl.gov) and ParaView
36 // (paraview.org) is also illustrated.
37 
38 #include "mfem.hpp"
39 #include <fstream>
40 #include <iostream>
41 #include <algorithm>
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  SparseMatrix &M, &K, A;
70  GMRESSolver linear_solver;
71  BlockILU prec;
72  double dt;
73 public:
74  DG_Solver(SparseMatrix &M_, SparseMatrix &K_, const FiniteElementSpace &fes)
75  : M(M_),
76  K(K_),
77  prec(fes.GetFE(0)->GetDof(),
78  BlockILU::Reordering::MINIMUM_DISCARDED_FILL),
79  dt(-1.0)
80  {
81  linear_solver.iterative_mode = false;
82  linear_solver.SetRelTol(1e-9);
83  linear_solver.SetAbsTol(0.0);
84  linear_solver.SetMaxIter(100);
85  linear_solver.SetPrintLevel(0);
86  linear_solver.SetPreconditioner(prec);
87  }
88 
89  void SetTimeStep(double dt_)
90  {
91  if (dt_ != dt)
92  {
93  dt = dt_;
94  // Form operator A = M - dt*K
95  A = K;
96  A *= -dt;
97  A += M;
98 
99  // this will also call SetOperator on the preconditioner
100  linear_solver.SetOperator(A);
101  }
102  }
103 
104  void SetOperator(const Operator &op)
105  {
106  linear_solver.SetOperator(op);
107  }
108 
109  virtual void Mult(const Vector &x, Vector &y) const
110  {
111  linear_solver.Mult(x, y);
112  }
113 };
114 
115 /** A time-dependent operator for the right-hand side of the ODE. The DG weak
116  form of du/dt = -v.grad(u) is M du/dt = K u + b, where M and K are the mass
117  and advection matrices, and b describes the flow on the boundary. This can
118  be written as a general ODE, du/dt = M^{-1} (K u + b), and this class is
119  used to evaluate the right-hand side. */
120 class FE_Evolution : public TimeDependentOperator
121 {
122 private:
123  BilinearForm &M, &K;
124  const Vector &b;
125  Solver *M_prec;
126  CGSolver M_solver;
127  DG_Solver *dg_solver;
128 
129  mutable Vector z;
130 
131 public:
132  FE_Evolution(BilinearForm &M_, BilinearForm &K_, const Vector &b_);
133 
134  virtual void Mult(const Vector &x, Vector &y) const;
135  virtual void ImplicitSolve(const double dt, const Vector &x, Vector &k);
136 
137  virtual ~FE_Evolution();
138 };
139 
140 
141 int main(int argc, char *argv[])
142 {
143  // 0. Initialize SUNDIALS.
144  Sundials::Init();
145 
146  // 1. Parse command-line options.
147  problem = 0;
148  const char *mesh_file = "../../data/periodic-hexagon.mesh";
149  int ref_levels = 2;
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 = 7;
156  double t_final = 10.0;
157  double dt = 0.01;
158  bool visualization = false;
159  bool visit = false;
160  bool paraview = false;
161  bool binary = false;
162  int vis_steps = 5;
163 
164  // Relative and absolute tolerances for CVODE and ARKODE.
165  const double reltol = 1e-2, abstol = 1e-2;
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(&ref_levels, "-r", "--refine",
176  "Number of times to refine the mesh uniformly.");
177  args.AddOption(&order, "-o", "--order",
178  "Order (degree) of the finite elements.");
179  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
180  "--no-partial-assembly", "Enable Partial Assembly.");
181  args.AddOption(&ea, "-ea", "--element-assembly", "-no-ea",
182  "--no-element-assembly", "Enable Element Assembly.");
183  args.AddOption(&fa, "-fa", "--full-assembly", "-no-fa",
184  "--no-full-assembly", "Enable Full Assembly.");
185  args.AddOption(&device_config, "-d", "--device",
186  "Device configuration string, see Device::Configure().");
187  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
188  "ODE solver:\n\t"
189  "1 - Forward Euler,\n\t"
190  "2 - RK2 SSP,\n\t"
191  "3 - RK3 SSP,\n\t"
192  "4 - RK4,\n\t"
193  "6 - RK6,\n\t"
194  "7 - CVODE (adaptive order implicit Adams),\n\t"
195  "8 - ARKODE default (4th order) explicit,\n\t"
196  "9 - ARKODE RK8.");
197  args.AddOption(&t_final, "-tf", "--t-final",
198  "Final time; start time is 0.");
199  args.AddOption(&dt, "-dt", "--time-step",
200  "Time step.");
201  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
202  "--no-visualization",
203  "Enable or disable GLVis visualization.");
204  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
205  "--no-visit-datafiles",
206  "Save data files for VisIt (visit.llnl.gov) visualization.");
207  args.AddOption(&paraview, "-paraview", "--paraview-datafiles", "-no-paraview",
208  "--no-paraview-datafiles",
209  "Save data files for ParaView (paraview.org) visualization.");
210  args.AddOption(&binary, "-binary", "--binary-datafiles", "-ascii",
211  "--ascii-datafiles",
212  "Use binary (Sidre) or ascii format for VisIt data files.");
213  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
214  "Visualize every n-th timestep.");
215  args.Parse();
216  if (!args.Good())
217  {
218  args.PrintUsage(cout);
219  return 1;
220  }
221  // check for valid ODE solver option
222  if (ode_solver_type < 1 || ode_solver_type > 9)
223  {
224  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
225  return 3;
226  }
227  args.PrintOptions(cout);
228 
229  Device device(device_config);
230  device.Print();
231 
232  // 2. Read the mesh from the given mesh file. We can handle geometrically
233  // periodic meshes in this code.
234  Mesh mesh(mesh_file, 1, 1);
235  int dim = mesh.Dimension();
236 
237  // 3. Refine the mesh to increase the resolution. In this example we do
238  // 'ref_levels' of uniform refinement, where 'ref_levels' is a
239  // command-line parameter. If the mesh is of NURBS type, we convert it to
240  // a (piecewise-polynomial) high-order mesh.
241  for (int lev = 0; lev < ref_levels; lev++)
242  {
243  mesh.UniformRefinement();
244  }
245  if (mesh.NURBSext)
246  {
247  mesh.SetCurvature(max(order, 1));
248  }
249  mesh.GetBoundingBox(bb_min, bb_max, max(order, 1));
250 
251  // 4. Define the discontinuous DG finite element space of the given
252  // polynomial order on the refined mesh.
253  DG_FECollection fec(order, dim, BasisType::GaussLobatto);
254  FiniteElementSpace fes(&mesh, &fec);
255 
256  cout << "Number of unknowns: " << fes.GetVSize() << endl;
257 
258  // 5. Set up and assemble the bilinear and linear forms corresponding to the
259  // DG discretization. The DGTraceIntegrator involves integrals over mesh
260  // interior faces.
264 
265  BilinearForm m(&fes);
266  BilinearForm k(&fes);
267  if (pa)
268  {
269  m.SetAssemblyLevel(AssemblyLevel::PARTIAL);
270  k.SetAssemblyLevel(AssemblyLevel::PARTIAL);
271  }
272  else if (ea)
273  {
274  m.SetAssemblyLevel(AssemblyLevel::ELEMENT);
275  k.SetAssemblyLevel(AssemblyLevel::ELEMENT);
276  }
277  else if (fa)
278  {
279  m.SetAssemblyLevel(AssemblyLevel::FULL);
280  k.SetAssemblyLevel(AssemblyLevel::FULL);
281  }
283  k.AddDomainIntegrator(new ConvectionIntegrator(velocity, -1.0));
285  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
287  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
288 
289  LinearForm b(&fes);
291  new BoundaryFlowIntegrator(inflow, velocity, -1.0, -0.5));
292 
293  m.Assemble();
294  int skip_zeros = 0;
295  k.Assemble(skip_zeros);
296  b.Assemble();
297  m.Finalize();
298  k.Finalize(skip_zeros);
299 
300  // 6. Define the initial conditions, save the corresponding grid function to
301  // a file and (optionally) save data in the VisIt format and initialize
302  // GLVis visualization.
303  GridFunction u(&fes);
304  u.ProjectCoefficient(u0);
305 
306  {
307  ofstream omesh("ex9.mesh");
308  omesh.precision(precision);
309  mesh.Print(omesh);
310  ofstream osol("ex9-init.gf");
311  osol.precision(precision);
312  u.Save(osol);
313  }
314 
315  // Create data collection for solution output: either VisItDataCollection for
316  // ascii data files, or SidreDataCollection for binary data files.
317  DataCollection *dc = NULL;
318  if (visit)
319  {
320  if (binary)
321  {
322 #ifdef MFEM_USE_SIDRE
323  dc = new SidreDataCollection("Example9", &mesh);
324 #else
325  MFEM_ABORT("Must build with MFEM_USE_SIDRE=YES for binary output.");
326 #endif
327  }
328  else
329  {
330  dc = new VisItDataCollection("Example9", &mesh);
331  dc->SetPrecision(precision);
332  }
333  dc->RegisterField("solution", &u);
334  dc->SetCycle(0);
335  dc->SetTime(0.0);
336  dc->Save();
337  }
338 
339  ParaViewDataCollection *pd = NULL;
340  if (paraview)
341  {
342  pd = new ParaViewDataCollection("Example9", &mesh);
343  pd->SetPrefixPath("ParaView");
344  pd->RegisterField("solution", &u);
345  pd->SetLevelsOfDetail(order);
346  pd->SetDataFormat(VTKFormat::BINARY);
347  pd->SetHighOrderOutput(true);
348  pd->SetCycle(0);
349  pd->SetTime(0.0);
350  pd->Save();
351  }
352 
353  socketstream sout;
354  if (visualization)
355  {
356  char vishost[] = "localhost";
357  int visport = 19916;
358  sout.open(vishost, visport);
359  if (!sout)
360  {
361  cout << "Unable to connect to GLVis server at "
362  << vishost << ':' << visport << endl;
363  visualization = false;
364  cout << "GLVis visualization disabled.\n";
365  }
366  else
367  {
368  sout.precision(precision);
369  sout << "solution\n" << mesh << u;
370  sout << "pause\n";
371  sout << flush;
372  cout << "GLVis visualization paused."
373  << " Press space (in the GLVis window) to resume it.\n";
374  }
375  }
376 
377  // 7. Define the time-dependent evolution operator describing the ODE
378  // right-hand side, and define the ODE solver used for time integration.
379  FE_Evolution adv(m, k, b);
380 
381  double t = 0.0;
382  adv.SetTime(t);
383 
384  // Create the time integrator
385  ODESolver *ode_solver = NULL;
386  CVODESolver *cvode = NULL;
387  ARKStepSolver *arkode = NULL;
388  switch (ode_solver_type)
389  {
390  case 1: ode_solver = new ForwardEulerSolver; break;
391  case 2: ode_solver = new RK2Solver(1.0); break;
392  case 3: ode_solver = new RK3SSPSolver; break;
393  case 4: ode_solver = new RK4Solver; break;
394  case 6: ode_solver = new RK6Solver; break;
395  case 7:
396  cvode = new CVODESolver(CV_ADAMS);
397  cvode->Init(adv);
398  cvode->SetSStolerances(reltol, abstol);
399  cvode->SetMaxStep(dt);
400  cvode->UseSundialsLinearSolver();
401  ode_solver = cvode; break;
402  case 8:
403  arkode = new ARKStepSolver(ARKStepSolver::EXPLICIT);
404  arkode->Init(adv);
405  arkode->SetSStolerances(reltol, abstol);
406  arkode->SetMaxStep(dt);
407  arkode->SetOrder(4);
408  ode_solver = arkode; break;
409  case 9:
410  arkode = new ARKStepSolver(ARKStepSolver::EXPLICIT);
411  arkode->Init(adv);
412  arkode->SetSStolerances(reltol, abstol);
413  arkode->SetMaxStep(dt);
415  ode_solver = arkode; break;
416  }
417 
418  // Initialize MFEM integrators, SUNDIALS integrators are initialized above
419  if (ode_solver_type < 7) { ode_solver->Init(adv); }
420 
421  // 8. Perform time-integration (looping over the time iterations, ti,
422  // with a time-step dt).
423  bool done = false;
424  for (int ti = 0; !done; )
425  {
426  double dt_real = min(dt, t_final - t);
427  ode_solver->Step(u, t, dt_real);
428  ti++;
429 
430  done = (t >= t_final - 1e-8*dt);
431 
432  if (done || ti % vis_steps == 0)
433  {
434  cout << "time step: " << ti << ", time: " << t << endl;
435  if (cvode) { cvode->PrintInfo(); }
436  if (arkode) { arkode->PrintInfo(); }
437 
438  if (visualization)
439  {
440  sout << "solution\n" << mesh << u << flush;
441  }
442 
443  if (visit)
444  {
445  dc->SetCycle(ti);
446  dc->SetTime(t);
447  dc->Save();
448  }
449 
450  if (paraview)
451  {
452  pd->SetCycle(ti);
453  pd->SetTime(t);
454  pd->Save();
455  }
456  }
457  }
458 
459  // 9. Save the final solution. This output can be viewed later using GLVis:
460  // "glvis -m ex9.mesh -g ex9-final.gf".
461  {
462  ofstream osol("ex9-final.gf");
463  osol.precision(precision);
464  u.Save(osol);
465  }
466 
467  // 10. Free the used memory.
468  delete ode_solver;
469  delete pd;
470  delete dc;
471 
472  return 0;
473 }
474 
475 
476 // Implementation of class FE_Evolution
478  : TimeDependentOperator(M_.Height()), M(M_), K(K_), b(b_), z(M_.Height())
479 {
480  Array<int> ess_tdof_list;
481  if (M.GetAssemblyLevel() == AssemblyLevel::LEGACY)
482  {
483  M_prec = new DSmoother(M.SpMat());
484  M_solver.SetOperator(M.SpMat());
485  dg_solver = new DG_Solver(M.SpMat(), K.SpMat(), *M.FESpace());
486  }
487  else
488  {
489  M_prec = new OperatorJacobiSmoother(M, ess_tdof_list);
490  M_solver.SetOperator(M);
491  dg_solver = NULL;
492  }
493  M_solver.SetPreconditioner(*M_prec);
494  M_solver.iterative_mode = false;
495  M_solver.SetRelTol(1e-9);
496  M_solver.SetAbsTol(0.0);
497  M_solver.SetMaxIter(100);
498  M_solver.SetPrintLevel(0);
499 }
500 
501 void FE_Evolution::Mult(const Vector &x, Vector &y) const
502 {
503  // y = M^{-1} (K x + b)
504  K.Mult(x, z);
505  z += b;
506  M_solver.Mult(z, y);
507 }
508 
509 void FE_Evolution::ImplicitSolve(const double dt, const Vector &x, Vector &k)
510 {
511  MFEM_VERIFY(dg_solver != NULL,
512  "Implicit time integration is not supported with partial assembly");
513  K.Mult(x, z);
514  z += b;
515  dg_solver->SetTimeStep(dt);
516  dg_solver->Mult(z, k);
517 }
518 
520 {
521  delete M_prec;
522  delete dg_solver;
523 }
524 
525 // Velocity coefficient
526 void velocity_function(const Vector &x, Vector &v)
527 {
528  int dim = x.Size();
529 
530  // map to the reference [-1,1] domain
531  Vector X(dim);
532  for (int i = 0; i < dim; i++)
533  {
534  double center = (bb_min[i] + bb_max[i]) * 0.5;
535  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
536  }
537 
538  switch (problem)
539  {
540  case 0:
541  {
542  // Translations in 1D, 2D, and 3D
543  switch (dim)
544  {
545  case 1: v(0) = 1.0; break;
546  case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.); break;
547  case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
548  break;
549  }
550  break;
551  }
552  case 1:
553  case 2:
554  {
555  // Clockwise rotation in 2D around the origin
556  const double w = M_PI/2;
557  switch (dim)
558  {
559  case 1: v(0) = 1.0; break;
560  case 2: v(0) = w*X(1); v(1) = -w*X(0); break;
561  case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0; break;
562  }
563  break;
564  }
565  case 3:
566  {
567  // Clockwise twisting rotation in 2D around the origin
568  const double w = M_PI/2;
569  double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
570  d = d*d;
571  switch (dim)
572  {
573  case 1: v(0) = 1.0; break;
574  case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0); break;
575  case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0; break;
576  }
577  break;
578  }
579  }
580 }
581 
582 // Initial condition
583 double u0_function(const Vector &x)
584 {
585  int dim = x.Size();
586 
587  // map to the reference [-1,1] domain
588  Vector X(dim);
589  for (int i = 0; i < dim; i++)
590  {
591  double center = (bb_min[i] + bb_max[i]) * 0.5;
592  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
593  }
594 
595  switch (problem)
596  {
597  case 0:
598  case 1:
599  {
600  switch (dim)
601  {
602  case 1:
603  return exp(-40.*pow(X(0)-0.5,2));
604  case 2:
605  case 3:
606  {
607  double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
608  if (dim == 3)
609  {
610  const double s = (1. + 0.25*cos(2*M_PI*X(2)));
611  rx *= s;
612  ry *= s;
613  }
614  return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
615  erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
616  }
617  }
618  }
619  case 2:
620  {
621  double x_ = X(0), y_ = X(1), rho, phi;
622  rho = hypot(x_, y_);
623  phi = atan2(y_, x_);
624  return pow(sin(M_PI*rho),2)*sin(3*phi);
625  }
626  case 3:
627  {
628  const double f = M_PI;
629  return sin(f*X(0))*sin(f*X(1));
630  }
631  }
632  return 0.0;
633 }
634 
635 // Inflow boundary condition (zero for the problems considered in this example)
636 double inflow_function(const Vector &x)
637 {
638  switch (problem)
639  {
640  case 0:
641  case 1:
642  case 2:
643  case 3: return 0.0;
644  }
645  return 0.0;
646 }
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:696
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:587
Conjugate gradient method.
Definition: solvers.hpp:465
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 SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
Definition: sundials.cpp:862
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:711
Helper class for ParaView visualization data.
void SetERKTableNum(ARKODE_ERKTableID table_id)
Choose a specific Butcher table for an explicit RK method.
Definition: sundials.cpp:1718
Base abstract class for first order time dependent operators.
Definition: operator.hpp:289
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:475
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
Definition: mesh.cpp:129
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:168
virtual void SetOperator(const Operator &a)
Set/update the solver for the given operator.
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:200
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition: sundials.cpp:1706
void SetOrder(int order)
Chooses integration order for all explicit / implicit / IMEX methods.
Definition: sundials.cpp:1712
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:333
Interface to ARKode&#39;s ARKStep module – additive Runge-Kutta methods.
Definition: sundials.hpp:661
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:3676
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
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:50
constexpr char vishost[]
virtual void Save()
Save the collection to disk.
Interface to the CVODE library – linear multi-step methods.
Definition: sundials.hpp:374
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:9816
constexpr int visport
Data collection with VisIt I/O routines.
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:5343
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition: sundials.cpp:880
void SetHighOrderOutput(bool high_order_output_)
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:85
int Dimension() const
Definition: mesh.hpp:1006
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)
Vector bb_min
Definition: ex9.cpp:67
A general vector function coefficient.
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:162
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 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:96
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
GMRES method.
Definition: solvers.hpp:497
constexpr ARKODE_ERKTableID ARKODE_FEHLBERG_13_7_8
Definition: sundials.hpp:54
virtual void Mult(const Vector &x, Vector &y) const
Perform the action of the operator: y = k = f(x, t), where k solves the algebraic equation F(x...
Definition: ex18.hpp:107
FE_Evolution(FiniteElementSpace &vfes_, Operator &A_, SparseMatrix &Aflux_)
Definition: ex18.hpp:81
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:272
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.
virtual void Print(std::ostream &os=mfem::out) const
Definition: mesh.hpp:1671
void PrintInfo() const
Print various ARKStep statistics.
Definition: sundials.cpp:1743
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:2781
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:324
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: gridfunc.cpp:2399
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_)
A general function coefficient.
Vector data type.
Definition: vector.hpp:60
Vector with associated FE space and LinearFormIntegrators.
Definition: linearform.hpp:24
virtual void Save() override
void PrintInfo() const
Print various CVODE statistics.
Definition: sundials.cpp:906
RefCoord s[3]
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
Base class for solvers.
Definition: operator.hpp:655
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
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults...
Definition: sundials.cpp:1461
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
Definition: sundials.cpp:1700
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:842
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:288
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)