43#ifndef MFEM_USE_SUNDIALS
44#error This example requires that MFEM is built with MFEM_USE_SUNDIALS=YES
73#if MFEM_HYPRE_VERSION >= 21800
78class AIR_prec :
public Solver
89 AIR_prec(
int blocksize_) : AIR_solver(NULL), blocksize(blocksize_) { }
97 MFEM_VERIFY(A != NULL,
"AIR_prec requires a HypreParMatrix.")
104 AIR_solver->SetAdvectiveOptions(1, "", "FA");
105 AIR_solver->SetPrintLevel(0);
106 AIR_solver->SetMaxLevels(50);
114 BlockInverseScaleJob::RHS_ONLY);
115 AIR_solver->Mult(z_s, y);
126class DG_Solver :
public Solver
141 linear_solver(M.GetComm()),
148 BlockILU::Reordering::MINIMUM_DISCARDED_FILL);
152#if MFEM_HYPRE_VERSION >= 21800
153 prec =
new AIR_prec(block_size);
155 MFEM_ABORT(
"Must have MFEM_HYPRE_VERSION >= 21800 to use AIR.\n");
168 void SetTimeStep(
double dt_)
175 A =
Add(-dt, K, 0.0, K);
178 A_diag.
Add(1.0, M_diag);
184 void SetOperator(
const Operator &op)
191 linear_solver.
Mult(x, y);
214 DG_Solver *dg_solver;
223 virtual void ImplicitSolve(
const double dt,
const Vector &x,
Vector &k);
225 virtual ~FE_Evolution();
229int main(
int argc,
char *argv[])
240 const char *mesh_file =
"../../data/periodic-hexagon.mesh";
241 int ser_ref_levels = 2;
242 int par_ref_levels = 0;
247 const char *device_config =
"cpu";
248 int ode_solver_type = 7;
249 double t_final = 10.0;
251 bool visualization =
false;
253 bool paraview =
false;
257#if MFEM_HYPRE_VERSION >= 21800
264 const double reltol = 1e-2, abstol = 1e-2;
267 cout.precision(precision);
270 args.
AddOption(&mesh_file,
"-m",
"--mesh",
271 "Mesh file to use.");
273 "Problem setup to use. See options in velocity_function().");
274 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
275 "Number of times to refine the mesh uniformly in serial.");
276 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
277 "Number of times to refine the mesh uniformly in parallel.");
279 "Order (degree) of the finite elements.");
280 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
281 "--no-partial-assembly",
"Enable Partial Assembly.");
282 args.
AddOption(&ea,
"-ea",
"--element-assembly",
"-no-ea",
283 "--no-element-assembly",
"Enable Element Assembly.");
284 args.
AddOption(&fa,
"-fa",
"--full-assembly",
"-no-fa",
285 "--no-full-assembly",
"Enable Full Assembly.");
286 args.
AddOption(&device_config,
"-d",
"--device",
287 "Device configuration string, see Device::Configure().");
288 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
290 "1 - Forward Euler,\n\t"
295 "7 - CVODE (adaptive order implicit Adams),\n\t"
296 "8 - ARKODE default (4th order) explicit,\n\t"
298 args.
AddOption(&t_final,
"-tf",
"--t-final",
299 "Final time; start time is 0.");
300 args.
AddOption(&dt,
"-dt",
"--time-step",
302 args.
AddOption((
int *)&prec_type,
"-pt",
"--prec-type",
"Preconditioner for "
303 "implicit solves. 0 for ILU, 1 for pAIR-AMG.");
304 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
305 "--no-visualization",
306 "Enable or disable GLVis visualization.");
307 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
308 "--no-visit-datafiles",
309 "Save data files for VisIt (visit.llnl.gov) visualization.");
310 args.
AddOption(¶view,
"-paraview",
"--paraview-datafiles",
"-no-paraview",
311 "--no-paraview-datafiles",
312 "Save data files for ParaView (paraview.org) visualization.");
313 args.
AddOption(&adios2,
"-adios2",
"--adios2-streams",
"-no-adios2",
314 "--no-adios2-streams",
315 "Save data using adios2 streams.");
316 args.
AddOption(&binary,
"-binary",
"--binary-datafiles",
"-ascii",
318 "Use binary (Sidre) or ascii format for VisIt data files.");
319 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
320 "Visualize every n-th timestep.");
336 if (ode_solver_type < 1 || ode_solver_type > 9)
340 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
345 Device device(device_config);
350 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
357 for (
int lev = 0; lev < ser_ref_levels; lev++)
372 for (
int lev = 0; lev < par_ref_levels; lev++)
385 cout <<
"Number of unknowns: " << global_vSize << endl;
414 constexpr double alpha = -1.0;
422 b->AddBdrFaceIntegrator(
438 u->ProjectCoefficient(u0);
442 ostringstream mesh_name, sol_name;
443 mesh_name <<
"ex9-mesh." << setfill(
'0') << setw(6) << myid;
444 sol_name <<
"ex9-init." << setfill(
'0') << setw(6) << myid;
445 ofstream omesh(mesh_name.str().c_str());
446 omesh.precision(precision);
448 ofstream osol(sol_name.str().c_str());
449 osol.precision(precision);
463 MFEM_ABORT(
"Must build with MFEM_USE_SIDRE=YES for binary output.");
495#ifdef MFEM_USE_ADIOS2
499 std::string postfix(mesh_file);
500 postfix.erase(0, std::string(
"../data/").size() );
501 postfix +=
"_o" + std::to_string(order);
502 const std::string collection_name =
"ex9-p-" + postfix +
".bp";
506 adios2_dc->
SetParameter(
"SubStreams", std::to_string(num_procs/2) );
524 cout <<
"Unable to connect to GLVis server at "
527 visualization =
false;
530 cout <<
"GLVis visualization disabled.\n";
535 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
536 sout.precision(precision);
537 sout <<
"solution\n" << *pmesh << *
u;
542 cout <<
"GLVis visualization paused."
543 <<
" Press space (in the GLVis window) to resume it.\n";
550 FE_Evolution adv(*m, *k, *B, prec_type);
559 switch (ode_solver_type)
562 case 2: ode_solver =
new RK2Solver(1.0);
break;
564 case 4: ode_solver =
new RK4Solver;
break;
565 case 6: ode_solver =
new RK6Solver;
break;
572 ode_solver = cvode;
break;
579 if (ode_solver_type == 9)
583 ode_solver = arkode;
break;
587 if (ode_solver_type < 7) { ode_solver->
Init(adv); }
592 for (
int ti = 0; !done; )
594 double dt_real = min(dt, t_final -
t);
595 ode_solver->
Step(*U,
t, dt_real);
598 done = (
t >= t_final - 1e-8*dt);
600 if (done || ti % vis_steps == 0)
604 cout <<
"time step: " << ti <<
", time: " <<
t << endl;
615 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
616 sout <<
"solution\n" << *pmesh << *
u << flush;
633#ifdef MFEM_USE_ADIOS2
649 ostringstream sol_name;
650 sol_name <<
"ex9-final." << setfill(
'0') << setw(6) << myid;
651 ofstream osol(sol_name.str().c_str());
652 osol.precision(precision);
667#ifdef MFEM_USE_ADIOS2
684 M_solver(M_.ParFESpace()->GetComm()),
698 M_solver.SetOperator(*M);
708 dg_solver =
new DG_Solver(M_mat, K_mat, *M_.
FESpace(), prec_type);
716 M_solver.SetPreconditioner(*M_prec);
717 M_solver.iterative_mode =
false;
718 M_solver.SetRelTol(1e-9);
719 M_solver.SetAbsTol(0.0);
720 M_solver.SetMaxIter(100);
721 M_solver.SetPrintLevel(0);
728void FE_Evolution::ImplicitSolve(
const double dt,
const Vector &x,
Vector &k)
732 dg_solver->SetTimeStep(dt);
733 dg_solver->Mult(z, k);
736void FE_Evolution::Mult(
const Vector &x,
Vector &y)
const
744FE_Evolution::~FE_Evolution()
758 for (
int i = 0; i <
dim; i++)
771 case 1: v(0) = 1.0;
break;
772 case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.);
break;
773 case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
782 const double w = M_PI/2;
785 case 1: v(0) = 1.0;
break;
786 case 2: v(0) = w*X(1); v(1) = -w*X(0);
break;
787 case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0;
break;
794 const double w = M_PI/2;
795 double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
799 case 1: v(0) = 1.0;
break;
800 case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0);
break;
801 case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0;
break;
815 for (
int i = 0; i <
dim; i++)
829 return exp(-40.*pow(X(0)-0.5,2));
833 double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
836 const double s = (1. + 0.25*cos(2*M_PI*X(2)));
840 return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
841 erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
847 double x_ = X(0), y_ = X(1), rho, phi;
850 return pow(sin(M_PI*rho),2)*sin(3*phi);
854 const double f = M_PI;
855 return sin(
f*X(0))*sin(
f*X(1));
void SetParameter(const std::string key, const std::string value) noexcept
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 Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults.
@ EXPLICIT
Explicit RK method.
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.
virtual void Mult(const Vector &b, Vector &x) const
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 SetMaxStep(double dt_max)
Set the maximum time step.
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
void PrintInfo() const
Print various CVODE statistics.
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
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...
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
The classical forward Euler method.
A general function coefficient.
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the GMRES method.
The BoomerAMG solver in hypre.
Wrapper for hypre's ParCSR matrix class.
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, real_t alpha=1.0, real_t beta=0.0) const
Computes y = alpha * A * x + beta * y.
Wrapper for hypre's parallel vector class.
Parallel smoothers in hypre.
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
virtual 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.
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 *)
static bool Root()
Return true if the rank in MPI_COMM_WORLD is zero.
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
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].
Pointer to an Operator of a specified type.
Jacobi smoothing for a given bilinear form (no matrix necessary).
int width
Dimension of the input / number of columns in the matrix.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
int height
Dimension of the output / number of rows in the matrix.
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
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.
Abstract parallel finite element space.
HYPRE_BigInt GlobalTrueVSize() const
Class for parallel grid function.
Class for parallel meshes.
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Helper class for ParaView visualization data.
void SetHighOrderOutput(bool high_order_output_)
void SetLevelsOfDetail(int levels_of_detail_)
virtual void Save() override
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.
void Add(const int i, const int j, const real_t val)
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)
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
void BlockInverseScale(const HypreParMatrix *A, HypreParMatrix *C, const Vector *b, HypreParVector *d, int blocksize, BlockInverseScaleJob job)
std::function< real_t(const Vector &)> f(real_t mass_coeff)
void Add(const DenseMatrix &A, const DenseMatrix &B, real_t alpha, DenseMatrix &C)
C = A + alpha*B.
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