32 #include "mfem-performance.hpp"
63 int main(
int argc,
char *argv[])
66 const char *mesh_file =
"../../data/fichera.mesh";
69 const char *basis_type =
"G";
70 bool static_cond =
false;
71 const char *pc =
"none";
73 bool matrix_free =
true;
74 bool visualization = 1;
77 args.
AddOption(&mesh_file,
"-m",
"--mesh",
79 args.
AddOption(&ref_levels,
"-r",
"--refine",
80 "Number of times to refine the mesh uniformly;"
81 " -1 = auto: <= 50,000 elements.");
83 "Finite element order (polynomial degree) or -1 for"
84 " isoparametric space.");
85 args.
AddOption(&basis_type,
"-b",
"--basis-type",
86 "Basis: G - Gauss-Lobatto, P - Positive, U - Uniform");
87 args.
AddOption(&perf,
"-perf",
"--hpc-version",
"-std",
"--standard-version",
88 "Enable high-performance, tensor-based, assembly/evaluation.");
89 args.
AddOption(&matrix_free,
"-mf",
"--matrix-free",
"-asm",
"--assembly",
90 "Use matrix-free evaluation or efficient matrix assembly in "
91 "the high-performance version.");
92 args.
AddOption(&pc,
"-pc",
"--preconditioner",
93 "Preconditioner: lor - low-order-refined (matrix-free) GS, "
94 "ho - high-order (assembled) GS, none.");
95 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
96 "--no-static-condensation",
"Enable static condensation.");
97 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
99 "Enable or disable GLVis visualization.");
106 if (static_cond && perf && matrix_free)
108 cout <<
"\nStatic condensation can not be used with matrix-free"
109 " evaluation!\n" << endl;
112 MFEM_VERIFY(perf || !matrix_free,
113 "--standard-version is not compatible with --matrix-free");
116 enum PCType {
NONE, LOR, HO };
118 if (!strcmp(pc,
"ho")) { pc_choice = HO; }
119 else if (!strcmp(pc,
"lor")) { pc_choice = LOR; }
120 else if (!strcmp(pc,
"none")) { pc_choice =
NONE; }
123 mfem_error(
"Invalid Preconditioner specified");
127 cout <<
"\nMFEM SIMD width: " << MFEM_SIMD_BYTES/
sizeof(double)
128 <<
" doubles\n" << endl;
131 int basis = BasisType::GetType(basis_type[0]);
132 cout <<
"Using " << BasisType::Name(basis) <<
" basis ..." << endl;
137 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
143 cout <<
"High-performance version using integration rule with "
144 << int_rule_t::qpts <<
" points ..." << endl;
145 if (!mesh_t::MatchesGeometry(*mesh))
147 cout <<
"The given mesh does not match the optimized 'geom' parameter.\n"
148 <<
"Recompile with suitable 'geom' value." << endl;
152 else if (!mesh_t::MatchesNodes(*mesh))
154 cout <<
"Switching the mesh curvature to match the "
155 <<
"optimized value (order " <<
mesh_p <<
") ..." << endl;
166 ref_levels = (ref_levels != -1) ? ref_levels :
167 (
int)floor(log(50000./mesh->
GetNE())/log(2.)/
dim);
168 for (
int l = 0; l < ref_levels; l++)
175 MFEM_VERIFY(pc_choice != LOR,
"triangle and tet meshes do not support"
176 " the LOR preconditioner yet");
190 cout <<
"Using isoparametric FEs: " << fec->
Name() << endl;
197 cout <<
"Number of finite element unknowns: "
202 Mesh *mesh_lor = NULL;
205 if (pc_choice == LOR)
207 int basis_lor = basis;
208 if (basis == BasisType::Positive) { basis_lor=BasisType::ClosedUniform; }
209 mesh_lor =
new Mesh(mesh, order, basis_lor);
215 if (perf && !sol_fes_t::Matches(*fespace))
217 cout <<
"The given order does not match the optimized parameter.\n"
218 <<
"Recompile with suitable 'sol_p' value." << endl;
256 if (pc_choice == LOR) { a_pc =
new BilinearForm(fespace_lor); }
257 if (pc_choice == HO) { a_pc =
new BilinearForm(fespace); }
266 MFEM_VERIFY(pc_choice != LOR,
267 "cannot use LOR preconditioner with static condensation");
270 cout <<
"Assembling the bilinear form ..." << flush;
307 if (perf && matrix_free)
310 cout <<
"Size of linear system: " << a_hpc->
Height() << endl;
315 cout <<
"Size of linear system: " << A.
Height() << endl;
320 cout <<
"Assembling the preconditioning matrix ..." << flush;
325 if (pc_choice == LOR)
333 else if (pc_choice == HO)
351 if (pc_choice !=
NONE)
354 PCG(*a_oper, M, B, X, 1, 500, 1e-12, 0.0);
358 CG(*a_oper, B, X, 1, 500, 1e-12, 0.0);
362 if (perf && matrix_free)
373 ofstream mesh_ofs(
"refined.mesh");
374 mesh_ofs.precision(8);
375 mesh->
Print(mesh_ofs);
376 ofstream sol_ofs(
"sol.gf");
377 sol_ofs.precision(8);
386 sol_sock.precision(8);
387 sol_sock <<
"solution\n" << *mesh << x << flush;
393 if (a_oper != &A) {
delete a_oper; }
400 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 Print(std::ostream &out=mfem::out) const
Class for grid function - Vector with associated FE space.
A coefficient that is constant across space and time.
void MakeRef(const SparseMatrix &master)
Clear the contents of the SparseMatrix and make it a reference to master.
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
int GetNE() const
Returns number of elements.
int main(int argc, char *argv[])
The Integrator class combines a kernel and a coefficient.
Data type for Gauss-Seidel smoother of sparse matrix.
void Stop()
Stop the stopwatch.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
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().
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 *)
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 CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Conjugate gradient method. (tolerances are squared)
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
void PrintUsage(std::ostream &out) const
Print the usage message.
void Start()
Clear the elapsed time and start the stopwatch.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
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...
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0)
Form a constrained linear system using a matrix-free approach.
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear sys...
virtual const char * Name() const
void PrintOptions(std::ostream &out) const
Print the options.
void GetNodes(Vector &node_coord) const
Arbitrary order H1-conforming (continuous) finite elements.
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.