77 #if MFEM_HYPRE_VERSION >= 21800
82 class 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);
130 class 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");
162 linear_solver.iterative_mode =
false;
163 linear_solver.SetRelTol(1e-9);
164 linear_solver.SetAbsTol(0.0);
165 linear_solver.SetMaxIter(100);
166 linear_solver.SetPrintLevel(0);
167 linear_solver.SetPreconditioner(*prec);
172 void SetTimeStep(
double dt_)
179 A =
Add(-dt, K, 0.0, K);
182 A_diag.
Add(1.0, M_diag);
184 linear_solver.SetOperator(*A);
188 void SetOperator(
const Operator &op)
190 linear_solver.SetOperator(op);
195 linear_solver.Mult(x, y);
218 DG_Solver *dg_solver;
227 virtual void ImplicitSolve(
const double dt,
const Vector &x,
Vector &k);
233 int main(
int argc,
char *argv[])
237 int num_procs = Mpi::WorldSize();
238 int myid = Mpi::WorldRank();
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;
252 double t_final = 10.0;
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);
332 if (Mpi::Root()) { device.
Print(); }
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;
428 constexpr
double alpha = -1.0;
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);
475 #ifdef MFEM_USE_SIDRE
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) );
535 sout.
open(vishost, visport);
539 cout <<
"Unable to connect to GLVis server at "
540 << vishost <<
':' << visport << endl;
541 visualization =
false;
544 cout <<
"GLVis visualization disabled.\n";
549 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
550 sout.precision(precision);
551 sout <<
"solution\n" << *pmesh << *
u;
555 cout <<
"GLVis visualization paused."
556 <<
" Press space (in the GLVis window) to resume it.\n";
567 ode_solver->
Init(adv);
570 for (
int ti = 0; !done; )
572 double dt_real = min(dt, t_final - t);
573 ode_solver->
Step(*U, t, dt_real);
576 done = (t >= t_final - 1e-8*dt);
578 if (done || ti % vis_steps == 0)
582 cout <<
"time step: " << ti <<
", time: " << t << endl;
591 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
592 sout <<
"solution\n" << *pmesh << *u << flush;
609 #ifdef MFEM_USE_ADIOS2
625 ostringstream sol_name;
626 sol_name <<
"ex9-final." << setfill(
'0') << setw(6) << myid;
627 ofstream osol(sol_name.str().c_str());
628 osol.precision(precision);
643 #ifdef MFEM_USE_ADIOS2
659 M_solver(M_.ParFESpace()->GetComm()),
683 dg_solver =
new DG_Solver(M_mat, K_mat, *M_.
FESpace(), prec_type);
707 dg_solver->SetTimeStep(dt);
708 dg_solver->Mult(z, k);
733 for (
int i = 0; i <
dim; i++)
746 case 1: v(0) = 1.0;
break;
747 case 2: v(0) =
sqrt(2./3.); v(1) =
sqrt(1./3.);
break;
748 case 3: v(0) =
sqrt(3./6.); v(1) =
sqrt(2./6.); v(2) =
sqrt(1./6.);
757 const double w = M_PI/2;
760 case 1: v(0) = 1.0;
break;
761 case 2: v(0) = w*X(1); v(1) = -w*X(0);
break;
762 case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0;
break;
769 const double w = M_PI/2;
770 double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
774 case 1: v(0) = 1.0;
break;
775 case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0);
break;
776 case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0;
break;
790 for (
int i = 0; i <
dim; i++)
804 return exp(-40.*
pow(X(0)-0.5,2));
808 double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
811 const double s = (1. + 0.25*
cos(2*M_PI*X(2)));
815 return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
816 erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
822 double x_ = X(0), y_ = X(1), rho, phi;
829 const double f = M_PI;
830 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.
Conjugate gradient method.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
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)
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
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 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.
FDualNumber< tbase > exp(const FDualNumber< tbase > &f)
exp([dual number])
void BlockInverseScale(const HypreParMatrix *A, HypreParMatrix *C, const Vector *b, HypreParVector *d, int blocksize, BlockInverseScaleJob job)
Backward Euler ODE solver. L-stable.
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.
The BoomerAMG solver in hypre.
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 ...
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
virtual void Save()
Save the collection to disk.
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 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)
FDualNumber< tbase > cos(const FDualNumber< tbase > &f)
cos([dual number])
FDualNumber< tbase > pow(const FDualNumber< tbase > &a, const FDualNumber< tbase > &b)
pow([dual number],[dual number])
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)
FDualNumber< tbase > sin(const FDualNumber< tbase > &f)
sin([dual number])
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...
Implicit midpoint method. A-stable, not L-stable.
FE_Evolution(FiniteElementSpace &vfes_, Operator &A_, SparseMatrix &Aflux_)
NURBSExtension * NURBSext
Optional NURBS mesh extension.
void GetDiag(Vector &d) const
Returns the Diagonal of A.
FDualNumber< tbase > sqrt(const FDualNumber< tbase > &f)
sqrt([dual number])
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 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.
Class for parallel meshes.
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)
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.