MFEM  v3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
miniapps/performance/ex1p.cpp
Go to the documentation of this file.
1 // MFEM Example 1 - Parallel High-Performance Version
2 //
3 // Compile with: make ex1p
4 //
5 // Sample runs: mpirun -np 4 ex1p -perf -m ../../data/fichera.mesh
6 // mpirun -np 4 ex1p -perf -m ../../data/amr-hex.mesh -sc
7 // mpirun -np 4 ex1p -perf -m ../../data/ball-nurbs.mesh -sc
8 // mpirun -np 4 ex1p -perf -m ../../data/pipe-nurbs.mesh
9 // mpirun -np 4 ex1p -std -m ../../data/square-disc.mesh
10 // mpirun -np 4 ex1p -std -m ../../data/star.mesh
11 // mpirun -np 4 ex1p -std -m ../../data/escher.mesh
12 // mpirun -np 4 ex1p -std -m ../../data/square-disc-p2.vtk -o 2
13 // mpirun -np 4 ex1p -std -m ../../data/square-disc-p3.mesh -o 3
14 // mpirun -np 4 ex1p -std -m ../../data/square-disc-nurbs.mesh -o -1
15 // mpirun -np 4 ex1p -std -m ../../data/disc-nurbs.mesh -o -1
16 // mpirun -np 4 ex1p -std -m ../../data/pipe-nurbs.mesh -o -1
17 // mpirun -np 4 ex1p -std -m ../../data/star-surf.mesh
18 // mpirun -np 4 ex1p -std -m ../../data/square-disc-surf.mesh
19 // mpirun -np 4 ex1p -std -m ../../data/inline-segment.mesh
20 // mpirun -np 4 ex1p -std -m ../../data/amr-quad.mesh
21 // mpirun -np 4 ex1p -std -m ../../data/amr-hex.mesh
22 // mpirun -np 4 ex1p -std -m ../../data/fichera-amr.mesh
23 // mpirun -np 4 ex1p -std -m ../../data/mobius-strip.mesh
24 // mpirun -np 4 ex1p -std -m ../../data/mobius-strip.mesh -o -1 -sc
25 //
26 // Description: This example code demonstrates the use of MFEM to define a
27 // simple finite element discretization of the Laplace problem
28 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
29 // Specifically, we discretize using a FE space of the specified
30 // order, or if order < 1 using an isoparametric/isogeometric
31 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
32 // NURBS mesh, etc.)
33 //
34 // The example highlights the use of mesh refinement, finite
35 // element grid functions, as well as linear and bilinear forms
36 // corresponding to the left-hand side and right-hand side of the
37 // discrete linear system. We also cover the explicit elimination
38 // of essential boundary conditions, static condensation, and the
39 // optional connection to the GLVis tool for visualization.
40 
41 #include "mfem-performance.hpp"
42 #include <fstream>
43 #include <iostream>
44 
45 using namespace std;
46 using namespace mfem;
47 
48 // Define template parameters for optimized build.
49 const Geometry::Type geom = Geometry::CUBE; // mesh elements (default: hex)
50 const int mesh_p = 3; // mesh curvature (default: 3)
51 const int sol_p = 3; // solution order (default: 3)
53 const int ir_order = 2*sol_p+rdim-1;
54 
55 // Static mesh type
59 
60 // Static solution finite element space type
63 
64 // Static quadrature, coefficient and integrator types
68 
69 // Static bilinear form type, combining the above types
71 
72 int main(int argc, char *argv[])
73 {
74  // 1. Initialize MPI.
75  int num_procs, myid;
76  MPI_Init(&argc, &argv);
77  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
78  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
79 
80  // 2. Parse command-line options.
81  const char *mesh_file = "../../data/fichera.mesh";
82  int order = sol_p;
83  bool static_cond = false;
84  bool visualization = 1;
85  bool perf = true;
86 
87  OptionsParser args(argc, argv);
88  args.AddOption(&mesh_file, "-m", "--mesh",
89  "Mesh file to use.");
90  args.AddOption(&order, "-o", "--order",
91  "Finite element order (polynomial degree) or -1 for"
92  " isoparametric space.");
93  args.AddOption(&perf, "-perf", "--hpc-version", "-std", "--standard-version",
94  "Enable high-performance, tensor-based, assembly.");
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",
98  "--no-visualization",
99  "Enable or disable GLVis visualization.");
100  args.Parse();
101  if (!args.Good())
102  {
103  if (myid == 0)
104  {
105  args.PrintUsage(cout);
106  }
107  MPI_Finalize();
108  return 1;
109  }
110  if (myid == 0)
111  {
112  args.PrintOptions(cout);
113  }
114 
115  // 3. Read the (serial) mesh from the given mesh file on all processors. We
116  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
117  // and volume meshes with the same code.
118  Mesh *mesh = new Mesh(mesh_file, 1, 1);
119  int dim = mesh->Dimension();
120 
121  // 4. Check if the optimized version matches the given mesh
122  if (perf)
123  {
124  if (myid == 0)
125  {
126  cout << "High-performance version using integration rule with "
127  << int_rule_t::qpts << " points ..." << endl;
128  }
129  if (!mesh_t::MatchesGeometry(*mesh))
130  {
131  if (myid == 0)
132  {
133  cout << "The given mesh does not match the optimized 'geom' parameter.\n"
134  << "Recompile with suitable 'geom' value." << endl;
135  }
136  delete mesh;
137  MPI_Finalize();
138  return 3;
139  }
140  else if (!mesh_t::MatchesNodes(*mesh))
141  {
142  if (myid == 0)
143  {
144  cout << "Switching the mesh curvature to match the "
145  << "optimized value (order " << mesh_p << ") ..." << endl;
146  }
147  mesh->SetCurvature(mesh_p, false, -1, Ordering::byNODES);
148  }
149  }
150 
151  // 5. Refine the serial mesh on all processors to increase the resolution. In
152  // this example we do 'ref_levels' of uniform refinement. We choose
153  // 'ref_levels' to be the largest number that gives a final mesh with no
154  // more than 10,000 elements.
155  {
156  int ref_levels =
157  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
158  for (int l = 0; l < ref_levels; l++)
159  {
160  mesh->UniformRefinement();
161  }
162  }
163 
164  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
165  // this mesh further in parallel to increase the resolution. Once the
166  // parallel mesh is defined, the serial mesh can be deleted.
167  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
168  delete mesh;
169  {
170  int par_ref_levels = 1;
171  for (int l = 0; l < par_ref_levels; l++)
172  {
173  pmesh->UniformRefinement();
174  }
175  }
176 
177  // 7. Define a parallel finite element space on the parallel mesh. Here we
178  // use continuous Lagrange finite elements of the specified order. If
179  // order < 1, we instead use an isoparametric/isogeometric space.
181  if (order > 0)
182  {
183  fec = new H1_FECollection(order, dim);
184  }
185  else if (pmesh->GetNodes())
186  {
187  fec = pmesh->GetNodes()->OwnFEC();
188  if (myid == 0)
189  {
190  cout << "Using isoparametric FEs: " << fec->Name() << endl;
191  }
192  }
193  else
194  {
195  fec = new H1_FECollection(order = 1, dim);
196  }
197  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
198  HYPRE_Int size = fespace->GlobalTrueVSize();
199  if (myid == 0)
200  {
201  cout << "Number of finite element unknowns: " << size << endl;
202  }
203 
204  // 8. Check if the optimized version matches the given space
205  if (perf && !sol_fes_t::Matches(*fespace))
206  {
207  if (myid == 0)
208  {
209  cout << "The given order does not match the optimized parameter.\n"
210  << "Recompile with suitable 'sol_p' value." << endl;
211  }
212  delete fespace;
213  delete fec;
214  delete mesh;
215  MPI_Finalize();
216  return 4;
217  }
218 
219  // 9. Determine the list of true (i.e. parallel conforming) essential
220  // boundary dofs. In this example, the boundary conditions are defined
221  // by marking all the boundary attributes from the mesh as essential
222  // (Dirichlet) and converting them to a list of true dofs.
223  Array<int> ess_tdof_list;
224  if (pmesh->bdr_attributes.Size())
225  {
226  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
227  ess_bdr = 1;
228  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
229  }
230 
231  // 10. Set up the parallel linear form b(.) which corresponds to the
232  // right-hand side of the FEM linear system, which in this case is
233  // (1,phi_i) where phi_i are the basis functions in fespace.
234  ParLinearForm *b = new ParLinearForm(fespace);
235  ConstantCoefficient one(1.0);
237  b->Assemble();
238 
239  // 11. Define the solution vector x as a parallel finite element grid
240  // function corresponding to fespace. Initialize x with initial guess of
241  // zero, which satisfies the boundary conditions.
242  ParGridFunction x(fespace);
243  x = 0.0;
244 
245  // 12. Set up the parallel bilinear form a(.,.) on the finite element space
246  // that will hold the matrix corresponding to the Laplacian operator.
247  ParBilinearForm *a = new ParBilinearForm(fespace);
248 
249  // 13. Assemble the parallel bilinear form and the corresponding linear
250  // system, applying any necessary transformations such as: parallel
251  // assembly, eliminating boundary conditions, applying conforming
252  // constraints for non-conforming AMR, static condensation, etc.
253  if (static_cond) { a->EnableStaticCondensation(); }
254 
255  if (myid == 0)
256  {
257  cout << "Assembling the matrix ..." << flush;
258  }
259  tic_toc.Clear();
260  tic_toc.Start();
261  // Pre-allocate sparsity assuming dense element matrices
263  if (!perf)
264  {
265  // Standard assembly using a diffusion domain integrator
267  a->Assemble();
268  }
269  else
270  {
271  // High-performance assembly using the templated operator type
272  oper_t a_oper(integ_t(coeff_t(1.0)), *fespace);
273  a_oper.AssembleBilinearForm(*a);
274  }
275  tic_toc.Stop();
276  if (myid == 0)
277  {
278  cout << " done, " << tic_toc.RealTime() << "s." << endl;
279  }
280 
281  HypreParMatrix A;
282  Vector B, X;
283  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
284 
285  if (myid == 0)
286  {
287  cout << "Size of linear system: " << A.GetGlobalNumRows() << endl;
288  }
289 
290  // 14. Define and apply a parallel PCG solver for AX=B with the BoomerAMG
291  // preconditioner from hypre.
292  HypreSolver *amg = new HypreBoomerAMG(A);
293  HyprePCG *pcg = new HyprePCG(A);
294  pcg->SetTol(1e-12);
295  pcg->SetMaxIter(200);
296  pcg->SetPrintLevel(2);
297  pcg->SetPreconditioner(*amg);
298  pcg->Mult(B, X);
299 
300  // 15. Recover the parallel grid function corresponding to X. This is the
301  // local finite element solution on each processor.
302  a->RecoverFEMSolution(X, *b, x);
303 
304  // 16. Save the refined mesh and the solution in parallel. This output can
305  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
306  {
307  ostringstream mesh_name, sol_name;
308  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
309  sol_name << "sol." << setfill('0') << setw(6) << myid;
310 
311  ofstream mesh_ofs(mesh_name.str().c_str());
312  mesh_ofs.precision(8);
313  pmesh->Print(mesh_ofs);
314 
315  ofstream sol_ofs(sol_name.str().c_str());
316  sol_ofs.precision(8);
317  x.Save(sol_ofs);
318  }
319 
320  // 17. Send the solution by socket to a GLVis server.
321  if (visualization)
322  {
323  char vishost[] = "localhost";
324  int visport = 19916;
325  socketstream sol_sock(vishost, visport);
326  sol_sock << "parallel " << num_procs << " " << myid << "\n";
327  sol_sock.precision(8);
328  sol_sock << "solution\n" << *pmesh << x << flush;
329  }
330 
331  // 18. Free the used memory.
332  delete pcg;
333  delete amg;
334  delete a;
335  delete b;
336  delete fespace;
337  if (order > 0) { delete fec; }
338  delete pmesh;
339 
340  MPI_Finalize();
341 
342  return 0;
343 }
void SetTol(double tol)
Definition: hypre.cpp:1964
int Size() const
Logical size of the array.
Definition: array.hpp:109
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:47
Subclass constant coefficient.
Definition: coefficient.hpp:57
const Geometry::Type geom
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:34
void AssembleBilinearForm(BilinearForm &a) const
const int rdim
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:496
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:312
StopWatch tic_toc
Definition: tic_toc.cpp:338
double RealTime()
Definition: tic_toc.cpp:317
Abstract parallel finite element space.
Definition: pfespace.hpp:28
int main(int argc, char *argv[])
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:1979
H1_FiniteElement< geom, sol_p > sol_fe_t
TConstantCoefficient coeff_t
The BoomerAMG solver in hypre.
Definition: hypre.hpp:705
Class for parallel linear form.
Definition: plinearform.hpp:26
HYPRE_Int GetGlobalNumRows() const
Definition: hypre.hpp:340
int dim
Definition: ex3.cpp:47
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6402
const int sol_p
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3020
T Max() const
Definition: array.cpp:90
void Assemble(int skip_zeros=1)
Assemble the local matrix.
int Dimension() const
Definition: mesh.hpp:523
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetMaxIter(int max_iter)
Definition: hypre.cpp:1969
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Array< int > bdr_attributes
Definition: mesh.hpp:141
H1_FiniteElementSpace< mesh_fe_t > mesh_fes_t
void UsePrecomputedSparsity(int ps=1)
H1_FiniteElementSpace< sol_fe_t > sol_fes_t
PCG solver in hypre.
Definition: hypre.hpp:565
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list)
Definition: pfespace.cpp:479
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
H1_FiniteElement< geom, mesh_p > mesh_fe_t
TMesh< mesh_fes_t > mesh_t
void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
virtual const char * Name() const
Definition: fe_coll.hpp:50
const int mesh_p
void FormLinearSystem(Array< int > &ess_tdof_list, Vector &x, Vector &b, HypreParMatrix &A, Vector &X, Vector &B, int copy_interior=0)
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:1984
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:529
Vector data type.
Definition: vector.hpp:33
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5114
TIntegrationRule< geom, ir_order > int_rule_t
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:77
TIntegrator< coeff_t, TDiffusionKernel > integ_t
Class for parallel grid function.
Definition: pgridfunc.hpp:31
TBilinearForm< mesh_t, sol_fes_t, int_rule_t, integ_t > oper_t
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:150
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2005
Class for parallel meshes.
Definition: pmesh.hpp:28
void EnableStaticCondensation()
const int ir_order
virtual void Print(std::ostream &out=std::cout) const
Definition: pmesh.cpp:2927
bool Good() const
Definition: optparser.hpp:120