43 #ifndef MFEM_USE_SUNDIALS
44 #error This example requires that MFEM is built with MFEM_USE_SUNDIALS=YES
66 class DG_Solver :
public Solver
78 BlockILU::Reordering::MINIMUM_DISCARDED_FILL),
81 linear_solver.iterative_mode =
false;
82 linear_solver.SetRelTol(1e-9);
83 linear_solver.SetAbsTol(0.0);
84 linear_solver.SetMaxIter(100);
85 linear_solver.SetPrintLevel(0);
86 linear_solver.SetPreconditioner(prec);
89 void SetTimeStep(
double dt_)
100 linear_solver.SetOperator(A);
104 void SetOperator(
const Operator &op)
106 linear_solver.SetOperator(op);
111 linear_solver.Mult(x, y);
127 DG_Solver *dg_solver;
135 virtual void ImplicitSolve(
const double dt,
const Vector &x,
Vector &k);
141 int main(
int argc,
char *argv[])
145 const char *mesh_file =
"../../data/periodic-hexagon.mesh";
151 const char *device_config =
"cpu";
152 int ode_solver_type = 7;
153 double t_final = 10.0;
155 bool visualization =
false;
157 bool paraview =
false;
162 const double reltol = 1e-2, abstol = 1e-2;
165 cout.precision(precision);
168 args.
AddOption(&mesh_file,
"-m",
"--mesh",
169 "Mesh file to use.");
171 "Problem setup to use. See options in velocity_function().");
172 args.
AddOption(&ref_levels,
"-r",
"--refine",
173 "Number of times to refine the mesh uniformly.");
175 "Order (degree) of the finite elements.");
176 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
177 "--no-partial-assembly",
"Enable Partial Assembly.");
178 args.
AddOption(&ea,
"-ea",
"--element-assembly",
"-no-ea",
179 "--no-element-assembly",
"Enable Element Assembly.");
180 args.
AddOption(&fa,
"-fa",
"--full-assembly",
"-no-fa",
181 "--no-full-assembly",
"Enable Full Assembly.");
182 args.
AddOption(&device_config,
"-d",
"--device",
183 "Device configuration string, see Device::Configure().");
184 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
186 "1 - Forward Euler,\n\t"
191 "7 - CVODE (adaptive order implicit Adams),\n\t"
192 "8 - ARKODE default (4th order) explicit,\n\t"
194 args.
AddOption(&t_final,
"-tf",
"--t-final",
195 "Final time; start time is 0.");
196 args.
AddOption(&dt,
"-dt",
"--time-step",
198 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
199 "--no-visualization",
200 "Enable or disable GLVis visualization.");
201 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
202 "--no-visit-datafiles",
203 "Save data files for VisIt (visit.llnl.gov) visualization.");
204 args.
AddOption(¶view,
"-paraview",
"--paraview-datafiles",
"-no-paraview",
205 "--no-paraview-datafiles",
206 "Save data files for ParaView (paraview.org) visualization.");
207 args.
AddOption(&binary,
"-binary",
"--binary-datafiles",
"-ascii",
209 "Use binary (Sidre) or ascii format for VisIt data files.");
210 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
211 "Visualize every n-th timestep.");
219 if (ode_solver_type < 1 || ode_solver_type > 9)
221 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
226 Device device(device_config);
231 Mesh mesh(mesh_file, 1, 1);
238 for (
int lev = 0; lev < ref_levels; lev++)
253 cout <<
"Number of unknowns: " << fes.
GetVSize() << endl;
304 ofstream omesh(
"ex9.mesh");
305 omesh.precision(precision);
307 ofstream osol(
"ex9-init.gf");
308 osol.precision(precision);
319 #ifdef MFEM_USE_SIDRE
322 MFEM_ABORT(
"Must build with MFEM_USE_SIDRE=YES for binary output.");
355 sout.
open(vishost, visport);
358 cout <<
"Unable to connect to GLVis server at "
359 << vishost <<
':' << visport << endl;
360 visualization =
false;
361 cout <<
"GLVis visualization disabled.\n";
365 sout.precision(precision);
366 sout <<
"solution\n" << mesh << u;
369 cout <<
"GLVis visualization paused."
370 <<
" Press space (in the GLVis window) to resume it.\n";
385 switch (ode_solver_type)
388 case 2: ode_solver =
new RK2Solver(1.0);
break;
390 case 4: ode_solver =
new RK4Solver;
break;
391 case 6: ode_solver =
new RK6Solver;
break;
398 ode_solver = cvode;
break;
405 ode_solver = arkode;
break;
412 ode_solver = arkode;
break;
416 if (ode_solver_type < 7) { ode_solver->
Init(adv); }
421 for (
int ti = 0; !done; )
423 double dt_real = min(dt, t_final - t);
424 ode_solver->
Step(u, t, dt_real);
427 done = (t >= t_final - 1e-8*dt);
429 if (done || ti % vis_steps == 0)
431 cout <<
"time step: " << ti <<
", time: " << t << endl;
437 sout <<
"solution\n" << mesh << u << flush;
459 ofstream osol(
"ex9-final.gf");
460 osol.precision(precision);
478 if (M.GetAssemblyLevel() == AssemblyLevel::LEGACYFULL)
482 dg_solver =
new DG_Solver(M.SpMat(), K.SpMat(), *M.FESpace());
487 M_solver.SetOperator(M);
490 M_solver.SetPreconditioner(*M_prec);
491 M_solver.iterative_mode =
false;
492 M_solver.SetRelTol(1e-9);
493 M_solver.SetAbsTol(0.0);
494 M_solver.SetMaxIter(100);
495 M_solver.SetPrintLevel(0);
508 MFEM_VERIFY(dg_solver != NULL,
509 "Implicit time integration is not supported with partial assembly");
512 dg_solver->SetTimeStep(dt);
513 dg_solver->Mult(z, k);
529 for (
int i = 0; i <
dim; i++)
542 case 1: v(0) = 1.0;
break;
543 case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.);
break;
544 case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
553 const double w = M_PI/2;
556 case 1: v(0) = 1.0;
break;
557 case 2: v(0) = w*X(1); v(1) = -w*X(0);
break;
558 case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0;
break;
565 const double w = M_PI/2;
566 double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
570 case 1: v(0) = 1.0;
break;
571 case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0);
break;
572 case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0;
break;
586 for (
int i = 0; i <
dim; i++)
600 return exp(-40.*pow(X(0)-0.5,2));
604 double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
607 const double s = (1. + 0.25*cos(2*M_PI*X(2)));
611 return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
612 erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
618 double x_ = X(0), y_ = X(1), rho, phi;
621 return pow(sin(M_PI*rho),2)*sin(3*phi);
625 const double f = M_PI;
626 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.
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
virtual void Print(std::ostream &out=mfem::out) const
Conjugate gradient method.
Class for grid function - Vector with associated FE space.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
Data type for scaled Jacobi-type smoother of sparse matrix.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
void SetDataFormat(VTKFormat fmt)
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Helper class for ParaView visualization data.
Base abstract class for first order 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)
Returns the minimum and maximum corners of the mesh bounding box.
virtual void SetOperator(const Operator &a)
Set/update the solver for the given operator.
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.
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)
void SetMaxStep(double dt_max)
Set the maximum time step.
void SetOrder(int order)
Chooses integration order for all explicit / implicit / IMEX methods.
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
int main(int argc, char *argv[])
FE_Evolution(FiniteElementSpace &_vfes, Operator &_A, SparseMatrix &_Aflux)
Interface to ARKode's ARKStep module – additive Runge-Kutta methods.
Data collection with Sidre routines following the Conduit mesh blueprint specification.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
virtual void Save()
Save the collection to disk.
Interface to the CVODE library – linear multi-step methods.
Jacobi smoothing for a given bilinear form (no matrix necessary).
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Data collection with VisIt I/O routines.
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
void SetMaxStep(double dt_max)
Set the maximum time step.
void SetHighOrderOutput(bool high_order_output_)
void PrintUsage(std::ostream &out) const
Print the usage message.
void SetTime(double t)
Set physical time (for time-dependent simulations)
A general vector function coefficient.
The classical explicit forth-order Runge-Kutta method, RK4.
virtual void ImplicitSolve(const double dt, const Vector &x, Vector &k)
Solve the equation: k = f(x + dt k, t), for the unknown k at the current time t.
void velocity_function(const Vector &x, Vector &v)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
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...
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 PrintInfo() const
Print various ARKStep statistics.
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
void PrintOptions(std::ostream &out) const
Print the options.
virtual void ProjectCoefficient(Coefficient &coeff)
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
double u0_function(const Vector &x)
void SetLevelsOfDetail(int levels_of_detail_)
A general function coefficient.
virtual void Save() override
void PrintInfo() const
Print various CVODE statistics.
void SetERKTableNum(int table_num)
Choose a specific Butcher table for an explicit RK method.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
The classical forward Euler method.
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults...
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
double inflow_function(const Vector &x)
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
Arbitrary order "L2-conforming" discontinuous finite elements.
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.