MFEM  v3.3
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
miniapps/performance/ex1.cpp
Go to the documentation of this file.
1 // MFEM Example 1 - High-Performance Version
2 //
3 // Compile with: make ex1
4 //
5 // Sample runs: ex1 -m ../../data/fichera.mesh -perf -mf -pc lor
6 // ex1 -m ../../data/fichera.mesh -perf -asm -pc ho
7 // ex1 -m ../../data/fichera.mesh -perf -asm -pc ho -sc
8 // ex1 -m ../../data/fichera.mesh -std -asm -pc ho
9 // ex1 -m ../../data/fichera.mesh -std -asm -pc ho -sc
10 // ex1 -m ../../data/amr-hex.mesh -perf -asm -pc ho -sc
11 // ex1 -m ../../data/amr-hex.mesh -std -asm -pc ho -sc
12 // ex1 -m ../../data/ball-nurbs.mesh -perf -asm -pc ho -sc
13 // ex1 -m ../../data/ball-nurbs.mesh -std -asm -pc ho -sc
14 // ex1 -m ../../data/pipe-nurbs.mesh -perf -mf -pc lor
15 // ex1 -m ../../data/pipe-nurbs.mesh -std -asm -pc ho -sc
16 //
17 // Description: This example code demonstrates the use of MFEM to define a
18 // simple finite element discretization of the Laplace problem
19 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
20 // Specifically, we discretize using a FE space of the specified
21 // order, or if order < 1 using an isoparametric/isogeometric
22 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
23 // NURBS mesh, etc.)
24 //
25 // The example highlights the use of mesh refinement, finite
26 // element grid functions, as well as linear and bilinear forms
27 // corresponding to the left-hand side and right-hand side of the
28 // discrete linear system. We also cover the explicit elimination
29 // of essential boundary conditions, static condensation, and the
30 // optional connection to the GLVis tool for visualization.
31 
32 #include "mfem-performance.hpp"
33 #include <fstream>
34 #include <iostream>
35 
36 using namespace std;
37 using namespace mfem;
38 
39 // Define template parameters for optimized build.
40 const Geometry::Type geom = Geometry::CUBE; // mesh elements (default: hex)
41 const int mesh_p = 3; // mesh curvature (default: 3)
42 const int sol_p = 3; // solution order (default: 3)
44 const int ir_order = 2*sol_p+rdim-1;
45 
46 // Static mesh type
50 
51 // Static solution finite element space type
54 
55 // Static quadrature, coefficient and integrator types
59 
60 // Static bilinear form type, combining the above types
62 
63 int main(int argc, char *argv[])
64 {
65  // 1. Parse command-line options.
66  const char *mesh_file = "../../data/fichera.mesh";
67  int order = sol_p;
68  const char *basis_type = "G"; // Gauss-Lobatto
69  bool static_cond = false;
70  const char *pc = "none";
71  bool perf = true;
72  bool matrix_free = true;
73  bool visualization = 1;
74 
75  OptionsParser args(argc, argv);
76  args.AddOption(&mesh_file, "-m", "--mesh",
77  "Mesh file to use.");
78  args.AddOption(&order, "-o", "--order",
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",
94  "--no-visualization",
95  "Enable or disable GLVis visualization.");
96  args.Parse();
97  if (!args.Good())
98  {
99  args.PrintUsage(cout);
100  return 1;
101  }
102  if (static_cond && perf && matrix_free)
103  {
104  cout << "\nStatic condensation can not be used with matrix-free"
105  " evaluation!\n" << endl;
106  return 2;
107  }
108  MFEM_VERIFY(perf || !matrix_free,
109  "--standard-version is not compatible with --matrix-free");
110  args.PrintOptions(cout);
111 
112  enum PCType { NONE, LOR, HO };
113  PCType pc_choice;
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; }
117  else
118  {
119  mfem_error("Invalid Preconditioner specified");
120  return 3;
121  }
122 
123  // See class BasisType in fem/fe_coll.hpp for available basis types
124  int basis = BasisType::GetType(basis_type[0]);
125  cout << "Using " << BasisType::Name(basis) << " basis ..." << endl;
126 
127  // 2. Read the mesh from the given mesh file. We can handle triangular,
128  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
129  // the same code.
130  Mesh *mesh = new Mesh(mesh_file, 1, 1);
131  int dim = mesh->Dimension();
132 
133  // 3. Check if the optimized version matches the given mesh
134  if (perf)
135  {
136  cout << "High-performance version using integration rule with "
137  << int_rule_t::qpts << " points ..." << endl;
138  if (!mesh_t::MatchesGeometry(*mesh))
139  {
140  cout << "The given mesh does not match the optimized 'geom' parameter.\n"
141  << "Recompile with suitable 'geom' value." << endl;
142  delete mesh;
143  return 4;
144  }
145  else if (!mesh_t::MatchesNodes(*mesh))
146  {
147  cout << "Switching the mesh curvature to match the "
148  << "optimized value (order " << mesh_p << ") ..." << endl;
149  mesh->SetCurvature(mesh_p, false, -1, Ordering::byNODES);
150  }
151  }
152 
153  // 4. Refine the mesh to increase the resolution. In this example we do
154  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
155  // largest number that gives a final mesh with no more than 50,000
156  // elements.
157  {
158  int ref_levels =
159  (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
160  for (int l = 0; l < ref_levels; l++)
161  {
162  mesh->UniformRefinement();
163  }
164  }
165  if (mesh->MeshGenerator() & 1) // simplex mesh
166  {
167  MFEM_VERIFY(pc_choice != LOR, "triangle and tet meshes do not support"
168  " the LOR preconditioner yet");
169  }
170 
171  // 5. Define a finite element space on the mesh. Here we use continuous
172  // Lagrange finite elements of the specified order. If order < 1, we
173  // instead use an isoparametric/isogeometric space.
175  if (order > 0)
176  {
177  fec = new H1_FECollection(order, dim, basis);
178  }
179  else if (mesh->GetNodes())
180  {
181  fec = mesh->GetNodes()->OwnFEC();
182  cout << "Using isoparametric FEs: " << fec->Name() << endl;
183  }
184  else
185  {
186  fec = new H1_FECollection(order = 1, dim, basis);
187  }
188  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
189  cout << "Number of finite element unknowns: "
190  << fespace->GetTrueVSize() << endl;
191 
192  // Create the LOR mesh and finite element space. In the settings of this
193  // example, we can transfer between HO and LOR with the identity operator.
194  Mesh *mesh_lor = NULL;
195  FiniteElementCollection *fec_lor = NULL;
196  FiniteElementSpace *fespace_lor = NULL;
197  if (pc_choice == LOR)
198  {
199  int basis_lor = basis;
200  if (basis == BasisType::Positive) { basis_lor=BasisType::ClosedUniform; }
201  mesh_lor = new Mesh(mesh, order, basis_lor);
202  fec_lor = new H1_FECollection(1, dim);
203  fespace_lor = new FiniteElementSpace(mesh_lor, fec_lor);
204  }
205 
206  // 6. Check if the optimized version matches the given space
207  if (perf && !sol_fes_t::Matches(*fespace))
208  {
209  cout << "The given order does not match the optimized parameter.\n"
210  << "Recompile with suitable 'sol_p' value." << endl;
211  delete fespace;
212  delete fec;
213  delete mesh;
214  return 5;
215  }
216 
217  // 7. Determine the list of true (i.e. conforming) essential boundary dofs.
218  // In this example, the boundary conditions are defined by marking all
219  // the boundary attributes from the mesh as essential (Dirichlet) and
220  // converting them to a list of true dofs.
221  Array<int> ess_tdof_list;
222  if (mesh->bdr_attributes.Size())
223  {
224  Array<int> ess_bdr(mesh->bdr_attributes.Max());
225  ess_bdr = 1;
226  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
227  }
228 
229  // 8. Set up the linear form b(.) which corresponds to the right-hand side of
230  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
231  // the basis functions in the finite element fespace.
232  LinearForm *b = new LinearForm(fespace);
233  ConstantCoefficient one(1.0);
235  b->Assemble();
236 
237  // 9. Define the solution vector x as a finite element grid function
238  // corresponding to fespace. Initialize x with initial guess of zero,
239  // which satisfies the boundary conditions.
240  GridFunction x(fespace);
241  x = 0.0;
242 
243  // 10. Set up the bilinear form a(.,.) on the finite element space that will
244  // hold the matrix corresponding to the Laplacian operator -Delta.
245  // Optionally setup a form to be assembled for preconditioning (a_pc).
246  BilinearForm *a = new BilinearForm(fespace);
247  BilinearForm *a_pc = NULL;
248  if (pc_choice == LOR) { a_pc = new BilinearForm(fespace_lor); }
249  if (pc_choice == HO) { a_pc = new BilinearForm(fespace); }
250 
251  // 11. Assemble the bilinear form and the corresponding linear system,
252  // applying any necessary transformations such as: eliminating boundary
253  // conditions, applying conforming constraints for non-conforming AMR,
254  // static condensation, etc.
255  if (static_cond)
256  {
258  MFEM_VERIFY(pc_choice != LOR,
259  "cannot use LOR preconditioner with static condensation");
260  }
261 
262  cout << "Assembling the bilinear form ..." << flush;
263  tic_toc.Clear();
264  tic_toc.Start();
265  // Pre-allocate sparsity assuming dense element matrices
267 
268  HPCBilinearForm *a_hpc = NULL;
269  Operator *a_oper = NULL;
270 
271  if (!perf)
272  {
273  // Standard assembly using a diffusion domain integrator
275  a->Assemble();
276  }
277  else
278  {
279  // High-performance assembly/evaluation using the templated operator type
280  a_hpc = new HPCBilinearForm(integ_t(coeff_t(1.0)), *fespace);
281  if (matrix_free)
282  {
283  a_hpc->Assemble(); // partial assembly
284  }
285  else
286  {
287  a_hpc->AssembleBilinearForm(*a); // full matrix assembly
288  }
289  }
290  tic_toc.Stop();
291  cout << " done, " << tic_toc.RealTime() << "s." << endl;
292 
293  // 12. Solve the system A X = B with CG. In the standard case, use a simple
294  // symmetric Gauss-Seidel preconditioner.
295 
296  // Setup the operator matrix (if applicable)
297  SparseMatrix A;
298  Vector B, X;
299  if (perf && matrix_free)
300  {
301  a_hpc->FormLinearSystem(ess_tdof_list, x, *b, a_oper, X, B);
302  cout << "Size of linear system: " << a_hpc->Height() << endl;
303  }
304  else
305  {
306  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
307  cout << "Size of linear system: " << A.Height() << endl;
308  a_oper = &A;
309  }
310 
311  // Setup the matrix used for preconditioning
312  cout << "Assembling the preconditioning matrix ..." << flush;
313  tic_toc.Clear();
314  tic_toc.Start();
315 
316  SparseMatrix A_pc;
317  if (pc_choice == LOR)
318  {
319  // TODO: assemble the LOR matrix using the performance code
321  a_pc->UsePrecomputedSparsity();
322  a_pc->Assemble();
323  a_pc->FormSystemMatrix(ess_tdof_list, A_pc);
324  }
325  else if (pc_choice == HO)
326  {
327  if (!matrix_free)
328  {
329  A_pc.MakeRef(A); // matrix already assembled, reuse it
330  }
331  else
332  {
333  a_pc->UsePrecomputedSparsity();
334  a_hpc->AssembleBilinearForm(*a_pc);
335  a_pc->FormSystemMatrix(ess_tdof_list, A_pc);
336  }
337  }
338 
339  tic_toc.Stop();
340  cout << " done, " << tic_toc.RealTime() << "s." << endl;
341 
342  // Solve with CG or PCG, depending if the matrix A_pc is available
343  if (pc_choice != NONE)
344  {
345  GSSmoother M(A_pc);
346  PCG(*a_oper, M, B, X, 1, 500, 1e-12, 0.0);
347  }
348  else
349  {
350  CG(*a_oper, B, X, 1, 500, 1e-12, 0.0);
351  }
352 
353  // 13. Recover the solution as a finite element grid function.
354  if (perf && matrix_free)
355  {
356  a_hpc->RecoverFEMSolution(X, *b, x);
357  }
358  else
359  {
360  a->RecoverFEMSolution(X, *b, x);
361  }
362 
363  // 14. Save the refined mesh and the solution. This output can be viewed later
364  // using GLVis: "glvis -m refined.mesh -g sol.gf".
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);
370  x.Save(sol_ofs);
371 
372  // 15. Send the solution by socket to a GLVis server.
373  if (visualization)
374  {
375  char vishost[] = "localhost";
376  int visport = 19916;
377  socketstream sol_sock(vishost, visport);
378  sol_sock.precision(8);
379  sol_sock << "solution\n" << *mesh << x << flush;
380  }
381 
382  // 16. Free the used memory.
383  delete a;
384  delete a_hpc;
385  if (a_oper != &A) { delete a_oper; }
386  delete a_pc;
387  delete b;
388  delete fespace;
389  delete fespace_lor;
390  delete fec_lor;
391  delete mesh_lor;
392  if (order > 0) { delete fec; }
393  delete mesh;
394 
395  return 0;
396 }
int Size() const
Logical size of the array.
Definition: array.hpp:109
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:47
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
Subclass constant coefficient.
Definition: coefficient.hpp:57
void MakeRef(const SparseMatrix &master)
Clear the contents of the SparseMatrix and make it a reference to master.
Definition: sparsemat.cpp:187
const Geometry::Type geom
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:42
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list)
Definition: fespace.cpp:295
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
void AssembleBilinearForm(BilinearForm &a) const
const int rdim
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:584
StopWatch tic_toc
Definition: tic_toc.cpp:338
void FormSystemMatrix(const Array< int > &ess_tdof_list, SparseMatrix &A)
Form the linear system matrix A, see FormLinearSystem for details.
double RealTime()
Definition: tic_toc.cpp:317
int main(int argc, char *argv[])
Data type for Gauss-Seidel smoother of sparse matrix.
H1_FiniteElement< geom, sol_p > sol_fe_t
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2185
TConstantCoefficient coeff_t
int dim
Definition: ex3.cpp:47
Data type sparse matrix.
Definition: sparsemat.hpp:38
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:36
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6684
const int sol_p
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3278
int MeshGenerator()
Get the mesh generator/type.
Definition: mesh.hpp:577
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:108
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)
Definition: solvers.cpp:443
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)
Definition: solvers.cpp:456
int Dimension() const
Definition: mesh.hpp:611
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:142
H1_FiniteElementSpace< mesh_fe_t > mesh_fes_t
void UsePrecomputedSparsity(int ps=1)
H1_FiniteElementSpace< sol_fe_t > sol_fes_t
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)
Definition: optparser.hpp:74
void mfem_error(const char *msg)
Definition: error.cpp:106
H1_FiniteElement< geom, mesh_p > mesh_fe_t
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.
Definition: operator.cpp:21
TMesh< mesh_fes_t > mesh_t
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...
Definition: operator.cpp:57
TBilinearForm< mesh_t, sol_fes_t, int_rule_t, integ_t > HPCBilinearForm
virtual int GetTrueVSize()
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:166
virtual const char * Name() const
Definition: fe_coll.hpp:117
const int mesh_p
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
Vector data type.
Definition: vector.hpp:36
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5368
TIntegrationRule< geom, ir_order > int_rule_t
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:146
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
TIntegrator< coeff_t, TDiffusionKernel > integ_t
Abstract operator.
Definition: operator.hpp:21
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, SparseMatrix &A, Vector &X, Vector &B, int copy_interior=0)
void EnableStaticCondensation()
const int ir_order
virtual void Print(std::ostream &out=std::cout) const
Definition: mesh.hpp:988
bool Good() const
Definition: optparser.hpp:120