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