41 #include "mfem-performance.hpp"
51 template <
int dim>
struct geom_t { };
53 struct geom_t<2> {
static const Geometry::Type value = Geometry::SQUARE; };
55 struct geom_t<3> {
static const Geometry::Type value = Geometry::CUBE; };
65 static const int ir_order = 2*
sol_p+rdim-1;
83 static int run(
Mesh *mesh,
int ser_ref_levels,
int par_ref_levels,
int order,
84 int basis,
bool static_cond,
PCType pc_choice,
bool perf,
85 bool matrix_free,
bool visualization);
88 int main(
int argc,
char *argv[])
91 Mpi::Init(argc, argv);
92 int myid = Mpi::WorldRank();
96 #ifdef MFEM_HPC_EX1_2D
97 const char *mesh_file =
"../../data/star.mesh";
99 const char *mesh_file =
"../../data/fichera.mesh";
101 int ser_ref_levels = -1;
102 int par_ref_levels = 1;
104 const char *basis_type =
"G";
105 bool static_cond =
false;
106 const char *pc =
"lor";
108 bool matrix_free =
true;
109 bool visualization = 1;
112 args.
AddOption(&mesh_file,
"-m",
"--mesh",
113 "Mesh file to use.");
114 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
115 "Number of times to refine the mesh uniformly in serial;"
116 " -1 = auto: <= 10,000 elements.");
117 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
118 "Number of times to refine the mesh uniformly in parallel.");
120 "Finite element order (polynomial degree) or -1 for"
121 " isoparametric space.");
122 args.
AddOption(&basis_type,
"-b",
"--basis-type",
123 "Basis: G - Gauss-Lobatto, P - Positive, U - Uniform");
124 args.
AddOption(&perf,
"-perf",
"--hpc-version",
"-std",
"--standard-version",
125 "Enable high-performance, tensor-based, assembly/evaluation.");
126 args.
AddOption(&matrix_free,
"-mf",
"--matrix-free",
"-asm",
"--assembly",
127 "Use matrix-free evaluation or efficient matrix assembly in "
128 "the high-performance version.");
129 args.
AddOption(&pc,
"-pc",
"--preconditioner",
130 "Preconditioner: lor - low-order-refined (matrix-free) AMG, "
131 "ho - high-order (assembled) AMG, none.");
132 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
133 "--no-static-condensation",
"Enable static condensation.");
134 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
135 "--no-visualization",
136 "Enable or disable GLVis visualization.");
146 if (static_cond && perf && matrix_free)
150 cout <<
"\nStatic condensation can not be used with matrix-free"
151 " evaluation!\n" << endl;
155 MFEM_VERIFY(perf || !matrix_free,
156 "--standard-version is not compatible with --matrix-free");
163 if (!strcmp(pc,
"ho")) { pc_choice =
PCType::HO; }
164 else if (!strcmp(pc,
"lor")) { pc_choice =
PCType::LOR; }
165 else if (!strcmp(pc,
"none")) { pc_choice =
PCType::NONE; }
168 mfem_error(
"Invalid Preconditioner specified");
174 cout <<
"\nMFEM SIMD width: " << MFEM_SIMD_BYTES/
sizeof(double)
175 <<
" doubles\n" << endl;
179 int basis = BasisType::GetType(basis_type[0]);
182 cout <<
"Using " << BasisType::Name(basis) <<
" basis ..." << endl;
188 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
193 return ex1_t<2>::run(mesh, ser_ref_levels, par_ref_levels, order, basis,
194 static_cond, pc_choice, perf, matrix_free,
199 return ex1_t<3>::run(mesh, ser_ref_levels, par_ref_levels, order,
200 basis, static_cond, pc_choice, perf, matrix_free,
205 MFEM_ABORT(
"Dimension must be 2 or 3.")
212 int ex1_t<dim>::run(
Mesh *mesh,
int ser_ref_levels,
int par_ref_levels,
213 int order,
int basis,
bool static_cond,
PCType pc_choice,
214 bool perf,
bool matrix_free,
bool visualization)
217 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
218 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
225 cout <<
"High-performance version using integration rule with "
226 << int_rule_t::qpts <<
" points ..." << endl;
228 if (!mesh_t::MatchesGeometry(*mesh))
232 cout <<
"The given mesh does not match the optimized 'geom' parameter.\n"
233 <<
"Recompile with suitable 'geom' value." << endl;
238 else if (!mesh_t::MatchesNodes(*mesh))
242 cout <<
"Switching the mesh curvature to match the "
243 <<
"optimized value (order " <<
mesh_p <<
") ..." << endl;
255 int ref_levels = (ser_ref_levels != -1) ? ser_ref_levels :
256 (
int)floor(log(10000./mesh->
GetNE())/log(2.)/
dim);
257 for (
int l = 0; l < ref_levels; l++)
267 cout <<
"NURBS mesh: switching the mesh curvature to be "
268 <<
"min(sol_p, mesh_p) = " << new_mesh_p <<
" ..." << endl;
270 mesh->
SetCurvature(new_mesh_p,
false, -1, Ordering::byNODES);
279 for (
int l = 0; l < par_ref_levels; l++)
286 MFEM_VERIFY(pc_choice !=
PCType::LOR,
"triangle and tet meshes do not "
287 "support the LOR preconditioner yet");
303 cout <<
"Using isoparametric FEs: " << fec->
Name() << endl;
314 cout <<
"Number of finite element unknowns: " << size << endl;
322 int basis_lor = basis;
323 if (basis == BasisType::Positive) { basis_lor=BasisType::ClosedUniform; }
324 pmesh_lor = ParMesh::MakeRefined(*pmesh, order, basis_lor);
330 if (perf && !sol_fes_t::Matches(*fespace))
334 cout <<
"The given order does not match the optimized parameter.\n"
335 <<
"Recompile with suitable 'sol_p' value." << endl;
384 "cannot use LOR preconditioner with static condensation");
389 cout <<
"Assembling the matrix ..." << flush;
396 HPCBilinearForm *a_hpc = NULL;
408 a_hpc =
new HPCBilinearForm(integ_t(coeff_t(1.0)), *fespace);
415 a_hpc->AssembleBilinearForm(*a);
430 if (perf && matrix_free)
432 a_hpc->FormLinearSystem(ess_tdof_list, x, *b, a_oper, X, B);
436 cout <<
"Size of linear system: " << glob_size << endl;
445 cout <<
"Size of linear system: " << glob_size << endl;
453 cout <<
"Assembling the preconditioning matrix ..." << flush;
476 a_hpc->AssembleBilinearForm(*a_pc);
514 cout <<
"Time per CG step: "
520 if (perf && matrix_free)
522 a_hpc->RecoverFEMSolution(X, *b, x);
532 ostringstream mesh_name, sol_name;
533 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
534 sol_name <<
"sol." << setfill(
'0') << setw(6) << myid;
536 ofstream mesh_ofs(mesh_name.str().c_str());
537 mesh_ofs.precision(8);
538 pmesh->
Print(mesh_ofs);
540 ofstream sol_ofs(sol_name.str().c_str());
541 sol_ofs.precision(8);
551 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
552 sol_sock.precision(8);
553 sol_sock <<
"solution\n" << *pmesh << x << flush;
559 if (a_oper != &A) {
delete a_oper; }
565 if (order > 0) {
delete fec; }
int Size() const
Return the logical size of the array.
Class for domain integration L(v) := (f, v)
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Conjugate gradient method.
int GetNumIterations() const
A coefficient that is constant across space and time.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
void MakeRef(const HypreParMatrix &master)
Make this HypreParMatrix a reference to 'master'.
HYPRE_BigInt GlobalTrueVSize() const
int GetNE() const
Returns number of elements.
Abstract parallel finite element space.
The Integrator class combines a kernel and a coefficient.
void Stop()
Stop the stopwatch.
The BoomerAMG solver in hypre.
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 ...
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void SetMaxIter(int max_it)
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
int MeshGenerator()
Get the mesh generator/type.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
void PrintUsage(std::ostream &out) const
Print the usage message.
void Start()
Start the stopwatch. The elapsed time is not cleared.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
void SetRelTol(double rtol)
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
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...
NURBSExtension * NURBSext
Optional NURBS mesh extension.
virtual const char * Name() const
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
void PrintOptions(std::ostream &out) const
Print the options.
Abstract class for hypre's solvers and preconditioners.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
void GetNodes(Vector &node_coord) const
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Arbitrary order H1-conforming (continuous) finite elements.
void Print(std::ostream &out=mfem::out) const override
Class for parallel grid function.
Wrapper for hypre's ParCSR matrix class.
Class for parallel meshes.
void Clear()
Clear the elapsed time on the stopwatch and restart it if it's running.
bool Good() const
Return true if the command line options were parsed successfully.