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 ref_levels,
int order,
int basis,
84 bool static_cond,
PCType pc_choice,
bool perf,
85 bool matrix_free,
bool visualization,
int visport = 19916);
88int 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;
100 bool visualization = 1;
103 args.
AddOption(&mesh_file,
"-m",
"--mesh",
104 "Mesh file to use.");
105 args.
AddOption(&ref_levels,
"-r",
"--refine",
106 "Number of times to refine the mesh uniformly;"
107 " -1 = auto: <= 50,000 elements.");
109 "Finite element order (polynomial degree) or -1 for"
110 " isoparametric space.");
111 args.
AddOption(&basis_type,
"-b",
"--basis-type",
112 "Basis: G - Gauss-Lobatto, P - Positive, U - Uniform");
113 args.
AddOption(&perf,
"-perf",
"--hpc-version",
"-std",
"--standard-version",
114 "Enable high-performance, tensor-based, assembly/evaluation.");
115 args.
AddOption(&matrix_free,
"-mf",
"--matrix-free",
"-asm",
"--assembly",
116 "Use matrix-free evaluation or efficient matrix assembly in "
117 "the high-performance version.");
118 args.
AddOption(&pc,
"-pc",
"--preconditioner",
119 "Preconditioner: lor - low-order-refined (matrix-free) GS, "
120 "ho - high-order (assembled) GS, none.");
121 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
122 "--no-static-condensation",
"Enable static condensation.");
123 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
124 "--no-visualization",
125 "Enable or disable GLVis visualization.");
126 args.
AddOption(&visport,
"-p",
"--send-port",
"Socket for GLVis.");
133 if (static_cond && perf && matrix_free)
135 cout <<
"\nStatic condensation can not be used with matrix-free"
136 " evaluation!\n" << endl;
139 MFEM_VERIFY(perf || !matrix_free,
140 "--standard-version is not compatible with --matrix-free");
144 if (!strcmp(pc,
"ho")) { pc_choice =
PCType::HO; }
145 else if (!strcmp(pc,
"lor")) { pc_choice =
PCType::LOR; }
146 else if (!strcmp(pc,
"none")) { pc_choice =
PCType::NONE; }
149 mfem_error(
"Invalid Preconditioner specified");
153 cout <<
"\nMFEM SIMD width: " << MFEM_SIMD_BYTES/
sizeof(double)
154 <<
" doubles\n" << endl;
163 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
168 return ex1_t<2>::run(mesh, ref_levels, order, basis, static_cond,
169 pc_choice, perf, matrix_free, visualization, visport);
173 return ex1_t<3>::run(mesh, ref_levels, order, basis, static_cond,
174 pc_choice, perf, matrix_free, visualization, visport);
178 MFEM_ABORT(
"Dimension must be 2 or 3.")
185int ex1_t<dim>::run(
Mesh *mesh,
int ref_levels,
int order,
int basis,
186 bool static_cond,
PCType pc_choice,
bool perf,
187 bool matrix_free,
bool visualization,
int visport)
192 cout <<
"High-performance version using integration rule with "
193 << int_rule_t::qpts <<
" points ..." << endl;
194 if (!mesh_t::MatchesGeometry(*mesh))
196 cout <<
"The given mesh does not match the optimized 'geom' parameter.\n"
197 <<
"Recompile with suitable 'geom' value." << endl;
201 else if (!mesh_t::MatchesNodes(*mesh))
203 cout <<
"Switching the mesh curvature to match the "
204 <<
"optimized value (order " <<
mesh_p <<
") ..." << endl;
215 ref_levels = (ref_levels != -1) ? ref_levels :
216 (int)floor(log(50000./mesh->GetNE())/log(2.)/
dim);
217 for (
int l = 0; l < ref_levels; l++)
224 MFEM_VERIFY(pc_choice !=
PCType::LOR,
"triangle and tet meshes do not "
225 " support the LOR preconditioner yet");
239 cout <<
"Using isoparametric FEs: " << fec->
Name() << endl;
246 cout <<
"Number of finite element unknowns: "
256 int basis_lor = basis;
264 if (perf && !sol_fes_t::Matches(*fespace))
266 cout <<
"The given order does not match the optimized parameter.\n"
267 <<
"Recompile with suitable 'sol_p' value." << endl;
314 a->EnableStaticCondensation();
316 "cannot use LOR preconditioner with static condensation");
319 cout <<
"Assembling the bilinear form ..." << flush;
323 a->UsePrecomputedSparsity();
325 HPCBilinearForm *a_hpc = NULL;
337 a_hpc =
new HPCBilinearForm(integ_t(coeff_t(1.0)), *fespace);
344 a_hpc->AssembleBilinearForm(*
a);
356 if (perf && matrix_free)
358 a_hpc->FormLinearSystem(ess_tdof_list, x, *
b, a_oper, X, B);
359 cout <<
"Size of linear system: " << a_hpc->Height() << endl;
363 a->FormLinearSystem(ess_tdof_list, x, *
b, A, X, B);
364 cout <<
"Size of linear system: " << A.
Height() << endl;
369 cout <<
"Assembling the preconditioning matrix ..." << flush;
391 a_hpc->AssembleBilinearForm(*a_pc);
403 PCG(*a_oper, M, B, X, 1, 500, 1e-12, 0.0);
407 CG(*a_oper, B, X, 1, 500, 1e-12, 0.0);
411 if (perf && matrix_free)
413 a_hpc->RecoverFEMSolution(X, *
b, x);
417 a->RecoverFEMSolution(X, *
b, x);
422 ofstream mesh_ofs(
"refined.mesh");
423 mesh_ofs.precision(8);
424 mesh->
Print(mesh_ofs);
425 ofstream sol_ofs(
"sol.gf");
426 sol_ofs.precision(8);
434 sol_sock.precision(8);
435 sol_sock <<
"solution\n" << *mesh << x << flush;
441 if (a_oper != &A) {
delete a_oper; }
447 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.
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
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Data type for Gauss-Seidel smoother of sparse matrix.
Class for grid function - Vector with associated FE space.
Arbitrary order H1-conforming (continuous) finite elements.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
int Dimension() const
Dimension of the reference space used within the elements.
static Mesh MakeRefined(Mesh &orig_mesh, int ref_factor, int ref_type)
Create a refined (by any factor) version of orig_mesh.
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 *)
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
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.
void MakeRef(const SparseMatrix &master)
Clear the contents of the SparseMatrix and make it a reference to master.
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)
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
void CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Conjugate gradient method. (tolerances are squared)