MFEM  v4.1.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 -pa -m ../data/periodic-cube.mesh
23 // ex9 -pa -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 using namespace std;
44 using namespace mfem;
45 
46 // Choice for the problem setup. The fluid velocity, initial condition and
47 // inflow boundary condition are chosen based on this parameter.
48 int problem;
49 
50 // Velocity coefficient
51 void velocity_function(const Vector &x, Vector &v);
52 
53 // Initial condition
54 double u0_function(const Vector &x);
55 
56 // Inflow boundary condition
57 double inflow_function(const Vector &x);
58 
59 // Mesh bounding box
61 
62 class DG_Solver : public Solver
63 {
64 private:
65  SparseMatrix &M, &K, A;
66  GMRESSolver linear_solver;
67  BlockILU prec;
68  double dt;
69 public:
70  DG_Solver(SparseMatrix &M_, SparseMatrix &K_, const FiniteElementSpace &fes)
71  : M(M_),
72  K(K_),
73  prec(fes.GetFE(0)->GetDof(),
74  BlockILU::Reordering::MINIMUM_DISCARDED_FILL),
75  dt(-1.0)
76  {
77  linear_solver.iterative_mode = false;
78  linear_solver.SetRelTol(1e-9);
79  linear_solver.SetAbsTol(0.0);
80  linear_solver.SetMaxIter(100);
81  linear_solver.SetPrintLevel(0);
82  linear_solver.SetPreconditioner(prec);
83  }
84 
85  void SetTimeStep(double dt_)
86  {
87  if (dt_ != dt)
88  {
89  dt = dt_;
90  // Form operator A = M - dt*K
91  A = K;
92  A *= -dt;
93  A += M;
94 
95  // this will also call SetOperator on the preconditioner
96  linear_solver.SetOperator(A);
97  }
98  }
99 
100  void SetOperator(const Operator &op)
101  {
102  linear_solver.SetOperator(op);
103  }
104 
105  virtual void Mult(const Vector &x, Vector &y) const
106  {
107  linear_solver.Mult(x, y);
108  }
109 };
110 
111 /** A time-dependent operator for the right-hand side of the ODE. The DG weak
112  form of du/dt = -v.grad(u) is M du/dt = K u + b, where M and K are the mass
113  and advection matrices, and b describes the flow on the boundary. This can
114  be written as a general ODE, du/dt = M^{-1} (K u + b), and this class is
115  used to evaluate the right-hand side. */
116 class FE_Evolution : public TimeDependentOperator
117 {
118 private:
119  BilinearForm &M, &K;
120  const Vector &b;
121  Solver *M_prec;
122  CGSolver M_solver;
123  DG_Solver *dg_solver;
124 
125  mutable Vector z;
126 
127 public:
128  FE_Evolution(BilinearForm &_M, BilinearForm &_K, const Vector &_b);
129 
130  virtual void Mult(const Vector &x, Vector &y) const;
131  virtual void ImplicitSolve(const double dt, const Vector &x, Vector &k);
132 
133  virtual ~FE_Evolution();
134 };
135 
136 
137 int main(int argc, char *argv[])
138 {
139  // 1. Parse command-line options.
140  problem = 0;
141  const char *mesh_file = "../data/periodic-hexagon.mesh";
142  int ref_levels = 2;
143  int order = 3;
144  bool pa = false;
145  const char *device_config = "cpu";
146  int ode_solver_type = 4;
147  double t_final = 10.0;
148  double dt = 0.01;
149  bool visualization = true;
150  bool visit = false;
151  bool paraview = false;
152  bool binary = false;
153  int vis_steps = 5;
154 
155  int precision = 8;
156  cout.precision(precision);
157 
158  OptionsParser args(argc, argv);
159  args.AddOption(&mesh_file, "-m", "--mesh",
160  "Mesh file to use.");
161  args.AddOption(&problem, "-p", "--problem",
162  "Problem setup to use. See options in velocity_function().");
163  args.AddOption(&ref_levels, "-r", "--refine",
164  "Number of times to refine the mesh uniformly.");
165  args.AddOption(&order, "-o", "--order",
166  "Order (degree) of the finite elements.");
167  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
168  "--no-partial-assembly", "Enable Partial Assembly.");
169  args.AddOption(&device_config, "-d", "--device",
170  "Device configuration string, see Device::Configure().");
171  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
172  "ODE solver: 1 - Forward Euler,\n\t"
173  " 2 - RK2 SSP, 3 - RK3 SSP, 4 - RK4, 6 - RK6,\n\t"
174  " 11 - Backward Euler,\n\t"
175  " 12 - SDIRK23 (L-stable), 13 - SDIRK33,\n\t"
176  " 22 - Implicit Midpoint Method,\n\t"
177  " 23 - SDIRK23 (A-stable), 24 - SDIRK34");
178  args.AddOption(&t_final, "-tf", "--t-final",
179  "Final time; start time is 0.");
180  args.AddOption(&dt, "-dt", "--time-step",
181  "Time step.");
182  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
183  "--no-visualization",
184  "Enable or disable GLVis visualization.");
185  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
186  "--no-visit-datafiles",
187  "Save data files for VisIt (visit.llnl.gov) visualization.");
188  args.AddOption(&paraview, "-paraview", "--paraview-datafiles", "-no-paraview",
189  "--no-paraview-datafiles",
190  "Save data files for ParaView (paraview.org) visualization.");
191  args.AddOption(&binary, "-binary", "--binary-datafiles", "-ascii",
192  "--ascii-datafiles",
193  "Use binary (Sidre) or ascii format for VisIt data files.");
194  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
195  "Visualize every n-th timestep.");
196  args.Parse();
197  if (!args.Good())
198  {
199  args.PrintUsage(cout);
200  return 1;
201  }
202  args.PrintOptions(cout);
203 
204  Device device(device_config);
205  device.Print();
206 
207  // 2. Read the mesh from the given mesh file. We can handle geometrically
208  // periodic meshes in this code.
209  Mesh mesh(mesh_file, 1, 1);
210  int dim = mesh.Dimension();
211 
212  // 3. Define the ODE solver used for time integration. Several explicit
213  // Runge-Kutta methods are available.
214  ODESolver *ode_solver = NULL;
215  switch (ode_solver_type)
216  {
217  // Explicit methods
218  case 1: ode_solver = new ForwardEulerSolver; break;
219  case 2: ode_solver = new RK2Solver(1.0); break;
220  case 3: ode_solver = new RK3SSPSolver; break;
221  case 4: ode_solver = new RK4Solver; break;
222  case 6: ode_solver = new RK6Solver; break;
223  // Implicit (L-stable) methods
224  case 11: ode_solver = new BackwardEulerSolver; break;
225  case 12: ode_solver = new SDIRK23Solver(2); break;
226  case 13: ode_solver = new SDIRK33Solver; break;
227  // Implicit A-stable methods (not L-stable)
228  case 22: ode_solver = new ImplicitMidpointSolver; break;
229  case 23: ode_solver = new SDIRK23Solver; break;
230  case 24: ode_solver = new SDIRK34Solver; break;
231 
232  default:
233  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
234  return 3;
235  }
236 
237  // 4. 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  // 5. 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  // 6. 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  }
273  k.AddDomainIntegrator(new ConvectionIntegrator(velocity, -1.0));
275  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
277  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
278 
279  LinearForm b(&fes);
281  new BoundaryFlowIntegrator(inflow, velocity, -1.0, -0.5));
282 
283  m.Assemble();
284  int skip_zeros = 0;
285  k.Assemble(skip_zeros);
286  b.Assemble();
287  m.Finalize();
288  k.Finalize(skip_zeros);
289 
290  // 7. Define the initial conditions, save the corresponding grid function to
291  // a file and (optionally) save data in the VisIt format and initialize
292  // GLVis visualization.
293  GridFunction u(&fes);
294  u.ProjectCoefficient(u0);
295 
296  {
297  ofstream omesh("ex9.mesh");
298  omesh.precision(precision);
299  mesh.Print(omesh);
300  ofstream osol("ex9-init.gf");
301  osol.precision(precision);
302  u.Save(osol);
303  }
304 
305  // Create data collection for solution output: either VisItDataCollection for
306  // ascii data files, or SidreDataCollection for binary data files.
307  DataCollection *dc = NULL;
308  if (visit)
309  {
310  if (binary)
311  {
312 #ifdef MFEM_USE_SIDRE
313  dc = new SidreDataCollection("Example9", &mesh);
314 #else
315  MFEM_ABORT("Must build with MFEM_USE_SIDRE=YES for binary output.");
316 #endif
317  }
318  else
319  {
320  dc = new VisItDataCollection("Example9", &mesh);
321  dc->SetPrecision(precision);
322  }
323  dc->RegisterField("solution", &u);
324  dc->SetCycle(0);
325  dc->SetTime(0.0);
326  dc->Save();
327  }
328 
329  ParaViewDataCollection *pd = NULL;
330  if (paraview)
331  {
332  pd = new ParaViewDataCollection("Example9", &mesh);
333  pd->SetPrefixPath("ParaView");
334  pd->RegisterField("solution", &u);
335  pd->SetLevelsOfDetail(order);
336  pd->SetDataFormat(VTKFormat::BINARY);
337  pd->SetHighOrderOutput(true);
338  pd->SetCycle(0);
339  pd->SetTime(0.0);
340  pd->Save();
341  }
342 
343  socketstream sout;
344  if (visualization)
345  {
346  char vishost[] = "localhost";
347  int visport = 19916;
348  sout.open(vishost, visport);
349  if (!sout)
350  {
351  cout << "Unable to connect to GLVis server at "
352  << vishost << ':' << visport << endl;
353  visualization = false;
354  cout << "GLVis visualization disabled.\n";
355  }
356  else
357  {
358  sout.precision(precision);
359  sout << "solution\n" << mesh << u;
360  sout << "pause\n";
361  sout << flush;
362  cout << "GLVis visualization paused."
363  << " Press space (in the GLVis window) to resume it.\n";
364  }
365  }
366 
367  // 8. Define the time-dependent evolution operator describing the ODE
368  // right-hand side, and perform time-integration (looping over the time
369  // iterations, ti, with a time-step dt).
370  FE_Evolution adv(m, k, b);
371 
372  double t = 0.0;
373  adv.SetTime(t);
374  ode_solver->Init(adv);
375 
376  bool done = false;
377  for (int ti = 0; !done; )
378  {
379  double dt_real = min(dt, t_final - t);
380  ode_solver->Step(u, t, dt_real);
381  ti++;
382 
383  done = (t >= t_final - 1e-8*dt);
384 
385  if (done || ti % vis_steps == 0)
386  {
387  cout << "time step: " << ti << ", time: " << t << endl;
388 
389  if (visualization)
390  {
391  sout << "solution\n" << mesh << u << flush;
392  }
393 
394  if (visit)
395  {
396  dc->SetCycle(ti);
397  dc->SetTime(t);
398  dc->Save();
399  }
400 
401  if (paraview)
402  {
403  pd->SetCycle(ti);
404  pd->SetTime(t);
405  pd->Save();
406  }
407  }
408  }
409 
410  // 9. Save the final solution. This output can be viewed later using GLVis:
411  // "glvis -m ex9.mesh -g ex9-final.gf".
412  {
413  ofstream osol("ex9-final.gf");
414  osol.precision(precision);
415  u.Save(osol);
416  }
417 
418  // 10. Free the used memory.
419  delete ode_solver;
420  delete pd;
421  delete dc;
422 
423  return 0;
424 }
425 
426 
427 // Implementation of class FE_Evolution
429  : TimeDependentOperator(_M.Height()), M(_M), K(_K), b(_b), z(_M.Height())
430 {
431  bool pa = M.GetAssemblyLevel() == AssemblyLevel::PARTIAL;
432  Array<int> ess_tdof_list;
433  if (pa)
434  {
435  M_prec = new OperatorJacobiSmoother(M, ess_tdof_list);
436  M_solver.SetOperator(M);
437  dg_solver = NULL;
438  }
439  else
440  {
441  M_prec = new DSmoother(M.SpMat());
442  dg_solver = new DG_Solver(M.SpMat(), K.SpMat(), *M.FESpace());
443  M_solver.SetOperator(M.SpMat());
444  }
445  M_solver.SetPreconditioner(*M_prec);
446  M_solver.iterative_mode = false;
447  M_solver.SetRelTol(1e-9);
448  M_solver.SetAbsTol(0.0);
449  M_solver.SetMaxIter(100);
450  M_solver.SetPrintLevel(0);
451 }
452 
453 void FE_Evolution::Mult(const Vector &x, Vector &y) const
454 {
455  // y = M^{-1} (K x + b)
456  K.Mult(x, z);
457  z += b;
458  M_solver.Mult(z, y);
459 }
460 
461 void FE_Evolution::ImplicitSolve(const double dt, const Vector &x, Vector &k)
462 {
463  MFEM_VERIFY(dg_solver != NULL,
464  "Implicit time integration is not supported with partial assembly");
465  K.Mult(x, z);
466  z += b;
467  dg_solver->SetTimeStep(dt);
468  dg_solver->Mult(z, k);
469 }
470 
472 {
473  delete M_prec;
474  delete dg_solver;
475 }
476 
477 // Velocity coefficient
478 void velocity_function(const Vector &x, Vector &v)
479 {
480  int dim = x.Size();
481 
482  // map to the reference [-1,1] domain
483  Vector X(dim);
484  for (int i = 0; i < dim; i++)
485  {
486  double center = (bb_min[i] + bb_max[i]) * 0.5;
487  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
488  }
489 
490  switch (problem)
491  {
492  case 0:
493  {
494  // Translations in 1D, 2D, and 3D
495  switch (dim)
496  {
497  case 1: v(0) = 1.0; break;
498  case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.); break;
499  case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
500  break;
501  }
502  break;
503  }
504  case 1:
505  case 2:
506  {
507  // Clockwise rotation in 2D around the origin
508  const double w = M_PI/2;
509  switch (dim)
510  {
511  case 1: v(0) = 1.0; break;
512  case 2: v(0) = w*X(1); v(1) = -w*X(0); break;
513  case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0; break;
514  }
515  break;
516  }
517  case 3:
518  {
519  // Clockwise twisting rotation in 2D around the origin
520  const double w = M_PI/2;
521  double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
522  d = d*d;
523  switch (dim)
524  {
525  case 1: v(0) = 1.0; break;
526  case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0); break;
527  case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0; break;
528  }
529  break;
530  }
531  }
532 }
533 
534 // Initial condition
535 double u0_function(const Vector &x)
536 {
537  int dim = x.Size();
538 
539  // map to the reference [-1,1] domain
540  Vector X(dim);
541  for (int i = 0; i < dim; i++)
542  {
543  double center = (bb_min[i] + bb_max[i]) * 0.5;
544  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
545  }
546 
547  switch (problem)
548  {
549  case 0:
550  case 1:
551  {
552  switch (dim)
553  {
554  case 1:
555  return exp(-40.*pow(X(0)-0.5,2));
556  case 2:
557  case 3:
558  {
559  double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
560  if (dim == 3)
561  {
562  const double s = (1. + 0.25*cos(2*M_PI*X(2)));
563  rx *= s;
564  ry *= s;
565  }
566  return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
567  erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
568  }
569  }
570  }
571  case 2:
572  {
573  double x_ = X(0), y_ = X(1), rho, phi;
574  rho = hypot(x_, y_);
575  phi = atan2(y_, x_);
576  return pow(sin(M_PI*rho),2)*sin(3*phi);
577  }
578  case 3:
579  {
580  const double f = M_PI;
581  return sin(f*X(0))*sin(f*X(1));
582  }
583  }
584  return 0.0;
585 }
586 
587 // Inflow boundary condition (zero for the problems considered in this example)
588 double inflow_function(const Vector &x)
589 {
590  switch (problem)
591  {
592  case 0:
593  case 1:
594  case 2:
595  case 3: return 0.0;
596  }
597  return 0.0;
598 }
Vector bb_max
Definition: ex9.cpp:60
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:381
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1188
Conjugate gradient method.
Definition: solvers.hpp:152
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
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:358
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:259
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:120
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 RegisterField(const std::string &field_name, mfem::GridFunction *gf) override
Add a grid function to the collection.
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:303
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
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:236
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::FULL.
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:504
int main(int argc, char *argv[])
Definition: ex1.cpp:62
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:360
FE_Evolution(FiniteElementSpace &_vfes, Operator &_A, SparseMatrix &_Aflux)
Definition: ex18.hpp:104
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:2731
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
Data type sparse matrix.
Definition: sparsemat.hpp:40
virtual void Save()
Save the collection to disk.
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:84
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7982
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
Definition: solvers.hpp:65
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3924
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:744
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetTime(double t)
Set physical time (for time-dependent simulations)
Vector bb_min
Definition: ex9.cpp:60
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:145
void SetAbsTol(double atol)
Definition: solvers.hpp:64
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:461
AssemblyLevel GetAssemblyLevel()
Get the assembly level.
void SetRelTol(double rtol)
Definition: solvers.hpp:63
void velocity_function(const Vector &x, Vector &v)
Definition: ex9.cpp:478
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:314
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)
Definition: optparser.hpp:76
int problem
Definition: ex15.cpp:54
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:132
GMRES method.
Definition: solvers.hpp:184
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:130
Implicit midpoint method. A-stable, not L-stable.
Definition: ode.hpp:373
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:191
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
int dim
Definition: ex24.cpp:43
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 associated with i&#39;th element.
Definition: fespace.cpp:1642
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1685
virtual ~FE_Evolution()
Definition: ex18.hpp:43
int open(const char hostname[], int port)
double u0_function(const Vector &x)
Definition: ex9.cpp:535
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:166
class for C-function coefficient
Vector data type.
Definition: vector.hpp:48
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:91
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
virtual void Save() override
const SparseMatrix & SpMat() const
Returns a reference to the sparse matrix.
Base class for solvers.
Definition: operator.hpp:500
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:114
The classical forward Euler method.
Definition: ode.hpp:99
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:588
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:143
bool Good() const
Definition: optparser.hpp:125
alpha (q . grad u, v)