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 ref_levels,
int order,
int basis,
84 bool static_cond,
PCType pc_choice,
bool perf,
85 bool matrix_free,
bool visualization);
88 int main(
int argc,
char *argv[])
91 const char *mesh_file =
"../../data/fichera.mesh";
94 const char *basis_type =
"G";
95 bool static_cond =
false;
96 const char *pc =
"none";
98 bool matrix_free =
true;
99 bool visualization = 1;
102 args.
AddOption(&mesh_file,
"-m",
"--mesh",
103 "Mesh file to use.");
104 args.
AddOption(&ref_levels,
"-r",
"--refine",
105 "Number of times to refine the mesh uniformly;"
106 " -1 = auto: <= 50,000 elements.");
108 "Finite element order (polynomial degree) or -1 for"
109 " isoparametric space.");
110 args.
AddOption(&basis_type,
"-b",
"--basis-type",
111 "Basis: G - Gauss-Lobatto, P - Positive, U - Uniform");
112 args.
AddOption(&perf,
"-perf",
"--hpc-version",
"-std",
"--standard-version",
113 "Enable high-performance, tensor-based, assembly/evaluation.");
114 args.
AddOption(&matrix_free,
"-mf",
"--matrix-free",
"-asm",
"--assembly",
115 "Use matrix-free evaluation or efficient matrix assembly in "
116 "the high-performance version.");
117 args.
AddOption(&pc,
"-pc",
"--preconditioner",
118 "Preconditioner: lor - low-order-refined (matrix-free) GS, "
119 "ho - high-order (assembled) GS, none.");
120 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
121 "--no-static-condensation",
"Enable static condensation.");
122 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
123 "--no-visualization",
124 "Enable or disable GLVis visualization.");
131 if (static_cond && perf && matrix_free)
133 cout <<
"\nStatic condensation can not be used with matrix-free"
134 " evaluation!\n" << endl;
137 MFEM_VERIFY(perf || !matrix_free,
138 "--standard-version is not compatible with --matrix-free");
142 if (!strcmp(pc,
"ho")) { pc_choice =
PCType::HO; }
143 else if (!strcmp(pc,
"lor")) { pc_choice =
PCType::LOR; }
144 else if (!strcmp(pc,
"none")) { pc_choice =
PCType::NONE; }
147 mfem_error(
"Invalid Preconditioner specified");
151 cout <<
"\nMFEM SIMD width: " << MFEM_SIMD_BYTES/
sizeof(double)
152 <<
" doubles\n" << endl;
155 int basis = BasisType::GetType(basis_type[0]);
156 cout <<
"Using " << BasisType::Name(basis) <<
" basis ..." << endl;
161 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
166 return ex1_t<2>::run(mesh, ref_levels, order, basis, static_cond,
167 pc_choice, perf, matrix_free, visualization);
171 return ex1_t<3>::run(mesh, ref_levels, order, basis, static_cond,
172 pc_choice, perf, matrix_free, visualization);
176 MFEM_ABORT(
"Dimension must be 2 or 3.")
183 int ex1_t<dim>::run(
Mesh *mesh,
int ref_levels,
int order,
int basis,
184 bool static_cond,
PCType pc_choice,
bool perf,
185 bool matrix_free,
bool visualization)
190 cout <<
"High-performance version using integration rule with "
191 << int_rule_t::qpts <<
" points ..." << endl;
192 if (!mesh_t::MatchesGeometry(*mesh))
194 cout <<
"The given mesh does not match the optimized 'geom' parameter.\n"
195 <<
"Recompile with suitable 'geom' value." << endl;
199 else if (!mesh_t::MatchesNodes(*mesh))
201 cout <<
"Switching the mesh curvature to match the "
202 <<
"optimized value (order " <<
mesh_p <<
") ..." << endl;
213 ref_levels = (ref_levels != -1) ? ref_levels :
215 for (
int l = 0; l < ref_levels; l++)
222 MFEM_VERIFY(pc_choice !=
PCType::LOR,
"triangle and tet meshes do not "
223 " support the LOR preconditioner yet");
237 cout <<
"Using isoparametric FEs: " << fec->
Name() << endl;
244 cout <<
"Number of finite element unknowns: "
254 int basis_lor = basis;
255 if (basis == BasisType::Positive) { basis_lor=BasisType::ClosedUniform; }
256 mesh_lor = Mesh::MakeRefined(*mesh, order, basis_lor);
262 if (perf && !sol_fes_t::Matches(*fespace))
264 cout <<
"The given order does not match the optimized parameter.\n"
265 <<
"Recompile with suitable 'sol_p' value." << endl;
314 "cannot use LOR preconditioner with static condensation");
317 cout <<
"Assembling the bilinear form ..." << flush;
323 HPCBilinearForm *a_hpc = NULL;
335 a_hpc =
new HPCBilinearForm(integ_t(coeff_t(1.0)), *fespace);
342 a_hpc->AssembleBilinearForm(*a);
354 if (perf && matrix_free)
356 a_hpc->FormLinearSystem(ess_tdof_list, x, *b, a_oper, X, B);
357 cout <<
"Size of linear system: " << a_hpc->Height() << endl;
362 cout <<
"Size of linear system: " << A.
Height() << endl;
367 cout <<
"Assembling the preconditioning matrix ..." << flush;
389 a_hpc->AssembleBilinearForm(*a_pc);
401 PCG(*a_oper, M, B, X, 1, 500, 1e-12, 0.0);
405 CG(*a_oper, B, X, 1, 500, 1e-12, 0.0);
409 if (perf && matrix_free)
411 a_hpc->RecoverFEMSolution(X, *b, x);
420 ofstream mesh_ofs(
"refined.mesh");
421 mesh_ofs.precision(8);
422 mesh->
Print(mesh_ofs);
423 ofstream sol_ofs(
"sol.gf");
424 sol_ofs.precision(8);
433 sol_sock.precision(8);
434 sol_sock <<
"solution\n" << *mesh << x << flush;
440 if (a_oper != &A) {
delete a_oper; }
446 if (order > 0) {
delete fec; }
int Size() const
Return the logical size of the array.
Class for domain integration L(v) := (f, v)
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.
The Integrator class combines a kernel and a coefficient.
Data type for Gauss-Seidel smoother of sparse matrix.
void Stop()
Stop the stopwatch.
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()
Start the stopwatch. The elapsed time is not cleared.
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...
FDualNumber< tbase > log(const FDualNumber< tbase > &f)
log([dual number])
virtual const char * Name() const
virtual void Print(std::ostream &os=mfem::out) 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.