MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex9p.cpp
Go to the documentation of this file.
1 // MFEM Example 9 - Parallel Version
2 // SUNDIALS Modification
3 //
4 // Compile with: make ex9p
5 //
6 // Sample runs:
7 // mpirun -np 4 ex9p -m ../../data/periodic-segment.mesh -p 1 -rp 1 -s 7 -dt 0.0025
8 // mpirun -np 4 ex9p -m ../../data/periodic-square.mesh -p 1 -rp 1 -s 8 -dt 0.0025 -tf 9
9 // mpirun -np 4 ex9p -m ../../data/periodic-hexagon.mesh -p 0 -rp 1 -s 7 -dt 0.0009 -vs 25
10 // mpirun -np 4 ex9p -m ../../data/periodic-hexagon.mesh -p 0 -rp 1 -s 9 -dt 0.005 -vs 15
11 // mpirun -np 4 ex9p -m ../../data/amr-quad.mesh -p 1 -rp 1 -s 9 -dt 0.001 -tf 9
12 // mpirun -np 4 ex9p -m ../../data/star-q3.mesh -p 1 -rp 1 -s 9 -dt 0.0025 -tf 9
13 // mpirun -np 4 ex9p -m ../../data/disc-nurbs.mesh -p 1 -rp 2 -s 7 -dt 0.0025 -tf 9
14 // mpirun -np 4 ex9p -m ../../data/periodic-cube.mesh -p 0 -rp 1 -s 8 -dt 0.01 -tf 8 -o 2
15 //
16 // Description: This example code solves the time-dependent advection equation
17 // du/dt + v.grad(u) = 0, where v is a given fluid velocity, and
18 // u0(x)=u(0,x) is a given initial condition.
19 //
20 // The example demonstrates the use of Discontinuous Galerkin (DG)
21 // bilinear forms in MFEM (face integrators), the use of explicit
22 // ODE time integrators, the definition of periodic boundary
23 // conditions through periodic meshes, as well as the use of GLVis
24 // for persistent visualization of a time-evolving solution. The
25 // saving of time-dependent data files for external visualization
26 // with VisIt (visit.llnl.gov) is also illustrated.
27 
28 #include "mfem.hpp"
29 #include <fstream>
30 #include <iostream>
31 
32 #ifndef MFEM_USE_SUNDIALS
33 #error This example requires that MFEM is built with MFEM_USE_SUNDIALS=YES
34 #endif
35 
36 using namespace std;
37 using namespace mfem;
38 
39 // Choice for the problem setup. The fluid velocity, initial condition and
40 // inflow boundary condition are chosen based on this parameter.
41 int problem;
42 
43 // Velocity coefficient
44 void velocity_function(const Vector &x, Vector &v);
45 
46 // Initial condition
47 double u0_function(const Vector &x);
48 
49 // Inflow boundary condition
50 double inflow_function(const Vector &x);
51 
52 // Mesh bounding box
54 
55 
56 /** A time-dependent operator for the right-hand side of the ODE. The DG weak
57  form of du/dt = -v.grad(u) is M du/dt = K u + b, where M and K are the mass
58  and advection matrices, and b describes the flow on the boundary. This can
59  be written as a general ODE, du/dt = M^{-1} (K u + b), and this class is
60  used to evaluate the right-hand side. */
62 {
63 private:
64  HypreParMatrix &M, &K;
65  const Vector &b;
66  HypreSmoother M_prec;
67  CGSolver M_solver;
68 
69  mutable Vector z;
70 
71 public:
72  FE_Evolution(HypreParMatrix &_M, HypreParMatrix &_K, const Vector &_b);
73 
74  virtual void Mult(const Vector &x, Vector &y) const;
75 
76  virtual ~FE_Evolution() { }
77 };
78 
79 
80 int main(int argc, char *argv[])
81 {
82  // 1. Initialize MPI.
83  int num_procs, myid;
84  MPI_Init(&argc, &argv);
85  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
86  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
87 
88  // 2. Parse command-line options.
89  problem = 0;
90  const char *mesh_file = "../../data/periodic-hexagon.mesh";
91  int ser_ref_levels = 2;
92  int par_ref_levels = 0;
93  int order = 3;
94  int ode_solver_type = 4;
95  double t_final = 10.0;
96  double dt = 0.01;
97  bool visualization = true;
98  bool visit = false;
99  bool binary = false;
100  int vis_steps = 5;
101 
102  // Relative and absolute tolerances for CVODE and ARKODE.
103  const double reltol = 1e-2, abstol = 1e-2;
104 
105  int precision = 8;
106  cout.precision(precision);
107 
108  OptionsParser args(argc, argv);
109  args.AddOption(&mesh_file, "-m", "--mesh",
110  "Mesh file to use.");
111  args.AddOption(&problem, "-p", "--problem",
112  "Problem setup to use. See options in velocity_function().");
113  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
114  "Number of times to refine the mesh uniformly in serial.");
115  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
116  "Number of times to refine the mesh uniformly in parallel.");
117  args.AddOption(&order, "-o", "--order",
118  "Order (degree) of the finite elements.");
119  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
120  "ODE solver:\n\t"
121  "1 - Forward Euler,\n\t"
122  "2 - RK2 SSP,\n\t"
123  "3 - RK3 SSP,\n\t"
124  "4 - RK4,\n\t"
125  "6 - RK6,\n\t"
126  "7 - CVODE (adaptive order implicit Adams),\n\t"
127  "8 - ARKODE default (4th order) explicit,\n\t"
128  "9 - ARKODE RK8.");
129  args.AddOption(&t_final, "-tf", "--t-final",
130  "Final time; start time is 0.");
131  args.AddOption(&dt, "-dt", "--time-step",
132  "Time step.");
133  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
134  "--no-visualization",
135  "Enable or disable GLVis visualization.");
136  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
137  "--no-visit-datafiles",
138  "Save data files for VisIt (visit.llnl.gov) visualization.");
139  args.AddOption(&binary, "-binary", "--binary-datafiles", "-ascii",
140  "--ascii-datafiles",
141  "Use binary (Sidre) or ascii format for VisIt data files.");
142  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
143  "Visualize every n-th timestep.");
144  args.Parse();
145  if (!args.Good())
146  {
147  if (myid == 0)
148  {
149  args.PrintUsage(cout);
150  }
151  MPI_Finalize();
152  return 1;
153  }
154  if (myid == 0)
155  {
156  args.PrintOptions(cout);
157  }
158  // check for vaild ODE solver option
159  if (ode_solver_type < 1 || ode_solver_type > 9)
160  {
161  if (myid == 0)
162  {
163  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
164  }
165  MPI_Finalize();
166  return 3;
167  }
168 
169  // 3. Read the serial mesh from the given mesh file on all processors. We can
170  // handle geometrically periodic meshes in this code.
171  Mesh *mesh = new Mesh(mesh_file, 1, 1);
172  int dim = mesh->Dimension();
173 
174  // 4. Refine the mesh in serial to increase the resolution. In this example
175  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
176  // a command-line parameter. If the mesh is of NURBS type, we convert it
177  // to a (piecewise-polynomial) high-order mesh.
178  for (int lev = 0; lev < ser_ref_levels; lev++)
179  {
180  mesh->UniformRefinement();
181  }
182  if (mesh->NURBSext)
183  {
184  mesh->SetCurvature(max(order, 1));
185  }
186  mesh->GetBoundingBox(bb_min, bb_max, max(order, 1));
187 
188  // 5. Define the parallel mesh by a partitioning of the serial mesh. Refine
189  // this mesh further in parallel to increase the resolution. Once the
190  // parallel mesh is defined, the serial mesh can be deleted.
191  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
192  delete mesh;
193  for (int lev = 0; lev < par_ref_levels; lev++)
194  {
195  pmesh->UniformRefinement();
196  }
197 
198  // 6. Define the parallel discontinuous DG finite element space on the
199  // parallel refined mesh of the given polynomial order.
200  DG_FECollection fec(order, dim);
201  ParFiniteElementSpace *fes = new ParFiniteElementSpace(pmesh, &fec);
202 
203  HYPRE_Int global_vSize = fes->GlobalTrueVSize();
204  if (myid == 0)
205  {
206  cout << "Number of unknowns: " << global_vSize << endl;
207  }
208 
209  // 7. Set up and assemble the parallel bilinear and linear forms (and the
210  // parallel hypre matrices) corresponding to the DG discretization. The
211  // DGTraceIntegrator involves integrals over mesh interior faces.
215 
216  ParBilinearForm *m = new ParBilinearForm(fes);
218  ParBilinearForm *k = new ParBilinearForm(fes);
219  k->AddDomainIntegrator(new ConvectionIntegrator(velocity, -1.0));
221  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
223  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
224 
225  ParLinearForm *b = new ParLinearForm(fes);
227  new BoundaryFlowIntegrator(inflow, velocity, -1.0, -0.5));
228 
229  m->Assemble();
230  m->Finalize();
231  int skip_zeros = 0;
232  k->Assemble(skip_zeros);
233  k->Finalize(skip_zeros);
234  b->Assemble();
235 
239 
240  // 8. Define the initial conditions, save the corresponding grid function to
241  // a file and (optionally) save data in the VisIt format and initialize
242  // GLVis visualization.
243  ParGridFunction *u = new ParGridFunction(fes);
244  u->ProjectCoefficient(u0);
245  HypreParVector *U = u->GetTrueDofs();
246 
247  {
248  ostringstream mesh_name, sol_name;
249  mesh_name << "ex9-mesh." << setfill('0') << setw(6) << myid;
250  sol_name << "ex9-init." << setfill('0') << setw(6) << myid;
251  ofstream omesh(mesh_name.str().c_str());
252  omesh.precision(precision);
253  pmesh->Print(omesh);
254  ofstream osol(sol_name.str().c_str());
255  osol.precision(precision);
256  u->Save(osol);
257  }
258 
259  // Create data collection for solution output: either VisItDataCollection for
260  // ascii data files, or SidreDataCollection for binary data files.
261  DataCollection *dc = NULL;
262  if (visit)
263  {
264  if (binary)
265  {
266 #ifdef MFEM_USE_SIDRE
267  dc = new SidreDataCollection("Example9-Parallel", pmesh);
268 #else
269  MFEM_ABORT("Must build with MFEM_USE_SIDRE=YES for binary output.");
270 #endif
271  }
272  else
273  {
274  dc = new VisItDataCollection("Example9-Parallel", pmesh);
275  dc->SetPrecision(precision);
276  }
277  dc->RegisterField("solution", u);
278  dc->SetCycle(0);
279  dc->SetTime(0.0);
280  dc->Save();
281  }
282 
283  socketstream sout;
284  if (visualization)
285  {
286  char vishost[] = "localhost";
287  int visport = 19916;
288  sout.open(vishost, visport);
289  if (!sout)
290  {
291  if (myid == 0)
292  cout << "Unable to connect to GLVis server at "
293  << vishost << ':' << visport << endl;
294  visualization = false;
295  if (myid == 0)
296  {
297  cout << "GLVis visualization disabled.\n";
298  }
299  }
300  else
301  {
302  sout << "parallel " << num_procs << " " << myid << "\n";
303  sout.precision(precision);
304  sout << "solution\n" << *pmesh << *u;
305  sout << "pause\n";
306  sout << flush;
307  if (myid == 0)
308  cout << "GLVis visualization paused."
309  << " Press space (in the GLVis window) to resume it.\n";
310  }
311  }
312 
313  // 9. Define the time-dependent evolution operator describing the ODE
314  // right-hand side, and define the ODE solver used for time integration.
315  FE_Evolution adv(*M, *K, *B);
316 
317  double t = 0.0;
318  adv.SetTime(t);
319 
320  // Create the time integrator
321  ODESolver *ode_solver = NULL;
322  CVODESolver *cvode = NULL;
323  ARKStepSolver *arkode = NULL;
324  switch (ode_solver_type)
325  {
326  case 1: ode_solver = new ForwardEulerSolver; break;
327  case 2: ode_solver = new RK2Solver(1.0); break;
328  case 3: ode_solver = new RK3SSPSolver; break;
329  case 4: ode_solver = new RK4Solver; break;
330  case 6: ode_solver = new RK6Solver; break;
331  case 7:
332  cvode = new CVODESolver(MPI_COMM_WORLD, CV_ADAMS);
333  cvode->Init(adv);
334  cvode->SetSStolerances(reltol, abstol);
335  cvode->SetMaxStep(dt);
336  cvode->UseSundialsLinearSolver();
337  ode_solver = cvode; break;
338  case 8:
339  case 9:
340  arkode = new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::EXPLICIT);
341  arkode->Init(adv);
342  arkode->SetSStolerances(reltol, abstol);
343  arkode->SetMaxStep(dt);
344  if (ode_solver_type == 9) { arkode->SetERKTableNum(FEHLBERG_13_7_8); }
345  ode_solver = arkode; break;
346  }
347 
348  // Initialize MFEM integrators, SUNDIALS integrators are initialized above
349  if (ode_solver_type < 7) { ode_solver->Init(adv); }
350 
351  // 10. Perform time-integration (looping over the time iterations, ti,
352  // with a time-step dt).
353  bool done = false;
354  for (int ti = 0; !done; )
355  {
356  double dt_real = min(dt, t_final - t);
357  ode_solver->Step(*U, t, dt_real);
358  ti++;
359 
360  done = (t >= t_final - 1e-8*dt);
361 
362  if (done || ti % vis_steps == 0)
363  {
364  if (myid == 0)
365  {
366  cout << "time step: " << ti << ", time: " << t << endl;
367  if (cvode) { cvode->PrintInfo(); }
368  if (arkode) { arkode->PrintInfo(); }
369  }
370 
371  // 11. Extract the parallel grid function corresponding to the finite
372  // element approximation U (the local solution on each processor).
373  *u = *U;
374 
375  if (visualization)
376  {
377  sout << "parallel " << num_procs << " " << myid << "\n";
378  sout << "solution\n" << *pmesh << *u << flush;
379  }
380 
381  if (visit)
382  {
383  dc->SetCycle(ti);
384  dc->SetTime(t);
385  dc->Save();
386  }
387  }
388  }
389 
390  // 12. Save the final solution in parallel. This output can be viewed later
391  // using GLVis: "glvis -np <np> -m ex9-mesh -g ex9-final".
392  {
393  *u = *U;
394  ostringstream sol_name;
395  sol_name << "ex9-final." << setfill('0') << setw(6) << myid;
396  ofstream osol(sol_name.str().c_str());
397  osol.precision(precision);
398  u->Save(osol);
399  }
400 
401  // 13. Free the used memory.
402  delete U;
403  delete u;
404  delete B;
405  delete b;
406  delete K;
407  delete k;
408  delete M;
409  delete m;
410  delete fes;
411  delete pmesh;
412  delete ode_solver;
413  delete dc;
414 
415  MPI_Finalize();
416  return 0;
417 }
418 
419 
420 // Implementation of class FE_Evolution
422  const Vector &_b)
423  : TimeDependentOperator(_M.Height()),
424  M(_M), K(_K), b(_b), M_solver(M.GetComm()), z(_M.Height())
425 {
426  M_prec.SetType(HypreSmoother::Jacobi);
427  M_solver.SetPreconditioner(M_prec);
428  M_solver.SetOperator(M);
429 
430  M_solver.iterative_mode = false;
431  M_solver.SetRelTol(1e-9);
432  M_solver.SetAbsTol(0.0);
433  M_solver.SetMaxIter(100);
434  M_solver.SetPrintLevel(0);
435 }
436 
437 void FE_Evolution::Mult(const Vector &x, Vector &y) const
438 {
439  // y = M^{-1} (K x + b)
440  K.Mult(x, z);
441  z += b;
442  M_solver.Mult(z, y);
443 }
444 
445 
446 // Velocity coefficient
447 void velocity_function(const Vector &x, Vector &v)
448 {
449  int dim = x.Size();
450 
451  // map to the reference [-1,1] domain
452  Vector X(dim);
453  for (int i = 0; i < dim; i++)
454  {
455  double center = (bb_min[i] + bb_max[i]) * 0.5;
456  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
457  }
458 
459  switch (problem)
460  {
461  case 0:
462  {
463  // Translations in 1D, 2D, and 3D
464  switch (dim)
465  {
466  case 1: v(0) = 1.0; break;
467  case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.); break;
468  case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
469  break;
470  }
471  break;
472  }
473  case 1:
474  case 2:
475  {
476  // Clockwise rotation in 2D around the origin
477  const double w = M_PI/2;
478  switch (dim)
479  {
480  case 1: v(0) = 1.0; break;
481  case 2: v(0) = w*X(1); v(1) = -w*X(0); break;
482  case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0; break;
483  }
484  break;
485  }
486  case 3:
487  {
488  // Clockwise twisting rotation in 2D around the origin
489  const double w = M_PI/2;
490  double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
491  d = d*d;
492  switch (dim)
493  {
494  case 1: v(0) = 1.0; break;
495  case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0); break;
496  case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0; break;
497  }
498  break;
499  }
500  }
501 }
502 
503 // Initial condition
504 double u0_function(const Vector &x)
505 {
506  int dim = x.Size();
507 
508  // map to the reference [-1,1] domain
509  Vector X(dim);
510  for (int i = 0; i < dim; i++)
511  {
512  double center = (bb_min[i] + bb_max[i]) * 0.5;
513  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
514  }
515 
516  switch (problem)
517  {
518  case 0:
519  case 1:
520  {
521  switch (dim)
522  {
523  case 1:
524  return exp(-40.*pow(X(0)-0.5,2));
525  case 2:
526  case 3:
527  {
528  double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
529  if (dim == 3)
530  {
531  const double s = (1. + 0.25*cos(2*M_PI*X(2)));
532  rx *= s;
533  ry *= s;
534  }
535  return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
536  erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
537  }
538  }
539  }
540  case 2:
541  {
542  double x_ = X(0), y_ = X(1), rho, phi;
543  rho = hypot(x_, y_);
544  phi = atan2(y_, x_);
545  return pow(sin(M_PI*rho),2)*sin(3*phi);
546  }
547  case 3:
548  {
549  const double f = M_PI;
550  return sin(f*X(0))*sin(f*X(1));
551  }
552  }
553  return 0.0;
554 }
555 
556 // Inflow boundary condition (zero for the problems considered in this example)
557 double inflow_function(const Vector &x)
558 {
559  switch (problem)
560  {
561  case 0:
562  case 1:
563  case 2:
564  case 3: return 0.0;
565  }
566  return 0.0;
567 }
Vector bb_max
Definition: ex9.cpp:60
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
Definition: sundials.cpp:151
Conjugate gradient method.
Definition: solvers.hpp:152
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
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:342
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:358
Base abstract class for first order time dependent operators.
Definition: operator.hpp:259
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
Definition: mesh.cpp:120
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
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
virtual void SetTime(const double _t)
Set the current time.
Definition: operator.hpp:303
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:485
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition: sundials.cpp:835
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:301
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:504
int main(int argc, char *argv[])
Definition: ex1.cpp:62
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
FE_Evolution(FiniteElementSpace &_vfes, Operator &_A, SparseMatrix &_Aflux)
Definition: ex18.hpp:104
Interface to ARKode&#39;s ARKStep module – additive Runge-Kutta methods.
Definition: sundials.hpp:197
Data collection with Sidre routines following the Conduit mesh blueprint specification.
Class for parallel linear form.
Definition: plinearform.hpp:26
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
virtual void Save()
Save the collection to disk.
Interface to the CVODE library – linear multi-step methods.
Definition: sundials.hpp:92
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7982
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
Definition: solvers.hpp:65
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3924
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:251
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition: sundials.cpp:348
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4102
Parallel smoothers in hypre.
Definition: hypre.hpp:581
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:66
int Dimension() const
Definition: mesh.hpp:744
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetTime(double t)
Set physical time (for time-dependent simulations)
Vector bb_min
Definition: ex9.cpp:60
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:73
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:145
void SetAbsTol(double atol)
Definition: solvers.hpp:64
void SetRelTol(double rtol)
Definition: solvers.hpp:63
void velocity_function(const Vector &x, Vector &v)
Definition: ex9.cpp:478
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)
Definition: optparser.hpp:76
int problem
Definition: ex15.cpp:54
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:132
virtual ~FE_Evolution()
Definition: ex9p.cpp:76
virtual void Mult(const Vector &x, Vector &y) const
Perform the action of the operator: y = k = f(x, t), where k solves the algebraic equation F(x...
Definition: ex18.hpp:130
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:191
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void PrintInfo() const
Print various ARKStep statistics.
Definition: sundials.cpp:871
int dim
Definition: ex24.cpp:43
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator. Assumes ownership of bfi.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
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)
double u0_function(const Vector &x)
Definition: ex9.cpp:535
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:166
class for C-function coefficient
Vector data type.
Definition: vector.hpp:48
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:91
void PrintInfo() const
Print various CVODE statistics.
Definition: sundials.cpp:360
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:847
The classical forward Euler method.
Definition: ode.hpp:99
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:187
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults...
Definition: sundials.cpp:565
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:829
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:34
double inflow_function(const Vector &x)
Definition: ex9.cpp:588
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
Definition: sundials.cpp:322
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator. Assumes ownership of bfi.
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:143
bool Good() const
Definition: optparser.hpp:125
alpha (q . grad u, v)