32 #include "mfem-performance.hpp"
63 int main(
int argc,
char *argv[])
66 const char *mesh_file =
"../../data/fichera.mesh";
68 const char *basis_type =
"G";
69 bool static_cond =
false;
70 const char *pc =
"none";
72 bool matrix_free =
true;
73 bool visualization = 1;
76 args.
AddOption(&mesh_file,
"-m",
"--mesh",
79 "Finite element order (polynomial degree) or -1 for"
80 " isoparametric space.");
81 args.
AddOption(&basis_type,
"-b",
"--basis-type",
82 "Basis: G - Gauss-Lobatto, P - Positive, U - Uniform");
83 args.
AddOption(&perf,
"-perf",
"--hpc-version",
"-std",
"--standard-version",
84 "Enable high-performance, tensor-based, assembly/evaluation.");
85 args.
AddOption(&matrix_free,
"-mf",
"--matrix-free",
"-asm",
"--assembly",
86 "Use matrix-free evaluation or efficient matrix assembly in "
87 "the high-performance version.");
88 args.
AddOption(&pc,
"-pc",
"--preconditioner",
89 "Preconditioner: lor - low-order-refined (matrix-free) GS, "
90 "ho - high-order (assembled) GS, none.");
91 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
92 "--no-static-condensation",
"Enable static condensation.");
93 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
95 "Enable or disable GLVis visualization.");
102 if (static_cond && perf && matrix_free)
104 cout <<
"\nStatic condensation can not be used with matrix-free"
105 " evaluation!\n" << endl;
108 MFEM_VERIFY(perf || !matrix_free,
109 "--standard-version is not compatible with --matrix-free");
112 enum PCType { NONE, LOR, HO };
114 if (!strcmp(pc,
"ho")) { pc_choice = HO; }
115 else if (!strcmp(pc,
"lor")) { pc_choice = LOR; }
116 else if (!strcmp(pc,
"none")) { pc_choice = NONE; }
119 mfem_error(
"Invalid Preconditioner specified");
124 int basis = BasisType::GetType(basis_type[0]);
125 cout <<
"Using " << BasisType::Name(basis) <<
" basis ..." << endl;
130 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
136 cout <<
"High-performance version using integration rule with "
137 << int_rule_t::qpts <<
" points ..." << endl;
138 if (!mesh_t::MatchesGeometry(*mesh))
140 cout <<
"The given mesh does not match the optimized 'geom' parameter.\n"
141 <<
"Recompile with suitable 'geom' value." << endl;
145 else if (!mesh_t::MatchesNodes(*mesh))
147 cout <<
"Switching the mesh curvature to match the "
148 <<
"optimized value (order " <<
mesh_p <<
") ..." << endl;
159 (int)floor(log(50000./mesh->
GetNE())/log(2.)/
dim);
160 for (
int l = 0; l < ref_levels; l++)
167 MFEM_VERIFY(pc_choice != LOR,
"triangle and tet meshes do not support"
168 " the LOR preconditioner yet");
182 cout <<
"Using isoparametric FEs: " << fec->
Name() << endl;
189 cout <<
"Number of finite element unknowns: "
194 Mesh *mesh_lor = NULL;
197 if (pc_choice == LOR)
199 int basis_lor = basis;
200 if (basis == BasisType::Positive) { basis_lor=BasisType::ClosedUniform; }
201 mesh_lor =
new Mesh(mesh, order, basis_lor);
207 if (perf && !sol_fes_t::Matches(*fespace))
209 cout <<
"The given order does not match the optimized parameter.\n"
210 <<
"Recompile with suitable 'sol_p' value." << endl;
248 if (pc_choice == LOR) { a_pc =
new BilinearForm(fespace_lor); }
249 if (pc_choice == HO) { a_pc =
new BilinearForm(fespace); }
258 MFEM_VERIFY(pc_choice != LOR,
259 "cannot use LOR preconditioner with static condensation");
262 cout <<
"Assembling the bilinear form ..." << flush;
299 if (perf && matrix_free)
302 cout <<
"Size of linear system: " << a_hpc->
Height() << endl;
307 cout <<
"Size of linear system: " << A.
Height() << endl;
312 cout <<
"Assembling the preconditioning matrix ..." << flush;
317 if (pc_choice == LOR)
325 else if (pc_choice == HO)
343 if (pc_choice != NONE)
346 PCG(*a_oper, M, B, X, 1, 500, 1e-12, 0.0);
350 CG(*a_oper, B, X, 1, 500, 1e-12, 0.0);
354 if (perf && matrix_free)
365 ofstream mesh_ofs(
"refined.mesh");
366 mesh_ofs.precision(8);
367 mesh->
Print(mesh_ofs);
368 ofstream sol_ofs(
"sol.gf");
369 sol_ofs.precision(8);
375 char vishost[] =
"localhost";
378 sol_sock.precision(8);
379 sol_sock <<
"solution\n" << *mesh << x << flush;
385 if (a_oper != &A) {
delete a_oper; }
392 if (order > 0) {
delete fec; }
int Size() const
Logical size of the array.
Class for domain integration L(v) := (f, v)
Class for grid function - Vector with associated FE space.
Subclass constant coefficient.
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 GetNE() const
Returns number of elements.
int main(int argc, char *argv[])
Data type for Gauss-Seidel smoother of sparse matrix.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
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)
void PrintUsage(std::ostream &out) const
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
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)
void mfem_error(const char *msg)
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 int GetTrueVSize()
Return the number of vector true (conforming) dofs.
virtual const char * Name() const
void PrintOptions(std::ostream &out) const
void GetNodes(Vector &node_coord) const
Arbitrary order H1-conforming (continuous) finite elements.
virtual void Print(std::ostream &out=std::cout) const