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
80 linear_solver(M.GetComm()),
82 BlockILU::Reordering::MINIMUM_DISCARDED_FILL),
85 linear_solver.iterative_mode =
false;
86 linear_solver.SetRelTol(1e-9);
87 linear_solver.SetAbsTol(0.0);
88 linear_solver.SetMaxIter(100);
89 linear_solver.SetPrintLevel(0);
90 linear_solver.SetPreconditioner(prec);
95 void SetTimeStep(
double dt_)
102 A =
Add(-dt, K, 0.0, K);
105 A_diag.
Add(1.0, M_diag);
107 linear_solver.SetOperator(*A);
111 void SetOperator(
const Operator &op)
113 linear_solver.SetOperator(op);
118 linear_solver.Mult(x, y);
139 DG_Solver *dg_solver;
147 virtual void ImplicitSolve(
const double dt,
const Vector &x,
Vector &k);
153 int main(
int argc,
char *argv[])
156 Mpi::Init(argc, argv);
157 int num_procs = Mpi::WorldSize();
158 int myid = Mpi::WorldRank();
164 const char *mesh_file =
"../../data/periodic-hexagon.mesh";
165 int ser_ref_levels = 2;
166 int par_ref_levels = 0;
171 const char *device_config =
"cpu";
172 int ode_solver_type = 7;
173 double t_final = 10.0;
175 bool visualization =
false;
177 bool paraview =
false;
183 const double reltol = 1e-2, abstol = 1e-2;
186 cout.precision(precision);
189 args.
AddOption(&mesh_file,
"-m",
"--mesh",
190 "Mesh file to use.");
192 "Problem setup to use. See options in velocity_function().");
193 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
194 "Number of times to refine the mesh uniformly in serial.");
195 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
196 "Number of times to refine the mesh uniformly in parallel.");
198 "Order (degree) of the finite elements.");
199 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
200 "--no-partial-assembly",
"Enable Partial Assembly.");
201 args.
AddOption(&ea,
"-ea",
"--element-assembly",
"-no-ea",
202 "--no-element-assembly",
"Enable Element Assembly.");
203 args.
AddOption(&fa,
"-fa",
"--full-assembly",
"-no-fa",
204 "--no-full-assembly",
"Enable Full Assembly.");
205 args.
AddOption(&device_config,
"-d",
"--device",
206 "Device configuration string, see Device::Configure().");
207 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
209 "1 - Forward Euler,\n\t"
214 "7 - CVODE (adaptive order implicit Adams),\n\t"
215 "8 - ARKODE default (4th order) explicit,\n\t"
217 args.
AddOption(&t_final,
"-tf",
"--t-final",
218 "Final time; start time is 0.");
219 args.
AddOption(&dt,
"-dt",
"--time-step",
221 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
222 "--no-visualization",
223 "Enable or disable GLVis visualization.");
224 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
225 "--no-visit-datafiles",
226 "Save data files for VisIt (visit.llnl.gov) visualization.");
227 args.
AddOption(¶view,
"-paraview",
"--paraview-datafiles",
"-no-paraview",
228 "--no-paraview-datafiles",
229 "Save data files for ParaView (paraview.org) visualization.");
230 args.
AddOption(&adios2,
"-adios2",
"--adios2-streams",
"-no-adios2",
231 "--no-adios2-streams",
232 "Save data using adios2 streams.");
233 args.
AddOption(&binary,
"-binary",
"--binary-datafiles",
"-ascii",
235 "Use binary (Sidre) or ascii format for VisIt data files.");
236 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
237 "Visualize every n-th timestep.");
253 if (ode_solver_type < 1 || ode_solver_type > 9)
257 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
262 Device device(device_config);
263 if (myid == 0) { device.
Print(); }
267 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
274 for (
int lev = 0; lev < ser_ref_levels; lev++)
289 for (
int lev = 0; lev < par_ref_levels; lev++)
302 cout <<
"Number of unknowns: " << global_vSize << endl;
358 ostringstream mesh_name, sol_name;
359 mesh_name <<
"ex9-mesh." << setfill(
'0') << setw(6) << myid;
360 sol_name <<
"ex9-init." << setfill(
'0') << setw(6) << myid;
361 ofstream omesh(mesh_name.str().c_str());
362 omesh.precision(precision);
364 ofstream osol(sol_name.str().c_str());
365 osol.precision(precision);
376 #ifdef MFEM_USE_SIDRE
379 MFEM_ABORT(
"Must build with MFEM_USE_SIDRE=YES for binary output.");
411 #ifdef MFEM_USE_ADIOS2
415 std::string postfix(mesh_file);
416 postfix.erase(0, std::string(
"../data/").size() );
417 postfix +=
"_o" + std::to_string(order);
418 const std::string collection_name =
"ex9-p-" + postfix +
".bp";
422 adios2_dc->
SetParameter(
"SubStreams", std::to_string(num_procs/2) );
435 sout.
open(vishost, visport);
439 cout <<
"Unable to connect to GLVis server at "
440 << vishost <<
':' << visport << endl;
441 visualization =
false;
444 cout <<
"GLVis visualization disabled.\n";
449 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
450 sout.precision(precision);
451 sout <<
"solution\n" << *pmesh << *
u;
455 cout <<
"GLVis visualization paused."
456 <<
" Press space (in the GLVis window) to resume it.\n";
471 switch (ode_solver_type)
474 case 2: ode_solver =
new RK2Solver(1.0);
break;
476 case 4: ode_solver =
new RK4Solver;
break;
477 case 6: ode_solver =
new RK6Solver;
break;
484 ode_solver = cvode;
break;
487 arkode =
new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::EXPLICIT);
491 if (ode_solver_type == 9)
495 ode_solver = arkode;
break;
499 if (ode_solver_type < 7) { ode_solver->
Init(adv); }
504 for (
int ti = 0; !done; )
506 double dt_real = min(dt, t_final - t);
507 ode_solver->
Step(*U, t, dt_real);
510 done = (t >= t_final - 1e-8*dt);
512 if (done || ti % vis_steps == 0)
516 cout <<
"time step: " << ti <<
", time: " << t << endl;
527 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
528 sout <<
"solution\n" << *pmesh << *u << flush;
545 #ifdef MFEM_USE_ADIOS2
561 ostringstream sol_name;
562 sol_name <<
"ex9-final." << setfill(
'0') << setw(6) << myid;
563 ofstream osol(sol_name.str().c_str());
564 osol.precision(precision);
579 #ifdef MFEM_USE_ADIOS2
596 M_solver(M_.ParFESpace()->GetComm()),
620 dg_solver =
new DG_Solver(M_mat, K_mat, *M_.
FESpace());
640 dg_solver->SetTimeStep(dt);
641 dg_solver->Mult(z, k);
666 for (
int i = 0; i <
dim; i++)
679 case 1: v(0) = 1.0;
break;
680 case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.);
break;
681 case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
690 const double w = M_PI/2;
693 case 1: v(0) = 1.0;
break;
694 case 2: v(0) = w*X(1); v(1) = -w*X(0);
break;
695 case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0;
break;
702 const double w = M_PI/2;
703 double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
707 case 1: v(0) = 1.0;
break;
708 case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0);
break;
709 case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0;
break;
723 for (
int i = 0; i <
dim; i++)
737 return exp(-40.*pow(X(0)-0.5,2));
741 double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
744 const double s = (1. + 0.25*cos(2*M_PI*X(2)));
748 return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
749 erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
755 double x_ = X(0), y_ = X(1), rho, phi;
758 return pow(sin(M_PI*rho),2)*sin(3*phi);
762 const double f = M_PI;
763 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.
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.
void SetDataFormat(VTKFormat fmt)
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Helper class for ParaView visualization data.
void SetERKTableNum(ARKODE_ERKTableID table_id)
Choose a specific Butcher table for an explicit RK method.
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.
Pointer to an Operator of a specified type.
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]...
HYPRE_BigInt GlobalTrueVSize() const
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
void SetMaxStep(double dt_max)
Set the maximum time step.
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
virtual void SetTime(const double t_)
Set the current time.
Abstract parallel finite element space.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Interface to ARKode's ARKStep module – additive Runge-Kutta methods.
Data collection with Sidre routines following the Conduit mesh blueprint specification.
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
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.
void SetMaxIter(int max_it)
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 SetParameter(const std::string key, const std::string value) noexcept
Parallel smoothers in hypre.
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)
void Add(const int i, const int j, const double val)
A general vector function coefficient.
Wrapper for hypre's parallel vector class.
The classical explicit forth-order Runge-Kutta method, RK4.
void SetAbsTol(double atol)
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 SetRelTol(double rtol)
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.
constexpr ARKODE_ERKTableID ARKODE_FEHLBERG_13_7_8
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...
FE_Evolution(FiniteElementSpace &vfes_, Operator &A_, SparseMatrix &Aflux_)
NURBSExtension * NURBSext
Optional NURBS mesh extension.
void GetDiag(Vector &d) const
Returns the Diagonal of A.
void PrintInfo() const
Print various ARKStep statistics.
virtual 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.
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
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_)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
A general function coefficient.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
virtual void Save() override
void PrintInfo() const
Print various CVODE statistics.
void Print(std::ostream &out=mfem::out) const override
double u(const Vector &xvec)
Class for parallel grid function.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
The classical forward Euler method.
Wrapper for hypre's ParCSR matrix class.
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults...
Class for parallel meshes.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
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.