MFEM  v3.2
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 -perf -m ../../data/fichera.mesh
6 // ex1 -perf -m ../../data/amr-hex.mesh -sc
7 // ex1 -perf -m ../../data/ball-nurbs.mesh -sc
8 // ex1 -perf -m ../../data/pipe-nurbs.mesh
9 // ex1 -std -m ../../data/square-disc.mesh
10 // ex1 -std -m ../../data/star.mesh
11 // ex1 -std -m ../../data/escher.mesh
12 // ex1 -std -m ../../data/square-disc-p2.vtk -o 2
13 // ex1 -std -m ../../data/square-disc-p3.mesh -o 3
14 // ex1 -std -m ../../data/square-disc-nurbs.mesh -o -1
15 // ex1 -std -m ../../data/disc-nurbs.mesh -o -1
16 // ex1 -std -m ../../data/pipe-nurbs.mesh -o -1
17 // ex1 -std -m ../../data/star-surf.mesh
18 // ex1 -std -m ../../data/square-disc-surf.mesh
19 // ex1 -std -m ../../data/inline-segment.mesh
20 // ex1 -std -m ../../data/amr-quad.mesh
21 // ex1 -std -m ../../data/amr-hex.mesh
22 // ex1 -std -m ../../data/fichera-amr.mesh
23 // ex1 -std -m ../../data/mobius-strip.mesh
24 // ex1 -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. Parse command-line options.
75  const char *mesh_file = "../../data/fichera.mesh";
76  int order = sol_p;
77  bool static_cond = false;
78  bool visualization = 1;
79  bool perf = true;
80 
81  OptionsParser args(argc, argv);
82  args.AddOption(&mesh_file, "-m", "--mesh",
83  "Mesh file to use.");
84  args.AddOption(&order, "-o", "--order",
85  "Finite element order (polynomial degree) or -1 for"
86  " isoparametric space.");
87  args.AddOption(&perf, "-perf", "--hpc-version", "-std", "--standard-version",
88  "Enable high-performance, tensor-based, assembly.");
89  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
90  "--no-static-condensation", "Enable static condensation.");
91  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
92  "--no-visualization",
93  "Enable or disable GLVis visualization.");
94  args.Parse();
95  if (!args.Good())
96  {
97  args.PrintUsage(cout);
98  return 1;
99  }
100  args.PrintOptions(cout);
101 
102  // 2. Read the mesh from the given mesh file. We can handle triangular,
103  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
104  // the same code.
105  Mesh *mesh = new Mesh(mesh_file, 1, 1);
106  int dim = mesh->Dimension();
107 
108  // 3. Check if the optimized version matches the given mesh
109  if (perf)
110  {
111  cout << "High-performance version using integration rule with "
112  << int_rule_t::qpts << " points ..." << endl;
113  if (!mesh_t::MatchesGeometry(*mesh))
114  {
115  cout << "The given mesh does not match the optimized 'geom' parameter.\n"
116  << "Recompile with suitable 'geom' value." << endl;
117  delete mesh;
118  return 3;
119  }
120  else if (!mesh_t::MatchesNodes(*mesh))
121  {
122  cout << "Switching the mesh curvature to match the "
123  << "optimized value (order " << mesh_p << ") ..." << endl;
124  mesh->SetCurvature(mesh_p, false, -1, Ordering::byNODES);
125  }
126  }
127 
128  // 4. Refine the mesh to increase the resolution. In this example we do
129  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
130  // largest number that gives a final mesh with no more than 50,000
131  // elements.
132  {
133  int ref_levels =
134  (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
135  for (int l = 0; l < ref_levels; l++)
136  {
137  mesh->UniformRefinement();
138  }
139  }
140 
141  // 5. Define a finite element space on the mesh. Here we use continuous
142  // Lagrange finite elements of the specified order. If order < 1, we
143  // instead use an isoparametric/isogeometric space.
145  if (order > 0)
146  {
147  fec = new H1_FECollection(order, dim);
148  }
149  else if (mesh->GetNodes())
150  {
151  fec = mesh->GetNodes()->OwnFEC();
152  cout << "Using isoparametric FEs: " << fec->Name() << endl;
153  }
154  else
155  {
156  fec = new H1_FECollection(order = 1, dim);
157  }
158  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
159  cout << "Number of finite element unknowns: "
160  << fespace->GetTrueVSize() << endl;
161 
162  // 6. Check if the optimized version matches the given space
163  if (perf && !sol_fes_t::Matches(*fespace))
164  {
165  cout << "The given order does not match the optimized parameter.\n"
166  << "Recompile with suitable 'sol_p' value." << endl;
167  delete fespace;
168  delete fec;
169  delete mesh;
170  return 4;
171  }
172 
173  // 7. Determine the list of true (i.e. conforming) essential boundary dofs.
174  // In this example, the boundary conditions are defined by marking all
175  // the boundary attributes from the mesh as essential (Dirichlet) and
176  // converting them to a list of true dofs.
177  Array<int> ess_tdof_list;
178  if (mesh->bdr_attributes.Size())
179  {
180  Array<int> ess_bdr(mesh->bdr_attributes.Max());
181  ess_bdr = 1;
182  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
183  }
184 
185  // 8. Set up the linear form b(.) which corresponds to the right-hand side of
186  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
187  // the basis functions in the finite element fespace.
188  LinearForm *b = new LinearForm(fespace);
189  ConstantCoefficient one(1.0);
191  b->Assemble();
192 
193  // 9. Define the solution vector x as a finite element grid function
194  // corresponding to fespace. Initialize x with initial guess of zero,
195  // which satisfies the boundary conditions.
196  GridFunction x(fespace);
197  x = 0.0;
198 
199  // 10. Set up the bilinear form a(.,.) on the finite element space that will
200  // hold the matrix corresponding to the Laplacian operator -Delta.
201  BilinearForm *a = new BilinearForm(fespace);
202 
203  // 11. Assemble the bilinear form and the corresponding linear system,
204  // applying any necessary transformations such as: eliminating boundary
205  // conditions, applying conforming constraints for non-conforming AMR,
206  // static condensation, etc.
207  if (static_cond) { a->EnableStaticCondensation(); }
208 
209  cout << "Assembling the matrix ..." << flush;
210  tic_toc.Clear();
211  tic_toc.Start();
212  // Pre-allocate sparsity assuming dense element matrices
214  if (!perf)
215  {
216  // Standard assembly using a diffusion domain integrator
218  a->Assemble();
219  }
220  else
221  {
222  // High-performance assembly using the templated operator type
223  oper_t a_oper(integ_t(coeff_t(1.0)), *fespace);
224  a_oper.AssembleBilinearForm(*a);
225  }
226  tic_toc.Stop();
227  cout << " done, " << tic_toc.RealTime() << "s." << endl;
228 
229  SparseMatrix A;
230  Vector B, X;
231  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
232 
233  cout << "Size of linear system: " << A.Height() << endl;
234 
235 #ifndef MFEM_USE_SUITESPARSE
236  // 12. Define a simple symmetric Gauss-Seidel preconditioner and use it to
237  // solve the system A X = B with PCG.
238  GSSmoother M(A);
239  PCG(A, M, B, X, 1, 500, 1e-12, 0.0);
240 #else
241  // 12. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
242  UMFPackSolver umf_solver;
243  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
244  umf_solver.SetOperator(A);
245  umf_solver.Mult(B, X);
246 #endif
247 
248  // 13. Recover the solution as a finite element grid function.
249  a->RecoverFEMSolution(X, *b, x);
250 
251  // 14. Save the refined mesh and the solution. This output can be viewed later
252  // using GLVis: "glvis -m refined.mesh -g sol.gf".
253  ofstream mesh_ofs("refined.mesh");
254  mesh_ofs.precision(8);
255  mesh->Print(mesh_ofs);
256  ofstream sol_ofs("sol.gf");
257  sol_ofs.precision(8);
258  x.Save(sol_ofs);
259 
260  // 15. Send the solution by socket to a GLVis server.
261  if (visualization)
262  {
263  char vishost[] = "localhost";
264  int visport = 19916;
265  socketstream sol_sock(vishost, visport);
266  sol_sock.precision(8);
267  sol_sock << "solution\n" << *mesh << x << flush;
268  }
269 
270  // 16. Free the used memory.
271  delete a;
272  delete b;
273  delete fespace;
274  if (order > 0) { delete fec; }
275  delete mesh;
276 
277  return 0;
278 }
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
const Geometry::Type geom
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:34
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list)
Definition: fespace.cpp:274
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:496
StopWatch tic_toc
Definition: tic_toc.cpp:338
double RealTime()
Definition: tic_toc.cpp:317
int main(int argc, char *argv[])
void FormLinearSystem(Array< int > &ess_tdof_list, Vector &x, Vector &b, SparseMatrix &A, Vector &X, Vector &B, int copy_interior=0)
Data type for Gauss-Seidel smoother of sparse matrix.
H1_FiniteElement< geom, sol_p > sol_fe_t
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:333
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2116
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:35
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 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:523
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
virtual void Print(std::ostream &out=std::cout) const
Print the mesh to the given stream using the default MFEM mesh format.
Definition: mesh.cpp:6736
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
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:344
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
H1_FiniteElement< geom, mesh_p > mesh_fe_t
TMesh< mesh_fes_t > mesh_t
virtual int GetTrueVSize()
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:164
virtual const char * Name() const
Definition: fe_coll.hpp:50
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: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
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
TIntegrator< coeff_t, TDiffusionKernel > integ_t
TBilinearForm< mesh_t, sol_fes_t, int_rule_t, integ_t > oper_t
virtual void Mult(const Vector &b, Vector &x) const
Operator application.
Definition: solvers.cpp:1741
void EnableStaticCondensation()
const int ir_order
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: solvers.cpp:1641
bool Good() const
Definition: optparser.hpp:120