MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
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
47using namespace std;
48using 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.
53
54// Velocity coefficient
55void velocity_function(const Vector &x, Vector &v);
56
57// Initial condition
58double u0_function(const Vector &x);
59
60// Inflow boundary condition
61double inflow_function(const Vector &x);
62
63// Mesh bounding box
65
66class DG_Solver : public Solver
67{
68private:
69 SparseMatrix &M, &K, A;
70 GMRESSolver linear_solver;
71 BlockILU prec;
72 double dt;
73public:
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. */
120class FE_Evolution : public TimeDependentOperator
121{
122private:
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
131public:
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
141int main(int argc, char *argv[])
142{
143 // 0. Initialize SUNDIALS.
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.
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 constexpr double alpha = -1.0;
289
290 LinearForm b(&fes);
291 b.AddBdrFaceIntegrator(
292 new BoundaryFlowIntegrator(inflow, velocity, alpha));
293
294 m.Assemble();
295 int skip_zeros = 0;
296 k.Assemble(skip_zeros);
297 b.Assemble();
298 m.Finalize();
299 k.Finalize(skip_zeros);
300
301 // 6. Define the initial conditions, save the corresponding grid function to
302 // a file and (optionally) save data in the VisIt format and initialize
303 // GLVis visualization.
304 GridFunction u(&fes);
305 u.ProjectCoefficient(u0);
306
307 {
308 ofstream omesh("ex9.mesh");
309 omesh.precision(precision);
310 mesh.Print(omesh);
311 ofstream osol("ex9-init.gf");
312 osol.precision(precision);
313 u.Save(osol);
314 }
315
316 // Create data collection for solution output: either VisItDataCollection for
317 // ascii data files, or SidreDataCollection for binary data files.
318 DataCollection *dc = NULL;
319 if (visit)
320 {
321 if (binary)
322 {
323#ifdef MFEM_USE_SIDRE
324 dc = new SidreDataCollection("Example9", &mesh);
325#else
326 MFEM_ABORT("Must build with MFEM_USE_SIDRE=YES for binary output.");
327#endif
328 }
329 else
330 {
331 dc = new VisItDataCollection("Example9", &mesh);
332 dc->SetPrecision(precision);
333 }
334 dc->RegisterField("solution", &u);
335 dc->SetCycle(0);
336 dc->SetTime(0.0);
337 dc->Save();
338 }
339
340 ParaViewDataCollection *pd = NULL;
341 if (paraview)
342 {
343 pd = new ParaViewDataCollection("Example9", &mesh);
344 pd->SetPrefixPath("ParaView");
345 pd->RegisterField("solution", &u);
346 pd->SetLevelsOfDetail(order);
347 pd->SetDataFormat(VTKFormat::BINARY);
348 pd->SetHighOrderOutput(true);
349 pd->SetCycle(0);
350 pd->SetTime(0.0);
351 pd->Save();
352 }
353
354 socketstream sout;
355 if (visualization)
356 {
357 char vishost[] = "localhost";
358 int visport = 19916;
359 sout.open(vishost, visport);
360 if (!sout)
361 {
362 cout << "Unable to connect to GLVis server at "
363 << vishost << ':' << visport << endl;
364 visualization = false;
365 cout << "GLVis visualization disabled.\n";
366 }
367 else
368 {
369 sout.precision(precision);
370 sout << "solution\n" << mesh << u;
371 sout << "pause\n";
372 sout << flush;
373 cout << "GLVis visualization paused."
374 << " Press space (in the GLVis window) to resume it.\n";
375 }
376 }
377
378 // 7. Define the time-dependent evolution operator describing the ODE
379 // right-hand side, and define the ODE solver used for time integration.
380 FE_Evolution adv(m, k, b);
381
382 double t = 0.0;
383 adv.SetTime(t);
384
385 // Create the time integrator
386 ODESolver *ode_solver = NULL;
387 CVODESolver *cvode = NULL;
388 ARKStepSolver *arkode = NULL;
389 switch (ode_solver_type)
390 {
391 case 1: ode_solver = new ForwardEulerSolver; break;
392 case 2: ode_solver = new RK2Solver(1.0); break;
393 case 3: ode_solver = new RK3SSPSolver; break;
394 case 4: ode_solver = new RK4Solver; break;
395 case 6: ode_solver = new RK6Solver; break;
396 case 7:
397 cvode = new CVODESolver(CV_ADAMS);
398 cvode->Init(adv);
399 cvode->SetSStolerances(reltol, abstol);
400 cvode->SetMaxStep(dt);
402 ode_solver = cvode; break;
403 case 8:
405 arkode->Init(adv);
406 arkode->SetSStolerances(reltol, abstol);
407 arkode->SetMaxStep(dt);
408 arkode->SetOrder(4);
409 ode_solver = arkode; break;
410 case 9:
412 arkode->Init(adv);
413 arkode->SetSStolerances(reltol, abstol);
414 arkode->SetMaxStep(dt);
416 ode_solver = arkode; break;
417 }
418
419 // Initialize MFEM integrators, SUNDIALS integrators are initialized above
420 if (ode_solver_type < 7) { ode_solver->Init(adv); }
421
422 // 8. Perform time-integration (looping over the time iterations, ti,
423 // with a time-step dt).
424 bool done = false;
425 for (int ti = 0; !done; )
426 {
427 double dt_real = min(dt, t_final - t);
428 ode_solver->Step(u, t, dt_real);
429 ti++;
430
431 done = (t >= t_final - 1e-8*dt);
432
433 if (done || ti % vis_steps == 0)
434 {
435 cout << "time step: " << ti << ", time: " << t << endl;
436 if (cvode) { cvode->PrintInfo(); }
437 if (arkode) { arkode->PrintInfo(); }
438
439 if (visualization)
440 {
441 sout << "solution\n" << mesh << u << flush;
442 }
443
444 if (visit)
445 {
446 dc->SetCycle(ti);
447 dc->SetTime(t);
448 dc->Save();
449 }
450
451 if (paraview)
452 {
453 pd->SetCycle(ti);
454 pd->SetTime(t);
455 pd->Save();
456 }
457 }
458 }
459
460 // 9. Save the final solution. This output can be viewed later using GLVis:
461 // "glvis -m ex9.mesh -g ex9-final.gf".
462 {
463 ofstream osol("ex9-final.gf");
464 osol.precision(precision);
465 u.Save(osol);
466 }
467
468 // 10. Free the used memory.
469 delete ode_solver;
470 delete pd;
471 delete dc;
472
473 return 0;
474}
475
476
477// Implementation of class FE_Evolution
478FE_Evolution::FE_Evolution(BilinearForm &M_, BilinearForm &K_, const Vector &b_)
479 : TimeDependentOperator(M_.FESpace()->GetTrueVSize()),
480 M(M_), K(K_), b(b_), z(height)
481{
482 Array<int> ess_tdof_list;
483 if (M.GetAssemblyLevel() == AssemblyLevel::LEGACY)
484 {
485 M_prec = new DSmoother(M.SpMat());
486 M_solver.SetOperator(M.SpMat());
487 dg_solver = new DG_Solver(M.SpMat(), K.SpMat(), *M.FESpace());
488 }
489 else
490 {
491 M_prec = new OperatorJacobiSmoother(M, ess_tdof_list);
492 M_solver.SetOperator(M);
493 dg_solver = NULL;
494 }
495 M_solver.SetPreconditioner(*M_prec);
496 M_solver.iterative_mode = false;
497 M_solver.SetRelTol(1e-9);
498 M_solver.SetAbsTol(0.0);
499 M_solver.SetMaxIter(100);
500 M_solver.SetPrintLevel(0);
501}
502
503void FE_Evolution::Mult(const Vector &x, Vector &y) const
504{
505 // y = M^{-1} (K x + b)
506 K.Mult(x, z);
507 z += b;
508 M_solver.Mult(z, y);
509}
510
511void FE_Evolution::ImplicitSolve(const double dt, const Vector &x, Vector &k)
512{
513 MFEM_VERIFY(dg_solver != NULL,
514 "Implicit time integration is not supported with partial assembly");
515 K.Mult(x, z);
516 z += b;
517 dg_solver->SetTimeStep(dt);
518 dg_solver->Mult(z, k);
519}
520
521FE_Evolution::~FE_Evolution()
522{
523 delete M_prec;
524 delete dg_solver;
525}
526
527// Velocity coefficient
529{
530 int dim = x.Size();
531
532 // map to the reference [-1,1] domain
533 Vector X(dim);
534 for (int i = 0; i < dim; i++)
535 {
536 double center = (bb_min[i] + bb_max[i]) * 0.5;
537 X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
538 }
539
540 switch (problem)
541 {
542 case 0:
543 {
544 // Translations in 1D, 2D, and 3D
545 switch (dim)
546 {
547 case 1: v(0) = 1.0; break;
548 case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.); break;
549 case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
550 break;
551 }
552 break;
553 }
554 case 1:
555 case 2:
556 {
557 // Clockwise rotation in 2D around the origin
558 const double w = M_PI/2;
559 switch (dim)
560 {
561 case 1: v(0) = 1.0; break;
562 case 2: v(0) = w*X(1); v(1) = -w*X(0); break;
563 case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0; break;
564 }
565 break;
566 }
567 case 3:
568 {
569 // Clockwise twisting rotation in 2D around the origin
570 const double w = M_PI/2;
571 double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
572 d = d*d;
573 switch (dim)
574 {
575 case 1: v(0) = 1.0; break;
576 case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0); break;
577 case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0; break;
578 }
579 break;
580 }
581 }
582}
583
584// Initial condition
585double u0_function(const Vector &x)
586{
587 int dim = x.Size();
588
589 // map to the reference [-1,1] domain
590 Vector X(dim);
591 for (int i = 0; i < dim; i++)
592 {
593 double center = (bb_min[i] + bb_max[i]) * 0.5;
594 X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
595 }
596
597 switch (problem)
598 {
599 case 0:
600 case 1:
601 {
602 switch (dim)
603 {
604 case 1:
605 return exp(-40.*pow(X(0)-0.5,2));
606 case 2:
607 case 3:
608 {
609 double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
610 if (dim == 3)
611 {
612 const double s = (1. + 0.25*cos(2*M_PI*X(2)));
613 rx *= s;
614 ry *= s;
615 }
616 return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
617 erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
618 }
619 }
620 }
621 case 2:
622 {
623 double x_ = X(0), y_ = X(1), rho, phi;
624 rho = hypot(x_, y_);
625 phi = atan2(y_, x_);
626 return pow(sin(M_PI*rho),2)*sin(3*phi);
627 }
628 case 3:
629 {
630 const double f = M_PI;
631 return sin(f*X(0))*sin(f*X(1));
632 }
633 }
634 return 0.0;
635}
636
637// Inflow boundary condition (zero for the problems considered in this example)
638double inflow_function(const Vector &x)
639{
640 switch (problem)
641 {
642 case 0:
643 case 1:
644 case 2:
645 case 3: return 0.0;
646 }
647 return 0.0;
648}
Interface to ARKode's ARKStep module – additive Runge-Kutta methods.
Definition sundials.hpp:673
void SetMaxStep(double dt_max)
Set the maximum time step.
void PrintInfo() const
Print various ARKStep statistics.
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults.
void SetOrder(int order)
Chooses integration order for all explicit / implicit / IMEX methods.
@ EXPLICIT
Explicit RK method.
Definition sundials.hpp:678
void SetERKTableNum(ARKODE_ERKTableID table_id)
Choose a specific Butcher table for an explicit RK method.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
@ GaussLobatto
Closed type.
Definition fe_base.hpp:32
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY....
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication: .
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator. Assumes ownership of bfi.
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator. Assumes ownership of bfi.
Conjugate gradient method.
Definition solvers.hpp:513
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the Conjugate Gradient method.
Definition solvers.cpp:718
Interface to the CVODE library – linear multi-step methods.
Definition sundials.hpp:386
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
Definition sundials.cpp:875
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition sundials.cpp:893
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
Definition sundials.cpp:709
void PrintInfo() const
Print various CVODE statistics.
Definition sundials.cpp:919
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
Definition sundials.cpp:855
Data type for scaled Jacobi-type smoother of sparse matrix.
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetTime(real_t t)
Set physical time (for time-dependent simulations)
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
virtual void Save()
Save the collection to disk.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition device.hpp:123
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition device.cpp:286
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition fespace.cpp:3168
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition fespace.hpp:713
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition fe_base.hpp:329
The classical forward Euler method.
Definition ode.hpp:118
A general function coefficient.
GMRES method.
Definition solvers.hpp:547
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the GMRES method.
Definition solvers.cpp:980
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
virtual void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
Definition solvers.cpp:179
void SetRelTol(real_t rtol)
Definition solvers.hpp:209
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition solvers.cpp:173
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition solvers.cpp:71
void SetMaxIter(int max_it)
Definition solvers.hpp:211
void SetAbsTol(real_t atol)
Definition solvers.hpp:210
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
Vector with associated FE space and LinearFormIntegrators.
Mesh data type.
Definition mesh.hpp:56
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition mesh.hpp:290
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Definition mesh.hpp:2288
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
Definition mesh.cpp:138
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
Definition mesh.cpp:6211
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition ode.hpp:24
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:18
virtual void Step(Vector &x, real_t &t, real_t &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition solvers.hpp:313
Abstract operator.
Definition operator.hpp:25
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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 'var' to receive the value. Enable/disable tags are used to set the bool...
Definition optparser.hpp:82
bool Good() const
Return true if the command line options were parsed successfully.
Helper class for ParaView visualization data.
void SetHighOrderOutput(bool high_order_output_)
void SetLevelsOfDetail(int levels_of_detail_)
virtual void Save() override
void SetDataFormat(VTKFormat fmt)
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition ode.hpp:151
The classical explicit forth-order Runge-Kutta method, RK4.
Definition ode.hpp:164
Data collection with Sidre routines following the Conduit mesh blueprint specification.
Base class for solvers.
Definition operator.hpp:683
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition operator.hpp:686
Data type sparse matrix.
Definition sparsemat.hpp:51
static void Init()
Definition sundials.cpp:164
Base abstract class for first order time dependent operators.
Definition operator.hpp:317
virtual void SetTime(const real_t t_)
Set the current time.
Definition operator.hpp:360
A general vector function coefficient.
Vector data type.
Definition vector.hpp:80
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
Data collection with VisIt I/O routines.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
const real_t alpha
Definition ex15.cpp:369
int dim
Definition ex24.cpp:53
int main()
real_t b
Definition lissajous.cpp:42
const int visport
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
const char vishost[]
RefCoord t[3]
RefCoord s[3]
void velocity_function(const Vector &x, Vector &v)
Definition ex9.cpp:528
double u0_function(const Vector &x)
Definition ex9.cpp:585
int problem
Definition ex9.cpp:52
double inflow_function(const Vector &x)
Definition ex9.cpp:638
Vector bb_min
Definition ex9.cpp:64
Vector bb_max
Definition ex9.cpp:64
constexpr ARKODE_ERKTableID ARKODE_FEHLBERG_13_7_8
Definition sundials.hpp:65