41#include "mfem-performance.hpp"
51template <
int dim>
struct geom_t { };
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,
int visport = 19916);
88int main(
int argc,
char *argv[])
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;
110 bool visualization = 1;
113 args.
AddOption(&mesh_file,
"-m",
"--mesh",
114 "Mesh file to use.");
115 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
116 "Number of times to refine the mesh uniformly in serial;"
117 " -1 = auto: <= 10,000 elements.");
118 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
119 "Number of times to refine the mesh uniformly in parallel.");
121 "Finite element order (polynomial degree) or -1 for"
122 " isoparametric space.");
123 args.
AddOption(&basis_type,
"-b",
"--basis-type",
124 "Basis: G - Gauss-Lobatto, P - Positive, U - Uniform");
125 args.
AddOption(&perf,
"-perf",
"--hpc-version",
"-std",
"--standard-version",
126 "Enable high-performance, tensor-based, assembly/evaluation.");
127 args.
AddOption(&matrix_free,
"-mf",
"--matrix-free",
"-asm",
"--assembly",
128 "Use matrix-free evaluation or efficient matrix assembly in "
129 "the high-performance version.");
130 args.
AddOption(&pc,
"-pc",
"--preconditioner",
131 "Preconditioner: lor - low-order-refined (matrix-free) AMG, "
132 "ho - high-order (assembled) AMG, none.");
133 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
134 "--no-static-condensation",
"Enable static condensation.");
135 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
136 "--no-visualization",
137 "Enable or disable GLVis visualization.");
138 args.
AddOption(&visport,
"-p",
"--send-port",
"Socket for GLVis.");
148 if (static_cond && perf && matrix_free)
152 cout <<
"\nStatic condensation can not be used with matrix-free"
153 " evaluation!\n" << endl;
157 MFEM_VERIFY(perf || !matrix_free,
158 "--standard-version is not compatible with --matrix-free");
165 if (!strcmp(pc,
"ho")) { pc_choice =
PCType::HO; }
166 else if (!strcmp(pc,
"lor")) { pc_choice =
PCType::LOR; }
167 else if (!strcmp(pc,
"none")) { pc_choice =
PCType::NONE; }
170 mfem_error(
"Invalid Preconditioner specified");
176 cout <<
"\nMFEM SIMD width: " << MFEM_SIMD_BYTES/
sizeof(double)
177 <<
" doubles\n" << endl;
190 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
195 return ex1_t<2>::run(mesh, ser_ref_levels, par_ref_levels, order, basis,
196 static_cond, pc_choice, perf, matrix_free,
197 visualization, visport);
201 return ex1_t<3>::run(mesh, ser_ref_levels, par_ref_levels, order,
202 basis, static_cond, pc_choice, perf, matrix_free,
203 visualization, visport);
207 MFEM_ABORT(
"Dimension must be 2 or 3.")
214int ex1_t<dim>::run(
Mesh *mesh,
int ser_ref_levels,
int par_ref_levels,
215 int order,
int basis,
bool static_cond,
PCType pc_choice,
216 bool perf,
bool matrix_free,
bool visualization,
int visport)
219 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
220 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
227 cout <<
"High-performance version using integration rule with "
228 << int_rule_t::qpts <<
" points ..." << endl;
230 if (!mesh_t::MatchesGeometry(*mesh))
234 cout <<
"The given mesh does not match the optimized 'geom' parameter.\n"
235 <<
"Recompile with suitable 'geom' value." << endl;
240 else if (!mesh_t::MatchesNodes(*mesh))
244 cout <<
"Switching the mesh curvature to match the "
245 <<
"optimized value (order " <<
mesh_p <<
") ..." << endl;
257 int ref_levels = (ser_ref_levels != -1) ? ser_ref_levels :
258 (int)floor(log(10000./mesh->GetNE())/log(2.)/
dim);
259 for (
int l = 0; l < ref_levels; l++)
269 cout <<
"NURBS mesh: switching the mesh curvature to be "
270 <<
"min(sol_p, mesh_p) = " << new_mesh_p <<
" ..." << endl;
281 for (
int l = 0; l < par_ref_levels; l++)
288 MFEM_VERIFY(pc_choice !=
PCType::LOR,
"triangle and tet meshes do not "
289 "support the LOR preconditioner yet");
305 cout <<
"Using isoparametric FEs: " << fec->
Name() << endl;
316 cout <<
"Number of finite element unknowns: " << size << endl;
324 int basis_lor = basis;
332 if (perf && !sol_fes_t::Matches(*fespace))
336 cout <<
"The given order does not match the optimized parameter.\n"
337 <<
"Recompile with suitable 'sol_p' value." << endl;
384 a->EnableStaticCondensation();
386 "cannot use LOR preconditioner with static condensation");
391 cout <<
"Assembling the matrix ..." << flush;
396 a->UsePrecomputedSparsity();
398 HPCBilinearForm *a_hpc = NULL;
410 a_hpc =
new HPCBilinearForm(integ_t(coeff_t(1.0)), *fespace);
417 a_hpc->AssembleBilinearForm(*
a);
432 if (perf && matrix_free)
434 a_hpc->FormLinearSystem(ess_tdof_list, x, *
b, a_oper, X, B);
438 cout <<
"Size of linear system: " << glob_size << endl;
443 a->FormLinearSystem(ess_tdof_list, x, *
b, A, X, B);
447 cout <<
"Size of linear system: " << glob_size << endl;
455 cout <<
"Assembling the preconditioning matrix ..." << flush;
478 a_hpc->AssembleBilinearForm(*a_pc);
516 cout <<
"Time per CG step: "
522 if (perf && matrix_free)
524 a_hpc->RecoverFEMSolution(X, *
b, x);
528 a->RecoverFEMSolution(X, *
b, x);
534 ostringstream mesh_name, sol_name;
535 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
536 sol_name <<
"sol." << setfill(
'0') << setw(6) << myid;
538 ofstream mesh_ofs(mesh_name.str().c_str());
539 mesh_ofs.precision(8);
540 pmesh->
Print(mesh_ofs);
542 ofstream sol_ofs(sol_name.str().c_str());
543 sol_ofs.precision(8);
552 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
553 sol_sock.precision(8);
554 sol_sock <<
"solution\n" << *pmesh << x << flush;
560 if (a_oper != &A) {
delete a_oper; }
566 if (order > 0) {
delete fec; }
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
int Size() const
Return the logical size of the array.
@ Positive
Bernstein polynomials.
@ ClosedUniform
Nodes: x_i = i/(n-1), i=0,...,n-1.
static const char * Name(int b_type)
Check and convert a BasisType constant to a string identifier.
static int GetType(char b_ident)
Convert char basis identifier to a BasisType constant.
Conjugate gradient method.
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the Conjugate Gradient method.
void SetOperator(const Operator &op) override
Set/update the solver for the given operator.
A coefficient that is constant across space and time.
Class for domain integration .
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
virtual const char * Name() const
Arbitrary order H1-conforming (continuous) finite elements.
The BoomerAMG solver in hypre.
Wrapper for hypre's ParCSR matrix class.
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
void MakeRef(const HypreParMatrix &master)
Make this HypreParMatrix a reference to 'master'.
Abstract class for hypre's solvers and preconditioners.
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
void SetRelTol(real_t rtol)
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void SetMaxIter(int max_it)
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
int Dimension() const
Dimension of the reference space used within the elements.
void GetNodes(Vector &node_coord) const
int MeshGenerator() const
Get the mesh generator/type.
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 int WorldRank()
Return the MPI rank in 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).
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.
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
HYPRE_BigInt GlobalTrueVSize() const
Class for parallel grid function.
Class for parallel meshes.
static ParMesh MakeRefined(ParMesh &orig_mesh, int ref_factor, int ref_type)
Create a uniformly refined (by any factor) version of orig_mesh.
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
double RealTime()
Return the number of real seconds elapsed since the stopwatch was started.
void Start()
Start the stopwatch. The elapsed time is not cleared.
void Stop()
Stop the stopwatch.
void Clear()
Clear the elapsed time on the stopwatch and restart it if it's running.
The Integrator class combines a kernel and a coefficient.
void mfem_error(const char *msg)