45#ifndef MFEM_USE_SUNDIALS
46#error This example requires that MFEM is built with MFEM_USE_SUNDIALS=YES
75#if MFEM_HYPRE_VERSION >= 21800
80class AIR_prec :
public Solver
91 AIR_prec(
int blocksize_) : AIR_solver(NULL), blocksize(blocksize_) { }
99 MFEM_VERIFY(A != NULL,
"AIR_prec requires a HypreParMatrix.")
106 AIR_solver->SetAdvectiveOptions(1, "", "FA");
107 AIR_solver->SetPrintLevel(0);
108 AIR_solver->SetMaxLevels(50);
116 BlockInverseScaleJob::RHS_ONLY);
117 AIR_solver->Mult(z_s, y);
128class DG_Solver :
public Solver
143 linear_solver(M.GetComm()),
150 BlockILU::Reordering::MINIMUM_DISCARDED_FILL);
154#if MFEM_HYPRE_VERSION >= 21800
155 prec =
new AIR_prec(block_size);
157 MFEM_ABORT(
"Must have MFEM_HYPRE_VERSION >= 21800 to use AIR.\n");
170 void SetTimeStep(
double dt_)
177 A =
Add(-dt, K, 0.0, K);
180 A_diag.
Add(1.0, M_diag);
186 void SetOperator(
const Operator &op)
193 linear_solver.
Mult(x, y);
216 DG_Solver *dg_solver;
225 virtual void ImplicitSolve(
const double dt,
const Vector &x,
Vector &k);
227 virtual ~FE_Evolution();
231int main(
int argc,
char *argv[])
242 const char *mesh_file =
"../../data/periodic-hexagon.mesh";
243 int ser_ref_levels = 2;
244 int par_ref_levels = 0;
249 const char *device_config =
"cpu";
250 int ode_solver_type = 7;
251 double t_final = 10.0;
253 bool visualization =
false;
255 bool paraview =
false;
259#if MFEM_HYPRE_VERSION >= 21800
266 const double reltol = 1e-2, abstol = 1e-2;
269 cout.precision(precision);
272 args.
AddOption(&mesh_file,
"-m",
"--mesh",
273 "Mesh file to use.");
275 "Problem setup to use. See options in velocity_function().");
276 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
277 "Number of times to refine the mesh uniformly in serial.");
278 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
279 "Number of times to refine the mesh uniformly in parallel.");
281 "Order (degree) of the finite elements.");
282 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
283 "--no-partial-assembly",
"Enable Partial Assembly.");
284 args.
AddOption(&ea,
"-ea",
"--element-assembly",
"-no-ea",
285 "--no-element-assembly",
"Enable Element Assembly.");
286 args.
AddOption(&fa,
"-fa",
"--full-assembly",
"-no-fa",
287 "--no-full-assembly",
"Enable Full Assembly.");
288 args.
AddOption(&device_config,
"-d",
"--device",
289 "Device configuration string, see Device::Configure().");
290 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
292 "1 - Forward Euler,\n\t"
297 "7 - CVODE (adaptive order implicit Adams),\n\t"
298 "8 - ARKODE default (4th order) explicit,\n\t"
300 args.
AddOption(&t_final,
"-tf",
"--t-final",
301 "Final time; start time is 0.");
302 args.
AddOption(&dt,
"-dt",
"--time-step",
304 args.
AddOption((
int *)&prec_type,
"-pt",
"--prec-type",
"Preconditioner for "
305 "implicit solves. 0 for ILU, 1 for pAIR-AMG.");
306 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
307 "--no-visualization",
308 "Enable or disable GLVis visualization.");
309 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
310 "--no-visit-datafiles",
311 "Save data files for VisIt (visit.llnl.gov) visualization.");
312 args.
AddOption(¶view,
"-paraview",
"--paraview-datafiles",
"-no-paraview",
313 "--no-paraview-datafiles",
314 "Save data files for ParaView (paraview.org) visualization.");
315 args.
AddOption(&adios2,
"-adios2",
"--adios2-streams",
"-no-adios2",
316 "--no-adios2-streams",
317 "Save data using adios2 streams.");
318 args.
AddOption(&binary,
"-binary",
"--binary-datafiles",
"-ascii",
320 "Use binary (Sidre) or ascii format for VisIt data files.");
321 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
322 "Visualize every n-th timestep.");
338 if (ode_solver_type < 1 || ode_solver_type > 9)
342 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
347 Device device(device_config);
352 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
359 for (
int lev = 0; lev < ser_ref_levels; lev++)
374 for (
int lev = 0; lev < par_ref_levels; lev++)
387 cout <<
"Number of unknowns: " << global_vSize << endl;
416 constexpr double alpha = -1.0;
424 b->AddBdrFaceIntegrator(
440 u->ProjectCoefficient(u0);
444 ostringstream mesh_name, sol_name;
445 mesh_name <<
"ex9-mesh." << setfill(
'0') << setw(6) << myid;
446 sol_name <<
"ex9-init." << setfill(
'0') << setw(6) << myid;
447 ofstream omesh(mesh_name.str().c_str());
448 omesh.precision(precision);
450 ofstream osol(sol_name.str().c_str());
451 osol.precision(precision);
465 MFEM_ABORT(
"Must build with MFEM_USE_SIDRE=YES for binary output.");
497#ifdef MFEM_USE_ADIOS2
501 std::string postfix(mesh_file);
502 postfix.erase(0, std::string(
"../data/").size() );
503 postfix +=
"_o" + std::to_string(order);
504 const std::string collection_name =
"ex9-p-" + postfix +
".bp";
508 adios2_dc->
SetParameter(
"SubStreams", std::to_string(num_procs/2) );
526 cout <<
"Unable to connect to GLVis server at "
527 <<
vishost <<
':' << visport << endl;
529 visualization =
false;
532 cout <<
"GLVis visualization disabled.\n";
537 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
538 sout.precision(precision);
539 sout <<
"solution\n" << *pmesh << *
u;
544 cout <<
"GLVis visualization paused."
545 <<
" Press space (in the GLVis window) to resume it.\n";
552 FE_Evolution adv(*m, *k, *B, prec_type);
561 switch (ode_solver_type)
564 case 2: ode_solver =
new RK2Solver(1.0);
break;
566 case 4: ode_solver =
new RK4Solver;
break;
567 case 6: ode_solver =
new RK6Solver;
break;
574 ode_solver = cvode;
break;
581 if (ode_solver_type == 9)
585 ode_solver = arkode;
break;
589 if (ode_solver_type < 7) { ode_solver->
Init(adv); }
594 for (
int ti = 0; !done; )
596 double dt_real = min(dt, t_final - t);
597 ode_solver->
Step(*U, t, dt_real);
600 done = (t >= t_final - 1e-8*dt);
602 if (done || ti % vis_steps == 0)
606 cout <<
"time step: " << ti <<
", time: " << t << endl;
617 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
618 sout <<
"solution\n" << *pmesh << *
u << flush;
635#ifdef MFEM_USE_ADIOS2
651 ostringstream sol_name;
652 sol_name <<
"ex9-final." << setfill(
'0') << setw(6) << myid;
653 ofstream osol(sol_name.str().c_str());
654 osol.precision(precision);
669#ifdef MFEM_USE_ADIOS2
686 M_solver(M_.ParFESpace()->GetComm()),
700 M_solver.SetOperator(*M);
710 dg_solver =
new DG_Solver(M_mat, K_mat, *M_.
FESpace(), prec_type);
718 M_solver.SetPreconditioner(*M_prec);
719 M_solver.iterative_mode =
false;
720 M_solver.SetRelTol(1e-9);
721 M_solver.SetAbsTol(0.0);
722 M_solver.SetMaxIter(100);
723 M_solver.SetPrintLevel(0);
730void FE_Evolution::ImplicitSolve(
const double dt,
const Vector &x,
Vector &k)
734 dg_solver->SetTimeStep(dt);
735 dg_solver->Mult(z, k);
738void FE_Evolution::Mult(
const Vector &x,
Vector &y)
const
746FE_Evolution::~FE_Evolution()
760 for (
int i = 0; i <
dim; i++)
773 case 1: v(0) = 1.0;
break;
774 case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.);
break;
775 case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
784 const double w = M_PI/2;
787 case 1: v(0) = 1.0;
break;
788 case 2: v(0) = w*X(1); v(1) = -w*X(0);
break;
789 case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0;
break;
796 const double w = M_PI/2;
797 double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
801 case 1: v(0) = 1.0;
break;
802 case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0);
break;
803 case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0;
break;
817 for (
int i = 0; i <
dim; i++)
831 return exp(-40.*pow(X(0)-0.5,2));
835 double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
838 const double s = (1. + 0.25*cos(2*M_PI*X(2)));
842 return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
843 erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
849 double x_ = X(0), y_ = X(1), rho, phi;
852 return pow(sin(M_PI*rho),2)*sin(3*phi);
856 const double f = M_PI;
857 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.
@ 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.
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...
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.
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(),...
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_)
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