MFEM  v3.3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
sundials/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 11 -dt 0.0025
8 // mpirun -np 4 ex9p -m ../../data/periodic-square.mesh -p 1 -rp 1 -s 12 -dt 0.0025 -tf 9
9 // mpirun -np 4 ex9p -m ../../data/periodic-hexagon.mesh -p 0 -rp 1 -s 11 -dt 0.0009 -vs 25
10 // mpirun -np 4 ex9p -m ../../data/periodic-hexagon.mesh -p 0 -rp 1 -s 13 -dt 0.005 -vs 15
11 // mpirun -np 4 ex9p -m ../../data/amr-quad.mesh -p 1 -rp 1 -s 13 -dt 0.001 -tf 9
12 // mpirun -np 4 ex9p -m ../../data/star-q3.mesh -p 1 -rp 1 -s 13 -dt 0.0025 -tf 9
13 // mpirun -np 4 ex9p -m ../../data/disc-nurbs.mesh -p 1 -rp 2 -s 11 -dt 0.0025 -tf 9
14 // mpirun -np 4 ex9p -m ../../data/periodic-cube.mesh -p 0 -rp 1 -s 12 -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: 1 - Forward Euler,\n\t"
121  " 2 - RK2 SSP, 3 - RK3 SSP, 4 - RK4, 6 - RK6,\n\t"
122  " 11 - CVODE (adaptive order) explicit,\n\t"
123  " 12 - ARKODE default (4th order) explicit,\n\t"
124  " 13 - ARKODE RK8.");
125  args.AddOption(&t_final, "-tf", "--t-final",
126  "Final time; start time is 0.");
127  args.AddOption(&dt, "-dt", "--time-step",
128  "Time step.");
129  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
130  "--no-visualization",
131  "Enable or disable GLVis visualization.");
132  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
133  "--no-visit-datafiles",
134  "Save data files for VisIt (visit.llnl.gov) visualization.");
135  args.AddOption(&binary, "-binary", "--binary-datafiles", "-ascii",
136  "--ascii-datafiles",
137  "Use binary (Sidre) or ascii format for VisIt data files.");
138  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
139  "Visualize every n-th timestep.");
140  args.Parse();
141  if (!args.Good())
142  {
143  if (myid == 0)
144  {
145  args.PrintUsage(cout);
146  }
147  MPI_Finalize();
148  return 1;
149  }
150  if (myid == 0)
151  {
152  args.PrintOptions(cout);
153  }
154 
155  // 3. Read the serial mesh from the given mesh file on all processors. We can
156  // handle geometrically periodic meshes in this code.
157  Mesh *mesh = new Mesh(mesh_file, 1, 1);
158  int dim = mesh->Dimension();
159 
160  // 4. Define the ODE solver used for time integration. Several explicit
161  // Runge-Kutta methods are available.
162  ODESolver *ode_solver = NULL;
163  CVODESolver *cvode = NULL;
164  ARKODESolver *arkode = NULL;
165  switch (ode_solver_type)
166  {
167  case 1: ode_solver = new ForwardEulerSolver; break;
168  case 2: ode_solver = new RK2Solver(1.0); break;
169  case 3: ode_solver = new RK3SSPSolver; break;
170  case 4: ode_solver = new RK4Solver; break;
171  case 6: ode_solver = new RK6Solver; break;
172  case 11:
173  cvode = new CVODESolver(MPI_COMM_WORLD, CV_ADAMS, CV_FUNCTIONAL);
174  cvode->SetSStolerances(reltol, abstol);
175  cvode->SetMaxStep(dt);
176  ode_solver = cvode; break;
177  case 12:
178  case 13:
179  arkode = new ARKODESolver(MPI_COMM_WORLD, ARKODESolver::EXPLICIT);
180  arkode->SetSStolerances(reltol, abstol);
181  arkode->SetMaxStep(dt);
182  if (ode_solver_type == 13) { arkode->SetERKTableNum(FEHLBERG_13_7_8); }
183  ode_solver = arkode; break;
184  default:
185  if (myid == 0)
186  {
187  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
188  }
189  MPI_Finalize();
190  return 3;
191  }
192 
193  // 5. Refine the mesh in serial to increase the resolution. In this example
194  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
195  // a command-line parameter. If the mesh is of NURBS type, we convert it
196  // to a (piecewise-polynomial) high-order mesh.
197  for (int lev = 0; lev < ser_ref_levels; lev++)
198  {
199  mesh->UniformRefinement();
200  }
201  if (mesh->NURBSext)
202  {
203  mesh->SetCurvature(max(order, 1));
204  }
205  mesh->GetBoundingBox(bb_min, bb_max, max(order, 1));
206 
207  // 6. Define the 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  // 7. Define the parallel discontinuous DG finite element space on the
218  // parallel refined mesh of the given polynomial order.
219  DG_FECollection fec(order, dim);
220  ParFiniteElementSpace *fes = new ParFiniteElementSpace(pmesh, &fec);
221 
222  HYPRE_Int global_vSize = fes->GlobalTrueVSize();
223  if (myid == 0)
224  {
225  cout << "Number of unknowns: " << global_vSize << endl;
226  }
227 
228  // 8. Set up and assemble the parallel bilinear and linear forms (and the
229  // parallel hypre matrices) corresponding to the DG discretization. The
230  // DGTraceIntegrator involves integrals over mesh interior faces.
234 
235  ParBilinearForm *m = new ParBilinearForm(fes);
237  ParBilinearForm *k = new ParBilinearForm(fes);
238  k->AddDomainIntegrator(new ConvectionIntegrator(velocity, -1.0));
240  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
242  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
243 
244  ParLinearForm *b = new ParLinearForm(fes);
246  new BoundaryFlowIntegrator(inflow, velocity, -1.0, -0.5));
247 
248  m->Assemble();
249  m->Finalize();
250  int skip_zeros = 0;
251  k->Assemble(skip_zeros);
252  k->Finalize(skip_zeros);
253  b->Assemble();
254 
258 
259  // 9. Define the initial conditions, save the corresponding grid function to
260  // a file and (optionally) save data in the VisIt format and initialize
261  // GLVis visualization.
262  ParGridFunction *u = new ParGridFunction(fes);
263  u->ProjectCoefficient(u0);
264  HypreParVector *U = u->GetTrueDofs();
265 
266  {
267  ostringstream mesh_name, sol_name;
268  mesh_name << "ex9-mesh." << setfill('0') << setw(6) << myid;
269  sol_name << "ex9-init." << setfill('0') << setw(6) << myid;
270  ofstream omesh(mesh_name.str().c_str());
271  omesh.precision(precision);
272  pmesh->Print(omesh);
273  ofstream osol(sol_name.str().c_str());
274  osol.precision(precision);
275  u->Save(osol);
276  }
277 
278  // Create data collection for solution output: either VisItDataCollection for
279  // ascii data files, or SidreDataCollection for binary data files.
280  DataCollection *dc = NULL;
281  if (visit)
282  {
283  if (binary)
284  {
285 #ifdef MFEM_USE_SIDRE
286  dc = new SidreDataCollection("Example9-Parallel", pmesh);
287 #else
288  MFEM_ABORT("Must build with MFEM_USE_SIDRE=YES for binary output.");
289 #endif
290  }
291  else
292  {
293  dc = new VisItDataCollection("Example9-Parallel", pmesh);
294  dc->SetPrecision(precision);
295  }
296  dc->RegisterField("solution", u);
297  dc->SetCycle(0);
298  dc->SetTime(0.0);
299  dc->Save();
300  }
301 
302  socketstream sout;
303  if (visualization)
304  {
305  char vishost[] = "localhost";
306  int visport = 19916;
307  sout.open(vishost, visport);
308  if (!sout)
309  {
310  if (myid == 0)
311  cout << "Unable to connect to GLVis server at "
312  << vishost << ':' << visport << endl;
313  visualization = false;
314  if (myid == 0)
315  {
316  cout << "GLVis visualization disabled.\n";
317  }
318  }
319  else
320  {
321  sout << "parallel " << num_procs << " " << myid << "\n";
322  sout.precision(precision);
323  sout << "solution\n" << *pmesh << *u;
324  sout << "pause\n";
325  sout << flush;
326  if (myid == 0)
327  cout << "GLVis visualization paused."
328  << " Press space (in the GLVis window) to resume it.\n";
329  }
330  }
331 
332  // 10. Define the time-dependent evolution operator describing the ODE
333  // right-hand side, and perform time-integration (looping over the time
334  // iterations, ti, with a time-step dt).
335  FE_Evolution adv(*M, *K, *B);
336 
337  double t = 0.0;
338  adv.SetTime(t);
339  ode_solver->Init(adv);
340 
341  bool done = false;
342  for (int ti = 0; !done; )
343  {
344  double dt_real = min(dt, t_final - t);
345  ode_solver->Step(*U, t, dt_real);
346  ti++;
347 
348  done = (t >= t_final - 1e-8*dt);
349 
350  if (done || ti % vis_steps == 0)
351  {
352  if (myid == 0)
353  {
354  cout << "time step: " << ti << ", time: " << t << endl;
355  if (cvode) { cvode->PrintInfo(); }
356  if (arkode) { arkode->PrintInfo(); }
357  }
358 
359  // 11. Extract the parallel grid function corresponding to the finite
360  // element approximation U (the local solution on each processor).
361  *u = *U;
362 
363  if (visualization)
364  {
365  sout << "parallel " << num_procs << " " << myid << "\n";
366  sout << "solution\n" << *pmesh << *u << flush;
367  }
368 
369  if (visit)
370  {
371  dc->SetCycle(ti);
372  dc->SetTime(t);
373  dc->Save();
374  }
375  }
376  }
377 
378  // 12. Save the final solution in parallel. This output can be viewed later
379  // using GLVis: "glvis -np <np> -m ex9-mesh -g ex9-final".
380  {
381  *u = *U;
382  ostringstream sol_name;
383  sol_name << "ex9-final." << setfill('0') << setw(6) << myid;
384  ofstream osol(sol_name.str().c_str());
385  osol.precision(precision);
386  u->Save(osol);
387  }
388 
389  // 13. Free the used memory.
390  delete U;
391  delete u;
392  delete B;
393  delete b;
394  delete K;
395  delete k;
396  delete M;
397  delete m;
398  delete fes;
399  delete pmesh;
400  delete ode_solver;
401  delete dc;
402 
403  MPI_Finalize();
404  return 0;
405 }
406 
407 
408 // Implementation of class FE_Evolution
410  const Vector &_b)
411  : TimeDependentOperator(_M.Height()),
412  M(_M), K(_K), b(_b), M_solver(M.GetComm()), z(_M.Height())
413 {
414  M_prec.SetType(HypreSmoother::Jacobi);
415  M_solver.SetPreconditioner(M_prec);
416  M_solver.SetOperator(M);
417 
418  M_solver.iterative_mode = false;
419  M_solver.SetRelTol(1e-9);
420  M_solver.SetAbsTol(0.0);
421  M_solver.SetMaxIter(100);
422  M_solver.SetPrintLevel(0);
423 }
424 
425 void FE_Evolution::Mult(const Vector &x, Vector &y) const
426 {
427  // y = M^{-1} (K x + b)
428  K.Mult(x, z);
429  z += b;
430  M_solver.Mult(z, y);
431 }
432 
433 
434 // Velocity coefficient
435 void velocity_function(const Vector &x, Vector &v)
436 {
437  int dim = x.Size();
438 
439  // map to the reference [-1,1] domain
440  Vector X(dim);
441  for (int i = 0; i < dim; i++)
442  {
443  double center = (bb_min[i] + bb_max[i]) * 0.5;
444  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
445  }
446 
447  switch (problem)
448  {
449  case 0:
450  {
451  // Translations in 1D, 2D, and 3D
452  switch (dim)
453  {
454  case 1: v(0) = 1.0; break;
455  case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.); break;
456  case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
457  break;
458  }
459  break;
460  }
461  case 1:
462  case 2:
463  {
464  // Clockwise rotation in 2D around the origin
465  const double w = M_PI/2;
466  switch (dim)
467  {
468  case 1: v(0) = 1.0; break;
469  case 2: v(0) = w*X(1); v(1) = -w*X(0); break;
470  case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0; break;
471  }
472  break;
473  }
474  case 3:
475  {
476  // Clockwise twisting rotation in 2D around the origin
477  const double w = M_PI/2;
478  double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
479  d = d*d;
480  switch (dim)
481  {
482  case 1: v(0) = 1.0; break;
483  case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0); break;
484  case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0; break;
485  }
486  break;
487  }
488  }
489 }
490 
491 // Initial condition
492 double u0_function(const Vector &x)
493 {
494  int dim = x.Size();
495 
496  // map to the reference [-1,1] domain
497  Vector X(dim);
498  for (int i = 0; i < dim; i++)
499  {
500  double center = (bb_min[i] + bb_max[i]) * 0.5;
501  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
502  }
503 
504  switch (problem)
505  {
506  case 0:
507  case 1:
508  {
509  switch (dim)
510  {
511  case 1:
512  return exp(-40.*pow(X(0)-0.5,2));
513  case 2:
514  case 3:
515  {
516  double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
517  if (dim == 3)
518  {
519  const double s = (1. + 0.25*cos(2*M_PI*X(2)));
520  rx *= s;
521  ry *= s;
522  }
523  return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
524  erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
525  }
526  }
527  }
528  case 2:
529  {
530  double x_ = X(0), y_ = X(1), rho, phi;
531  rho = hypot(x_, y_);
532  phi = atan2(y_, x_);
533  return pow(sin(M_PI*rho),2)*sin(3*phi);
534  }
535  case 3:
536  {
537  const double f = M_PI;
538  return sin(f*X(0))*sin(f*X(1));
539  }
540  }
541  return 0.0;
542 }
543 
544 // Inflow boundary condition (zero for the problems considered in this example)
545 double inflow_function(const Vector &x)
546 {
547  switch (problem)
548  {
549  case 0:
550  case 1:
551  case 2:
552  case 3: return 0.0;
553  }
554  return 0.0;
555 }
Vector bb_max
Definition: ex9.cpp:52
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
Conjugate gradient method.
Definition: solvers.hpp:111
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:224
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:295
virtual void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.hpp:37
Base abstract class for time dependent operators.
Definition: operator.hpp:142
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:468
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Definition: mesh.cpp:95
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:51
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 SetTime(const double _t)
Set the current time.
Definition: operator.hpp:172
int Size() const
Returns the size of the vector.
Definition: vector.hpp:113
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:361
Abstract parallel finite element space.
Definition: pfespace.hpp:31
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:275
void SetSStolerances(double reltol, double abstol)
Specify the scalar relative and scalar absolute tolerances.
Definition: sundials.cpp:477
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void SetERKTableNum(int table_num)
Choose a specific Butcher table for explicit RK method.
Definition: sundials.cpp:528
FE_Evolution(FiniteElementSpace &_vfes, Operator &_A, SparseMatrix &_Aflux)
Definition: ex18.hpp:104
Data collection with Sidre routines following the Conduit mesh blueprint specification.
Class for parallel linear form.
Definition: plinearform.hpp:26
int dim
Definition: ex3.cpp:47
Wrapper for SUNDIALS&#39; CVODE library – Multi-step time integration.
Definition: sundials.hpp:133
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6718
Data collection with VisIt I/O routines.
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3315
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:178
void SetMaxStep(double dt_max)
Set the maximum time step of the linear multistep method.
Definition: sundials.hpp:183
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:3399
Parallel smoothers in hypre.
Definition: hypre.hpp:504
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator.
Definition: linearform.cpp:38
int Dimension() const
Definition: mesh.hpp:641
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetTime(double t)
Set physical time (for time-dependent simulations)
Wrapper for SUNDIALS&#39; ARKODE library – Runge-Kutta time integration.
Definition: sundials.hpp:220
Vector bb_min
Definition: ex9.cpp:52
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:73
int main()
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:147
void velocity_function(const Vector &x, Vector &v)
Definition: ex9.cpp:339
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:475
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:74
int problem
Definition: ex15.cpp:57
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:134
virtual ~FE_Evolution()
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:174
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:126
Class for parallel bilinear form.
int open(const char hostname[], int port)
double u0_function(const Vector &x)
Definition: ex9.cpp:396
class for C-function coefficient
Vector data type.
Definition: vector.hpp:41
void PrintInfo() const
Print CVODE statistics.
Definition: sundials.cpp:400
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
Class for parallel grid function.
Definition: pgridfunc.hpp:32
The classical forward Euler method.
Definition: ode.hpp:101
void SetMaxStep(double dt_max)
Set the maximum time step of the Runge-Kutta method.
Definition: sundials.hpp:276
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:175
Class for parallel meshes.
Definition: pmesh.hpp:29
void PrintInfo() const
Print ARKODE statistics.
Definition: sundials.cpp:691
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:449
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator.
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:195
bool Good() const
Definition: optparser.hpp:120
alpha (q . grad u, v)