MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
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 
38 using namespace std;
39 using 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  */
51 class ConductionOperator : public TimeDependentOperator
52 {
53 protected:
54  ParFiniteElementSpace &fespace;
55  Array<int> ess_tdof_list; // this list remains empty for pure Neumann b.c.
56 
57  ParBilinearForm *M;
58  ParBilinearForm *K;
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 
75 public:
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 
100 double InitialTemperature(const Vector &x);
101 
102 int main(int argc, char *argv[])
103 {
104  // 1. Initialize MPI.
105  int num_procs, myid;
106  MPI_Init(&argc, &argv);
107  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
108  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
109 
110  // 2. Parse command-line options.
111  const char *mesh_file = "../../data/star.mesh";
112  int ser_ref_levels = 2;
113  int par_ref_levels = 1;
114  int order = 2;
115  int ode_solver_type = 9; // CVODE implicit BDF
116  double t_final = 0.5;
117  double dt = 1.0e-2;
118  double alpha = 1.0e-2;
119  double kappa = 0.5;
120  bool visualization = true;
121  bool visit = false;
122  int vis_steps = 5;
123 
124  // Relative and absolute tolerances for CVODE and ARKODE.
125  const double reltol = 1e-4, abstol = 1e-4;
126 
127  int precision = 8;
128  cout.precision(precision);
129 
130  OptionsParser args(argc, argv);
131  args.AddOption(&mesh_file, "-m", "--mesh",
132  "Mesh file to use.");
133  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
134  "Number of times to refine the mesh uniformly in serial.");
135  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
136  "Number of times to refine the mesh uniformly in parallel.");
137  args.AddOption(&order, "-o", "--order",
138  "Order (degree) of the finite elements.");
139  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
140  "ODE solver:\n\t"
141  "1 - Forward Euler,\n\t"
142  "2 - RK2,\n\t"
143  "3 - RK3 SSP,\n\t"
144  "4 - RK4,\n\t"
145  "5 - Backward Euler,\n\t"
146  "6 - SDIRK 2,\n\t"
147  "7 - SDIRK 3,\n\t"
148  "8 - CVODE (implicit Adams),\n\t"
149  "9 - CVODE (implicit BDF),\n\t"
150  "10 - ARKODE (default explicit),\n\t"
151  "11 - ARKODE (explicit Fehlberg-6-4-5),\n\t"
152  "12 - ARKODE (default impicit).");
153  args.AddOption(&t_final, "-tf", "--t-final",
154  "Final time; start time is 0.");
155  args.AddOption(&dt, "-dt", "--time-step",
156  "Time step.");
157  args.AddOption(&alpha, "-a", "--alpha",
158  "Alpha coefficient.");
159  args.AddOption(&kappa, "-k", "--kappa",
160  "Kappa coefficient offset.");
161  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
162  "--no-visualization",
163  "Enable or disable GLVis visualization.");
164  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
165  "--no-visit-datafiles",
166  "Save data files for VisIt (visit.llnl.gov) visualization.");
167  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
168  "Visualize every n-th timestep.");
169  args.Parse();
170  if (!args.Good())
171  {
172  args.PrintUsage(cout);
173  MPI_Finalize();
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  MPI_Finalize();
190  return 1;
191  }
192 
193  // 3. Read the serial mesh from the given mesh file on all processors. We can
194  // handle triangular, quadrilateral, tetrahedral and hexahedral meshes
195  // with the same code.
196  Mesh *mesh = new Mesh(mesh_file, 1, 1);
197  int dim = mesh->Dimension();
198 
199  // 4. Refine the mesh in serial to increase the resolution. In this example
200  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
201  // a command-line parameter.
202  for (int lev = 0; lev < ser_ref_levels; lev++)
203  {
204  mesh->UniformRefinement();
205  }
206 
207  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
208  // this mesh further in parallel to increase the resolution. Once the
209  // parallel mesh is defined, the serial mesh can be deleted.
210  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
211  delete mesh;
212  for (int lev = 0; lev < par_ref_levels; lev++)
213  {
214  pmesh->UniformRefinement();
215  }
216 
217  // 6. Define the vector finite element space representing the current and the
218  // initial temperature, u_ref.
219  H1_FECollection fe_coll(order, dim);
220  ParFiniteElementSpace fespace(pmesh, &fe_coll);
221 
222  int fe_size = fespace.GlobalTrueVSize();
223  if (myid == 0)
224  {
225  cout << "Number of temperature unknowns: " << fe_size << endl;
226  }
227 
228  ParGridFunction u_gf(&fespace);
229 
230  // 7. Set the initial conditions for u. All boundaries are considered
231  // natural.
233  u_gf.ProjectCoefficient(u_0);
234  Vector u;
235  u_gf.GetTrueDofs(u);
236 
237  // 8. Initialize the conduction operator and the VisIt visualization.
238  ConductionOperator oper(fespace, alpha, kappa, u);
239 
240  u_gf.SetFromTrueDofs(u);
241  {
242  ostringstream mesh_name, sol_name;
243  mesh_name << "ex16-mesh." << setfill('0') << setw(6) << myid;
244  sol_name << "ex16-init." << setfill('0') << setw(6) << myid;
245  ofstream omesh(mesh_name.str().c_str());
246  omesh.precision(precision);
247  pmesh->Print(omesh);
248  ofstream osol(sol_name.str().c_str());
249  osol.precision(precision);
250  u_gf.Save(osol);
251  }
252 
253  VisItDataCollection visit_dc("Example16-Parallel", pmesh);
254  visit_dc.RegisterField("temperature", &u_gf);
255  if (visit)
256  {
257  visit_dc.SetCycle(0);
258  visit_dc.SetTime(0.0);
259  visit_dc.Save();
260  }
261 
262  socketstream sout;
263  if (visualization)
264  {
265  char vishost[] = "localhost";
266  int visport = 19916;
267  sout.open(vishost, visport);
268  sout << "parallel " << num_procs << " " << myid << endl;
269  int good = sout.good(), all_good;
270  MPI_Allreduce(&good, &all_good, 1, MPI_INT, MPI_MIN, pmesh->GetComm());
271  if (!all_good)
272  {
273  sout.close();
274  visualization = false;
275  if (myid == 0)
276  {
277  cout << "Unable to connect to GLVis server at "
278  << vishost << ':' << visport << endl;
279  cout << "GLVis visualization disabled.\n";
280  }
281  }
282  else
283  {
284  sout.precision(precision);
285  sout << "solution\n" << *pmesh << u_gf;
286  sout << "pause\n";
287  sout << flush;
288  if (myid == 0)
289  {
290  cout << "GLVis visualization paused."
291  << " Press space (in the GLVis window) to resume it.\n";
292  }
293  }
294  }
295 
296  // 9. Define the ODE solver used for time integration.
297  double t = 0.0;
298  ODESolver *ode_solver = NULL;
299  CVODESolver *cvode = NULL;
300  ARKStepSolver *arkode = NULL;
301  switch (ode_solver_type)
302  {
303  // MFEM explicit methods
304  case 1: ode_solver = new ForwardEulerSolver; break;
305  case 2: ode_solver = new RK2Solver(0.5); break; // midpoint method
306  case 3: ode_solver = new RK3SSPSolver; break;
307  case 4: ode_solver = new RK4Solver; break;
308  // MFEM implicit L-stable methods
309  case 5: ode_solver = new BackwardEulerSolver; break;
310  case 6: ode_solver = new SDIRK23Solver(2); break;
311  case 7: ode_solver = new SDIRK33Solver; break;
312  // CVODE
313  case 8:
314  cvode = new CVODESolver(MPI_COMM_WORLD, CV_ADAMS);
315  cvode->Init(oper);
316  cvode->SetSStolerances(reltol, abstol);
317  cvode->SetMaxStep(dt);
318  ode_solver = cvode; break;
319  case 9:
320  cvode = new CVODESolver(MPI_COMM_WORLD, CV_BDF);
321  cvode->Init(oper);
322  cvode->SetSStolerances(reltol, abstol);
323  cvode->SetMaxStep(dt);
324  ode_solver = cvode; break;
325  // ARKODE
326  case 10:
327  case 11:
328  arkode = new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::EXPLICIT);
329  arkode->Init(oper);
330  arkode->SetSStolerances(reltol, abstol);
331  arkode->SetMaxStep(dt);
332  if (ode_solver_type == 11) { arkode->SetERKTableNum(FEHLBERG_13_7_8); }
333  ode_solver = arkode; break;
334  case 12:
335  arkode = new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::IMPLICIT);
336  arkode->Init(oper);
337  arkode->SetSStolerances(reltol, abstol);
338  arkode->SetMaxStep(dt);
339  ode_solver = arkode; break;
340  }
341 
342  // Initialize MFEM integrators, SUNDIALS integrators are initialized above
343  if (ode_solver_type < 8) { ode_solver->Init(oper); }
344 
345  // Since we want to update the diffusion coefficient after every time step,
346  // we need to use the "one-step" mode of the SUNDIALS solvers.
347  if (cvode) { cvode->SetStepMode(CV_ONE_STEP); }
348  if (arkode) { arkode->SetStepMode(ARK_ONE_STEP); }
349 
350  // 10. Perform time-integration (looping over the time iterations, ti, with a
351  // time-step dt).
352  if (myid == 0)
353  {
354  cout << "Integrating the ODE ..." << endl;
355  }
356  tic_toc.Clear();
357  tic_toc.Start();
358 
359  bool last_step = false;
360  for (int ti = 1; !last_step; ti++)
361  {
362  double dt_real = min(dt, t_final - t);
363 
364  // Note that since we are using the "one-step" mode of the SUNDIALS
365  // solvers, they will, generally, step over the final time and will not
366  // explicitly perform the interpolation to t_final as they do in the
367  // "normal" step mode.
368 
369  ode_solver->Step(u, t, dt_real);
370 
371  last_step = (t >= t_final - 1e-8*dt);
372 
373  if (last_step || (ti % vis_steps) == 0)
374  {
375  if (myid == 0)
376  {
377  cout << "step " << ti << ", t = " << t << endl;
378  if (cvode) { cvode->PrintInfo(); }
379  if (arkode) { arkode->PrintInfo(); }
380  }
381 
382  u_gf.SetFromTrueDofs(u);
383  if (visualization)
384  {
385  sout << "parallel " << num_procs << " " << myid << "\n";
386  sout << "solution\n" << *pmesh << u_gf << flush;
387  }
388 
389  if (visit)
390  {
391  visit_dc.SetCycle(ti);
392  visit_dc.SetTime(t);
393  visit_dc.Save();
394  }
395  }
396  oper.SetParameters(u);
397  }
398  tic_toc.Stop();
399  if (myid == 0)
400  {
401  cout << "Done, " << tic_toc.RealTime() << "s." << endl;
402  }
403 
404  // 11. Save the final solution in parallel. This output can be viewed later
405  // using GLVis: "glvis -np <np> -m ex16-mesh -g ex16-final".
406  {
407  ostringstream sol_name;
408  sol_name << "ex16-final." << setfill('0') << setw(6) << myid;
409  ofstream osol(sol_name.str().c_str());
410  osol.precision(precision);
411  u_gf.Save(osol);
412  }
413 
414  // 12. Free the used memory.
415  delete ode_solver;
416  delete pmesh;
417 
418  MPI_Finalize();
419 
420  return 0;
421 }
422 
423 ConductionOperator::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 
458 void 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 
468 void 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 
482 int 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 
494 int 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 
502 void 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 
516  K->AddDomainIntegrator(new DiffusionIntegrator(u_coeff));
517  K->Assemble(0); // keep sparsity pattern of M and K the same
518  K->FormSystemMatrix(ess_tdof_list, Kmat);
519 }
520 
521 ConductionOperator::~ConductionOperator()
522 {
523  delete T;
524  delete M;
525  delete K;
526 }
527 
528 double InitialTemperature(const Vector &x)
529 {
530  if (x.Norml2() < 0.5)
531  {
532  return 2.0;
533  }
534  else
535  {
536  return 1.0;
537  }
538 }
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
Definition: sundials.cpp:532
Conjugate gradient method.
Definition: solvers.hpp:258
double InitialTemperature(const Vector &x)
Definition: ex16.cpp:382
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
Definition: sundials.cpp:698
void SetStepMode(int itask)
Select the ARKode step mode: ARK_NORMAL (default) or ARK_ONE_STEP.
Definition: sundials.cpp:1529
Base abstract class for first order time dependent operators.
Definition: operator.hpp:267
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:744
virtual void Step(Vector &x, double &t, double &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:18
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
virtual void Save()
Save the collection and a VisIt root file.
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:819
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition: sundials.cpp:1540
StopWatch tic_toc
Definition: tic_toc.cpp:447
double RealTime()
Definition: tic_toc.cpp:426
double kappa
Definition: ex24.cpp:54
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:474
int main(int argc, char *argv[])
Definition: ex1.cpp:66
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:389
int close()
Close the socketstream.
void Stop()
Stop the stopwatch.
Definition: tic_toc.cpp:416
Interface to ARKode&#39;s ARKStep module – additive Runge-Kutta methods.
Definition: sundials.hpp:539
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:1928
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: pgridfunc.hpp:164
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
constexpr char vishost[]
Interface to the CVODE library – linear multi-step methods.
Definition: sundials.hpp:252
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:8382
constexpr int visport
Data collection with VisIt I/O routines.
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:251
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition: sundials.cpp:716
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4169
Parallel smoothers in hypre.
Definition: hypre.hpp:592
int Dimension() const
Definition: mesh.hpp:788
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:434
void SetTime(double t)
Set physical time (for time-dependent simulations)
void Start()
Clear the elapsed time and start the stopwatch.
Definition: tic_toc.cpp:411
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:162
MPI_Comm GetComm() const
Definition: pmesh.hpp:236
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:149
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
void PrintInfo() const
Print various ARKStep statistics.
Definition: sundials.cpp:1576
int dim
Definition: ex24.cpp:53
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:304
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:143
Class for parallel bilinear form.
int open(const char hostname[], int port)
Open the socket stream on &#39;port&#39; at &#39;hostname&#39;.
const double alpha
Definition: ex15.cpp:336
A general function coefficient.
Vector data type.
Definition: vector.hpp:51
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:159
void PrintInfo() const
Print various CVODE statistics.
Definition: sundials.cpp:742
Class for parallel grid function.
Definition: pgridfunc.hpp:32
void SetERKTableNum(int table_num)
Choose a specific Butcher table for an explicit RK method.
Definition: sundials.cpp:1552
The classical forward Euler method.
Definition: ode.hpp:116
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults...
Definition: sundials.cpp:1297
Class for parallel meshes.
Definition: pmesh.hpp:32
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
Definition: sundials.cpp:1534
void Clear()
Clear the elapsed time on the stopwatch and restart it if it&#39;s running.
Definition: tic_toc.cpp:406
double f(const Vector &p)
void SetStepMode(int itask)
Select the CVODE step mode: CV_NORMAL (default) or CV_ONE_STEP.
Definition: sundials.cpp:693
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:145