32 #ifndef MFEM_USE_SUNDIALS
33 #error This example requires that MFEM is built with MFEM_USE_SUNDIALS=YES
80 int main(
int argc,
char *argv[])
84 MPI_Init(&argc, &argv);
85 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
86 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
90 const char *mesh_file =
"../../data/periodic-hexagon.mesh";
91 int ser_ref_levels = 2;
92 int par_ref_levels = 0;
94 int ode_solver_type = 4;
95 double t_final = 10.0;
97 bool visualization =
true;
103 const double reltol = 1e-2, abstol = 1e-2;
106 cout.precision(precision);
109 args.
AddOption(&mesh_file,
"-m",
"--mesh",
110 "Mesh file to use.");
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.");
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",
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",
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.");
157 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
165 switch (ode_solver_type)
168 case 2: ode_solver =
new RK2Solver(1.0);
break;
170 case 4: ode_solver =
new RK4Solver;
break;
171 case 6: ode_solver =
new RK6Solver;
break;
173 cvode =
new CVODESolver(MPI_COMM_WORLD, CV_ADAMS, CV_FUNCTIONAL);
176 ode_solver = cvode;
break;
179 arkode =
new ARKODESolver(MPI_COMM_WORLD, ARKODESolver::EXPLICIT);
182 if (ode_solver_type == 13) { arkode->
SetERKTableNum(FEHLBERG_13_7_8); }
183 ode_solver = arkode;
break;
187 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
197 for (
int lev = 0; lev < ser_ref_levels; lev++)
212 for (
int lev = 0; lev < par_ref_levels; lev++)
225 cout <<
"Number of unknowns: " << global_vSize << endl;
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);
273 ofstream osol(sol_name.str().c_str());
274 osol.precision(precision);
285 #ifdef MFEM_USE_SIDRE
288 MFEM_ABORT(
"Must build with MFEM_USE_SIDRE=YES for binary output.");
305 char vishost[] =
"localhost";
307 sout.
open(vishost, visport);
311 cout <<
"Unable to connect to GLVis server at "
312 << vishost <<
':' << visport << endl;
313 visualization =
false;
316 cout <<
"GLVis visualization disabled.\n";
321 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
322 sout.precision(precision);
323 sout <<
"solution\n" << *pmesh << *u;
327 cout <<
"GLVis visualization paused."
328 <<
" Press space (in the GLVis window) to resume it.\n";
339 ode_solver->
Init(adv);
342 for (
int ti = 0; !done; )
344 double dt_real = min(dt, t_final - t);
345 ode_solver->
Step(*U, t, dt_real);
348 done = (t >= t_final - 1e-8*dt);
350 if (done || ti % vis_steps == 0)
354 cout <<
"time step: " << ti <<
", time: " << t << endl;
365 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
366 sout <<
"solution\n" << *pmesh << *u << flush;
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);
412 M(_M), K(_K), b(_b), M_solver(M.GetComm()), z(_M.Height())
414 M_prec.SetType(HypreSmoother::Jacobi);
415 M_solver.SetPreconditioner(M_prec);
416 M_solver.SetOperator(M);
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);
441 for (
int i = 0; i <
dim; i++)
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.);
465 const double w = M_PI/2;
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;
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.);
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;
498 for (
int i = 0; i <
dim; i++)
512 return exp(-40.*pow(X(0)-0.5,2));
516 double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
519 const double s = (1. + 0.25*cos(2*M_PI*X(2)));
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;
530 double x_ = X(0), y_ = X(1), rho, phi;
533 return pow(sin(M_PI*rho),2)*sin(3*phi);
537 const double f = M_PI;
538 return sin(f*X(0))*sin(f*X(1));
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
Conjugate gradient method.
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.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
virtual void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Base abstract class for time dependent operators.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
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.
int Size() const
Returns the size of the vector.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
virtual void Save(std::ostream &out) const
Abstract parallel finite element space.
virtual void ProjectCoefficient(Coefficient &coeff)
void SetSStolerances(double reltol, double abstol)
Specify the scalar relative and scalar absolute tolerances.
void SetERKTableNum(int table_num)
Choose a specific Butcher table for explicit RK method.
FE_Evolution(FiniteElementSpace &_vfes, Operator &_A, SparseMatrix &_Aflux)
Data collection with Sidre routines following the Conduit mesh blueprint specification.
Wrapper for SUNDIALS' CVODE library – Multi-step time integration.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Data collection with VisIt I/O routines.
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
HYPRE_Int GlobalTrueVSize() const
void SetMaxStep(double dt_max)
Set the maximum time step of the linear multistep method.
virtual void Print(std::ostream &out=mfem::out) const
Parallel smoothers in hypre.
void PrintUsage(std::ostream &out) const
void SetTime(double t)
Set physical time (for time-dependent simulations)
Wrapper for SUNDIALS' ARKODE library – Runge-Kutta time integration.
Wrapper for hypre's parallel vector class.
The classical explicit forth-order Runge-Kutta method, RK4.
void velocity_function(const Vector &x, Vector &v)
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
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)
Third-order, strong stability preserving (SSP) Runge-Kutta method.
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...
NURBSExtension * NURBSext
Optional NURBS mesh extension.
void PrintOptions(std::ostream &out) const
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
int open(const char hostname[], int port)
double u0_function(const Vector &x)
class for C-function coefficient
void PrintInfo() const
Print CVODE statistics.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
Class for parallel grid function.
The classical forward Euler method.
void SetMaxStep(double dt_max)
Set the maximum time step of the Runge-Kutta method.
Wrapper for hypre's ParCSR matrix class.
Class for parallel meshes.
void PrintInfo() const
Print ARKODE statistics.
double inflow_function(const Vector &x)
Arbitrary order "L2-conforming" discontinuous finite elements.