MFEM  v4.4.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 and HYPRE.
105  Mpi::Init(argc, argv);
106  int num_procs = Mpi::WorldSize();
107  int myid = Mpi::WorldRank();
108  Hypre::Init();
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  return 1;
174  }
175 
176  if (myid == 0)
177  {
178  args.PrintOptions(cout);
179  }
180 
181  // check for valid ODE solver option
182  if (ode_solver_type < 1 || ode_solver_type > 12)
183  {
184  if (myid == 0)
185  {
186  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
187  }
188  return 1;
189  }
190 
191  // 3. Read the serial mesh from the given mesh file on all processors. We can
192  // handle triangular, quadrilateral, tetrahedral and hexahedral meshes
193  // with the same code.
194  Mesh *mesh = new Mesh(mesh_file, 1, 1);
195  int dim = mesh->Dimension();
196 
197  // 4. Refine the mesh in serial to increase the resolution. In this example
198  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
199  // a command-line parameter.
200  for (int lev = 0; lev < ser_ref_levels; lev++)
201  {
202  mesh->UniformRefinement();
203  }
204 
205  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
206  // this mesh further in parallel to increase the resolution. Once the
207  // parallel mesh is defined, the serial mesh can be deleted.
208  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
209  delete mesh;
210  for (int lev = 0; lev < par_ref_levels; lev++)
211  {
212  pmesh->UniformRefinement();
213  }
214 
215  // 6. Define the vector finite element space representing the current and the
216  // initial temperature, u_ref.
217  H1_FECollection fe_coll(order, dim);
218  ParFiniteElementSpace fespace(pmesh, &fe_coll);
219 
220  int fe_size = fespace.GlobalTrueVSize();
221  if (myid == 0)
222  {
223  cout << "Number of temperature unknowns: " << fe_size << endl;
224  }
225 
226  ParGridFunction u_gf(&fespace);
227 
228  // 7. Set the initial conditions for u. All boundaries are considered
229  // natural.
231  u_gf.ProjectCoefficient(u_0);
232  Vector u;
233  u_gf.GetTrueDofs(u);
234 
235  // 8. Initialize the conduction operator and the VisIt visualization.
236  ConductionOperator oper(fespace, alpha, kappa, u);
237 
238  u_gf.SetFromTrueDofs(u);
239  {
240  ostringstream mesh_name, sol_name;
241  mesh_name << "ex16-mesh." << setfill('0') << setw(6) << myid;
242  sol_name << "ex16-init." << setfill('0') << setw(6) << myid;
243  ofstream omesh(mesh_name.str().c_str());
244  omesh.precision(precision);
245  pmesh->Print(omesh);
246  ofstream osol(sol_name.str().c_str());
247  osol.precision(precision);
248  u_gf.Save(osol);
249  }
250 
251  VisItDataCollection visit_dc("Example16-Parallel", pmesh);
252  visit_dc.RegisterField("temperature", &u_gf);
253  if (visit)
254  {
255  visit_dc.SetCycle(0);
256  visit_dc.SetTime(0.0);
257  visit_dc.Save();
258  }
259 
260  socketstream sout;
261  if (visualization)
262  {
263  char vishost[] = "localhost";
264  int visport = 19916;
265  sout.open(vishost, visport);
266  sout << "parallel " << num_procs << " " << myid << endl;
267  int good = sout.good(), all_good;
268  MPI_Allreduce(&good, &all_good, 1, MPI_INT, MPI_MIN, pmesh->GetComm());
269  if (!all_good)
270  {
271  sout.close();
272  visualization = false;
273  if (myid == 0)
274  {
275  cout << "Unable to connect to GLVis server at "
276  << vishost << ':' << visport << endl;
277  cout << "GLVis visualization disabled.\n";
278  }
279  }
280  else
281  {
282  sout.precision(precision);
283  sout << "solution\n" << *pmesh << u_gf;
284  sout << "pause\n";
285  sout << flush;
286  if (myid == 0)
287  {
288  cout << "GLVis visualization paused."
289  << " Press space (in the GLVis window) to resume it.\n";
290  }
291  }
292  }
293 
294  // 9. Define the ODE solver used for time integration.
295  double t = 0.0;
296  ODESolver *ode_solver = NULL;
297  CVODESolver *cvode = NULL;
298  ARKStepSolver *arkode = NULL;
299  switch (ode_solver_type)
300  {
301  // MFEM explicit methods
302  case 1: ode_solver = new ForwardEulerSolver; break;
303  case 2: ode_solver = new RK2Solver(0.5); break; // midpoint method
304  case 3: ode_solver = new RK3SSPSolver; break;
305  case 4: ode_solver = new RK4Solver; break;
306  // MFEM implicit L-stable methods
307  case 5: ode_solver = new BackwardEulerSolver; break;
308  case 6: ode_solver = new SDIRK23Solver(2); break;
309  case 7: ode_solver = new SDIRK33Solver; break;
310  // CVODE
311  case 8:
312  cvode = new CVODESolver(MPI_COMM_WORLD, CV_ADAMS);
313  cvode->Init(oper);
314  cvode->SetSStolerances(reltol, abstol);
315  cvode->SetMaxStep(dt);
316  ode_solver = cvode; break;
317  case 9:
318  cvode = new CVODESolver(MPI_COMM_WORLD, CV_BDF);
319  cvode->Init(oper);
320  cvode->SetSStolerances(reltol, abstol);
321  cvode->SetMaxStep(dt);
322  ode_solver = cvode; break;
323  // ARKODE
324  case 10:
325  case 11:
326  arkode = new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::EXPLICIT);
327  arkode->Init(oper);
328  arkode->SetSStolerances(reltol, abstol);
329  arkode->SetMaxStep(dt);
330  if (ode_solver_type == 11) { arkode->SetERKTableNum(FEHLBERG_13_7_8); }
331  ode_solver = arkode; break;
332  case 12:
333  arkode = new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::IMPLICIT);
334  arkode->Init(oper);
335  arkode->SetSStolerances(reltol, abstol);
336  arkode->SetMaxStep(dt);
337  ode_solver = arkode; break;
338  }
339 
340  // Initialize MFEM integrators, SUNDIALS integrators are initialized above
341  if (ode_solver_type < 8) { ode_solver->Init(oper); }
342 
343  // Since we want to update the diffusion coefficient after every time step,
344  // we need to use the "one-step" mode of the SUNDIALS solvers.
345  if (cvode) { cvode->SetStepMode(CV_ONE_STEP); }
346  if (arkode) { arkode->SetStepMode(ARK_ONE_STEP); }
347 
348  // 10. Perform time-integration (looping over the time iterations, ti, with a
349  // time-step dt).
350  if (myid == 0)
351  {
352  cout << "Integrating the ODE ..." << endl;
353  }
354  tic_toc.Clear();
355  tic_toc.Start();
356 
357  bool last_step = false;
358  for (int ti = 1; !last_step; ti++)
359  {
360  double dt_real = min(dt, t_final - t);
361 
362  // Note that since we are using the "one-step" mode of the SUNDIALS
363  // solvers, they will, generally, step over the final time and will not
364  // explicitly perform the interpolation to t_final as they do in the
365  // "normal" step mode.
366 
367  ode_solver->Step(u, t, dt_real);
368 
369  last_step = (t >= t_final - 1e-8*dt);
370 
371  if (last_step || (ti % vis_steps) == 0)
372  {
373  if (myid == 0)
374  {
375  cout << "step " << ti << ", t = " << t << endl;
376  if (cvode) { cvode->PrintInfo(); }
377  if (arkode) { arkode->PrintInfo(); }
378  }
379 
380  u_gf.SetFromTrueDofs(u);
381  if (visualization)
382  {
383  sout << "parallel " << num_procs << " " << myid << "\n";
384  sout << "solution\n" << *pmesh << u_gf << flush;
385  }
386 
387  if (visit)
388  {
389  visit_dc.SetCycle(ti);
390  visit_dc.SetTime(t);
391  visit_dc.Save();
392  }
393  }
394  oper.SetParameters(u);
395  }
396  tic_toc.Stop();
397  if (myid == 0)
398  {
399  cout << "Done, " << tic_toc.RealTime() << "s." << endl;
400  }
401 
402  // 11. Save the final solution in parallel. This output can be viewed later
403  // using GLVis: "glvis -np <np> -m ex16-mesh -g ex16-final".
404  {
405  ostringstream sol_name;
406  sol_name << "ex16-final." << setfill('0') << setw(6) << myid;
407  ofstream osol(sol_name.str().c_str());
408  osol.precision(precision);
409  u_gf.Save(osol);
410  }
411 
412  // 12. Free the used memory.
413  delete ode_solver;
414  delete pmesh;
415 
416  return 0;
417 }
418 
419 ConductionOperator::ConductionOperator(ParFiniteElementSpace &f, double al,
420  double kap, const Vector &u)
421  : TimeDependentOperator(f.GetTrueVSize(), 0.0), fespace(f), M(NULL), K(NULL),
422  T(NULL),
423  M_solver(f.GetComm()), T_solver(f.GetComm()), z(height)
424 {
425  const double rel_tol = 1e-8;
426 
427  M = new ParBilinearForm(&fespace);
428  M->AddDomainIntegrator(new MassIntegrator());
429  M->Assemble(0); // keep sparsity pattern of M and K the same
430  M->FormSystemMatrix(ess_tdof_list, Mmat);
431 
432  M_solver.iterative_mode = false;
433  M_solver.SetRelTol(rel_tol);
434  M_solver.SetAbsTol(0.0);
435  M_solver.SetMaxIter(100);
436  M_solver.SetPrintLevel(0);
437  M_prec.SetType(HypreSmoother::Jacobi);
438  M_solver.SetPreconditioner(M_prec);
439  M_solver.SetOperator(Mmat);
440 
441  alpha = al;
442  kappa = kap;
443 
444  T_solver.iterative_mode = false;
445  T_solver.SetRelTol(rel_tol);
446  T_solver.SetAbsTol(0.0);
447  T_solver.SetMaxIter(100);
448  T_solver.SetPrintLevel(0);
449  T_solver.SetPreconditioner(T_prec);
450 
451  SetParameters(u);
452 }
453 
454 void ConductionOperator::Mult(const Vector &u, Vector &du_dt) const
455 {
456  // Compute:
457  // du_dt = M^{-1}*-K(u)
458  // for du_dt
459  Kmat.Mult(u, z);
460  z.Neg(); // z = -z
461  M_solver.Mult(z, du_dt);
462 }
463 
464 void ConductionOperator::ImplicitSolve(const double dt,
465  const Vector &u, Vector &du_dt)
466 {
467  // Solve the equation:
468  // du_dt = M^{-1}*[-K(u + dt*du_dt)]
469  // for du_dt
470  if (T) { delete T; }
471  T = Add(1.0, Mmat, dt, Kmat);
472  T_solver.SetOperator(*T);
473  Kmat.Mult(u, z);
474  z.Neg();
475  T_solver.Mult(z, du_dt);
476 }
477 
478 int ConductionOperator::SUNImplicitSetup(const Vector &x,
479  const Vector &fx, int jok, int *jcur,
480  double gamma)
481 {
482  // Setup the ODE Jacobian T = M + gamma K.
483  if (T) { delete T; }
484  T = Add(1.0, Mmat, gamma, Kmat);
485  T_solver.SetOperator(*T);
486  *jcur = 1;
487  return (0);
488 }
489 
490 int ConductionOperator::SUNImplicitSolve(const Vector &b, Vector &x, double tol)
491 {
492  // Solve the system A x = z => (M - gamma K) x = M b.
493  Mmat.Mult(b, z);
494  T_solver.Mult(z, x);
495  return (0);
496 }
497 
498 void ConductionOperator::SetParameters(const Vector &u)
499 {
500  ParGridFunction u_alpha_gf(&fespace);
501  u_alpha_gf.SetFromTrueDofs(u);
502  for (int i = 0; i < u_alpha_gf.Size(); i++)
503  {
504  u_alpha_gf(i) = kappa + alpha*u_alpha_gf(i);
505  }
506 
507  delete K;
508  K = new ParBilinearForm(&fespace);
509 
510  GridFunctionCoefficient u_coeff(&u_alpha_gf);
511 
512  K->AddDomainIntegrator(new DiffusionIntegrator(u_coeff));
513  K->Assemble(0); // keep sparsity pattern of M and K the same
514  K->FormSystemMatrix(ess_tdof_list, Kmat);
515 }
516 
517 ConductionOperator::~ConductionOperator()
518 {
519  delete T;
520  delete M;
521  delete K;
522 }
523 
524 double InitialTemperature(const Vector &x)
525 {
526  if (x.Norml2() < 0.5)
527  {
528  return 2.0;
529  }
530  else
531  {
532  return 1.0;
533  }
534 }
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:465
double InitialTemperature(const Vector &x)
Definition: ex16.cpp:385
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:285
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:472
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:792
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]...
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:285
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:873
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)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:525
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:391
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:1916
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: pgridfunc.hpp:169
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
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:9737
constexpr int visport
Data collection with VisIt I/O routines.
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition: sundials.cpp:716
Parallel smoothers in hypre.
Definition: hypre.hpp:896
int Dimension() const
Definition: mesh.hpp:999
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
void SetTime(double t)
Set physical time (for time-dependent simulations)
void Start()
Start the stopwatch. The elapsed time is not cleared.
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:288
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:324
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;.
RefCoord t[3]
const double alpha
Definition: ex15.cpp:369
A general function coefficient.
Vector data type.
Definition: vector.hpp:60
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
void PrintInfo() const
Print various CVODE statistics.
Definition: sundials.cpp:742
void Print(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4684
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
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:337
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
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:18
int main()
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:150