MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ex16p.cpp
Go to the documentation of this file.
1// MFEM Example 16 - Parallel Version
2// SUNDIALS Modification
3//
4// Compile with: make ex16p
5//
6// Sample runs:
7// mpirun -np 4 ex16p
8// mpirun -np 4 ex16p -m ../../data/inline-tri.mesh
9// mpirun -np 4 ex16p -m ../../data/disc-nurbs.mesh -tf 2
10// mpirun -np 4 ex16p -s 12 -a 0.0 -k 1.0
11// mpirun -np 4 ex16p -s 8 -a 1.0 -k 0.0 -dt 4e-6 -tf 2e-2 -vs 50
12// mpirun -np 8 ex16p -s 9 -a 0.5 -k 0.5 -o 4 -dt 8e-6 -tf 2e-2 -vs 50
13// mpirun -np 4 ex16p -s 10 -dt 2.0e-4 -tf 4.0e-2
14// mpirun -np 16 ex16p -m ../../data/fichera-q2.mesh
15// mpirun -np 16 ex16p -m ../../data/escher-p2.mesh
16// mpirun -np 8 ex16p -m ../../data/beam-tet.mesh -tf 10 -dt 0.1
17// mpirun -np 4 ex16p -m ../../data/amr-quad.mesh -o 4 -rs 0 -rp 0
18// mpirun -np 4 ex16p -m ../../data/amr-hex.mesh -o 2 -rs 0 -rp 0
19//
20// Description: This example solves a time dependent nonlinear heat equation
21// problem of the form du/dt = C(u), with a non-linear diffusion
22// operator C(u) = \nabla \cdot (\kappa + \alpha u) \nabla u.
23//
24// The example demonstrates the use of nonlinear operators (the
25// class ConductionOperator defining C(u)), as well as their
26// implicit time integration. Note that implementing the method
27// ConductionOperator::ImplicitSolve is the only requirement for
28// high-order implicit (SDIRK) time integration. By default, this
29// example uses the SUNDIALS ODE solvers from CVODE and ARKODE.
30//
31// We recommend viewing examples 2, 9 and 10 before viewing this
32// example.
33
34#include "mfem.hpp"
35#include <fstream>
36#include <iostream>
37
38using namespace std;
39using namespace mfem;
40
41/** After spatial discretization, the conduction model can be written as:
42 *
43 * du/dt = M^{-1}(-Ku)
44 *
45 * where u is the vector representing the temperature, M is the mass matrix,
46 * and K is the diffusion operator with diffusivity depending on u:
47 * (\kappa + \alpha u).
48 *
49 * Class ConductionOperator represents the right-hand side of the above ODE.
50 */
51class ConductionOperator : public TimeDependentOperator
52{
53protected:
54 ParFiniteElementSpace &fespace;
55 Array<int> ess_tdof_list; // this list remains empty for pure Neumann b.c.
56
59
60 HypreParMatrix Mmat;
61 HypreParMatrix Kmat;
62 HypreParMatrix *T; // T = M + dt K
63 double current_dt;
64
65 CGSolver M_solver; // Krylov solver for inverting the mass matrix M
66 HypreSmoother M_prec; // Preconditioner for the mass matrix M
67
68 CGSolver T_solver; // Implicit solver for T = M + dt K
69 HypreSmoother T_prec; // Preconditioner for the implicit solver
70
71 double alpha, kappa;
72
73 mutable Vector z; // auxiliary vector
74
75public:
76 ConductionOperator(ParFiniteElementSpace &f, double alpha, double kappa,
77 const Vector &u);
78
79 virtual void Mult(const Vector &u, Vector &du_dt) const;
80
81 /** Solve the Backward-Euler equation: k = f(u + dt*k, t), for the unknown k.
82 This is the only requirement for high-order SDIRK implicit integration.*/
83 virtual void ImplicitSolve(const double dt, const Vector &u, Vector &k);
84
85 /** Setup the system (M + dt K) x = M b. This method is used by the implicit
86 SUNDIALS solvers. */
87 virtual int SUNImplicitSetup(const Vector &x, const Vector &fx,
88 int jok, int *jcur, double gamma);
89
90 /** Solve the system (M + dt K) x = M b. This method is used by the implicit
91 SUNDIALS solvers. */
92 virtual int SUNImplicitSolve(const Vector &b, Vector &x, double tol);
93
94 /// Update the diffusion BilinearForm K using the given true-dof vector `u`.
95 void SetParameters(const Vector &u);
96
97 virtual ~ConductionOperator();
98};
99
100double InitialTemperature(const Vector &x);
101
102int main(int argc, char *argv[])
103{
104 // 1. Initialize MPI, HYPRE, and SUNDIALS.
105 Mpi::Init(argc, argv);
106 int num_procs = Mpi::WorldSize();
107 int myid = Mpi::WorldRank();
108 Hypre::Init();
110
111 // 2. Parse command-line options.
112 const char *mesh_file = "../../data/star.mesh";
113 int ser_ref_levels = 2;
114 int par_ref_levels = 1;
115 int order = 2;
116 int ode_solver_type = 9; // CVODE implicit BDF
117 double t_final = 0.5;
118 double dt = 1.0e-2;
119 double alpha = 1.0e-2;
120 double kappa = 0.5;
121 bool visualization = true;
122 bool visit = false;
123 int vis_steps = 5;
124
125 // Relative and absolute tolerances for CVODE and ARKODE.
126 const double reltol = 1e-4, abstol = 1e-4;
127
128 int precision = 8;
129 cout.precision(precision);
130
131 OptionsParser args(argc, argv);
132 args.AddOption(&mesh_file, "-m", "--mesh",
133 "Mesh file to use.");
134 args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
135 "Number of times to refine the mesh uniformly in serial.");
136 args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
137 "Number of times to refine the mesh uniformly in parallel.");
138 args.AddOption(&order, "-o", "--order",
139 "Order (degree) of the finite elements.");
140 args.AddOption(&ode_solver_type, "-s", "--ode-solver",
141 "ODE solver:\n\t"
142 "1 - Forward Euler,\n\t"
143 "2 - RK2,\n\t"
144 "3 - RK3 SSP,\n\t"
145 "4 - RK4,\n\t"
146 "5 - Backward Euler,\n\t"
147 "6 - SDIRK 2,\n\t"
148 "7 - SDIRK 3,\n\t"
149 "8 - CVODE (implicit Adams),\n\t"
150 "9 - CVODE (implicit BDF),\n\t"
151 "10 - ARKODE (default explicit),\n\t"
152 "11 - ARKODE (explicit Fehlberg-6-4-5),\n\t"
153 "12 - ARKODE (default impicit).");
154 args.AddOption(&t_final, "-tf", "--t-final",
155 "Final time; start time is 0.");
156 args.AddOption(&dt, "-dt", "--time-step",
157 "Time step.");
158 args.AddOption(&alpha, "-a", "--alpha",
159 "Alpha coefficient.");
160 args.AddOption(&kappa, "-k", "--kappa",
161 "Kappa coefficient offset.");
162 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
163 "--no-visualization",
164 "Enable or disable GLVis visualization.");
165 args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
166 "--no-visit-datafiles",
167 "Save data files for VisIt (visit.llnl.gov) visualization.");
168 args.AddOption(&vis_steps, "-vs", "--visualization-steps",
169 "Visualize every n-th timestep.");
170 args.Parse();
171 if (!args.Good())
172 {
173 args.PrintUsage(cout);
174 return 1;
175 }
176
177 if (myid == 0)
178 {
179 args.PrintOptions(cout);
180 }
181
182 // check for valid ODE solver option
183 if (ode_solver_type < 1 || ode_solver_type > 12)
184 {
185 if (myid == 0)
186 {
187 cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
188 }
189 return 1;
190 }
191
192 // 3. Read the serial mesh from the given mesh file on all processors. We can
193 // handle triangular, quadrilateral, tetrahedral and hexahedral meshes
194 // with the same code.
195 Mesh *mesh = new Mesh(mesh_file, 1, 1);
196 int dim = mesh->Dimension();
197
198 // 4. Refine the mesh in serial to increase the resolution. In this example
199 // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
200 // a command-line parameter.
201 for (int lev = 0; lev < ser_ref_levels; lev++)
202 {
203 mesh->UniformRefinement();
204 }
205
206 // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
207 // this mesh further in parallel to increase the resolution. Once the
208 // parallel mesh is defined, the serial mesh can be deleted.
209 ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
210 delete mesh;
211 for (int lev = 0; lev < par_ref_levels; lev++)
212 {
213 pmesh->UniformRefinement();
214 }
215
216 // 6. Define the vector finite element space representing the current and the
217 // initial temperature, u_ref.
218 H1_FECollection fe_coll(order, dim);
219 ParFiniteElementSpace fespace(pmesh, &fe_coll);
220
221 int fe_size = fespace.GlobalTrueVSize();
222 if (myid == 0)
223 {
224 cout << "Number of temperature unknowns: " << fe_size << endl;
225 }
226
227 ParGridFunction u_gf(&fespace);
228
229 // 7. Set the initial conditions for u. All boundaries are considered
230 // natural.
232 u_gf.ProjectCoefficient(u_0);
233 Vector u;
234 u_gf.GetTrueDofs(u);
235
236 // 8. Initialize the conduction operator and the VisIt visualization.
237 ConductionOperator oper(fespace, alpha, kappa, u);
238
239 u_gf.SetFromTrueDofs(u);
240 {
241 ostringstream mesh_name, sol_name;
242 mesh_name << "ex16-mesh." << setfill('0') << setw(6) << myid;
243 sol_name << "ex16-init." << setfill('0') << setw(6) << myid;
244 ofstream omesh(mesh_name.str().c_str());
245 omesh.precision(precision);
246 pmesh->Print(omesh);
247 ofstream osol(sol_name.str().c_str());
248 osol.precision(precision);
249 u_gf.Save(osol);
250 }
251
252 VisItDataCollection visit_dc("Example16-Parallel", pmesh);
253 visit_dc.RegisterField("temperature", &u_gf);
254 if (visit)
255 {
256 visit_dc.SetCycle(0);
257 visit_dc.SetTime(0.0);
258 visit_dc.Save();
259 }
260
261 socketstream sout;
262 if (visualization)
263 {
264 char vishost[] = "localhost";
265 int visport = 19916;
266 sout.open(vishost, visport);
267 sout << "parallel " << num_procs << " " << myid << endl;
268 int good = sout.good(), all_good;
269 MPI_Allreduce(&good, &all_good, 1, MPI_INT, MPI_MIN, pmesh->GetComm());
270 if (!all_good)
271 {
272 sout.close();
273 visualization = false;
274 if (myid == 0)
275 {
276 cout << "Unable to connect to GLVis server at "
277 << vishost << ':' << visport << endl;
278 cout << "GLVis visualization disabled.\n";
279 }
280 }
281 else
282 {
283 sout.precision(precision);
284 sout << "solution\n" << *pmesh << u_gf;
285 sout << "pause\n";
286 sout << flush;
287 if (myid == 0)
288 {
289 cout << "GLVis visualization paused."
290 << " Press space (in the GLVis window) to resume it.\n";
291 }
292 }
293 }
294
295 // 9. Define the ODE solver used for time integration.
296 double t = 0.0;
297 ODESolver *ode_solver = NULL;
298 CVODESolver *cvode = NULL;
299 ARKStepSolver *arkode = NULL;
300 switch (ode_solver_type)
301 {
302 // MFEM explicit methods
303 case 1: ode_solver = new ForwardEulerSolver; break;
304 case 2: ode_solver = new RK2Solver(0.5); break; // midpoint method
305 case 3: ode_solver = new RK3SSPSolver; break;
306 case 4: ode_solver = new RK4Solver; break;
307 // MFEM implicit L-stable methods
308 case 5: ode_solver = new BackwardEulerSolver; break;
309 case 6: ode_solver = new SDIRK23Solver(2); break;
310 case 7: ode_solver = new SDIRK33Solver; break;
311 // CVODE
312 case 8:
313 cvode = new CVODESolver(MPI_COMM_WORLD, CV_ADAMS);
314 cvode->Init(oper);
315 cvode->SetSStolerances(reltol, abstol);
316 cvode->SetMaxStep(dt);
317 ode_solver = cvode; break;
318 case 9:
319 cvode = new CVODESolver(MPI_COMM_WORLD, CV_BDF);
320 cvode->Init(oper);
321 cvode->SetSStolerances(reltol, abstol);
322 cvode->SetMaxStep(dt);
323 ode_solver = cvode; break;
324 // ARKODE
325 case 10:
326 case 11:
327 arkode = new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::EXPLICIT);
328 arkode->Init(oper);
329 arkode->SetSStolerances(reltol, abstol);
330 arkode->SetMaxStep(dt);
331 if (ode_solver_type == 11)
332 {
334 }
335 ode_solver = arkode; break;
336 case 12:
337 arkode = new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::IMPLICIT);
338 arkode->Init(oper);
339 arkode->SetSStolerances(reltol, abstol);
340 arkode->SetMaxStep(dt);
341 ode_solver = arkode; break;
342 }
343
344 // Initialize MFEM integrators, SUNDIALS integrators are initialized above
345 if (ode_solver_type < 8) { ode_solver->Init(oper); }
346
347 // Since we want to update the diffusion coefficient after every time step,
348 // we need to use the "one-step" mode of the SUNDIALS solvers.
349 if (cvode) { cvode->SetStepMode(CV_ONE_STEP); }
350 if (arkode) { arkode->SetStepMode(ARK_ONE_STEP); }
351
352 // 10. Perform time-integration (looping over the time iterations, ti, with a
353 // time-step dt).
354 if (myid == 0)
355 {
356 cout << "Integrating the ODE ..." << endl;
357 }
358 tic_toc.Clear();
359 tic_toc.Start();
360
361 bool last_step = false;
362 for (int ti = 1; !last_step; ti++)
363 {
364 double dt_real = min(dt, t_final - t);
365
366 // Note that since we are using the "one-step" mode of the SUNDIALS
367 // solvers, they will, generally, step over the final time and will not
368 // explicitly perform the interpolation to t_final as they do in the
369 // "normal" step mode.
370
371 ode_solver->Step(u, t, dt_real);
372
373 last_step = (t >= t_final - 1e-8*dt);
374
375 if (last_step || (ti % vis_steps) == 0)
376 {
377 if (myid == 0)
378 {
379 cout << "step " << ti << ", t = " << t << endl;
380 if (cvode) { cvode->PrintInfo(); }
381 if (arkode) { arkode->PrintInfo(); }
382 }
383
384 u_gf.SetFromTrueDofs(u);
385 if (visualization)
386 {
387 sout << "parallel " << num_procs << " " << myid << "\n";
388 sout << "solution\n" << *pmesh << u_gf << flush;
389 }
390
391 if (visit)
392 {
393 visit_dc.SetCycle(ti);
394 visit_dc.SetTime(t);
395 visit_dc.Save();
396 }
397 }
398 oper.SetParameters(u);
399 }
400 tic_toc.Stop();
401 if (myid == 0)
402 {
403 cout << "Done, " << tic_toc.RealTime() << "s." << endl;
404 }
405
406 // 11. Save the final solution in parallel. This output can be viewed later
407 // using GLVis: "glvis -np <np> -m ex16-mesh -g ex16-final".
408 {
409 ostringstream sol_name;
410 sol_name << "ex16-final." << setfill('0') << setw(6) << myid;
411 ofstream osol(sol_name.str().c_str());
412 osol.precision(precision);
413 u_gf.Save(osol);
414 }
415
416 // 12. Free the used memory.
417 delete ode_solver;
418 delete pmesh;
419
420 return 0;
421}
422
423ConductionOperator::ConductionOperator(ParFiniteElementSpace &f, double al,
424 double kap, const Vector &u)
425 : TimeDependentOperator(f.GetTrueVSize(), 0.0), fespace(f), M(NULL), K(NULL),
426 T(NULL),
427 M_solver(f.GetComm()), T_solver(f.GetComm()), z(height)
428{
429 const double rel_tol = 1e-8;
430
431 M = new ParBilinearForm(&fespace);
432 M->AddDomainIntegrator(new MassIntegrator());
433 M->Assemble(0); // keep sparsity pattern of M and K the same
434 M->FormSystemMatrix(ess_tdof_list, Mmat);
435
436 M_solver.iterative_mode = false;
437 M_solver.SetRelTol(rel_tol);
438 M_solver.SetAbsTol(0.0);
439 M_solver.SetMaxIter(100);
440 M_solver.SetPrintLevel(0);
441 M_prec.SetType(HypreSmoother::Jacobi);
442 M_solver.SetPreconditioner(M_prec);
443 M_solver.SetOperator(Mmat);
444
445 alpha = al;
446 kappa = kap;
447
448 T_solver.iterative_mode = false;
449 T_solver.SetRelTol(rel_tol);
450 T_solver.SetAbsTol(0.0);
451 T_solver.SetMaxIter(100);
452 T_solver.SetPrintLevel(0);
453 T_solver.SetPreconditioner(T_prec);
454
455 SetParameters(u);
456}
457
458void ConductionOperator::Mult(const Vector &u, Vector &du_dt) const
459{
460 // Compute:
461 // du_dt = M^{-1}*-K(u)
462 // for du_dt
463 Kmat.Mult(u, z);
464 z.Neg(); // z = -z
465 M_solver.Mult(z, du_dt);
466}
467
468void ConductionOperator::ImplicitSolve(const double dt,
469 const Vector &u, Vector &du_dt)
470{
471 // Solve the equation:
472 // du_dt = M^{-1}*[-K(u + dt*du_dt)]
473 // for du_dt
474 if (T) { delete T; }
475 T = Add(1.0, Mmat, dt, Kmat);
476 T_solver.SetOperator(*T);
477 Kmat.Mult(u, z);
478 z.Neg();
479 T_solver.Mult(z, du_dt);
480}
481
482int ConductionOperator::SUNImplicitSetup(const Vector &x,
483 const Vector &fx, int jok, int *jcur,
484 double gamma)
485{
486 // Setup the ODE Jacobian T = M + gamma K.
487 if (T) { delete T; }
488 T = Add(1.0, Mmat, gamma, Kmat);
489 T_solver.SetOperator(*T);
490 *jcur = 1;
491 return (0);
492}
493
494int ConductionOperator::SUNImplicitSolve(const Vector &b, Vector &x, double tol)
495{
496 // Solve the system A x = z => (M - gamma K) x = M b.
497 Mmat.Mult(b, z);
498 T_solver.Mult(z, x);
499 return (0);
500}
501
502void ConductionOperator::SetParameters(const Vector &u)
503{
504 ParGridFunction u_alpha_gf(&fespace);
505 u_alpha_gf.SetFromTrueDofs(u);
506 for (int i = 0; i < u_alpha_gf.Size(); i++)
507 {
508 u_alpha_gf(i) = kappa + alpha*u_alpha_gf(i);
509 }
510
511 delete K;
512 K = new ParBilinearForm(&fespace);
513
514 GridFunctionCoefficient u_coeff(&u_alpha_gf);
515
517 K->Assemble(0); // keep sparsity pattern of M and K the same
518 K->FormSystemMatrix(ess_tdof_list, Kmat);
519}
520
521ConductionOperator::~ConductionOperator()
522{
523 delete T;
524 delete M;
525 delete K;
526}
527
529{
530 if (x.Norml2() < 0.5)
531 {
532 return 2.0;
533 }
534 else
535 {
536 return 1.0;
537 }
538}
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 SetStepMode(int itask)
Select the ARKode step mode: ARK_NORMAL (default) or ARK_ONE_STEP.
@ IMPLICIT
Implicit RK method.
Definition sundials.hpp:679
@ 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.
Backward Euler ODE solver. L-stable.
Definition ode.hpp:412
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.
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
Conjugate gradient method.
Definition solvers.hpp:513
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition solvers.hpp:526
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 SetStepMode(int itask)
Select the CVODE step mode: CV_NORMAL (default) or CV_ONE_STEP.
Definition sundials.cpp:870
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 SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetTime(real_t t)
Set physical time (for time-dependent simulations)
The classical forward Euler method.
Definition ode.hpp:118
A general function coefficient.
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
Parallel smoothers in hypre.
Definition hypre.hpp:1020
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.hpp:74
Mesh data type.
Definition mesh.hpp:56
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
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].
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.
Class for parallel bilinear form.
Abstract parallel finite element space.
Definition pfespace.hpp:29
HYPRE_BigInt GlobalTrueVSize() const
Definition pfespace.hpp:285
Class for parallel grid function.
Definition pgridfunc.hpp:33
void Save(std::ostream &out) const override
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
void SetFromTrueDofs(const Vector &tv) override
Set the GridFunction from the given true-dof vector.
void ProjectCoefficient(Coefficient &coeff) override
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Class for parallel meshes.
Definition pmesh.hpp:34
MPI_Comm GetComm() const
Definition pmesh.hpp:402
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Definition pmesh.cpp:4801
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
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
double RealTime()
Return the number of real seconds elapsed since the stopwatch was started.
Definition tic_toc.cpp:432
void Start()
Start the stopwatch. The elapsed time is not cleared.
Definition tic_toc.cpp:411
void Stop()
Stop the stopwatch.
Definition tic_toc.cpp:422
void Clear()
Clear the elapsed time on the stopwatch and restart it if it's running.
Definition tic_toc.cpp:406
static void Init()
Definition sundials.cpp:164
Base abstract class for first order time dependent operators.
Definition operator.hpp:317
Vector data type.
Definition vector.hpp:80
void Neg()
(*this) = -(*this)
Definition vector.cpp:300
real_t Norml2() const
Returns the l2 norm of the vector.
Definition vector.cpp:832
Data collection with VisIt I/O routines.
virtual void Save() override
Save the collection and a VisIt root file.
virtual void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
int close()
Close the socketstream.
const real_t alpha
Definition ex15.cpp:369
real_t kappa
Definition ex24.cpp:54
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
StopWatch tic_toc
Definition tic_toc.cpp:450
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
void Add(const DenseMatrix &A, const DenseMatrix &B, real_t alpha, DenseMatrix &C)
C = A + alpha*B.
const char vishost[]
RefCoord t[3]
double InitialTemperature(const Vector &x)
Definition ex16p.cpp:528
constexpr ARKODE_ERKTableID ARKODE_FEHLBERG_13_7_8
Definition sundials.hpp:65