45#ifndef MFEM_USE_SUNDIALS
46#error This example requires that MFEM is built with MFEM_USE_SUNDIALS=YES
68class DG_Solver :
public Solver
80 BlockILU::Reordering::MINIMUM_DISCARDED_FILL),
91 void SetTimeStep(
double dt_)
106 void SetOperator(
const Operator &op)
113 linear_solver.
Mult(x, y);
129 DG_Solver *dg_solver;
137 virtual void ImplicitSolve(
const double dt,
const Vector &x,
Vector &k);
139 virtual ~FE_Evolution();
143int main(
int argc,
char *argv[])
150 const char *mesh_file =
"../../data/periodic-hexagon.mesh";
156 const char *device_config =
"cpu";
157 int ode_solver_type = 7;
158 double t_final = 10.0;
160 bool visualization =
false;
162 bool paraview =
false;
167 const double reltol = 1e-2, abstol = 1e-2;
170 cout.precision(precision);
173 args.
AddOption(&mesh_file,
"-m",
"--mesh",
174 "Mesh file to use.");
176 "Problem setup to use. See options in velocity_function().");
177 args.
AddOption(&ref_levels,
"-r",
"--refine",
178 "Number of times to refine the mesh uniformly.");
180 "Order (degree) of the finite elements.");
181 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
182 "--no-partial-assembly",
"Enable Partial Assembly.");
183 args.
AddOption(&ea,
"-ea",
"--element-assembly",
"-no-ea",
184 "--no-element-assembly",
"Enable Element Assembly.");
185 args.
AddOption(&fa,
"-fa",
"--full-assembly",
"-no-fa",
186 "--no-full-assembly",
"Enable Full Assembly.");
187 args.
AddOption(&device_config,
"-d",
"--device",
188 "Device configuration string, see Device::Configure().");
189 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
191 "1 - Forward Euler,\n\t"
196 "7 - CVODE (adaptive order implicit Adams),\n\t"
197 "8 - ARKODE default (4th order) explicit,\n\t"
199 args.
AddOption(&t_final,
"-tf",
"--t-final",
200 "Final time; start time is 0.");
201 args.
AddOption(&dt,
"-dt",
"--time-step",
203 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
204 "--no-visualization",
205 "Enable or disable GLVis visualization.");
206 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
207 "--no-visit-datafiles",
208 "Save data files for VisIt (visit.llnl.gov) visualization.");
209 args.
AddOption(¶view,
"-paraview",
"--paraview-datafiles",
"-no-paraview",
210 "--no-paraview-datafiles",
211 "Save data files for ParaView (paraview.org) visualization.");
212 args.
AddOption(&binary,
"-binary",
"--binary-datafiles",
"-ascii",
214 "Use binary (Sidre) or ascii format for VisIt data files.");
215 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
216 "Visualize every n-th timestep.");
224 if (ode_solver_type < 1 || ode_solver_type > 9)
226 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
231 Device device(device_config);
236 Mesh mesh(mesh_file, 1, 1);
243 for (
int lev = 0; lev < ref_levels; lev++)
258 cout <<
"Number of unknowns: " << fes.
GetVSize() << endl;
285 constexpr double alpha = -1.0;
293 b.AddBdrFaceIntegrator(
307 u.ProjectCoefficient(u0);
310 ofstream omesh(
"ex9.mesh");
311 omesh.precision(precision);
313 ofstream osol(
"ex9-init.gf");
314 osol.precision(precision);
328 MFEM_ABORT(
"Must build with MFEM_USE_SIDRE=YES for binary output.");
364 cout <<
"Unable to connect to GLVis server at "
365 <<
vishost <<
':' << visport << endl;
366 visualization =
false;
367 cout <<
"GLVis visualization disabled.\n";
371 sout.precision(precision);
372 sout <<
"solution\n" << mesh <<
u;
375 cout <<
"GLVis visualization paused."
376 <<
" Press space (in the GLVis window) to resume it.\n";
382 FE_Evolution adv(m, k,
b);
391 switch (ode_solver_type)
394 case 2: ode_solver =
new RK2Solver(1.0);
break;
396 case 4: ode_solver =
new RK4Solver;
break;
397 case 6: ode_solver =
new RK6Solver;
break;
404 ode_solver = cvode;
break;
411 ode_solver = arkode;
break;
418 ode_solver = arkode;
break;
422 if (ode_solver_type < 7) { ode_solver->
Init(adv); }
427 for (
int ti = 0; !done; )
429 double dt_real = min(dt, t_final - t);
430 ode_solver->
Step(
u, t, dt_real);
433 done = (t >= t_final - 1e-8*dt);
435 if (done || ti % vis_steps == 0)
437 cout <<
"time step: " << ti <<
", time: " << t << endl;
443 sout <<
"solution\n" << mesh <<
u << flush;
465 ofstream osol(
"ex9-final.gf");
466 osol.precision(precision);
482 M(M_), K(K_),
b(b_), z(height)
485 if (M.GetAssemblyLevel() == AssemblyLevel::LEGACY)
488 M_solver.SetOperator(M.SpMat());
489 dg_solver =
new DG_Solver(M.SpMat(), K.SpMat(), *M.FESpace());
494 M_solver.SetOperator(M);
497 M_solver.SetPreconditioner(*M_prec);
498 M_solver.iterative_mode =
false;
499 M_solver.SetRelTol(1e-9);
500 M_solver.SetAbsTol(0.0);
501 M_solver.SetMaxIter(100);
502 M_solver.SetPrintLevel(0);
505void FE_Evolution::Mult(
const Vector &x,
Vector &y)
const
513void FE_Evolution::ImplicitSolve(
const double dt,
const Vector &x,
Vector &k)
515 MFEM_VERIFY(dg_solver != NULL,
516 "Implicit time integration is not supported with partial assembly");
519 dg_solver->SetTimeStep(dt);
520 dg_solver->Mult(z, k);
523FE_Evolution::~FE_Evolution()
536 for (
int i = 0; i <
dim; i++)
549 case 1: v(0) = 1.0;
break;
550 case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.);
break;
551 case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
560 const double w = M_PI/2;
563 case 1: v(0) = 1.0;
break;
564 case 2: v(0) = w*X(1); v(1) = -w*X(0);
break;
565 case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0;
break;
572 const double w = M_PI/2;
573 double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
577 case 1: v(0) = 1.0;
break;
578 case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0);
break;
579 case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0;
break;
593 for (
int i = 0; i <
dim; i++)
607 return exp(-40.*pow(X(0)-0.5,2));
611 double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
614 const double s = (1. + 0.25*cos(2*M_PI*X(2)));
618 return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
619 erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
625 double x_ = X(0), y_ = X(1), rho, phi;
628 return pow(sin(M_PI*rho),2)*sin(3*phi);
632 const double f = M_PI;
633 return sin(
f*X(0))*sin(
f*X(1));
Interface to ARKode's ARKStep module – additive Runge-Kutta methods.
void SetMaxStep(double dt_max)
Set the maximum time step.
void PrintInfo() const
Print various ARKStep statistics.
void SetOrder(int order)
Chooses integration order for all explicit / implicit / IMEX methods.
@ EXPLICIT
Explicit RK method.
void Init(TimeDependentOperator &f_) override
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults.
void SetERKTableNum(ARKODE_ERKTableID table_id)
Choose a specific Butcher table for an explicit RK method.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
@ GaussLobatto
Closed type.
Conjugate gradient method.
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the Conjugate Gradient method.
Interface to the CVODE library – linear multi-step methods.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
void Init(TimeDependentOperator &f_) override
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
void SetMaxStep(double dt_max)
Set the maximum time step.
void PrintInfo() const
Print various CVODE statistics.
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
Data type for scaled Jacobi-type smoother of sparse matrix.
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetTime(real_t t)
Set physical time (for time-dependent simulations)
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
virtual void Save()
Save the collection to disk.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
const FiniteElement * GetTypicalFE() const
Return GetFE(0) if the local mesh is not empty; otherwise return a typical FE based on the Geometry t...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
The classical forward Euler method.
A general function coefficient.
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the GMRES method.
Class for grid function - Vector with associated FE space.
void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
void SetRelTol(real_t rtol)
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void SetMaxIter(int max_it)
void SetAbsTol(real_t atol)
Arbitrary order "L2-conforming" discontinuous finite elements.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
int Dimension() const
Dimension of the reference space used within the elements.
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
virtual void Step(Vector &x, real_t &t, real_t &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Jacobi smoothing for a given bilinear form (no matrix necessary).
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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 'var' to receive the value. Enable/disable tags are used to set the bool...
bool Good() const
Return true if the command line options were parsed successfully.
Helper class for ParaView visualization data.
void SetHighOrderOutput(bool high_order_output_)
void SetLevelsOfDetail(int levels_of_detail_)
void SetDataFormat(VTKFormat fmt)
Third-order, strong stability preserving (SSP) Runge-Kutta method.
The classical explicit forth-order Runge-Kutta method, RK4.
Data collection with Sidre routines following the Conduit mesh blueprint specification.
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Base abstract class for first order time dependent operators.
virtual void SetTime(const real_t t_)
Set the current time.
A general vector function coefficient.
int Size() const
Returns the size of the vector.
Data collection with VisIt I/O routines.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
real_t u(const Vector &xvec)
std::function< real_t(const Vector &)> f(real_t mass_coeff)
void velocity_function(const Vector &x, Vector &v)
double u0_function(const Vector &x)
double inflow_function(const Vector &x)
constexpr ARKODE_ERKTableID ARKODE_FEHLBERG_13_7_8