MFEM  v4.2.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 //
20 // Device sample runs:
21 // ex9 -pa
22 // ex9 -ea
23 // ex9 -fa
24 // ex9 -pa -m ../data/periodic-cube.mesh
25 // ex9 -pa -m ../data/periodic-cube.mesh -d cuda
26 // ex9 -ea -m ../data/periodic-cube.mesh -d cuda
27 // ex9 -fa -m ../data/periodic-cube.mesh -d cuda
28 //
29 // Description: This example code solves the time-dependent advection equation
30 // du/dt + v.grad(u) = 0, where v is a given fluid velocity, and
31 // u0(x)=u(0,x) is a given initial condition.
32 //
33 // The example demonstrates the use of Discontinuous Galerkin (DG)
34 // bilinear forms in MFEM (face integrators), the use of implicit
35 // and explicit ODE time integrators, the definition of periodic
36 // boundary conditions through periodic meshes, as well as the use
37 // of GLVis for persistent visualization of a time-evolving
38 // solution. The saving of time-dependent data files for external
39 // visualization with VisIt (visit.llnl.gov) and ParaView
40 // (paraview.org) is also illustrated.
41 
42 #include "mfem.hpp"
43 #include <fstream>
44 #include <iostream>
45 #include <algorithm>
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  // 1. Parse command-line options.
144  problem = 0;
145  const char *mesh_file = "../data/periodic-hexagon.mesh";
146  int ref_levels = 2;
147  int order = 3;
148  bool pa = false;
149  bool ea = false;
150  bool fa = false;
151  const char *device_config = "cpu";
152  int ode_solver_type = 4;
153  double t_final = 10.0;
154  double dt = 0.01;
155  bool visualization = true;
156  bool visit = false;
157  bool paraview = false;
158  bool binary = false;
159  int vis_steps = 5;
160 
161  int precision = 8;
162  cout.precision(precision);
163 
164  OptionsParser args(argc, argv);
165  args.AddOption(&mesh_file, "-m", "--mesh",
166  "Mesh file to use.");
167  args.AddOption(&problem, "-p", "--problem",
168  "Problem setup to use. See options in velocity_function().");
169  args.AddOption(&ref_levels, "-r", "--refine",
170  "Number of times to refine the mesh uniformly.");
171  args.AddOption(&order, "-o", "--order",
172  "Order (degree) of the finite elements.");
173  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
174  "--no-partial-assembly", "Enable Partial Assembly.");
175  args.AddOption(&ea, "-ea", "--element-assembly", "-no-ea",
176  "--no-element-assembly", "Enable Element Assembly.");
177  args.AddOption(&fa, "-fa", "--full-assembly", "-no-fa",
178  "--no-full-assembly", "Enable Full Assembly.");
179  args.AddOption(&device_config, "-d", "--device",
180  "Device configuration string, see Device::Configure().");
181  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
182  "ODE solver: 1 - Forward Euler,\n\t"
183  " 2 - RK2 SSP, 3 - RK3 SSP, 4 - RK4, 6 - RK6,\n\t"
184  " 11 - Backward Euler,\n\t"
185  " 12 - SDIRK23 (L-stable), 13 - SDIRK33,\n\t"
186  " 22 - Implicit Midpoint Method,\n\t"
187  " 23 - SDIRK23 (A-stable), 24 - SDIRK34");
188  args.AddOption(&t_final, "-tf", "--t-final",
189  "Final time; start time is 0.");
190  args.AddOption(&dt, "-dt", "--time-step",
191  "Time step.");
192  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
193  "--no-visualization",
194  "Enable or disable GLVis visualization.");
195  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
196  "--no-visit-datafiles",
197  "Save data files for VisIt (visit.llnl.gov) visualization.");
198  args.AddOption(&paraview, "-paraview", "--paraview-datafiles", "-no-paraview",
199  "--no-paraview-datafiles",
200  "Save data files for ParaView (paraview.org) visualization.");
201  args.AddOption(&binary, "-binary", "--binary-datafiles", "-ascii",
202  "--ascii-datafiles",
203  "Use binary (Sidre) or ascii format for VisIt data files.");
204  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
205  "Visualize every n-th timestep.");
206  args.Parse();
207  if (!args.Good())
208  {
209  args.PrintUsage(cout);
210  return 1;
211  }
212  args.PrintOptions(cout);
213 
214  Device device(device_config);
215  device.Print();
216 
217  // 2. Read the mesh from the given mesh file. We can handle geometrically
218  // periodic meshes in this code.
219  Mesh mesh(mesh_file, 1, 1);
220  int dim = mesh.Dimension();
221 
222  // 3. Define the ODE solver used for time integration. Several explicit
223  // Runge-Kutta methods are available.
224  ODESolver *ode_solver = NULL;
225  switch (ode_solver_type)
226  {
227  // Explicit methods
228  case 1: ode_solver = new ForwardEulerSolver; break;
229  case 2: ode_solver = new RK2Solver(1.0); break;
230  case 3: ode_solver = new RK3SSPSolver; break;
231  case 4: ode_solver = new RK4Solver; break;
232  case 6: ode_solver = new RK6Solver; break;
233  // Implicit (L-stable) methods
234  case 11: ode_solver = new BackwardEulerSolver; break;
235  case 12: ode_solver = new SDIRK23Solver(2); break;
236  case 13: ode_solver = new SDIRK33Solver; break;
237  // Implicit A-stable methods (not L-stable)
238  case 22: ode_solver = new ImplicitMidpointSolver; break;
239  case 23: ode_solver = new SDIRK23Solver; break;
240  case 24: ode_solver = new SDIRK34Solver; break;
241 
242  default:
243  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
244  return 3;
245  }
246 
247  // 4. Refine the mesh to increase the resolution. In this example we do
248  // 'ref_levels' of uniform refinement, where 'ref_levels' is a
249  // command-line parameter. If the mesh is of NURBS type, we convert it to
250  // a (piecewise-polynomial) high-order mesh.
251  for (int lev = 0; lev < ref_levels; lev++)
252  {
253  mesh.UniformRefinement();
254  }
255  if (mesh.NURBSext)
256  {
257  mesh.SetCurvature(max(order, 1));
258  }
259  mesh.GetBoundingBox(bb_min, bb_max, max(order, 1));
260 
261  // 5. Define the discontinuous DG finite element space of the given
262  // polynomial order on the refined mesh.
263  DG_FECollection fec(order, dim, BasisType::GaussLobatto);
264  FiniteElementSpace fes(&mesh, &fec);
265 
266  cout << "Number of unknowns: " << fes.GetVSize() << endl;
267 
268  // 6. Set up and assemble the bilinear and linear forms corresponding to the
269  // DG discretization. The DGTraceIntegrator involves integrals over mesh
270  // interior faces.
274 
275  BilinearForm m(&fes);
276  BilinearForm k(&fes);
277  if (pa)
278  {
279  m.SetAssemblyLevel(AssemblyLevel::PARTIAL);
280  k.SetAssemblyLevel(AssemblyLevel::PARTIAL);
281  }
282  else if (ea)
283  {
284  m.SetAssemblyLevel(AssemblyLevel::ELEMENT);
285  k.SetAssemblyLevel(AssemblyLevel::ELEMENT);
286  }
287  else if (fa)
288  {
289  m.SetAssemblyLevel(AssemblyLevel::FULL);
290  k.SetAssemblyLevel(AssemblyLevel::FULL);
291  }
293  k.AddDomainIntegrator(new ConvectionIntegrator(velocity, -1.0));
295  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
297  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
298 
299  LinearForm b(&fes);
301  new BoundaryFlowIntegrator(inflow, velocity, -1.0, -0.5));
302 
303  m.Assemble();
304  int skip_zeros = 0;
305  k.Assemble(skip_zeros);
306  b.Assemble();
307  m.Finalize();
308  k.Finalize(skip_zeros);
309 
310  // 7. Define the initial conditions, save the corresponding grid function to
311  // a file and (optionally) save data in the VisIt format and initialize
312  // GLVis visualization.
313  GridFunction u(&fes);
314  u.ProjectCoefficient(u0);
315 
316  {
317  ofstream omesh("ex9.mesh");
318  omesh.precision(precision);
319  mesh.Print(omesh);
320  ofstream osol("ex9-init.gf");
321  osol.precision(precision);
322  u.Save(osol);
323  }
324 
325  // Create data collection for solution output: either VisItDataCollection for
326  // ascii data files, or SidreDataCollection for binary data files.
327  DataCollection *dc = NULL;
328  if (visit)
329  {
330  if (binary)
331  {
332 #ifdef MFEM_USE_SIDRE
333  dc = new SidreDataCollection("Example9", &mesh);
334 #else
335  MFEM_ABORT("Must build with MFEM_USE_SIDRE=YES for binary output.");
336 #endif
337  }
338  else
339  {
340  dc = new VisItDataCollection("Example9", &mesh);
341  dc->SetPrecision(precision);
342  }
343  dc->RegisterField("solution", &u);
344  dc->SetCycle(0);
345  dc->SetTime(0.0);
346  dc->Save();
347  }
348 
349  ParaViewDataCollection *pd = NULL;
350  if (paraview)
351  {
352  pd = new ParaViewDataCollection("Example9", &mesh);
353  pd->SetPrefixPath("ParaView");
354  pd->RegisterField("solution", &u);
355  pd->SetLevelsOfDetail(order);
356  pd->SetDataFormat(VTKFormat::BINARY);
357  pd->SetHighOrderOutput(true);
358  pd->SetCycle(0);
359  pd->SetTime(0.0);
360  pd->Save();
361  }
362 
363  socketstream sout;
364  if (visualization)
365  {
366  char vishost[] = "localhost";
367  int visport = 19916;
368  sout.open(vishost, visport);
369  if (!sout)
370  {
371  cout << "Unable to connect to GLVis server at "
372  << vishost << ':' << visport << endl;
373  visualization = false;
374  cout << "GLVis visualization disabled.\n";
375  }
376  else
377  {
378  sout.precision(precision);
379  sout << "solution\n" << mesh << u;
380  sout << "pause\n";
381  sout << flush;
382  cout << "GLVis visualization paused."
383  << " Press space (in the GLVis window) to resume it.\n";
384  }
385  }
386 
387  // 8. Define the time-dependent evolution operator describing the ODE
388  // right-hand side, and perform time-integration (looping over the time
389  // iterations, ti, with a time-step dt).
390  FE_Evolution adv(m, k, b);
391 
392  double t = 0.0;
393  adv.SetTime(t);
394  ode_solver->Init(adv);
395 
396  bool done = false;
397  for (int ti = 0; !done; )
398  {
399  double dt_real = min(dt, t_final - t);
400  ode_solver->Step(u, t, dt_real);
401  ti++;
402 
403  done = (t >= t_final - 1e-8*dt);
404 
405  if (done || ti % vis_steps == 0)
406  {
407  cout << "time step: " << ti << ", time: " << t << endl;
408 
409  if (visualization)
410  {
411  sout << "solution\n" << mesh << u << flush;
412  }
413 
414  if (visit)
415  {
416  dc->SetCycle(ti);
417  dc->SetTime(t);
418  dc->Save();
419  }
420 
421  if (paraview)
422  {
423  pd->SetCycle(ti);
424  pd->SetTime(t);
425  pd->Save();
426  }
427  }
428  }
429 
430  // 9. Save the final solution. This output can be viewed later using GLVis:
431  // "glvis -m ex9.mesh -g ex9-final.gf".
432  {
433  ofstream osol("ex9-final.gf");
434  osol.precision(precision);
435  u.Save(osol);
436  }
437 
438  // 10. Free the used memory.
439  delete ode_solver;
440  delete pd;
441  delete dc;
442 
443  return 0;
444 }
445 
446 
447 // Implementation of class FE_Evolution
449  : TimeDependentOperator(_M.Height()), M(_M), K(_K), b(_b), z(_M.Height())
450 {
451  Array<int> ess_tdof_list;
452  if (M.GetAssemblyLevel() == AssemblyLevel::LEGACYFULL)
453  {
454  M_prec = new DSmoother(M.SpMat());
455  M_solver.SetOperator(M.SpMat());
456  dg_solver = new DG_Solver(M.SpMat(), K.SpMat(), *M.FESpace());
457  }
458  else
459  {
460  M_prec = new OperatorJacobiSmoother(M, ess_tdof_list);
461  M_solver.SetOperator(M);
462  dg_solver = NULL;
463  }
464  M_solver.SetPreconditioner(*M_prec);
465  M_solver.iterative_mode = false;
466  M_solver.SetRelTol(1e-9);
467  M_solver.SetAbsTol(0.0);
468  M_solver.SetMaxIter(100);
469  M_solver.SetPrintLevel(0);
470 }
471 
472 void FE_Evolution::Mult(const Vector &x, Vector &y) const
473 {
474  // y = M^{-1} (K x + b)
475  K.Mult(x, z);
476  z += b;
477  M_solver.Mult(z, y);
478 }
479 
480 void FE_Evolution::ImplicitSolve(const double dt, const Vector &x, Vector &k)
481 {
482  MFEM_VERIFY(dg_solver != NULL,
483  "Implicit time integration is not supported with partial assembly");
484  K.Mult(x, z);
485  z += b;
486  dg_solver->SetTimeStep(dt);
487  dg_solver->Mult(z, k);
488 }
489 
491 {
492  delete M_prec;
493  delete dg_solver;
494 }
495 
496 // Velocity coefficient
497 void velocity_function(const Vector &x, Vector &v)
498 {
499  int dim = x.Size();
500 
501  // map to the reference [-1,1] domain
502  Vector X(dim);
503  for (int i = 0; i < dim; i++)
504  {
505  double center = (bb_min[i] + bb_max[i]) * 0.5;
506  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
507  }
508 
509  switch (problem)
510  {
511  case 0:
512  {
513  // Translations in 1D, 2D, and 3D
514  switch (dim)
515  {
516  case 1: v(0) = 1.0; break;
517  case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.); break;
518  case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
519  break;
520  }
521  break;
522  }
523  case 1:
524  case 2:
525  {
526  // Clockwise rotation in 2D around the origin
527  const double w = M_PI/2;
528  switch (dim)
529  {
530  case 1: v(0) = 1.0; break;
531  case 2: v(0) = w*X(1); v(1) = -w*X(0); break;
532  case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0; break;
533  }
534  break;
535  }
536  case 3:
537  {
538  // Clockwise twisting rotation in 2D around the origin
539  const double w = M_PI/2;
540  double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
541  d = d*d;
542  switch (dim)
543  {
544  case 1: v(0) = 1.0; break;
545  case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0); break;
546  case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0; break;
547  }
548  break;
549  }
550  }
551 }
552 
553 // Initial condition
554 double u0_function(const Vector &x)
555 {
556  int dim = x.Size();
557 
558  // map to the reference [-1,1] domain
559  Vector X(dim);
560  for (int i = 0; i < dim; i++)
561  {
562  double center = (bb_min[i] + bb_max[i]) * 0.5;
563  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
564  }
565 
566  switch (problem)
567  {
568  case 0:
569  case 1:
570  {
571  switch (dim)
572  {
573  case 1:
574  return exp(-40.*pow(X(0)-0.5,2));
575  case 2:
576  case 3:
577  {
578  double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
579  if (dim == 3)
580  {
581  const double s = (1. + 0.25*cos(2*M_PI*X(2)));
582  rx *= s;
583  ry *= s;
584  }
585  return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
586  erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
587  }
588  }
589  }
590  case 2:
591  {
592  double x_ = X(0), y_ = X(1), rho, phi;
593  rho = hypot(x_, y_);
594  phi = atan2(y_, x_);
595  return pow(sin(M_PI*rho),2)*sin(3*phi);
596  }
597  case 3:
598  {
599  const double f = M_PI;
600  return sin(f*X(0))*sin(f*X(1));
601  }
602  }
603  return 0.0;
604 }
605 
606 // Inflow boundary condition (zero for the problems considered in this example)
607 double inflow_function(const Vector &x)
608 {
609  switch (problem)
610  {
611  case 0:
612  case 1:
613  case 2:
614  case 3: return 0.0;
615  }
616  return 0.0;
617 }
Vector bb_max
Definition: ex9.cpp:64
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:400
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1234
Conjugate gradient method.
Definition: solvers.hpp:258
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:535
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:267
AssemblyLevel GetAssemblyLevel() const
Returns the assembly level.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
Definition: mesh.cpp:118
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
virtual void Step(Vector &x, double &t, double &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:18
virtual void SetTime(const double _t)
Set the current time.
Definition: operator.hpp:311
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:261
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:638
int main(int argc, char *argv[])
Definition: ex1.cpp:66
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:389
FE_Evolution(FiniteElementSpace &_vfes, Operator &_A, SparseMatrix &_Aflux)
Definition: ex18.hpp:102
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:3417
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
Data type sparse matrix.
Definition: sparsemat.hpp:46
constexpr char vishost[]
virtual void Save()
Save the collection to disk.
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:128
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:8382
constexpr int visport
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
Definition: solvers.hpp:98
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:4165
void SetHighOrderOutput(bool high_order_output_)
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:66
int Dimension() const
Definition: mesh.hpp:788
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:434
void SetTime(double t)
Set physical time (for time-dependent simulations)
Vector bb_min
Definition: ex9.cpp:64
A general vector function coefficient.
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:162
void SetAbsTol(double atol)
Definition: solvers.hpp:97
virtual void ImplicitSolve(const double dt, const Vector &x, Vector &k)
Solve the equation: k = f(x + dt k, t), for the unknown k at the current time t.
Definition: ex9.cpp:480
void SetRelTol(double rtol)
Definition: solvers.hpp:96
void velocity_function(const Vector &x, Vector &v)
Definition: ex9.cpp:497
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:315
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
int problem
Definition: ex15.cpp:54
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:149
GMRES method.
Definition: solvers.hpp:290
virtual void Mult(const Vector &x, Vector &y) const
Perform the action of the operator: y = k = f(x, t), where k solves the algebraic equation F(x...
Definition: ex18.hpp:128
Implicit midpoint method. A-stable, not L-stable.
Definition: ode.hpp:402
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:203
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.
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:1798
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:304
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:2252
virtual ~FE_Evolution()
Definition: ex18.hpp:43
int open(const char hostname[], int port)
Open the socket stream on &#39;port&#39; at &#39;hostname&#39;.
double u0_function(const Vector &x)
Definition: ex9.cpp:554
void SetLevelsOfDetail(int levels_of_detail_)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:272
A general function coefficient.
Vector data type.
Definition: vector.hpp:51
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:91
Vector with associated FE space and LinearFormIntegrators.
Definition: linearform.hpp:23
virtual void Save() override
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix.
Base class for solvers.
Definition: operator.hpp:634
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:118
The classical forward Euler method.
Definition: ode.hpp:116
Abstract operator.
Definition: operator.hpp:24
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
double inflow_function(const Vector &x)
Definition: ex9.cpp:607
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator. Assumes ownership of bfi.
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:221
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:145
alpha (q . grad u, v)