77#if MFEM_HYPRE_VERSION >= 21800
82class AIR_prec :
public Solver
93 AIR_prec(
int blocksize_) : AIR_solver(NULL), blocksize(blocksize_) { }
101 MFEM_VERIFY(A != NULL,
"AIR_prec requires a HypreParMatrix.")
108 AIR_solver->SetAdvectiveOptions(1, "", "FA");
109 AIR_solver->SetPrintLevel(0);
110 AIR_solver->SetMaxLevels(50);
118 BlockInverseScaleJob::RHS_ONLY);
119 AIR_solver->Mult(z_s, y);
130class DG_Solver :
public Solver
145 linear_solver(M.GetComm()),
152 BlockILU::Reordering::MINIMUM_DISCARDED_FILL);
156#if MFEM_HYPRE_VERSION >= 21800
157 prec =
new AIR_prec(block_size);
159 MFEM_ABORT(
"Must have MFEM_HYPRE_VERSION >= 21800 to use AIR.\n");
172 void SetTimeStep(
real_t dt_)
179 A =
Add(-dt, K, 0.0, K);
182 A_diag.
Add(1.0, M_diag);
188 void SetOperator(
const Operator &op)
195 linear_solver.
Mult(x, y);
218 DG_Solver *dg_solver;
229 virtual ~FE_Evolution();
233int main(
int argc,
char *argv[])
243 const char *mesh_file =
"../data/periodic-hexagon.mesh";
244 int ser_ref_levels = 2;
245 int par_ref_levels = 0;
250 const char *device_config =
"cpu";
251 int ode_solver_type = 4;
254 bool visualization =
true;
256 bool paraview =
false;
260#if MFEM_HYPRE_VERSION >= 21800
266 cout.precision(precision);
269 args.
AddOption(&mesh_file,
"-m",
"--mesh",
270 "Mesh file to use.");
272 "Problem setup to use. See options in velocity_function().");
273 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
274 "Number of times to refine the mesh uniformly in serial.");
275 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
276 "Number of times to refine the mesh uniformly in parallel.");
278 "Order (degree) of the finite elements.");
279 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
280 "--no-partial-assembly",
"Enable Partial Assembly.");
281 args.
AddOption(&ea,
"-ea",
"--element-assembly",
"-no-ea",
282 "--no-element-assembly",
"Enable Element Assembly.");
283 args.
AddOption(&fa,
"-fa",
"--full-assembly",
"-no-fa",
284 "--no-full-assembly",
"Enable Full Assembly.");
285 args.
AddOption(&device_config,
"-d",
"--device",
286 "Device configuration string, see Device::Configure().");
287 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
288 "ODE solver: 1 - Forward Euler,\n\t"
289 " 2 - RK2 SSP, 3 - RK3 SSP, 4 - RK4, 6 - RK6,\n\t"
290 " 11 - Backward Euler,\n\t"
291 " 12 - SDIRK23 (L-stable), 13 - SDIRK33,\n\t"
292 " 22 - Implicit Midpoint Method,\n\t"
293 " 23 - SDIRK23 (A-stable), 24 - SDIRK34");
294 args.
AddOption(&t_final,
"-tf",
"--t-final",
295 "Final time; start time is 0.");
296 args.
AddOption(&dt,
"-dt",
"--time-step",
298 args.
AddOption((
int *)&prec_type,
"-pt",
"--prec-type",
"Preconditioner for "
299 "implicit solves. 0 for ILU, 1 for pAIR-AMG.");
300 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
301 "--no-visualization",
302 "Enable or disable GLVis visualization.");
303 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
304 "--no-visit-datafiles",
305 "Save data files for VisIt (visit.llnl.gov) visualization.");
306 args.
AddOption(¶view,
"-paraview",
"--paraview-datafiles",
"-no-paraview",
307 "--no-paraview-datafiles",
308 "Save data files for ParaView (paraview.org) visualization.");
309 args.
AddOption(&adios2,
"-adios2",
"--adios2-streams",
"-no-adios2",
310 "--no-adios2-streams",
311 "Save data using adios2 streams.");
312 args.
AddOption(&binary,
"-binary",
"--binary-datafiles",
"-ascii",
314 "Use binary (Sidre) or ascii format for VisIt data files.");
315 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
316 "Visualize every n-th timestep.");
331 Device device(device_config);
336 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
342 switch (ode_solver_type)
346 case 2: ode_solver =
new RK2Solver(1.0);
break;
348 case 4: ode_solver =
new RK4Solver;
break;
349 case 6: ode_solver =
new RK6Solver;
break;
361 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
371 for (
int lev = 0; lev < ser_ref_levels; lev++)
386 for (
int lev = 0; lev < par_ref_levels; lev++)
399 cout <<
"Number of unknowns: " << global_vSize << endl;
436 b->AddBdrFaceIntegrator(
453 u->ProjectCoefficient(u0);
457 ostringstream mesh_name, sol_name;
458 mesh_name <<
"ex9-mesh." << setfill(
'0') << setw(6) << myid;
459 sol_name <<
"ex9-init." << setfill(
'0') << setw(6) << myid;
460 ofstream omesh(mesh_name.str().c_str());
461 omesh.precision(precision);
463 ofstream osol(sol_name.str().c_str());
464 osol.precision(precision);
478 MFEM_ABORT(
"Must build with MFEM_USE_SIDRE=YES for binary output.");
510#ifdef MFEM_USE_ADIOS2
514 std::string postfix(mesh_file);
515 postfix.erase(0, std::string(
"../data/").size() );
516 postfix +=
"_o" + std::to_string(order);
517 const std::string collection_name =
"ex9-p-" + postfix +
".bp";
521 adios2_dc->
SetParameter(
"SubStreams", std::to_string(num_procs/2) );
540 cout <<
"Unable to connect to GLVis server at "
543 visualization =
false;
546 cout <<
"GLVis visualization disabled.\n";
551 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
552 sout.precision(precision);
553 sout <<
"solution\n" << *pmesh << *
u;
558 cout <<
"GLVis visualization paused."
559 <<
" Press space (in the GLVis window) to resume it.\n";
567 FE_Evolution adv(*m, *k, *B, prec_type);
571 ode_solver->
Init(adv);
574 for (
int ti = 0; !done; )
576 real_t dt_real = min(dt, t_final -
t);
577 ode_solver->
Step(*U,
t, dt_real);
580 done = (
t >= t_final - 1e-8*dt);
582 if (done || ti % vis_steps == 0)
586 cout <<
"time step: " << ti <<
", time: " <<
t << endl;
595 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
596 sout <<
"solution\n" << *pmesh << *
u << flush;
613#ifdef MFEM_USE_ADIOS2
629 ostringstream sol_name;
630 sol_name <<
"ex9-final." << setfill(
'0') << setw(6) << myid;
631 ofstream osol(sol_name.str().c_str());
632 osol.precision(precision);
647#ifdef MFEM_USE_ADIOS2
663 M_solver(M_.ParFESpace()->GetComm()),
677 M_solver.SetOperator(*M);
687 dg_solver =
new DG_Solver(M_mat, K_mat, *M_.
FESpace(), prec_type);
695 M_solver.SetPreconditioner(*M_prec);
696 M_solver.iterative_mode =
false;
697 M_solver.SetRelTol(1e-9);
698 M_solver.SetAbsTol(0.0);
699 M_solver.SetMaxIter(100);
700 M_solver.SetPrintLevel(0);
711 dg_solver->SetTimeStep(dt);
712 dg_solver->Mult(z, k);
715void FE_Evolution::Mult(
const Vector &x,
Vector &y)
const
723FE_Evolution::~FE_Evolution()
737 for (
int i = 0; i <
dim; i++)
750 case 1: v(0) = 1.0;
break;
751 case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.);
break;
752 case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
764 case 1: v(0) = 1.0;
break;
765 case 2: v(0) = w*X(1); v(1) = -w*X(0);
break;
766 case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0;
break;
774 real_t d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
778 case 1: v(0) = 1.0;
break;
779 case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0);
break;
780 case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0;
break;
794 for (
int i = 0; i <
dim; i++)
808 return exp(-40.*pow(X(0)-0.5,2));
812 real_t rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
815 const real_t s = (1. + 0.25*cos(2*M_PI*X(2)));
819 return ( std::erfc(w*(X(0)-cx-rx))*std::erfc(-w*(X(0)-cx+rx)) *
820 std::erfc(w*(X(1)-cy-ry))*std::erfc(-w*(X(1)-cy+ry)) )/16;
826 real_t x_ = X(0), y_ = X(1), rho, phi;
827 rho = std::hypot(x_, y_);
829 return pow(sin(M_PI*rho),2)*sin(3*phi);
834 return sin(
f*X(0))*sin(
f*X(1));
void SetParameter(const std::string key, const std::string value) noexcept
Backward Euler ODE solver. L-stable.
@ 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.
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(),...
Implicit midpoint method. A-stable, not L-stable.
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'.
void velocity_function(const Vector &x, Vector &v)
real_t inflow_function(const Vector &x)
real_t u0_function(const Vector &x)
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)