MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex2p.cpp
Go to the documentation of this file.
1 // MFEM Example 2 - Parallel Version
2 //
3 // Compile with: make ex2p
4 //
5 // Sample runs: mpirun -np 4 ex2p -m ../data/beam-tri.mesh
6 // mpirun -np 4 ex2p -m ../data/beam-quad.mesh
7 // mpirun -np 4 ex2p -m ../data/beam-tet.mesh
8 // mpirun -np 4 ex2p -m ../data/beam-hex.mesh
9 // mpirun -np 4 ex2p -m ../data/beam-wedge.mesh
10 // mpirun -np 4 ex2p -m ../data/beam-tri.mesh -o 2 -sys
11 // mpirun -np 4 ex2p -m ../data/beam-quad.mesh -o 3 -elast
12 // mpirun -np 4 ex2p -m ../data/beam-quad.mesh -o 3 -sc
13 // mpirun -np 4 ex2p -m ../data/beam-quad-nurbs.mesh
14 // mpirun -np 4 ex2p -m ../data/beam-hex-nurbs.mesh
15 //
16 // Description: This example code solves a simple linear elasticity problem
17 // describing a multi-material cantilever beam.
18 //
19 // Specifically, we approximate the weak form of -div(sigma(u))=0
20 // where sigma(u)=lambda*div(u)*I+mu*(grad*u+u*grad) is the stress
21 // tensor corresponding to displacement field u, and lambda and mu
22 // are the material Lame constants. The boundary conditions are
23 // u=0 on the fixed part of the boundary with attribute 1, and
24 // sigma(u).n=f on the remainder with f being a constant pull down
25 // vector on boundary elements with attribute 2, and zero
26 // otherwise. The geometry of the domain is assumed to be as
27 // follows:
28 //
29 // +----------+----------+
30 // boundary --->| material | material |<--- boundary
31 // attribute 1 | 1 | 2 | attribute 2
32 // (fixed) +----------+----------+ (pull down)
33 //
34 // The example demonstrates the use of high-order and NURBS vector
35 // finite element spaces with the linear elasticity bilinear form,
36 // meshes with curved elements, and the definition of piece-wise
37 // constant and vector coefficient objects. Static condensation is
38 // also illustrated.
39 //
40 // We recommend viewing Example 1 before viewing this example.
41 
42 #include "mfem.hpp"
43 #include <fstream>
44 #include <iostream>
45 
46 using namespace std;
47 using namespace mfem;
48 
49 int main(int argc, char *argv[])
50 {
51  // 1. Initialize MPI.
52  int num_procs, myid;
53  MPI_Init(&argc, &argv);
54  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
55  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
56 
57  // 2. Parse command-line options.
58  const char *mesh_file = "../data/beam-tri.mesh";
59  int order = 1;
60  bool static_cond = false;
61  bool visualization = 1;
62  bool amg_elast = 0;
63  bool reorder_space = false;
64  const char *device_config = "cpu";
65 
66  OptionsParser args(argc, argv);
67  args.AddOption(&mesh_file, "-m", "--mesh",
68  "Mesh file to use.");
69  args.AddOption(&order, "-o", "--order",
70  "Finite element order (polynomial degree).");
71  args.AddOption(&amg_elast, "-elast", "--amg-for-elasticity", "-sys",
72  "--amg-for-systems",
73  "Use the special AMG elasticity solver (GM/LN approaches), "
74  "or standard AMG for systems (unknown approach).");
75  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
76  "--no-static-condensation", "Enable static condensation.");
77  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
78  "--no-visualization",
79  "Enable or disable GLVis visualization.");
80  args.AddOption(&reorder_space, "-nodes", "--by-nodes", "-vdim", "--by-vdim",
81  "Use byNODES ordering of vector space instead of byVDIM");
82  args.AddOption(&device_config, "-d", "--device",
83  "Device configuration string, see Device::Configure().");
84  args.Parse();
85  if (!args.Good())
86  {
87  if (myid == 0)
88  {
89  args.PrintUsage(cout);
90  }
91  MPI_Finalize();
92  return 1;
93  }
94  if (myid == 0)
95  {
96  args.PrintOptions(cout);
97  }
98 
99  // 3. Enable hardware devices such as GPUs, and programming models such as
100  // CUDA, OCCA, RAJA and OpenMP based on command line options.
101  Device device(device_config);
102  if (myid == 0) { device.Print(); }
103 
104  // 4. Read the (serial) mesh from the given mesh file on all processors. We
105  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
106  // and volume meshes with the same code.
107  Mesh *mesh = new Mesh(mesh_file, 1, 1);
108  int dim = mesh->Dimension();
109 
110  if (mesh->attributes.Max() < 2 || mesh->bdr_attributes.Max() < 2)
111  {
112  if (myid == 0)
113  cerr << "\nInput mesh should have at least two materials and "
114  << "two boundary attributes! (See schematic in ex2.cpp)\n"
115  << endl;
116  MPI_Finalize();
117  return 3;
118  }
119 
120  // 5. Select the order of the finite element discretization space. For NURBS
121  // meshes, we increase the order by degree elevation.
122  if (mesh->NURBSext)
123  {
124  mesh->DegreeElevate(order, order);
125  }
126 
127  // 6. Refine the serial mesh on all processors to increase the resolution. In
128  // this example we do 'ref_levels' of uniform refinement. We choose
129  // 'ref_levels' to be the largest number that gives a final mesh with no
130  // more than 1,000 elements.
131  {
132  int ref_levels =
133  (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
134  for (int l = 0; l < ref_levels; l++)
135  {
136  mesh->UniformRefinement();
137  }
138  }
139 
140  // 7. Define a parallel mesh by a partitioning of the serial mesh. Refine
141  // this mesh further in parallel to increase the resolution. Once the
142  // parallel mesh is defined, the serial mesh can be deleted.
143  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
144  delete mesh;
145  {
146  int par_ref_levels = 1;
147  for (int l = 0; l < par_ref_levels; l++)
148  {
149  pmesh->UniformRefinement();
150  }
151  }
152 
153  // 8. Define a parallel finite element space on the parallel mesh. Here we
154  // use vector finite elements, i.e. dim copies of a scalar finite element
155  // space. We use the ordering by vector dimension (the last argument of
156  // the FiniteElementSpace constructor) which is expected in the systems
157  // version of BoomerAMG preconditioner. For NURBS meshes, we use the
158  // (degree elevated) NURBS space associated with the mesh nodes.
160  ParFiniteElementSpace *fespace;
161  const bool use_nodal_fespace = pmesh->NURBSext && !amg_elast;
162  if (use_nodal_fespace)
163  {
164  fec = NULL;
165  fespace = (ParFiniteElementSpace *)pmesh->GetNodes()->FESpace();
166  }
167  else
168  {
169  fec = new H1_FECollection(order, dim);
170  if (reorder_space)
171  {
172  fespace = new ParFiniteElementSpace(pmesh, fec, dim, Ordering::byNODES);
173  }
174  else
175  {
176  fespace = new ParFiniteElementSpace(pmesh, fec, dim, Ordering::byVDIM);
177  }
178  }
179  HYPRE_BigInt size = fespace->GlobalTrueVSize();
180  if (myid == 0)
181  {
182  cout << "Number of finite element unknowns: " << size << endl
183  << "Assembling: " << flush;
184  }
185 
186  // 9. Determine the list of true (i.e. parallel conforming) essential
187  // boundary dofs. In this example, the boundary conditions are defined by
188  // marking only boundary attribute 1 from the mesh as essential and
189  // converting it to a list of true dofs.
190  Array<int> ess_tdof_list, ess_bdr(pmesh->bdr_attributes.Max());
191  ess_bdr = 0;
192  ess_bdr[0] = 1;
193  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
194 
195  // 10. Set up the parallel linear form b(.) which corresponds to the
196  // right-hand side of the FEM linear system. In this case, b_i equals the
197  // boundary integral of f*phi_i where f represents a "pull down" force on
198  // the Neumann part of the boundary and phi_i are the basis functions in
199  // the finite element fespace. The force is defined by the object f, which
200  // is a vector of Coefficient objects. The fact that f is non-zero on
201  // boundary attribute 2 is indicated by the use of piece-wise constants
202  // coefficient for its last component.
204  for (int i = 0; i < dim-1; i++)
205  {
206  f.Set(i, new ConstantCoefficient(0.0));
207  }
208  {
209  Vector pull_force(pmesh->bdr_attributes.Max());
210  pull_force = 0.0;
211  pull_force(1) = -1.0e-2;
212  f.Set(dim-1, new PWConstCoefficient(pull_force));
213  }
214 
215  ParLinearForm *b = new ParLinearForm(fespace);
217  if (myid == 0)
218  {
219  cout << "r.h.s. ... " << flush;
220  }
221  b->Assemble();
222 
223  // 11. Define the solution vector x as a parallel finite element grid
224  // function corresponding to fespace. Initialize x with initial guess of
225  // zero, which satisfies the boundary conditions.
226  ParGridFunction x(fespace);
227  x = 0.0;
228 
229  // 12. Set up the parallel bilinear form a(.,.) on the finite element space
230  // corresponding to the linear elasticity integrator with piece-wise
231  // constants coefficient lambda and mu.
232  Vector lambda(pmesh->attributes.Max());
233  lambda = 1.0;
234  lambda(0) = lambda(1)*50;
235  PWConstCoefficient lambda_func(lambda);
236  Vector mu(pmesh->attributes.Max());
237  mu = 1.0;
238  mu(0) = mu(1)*50;
239  PWConstCoefficient mu_func(mu);
240 
241  ParBilinearForm *a = new ParBilinearForm(fespace);
242  a->AddDomainIntegrator(new ElasticityIntegrator(lambda_func, mu_func));
243 
244  // 13. Assemble the parallel bilinear form and the corresponding linear
245  // system, applying any necessary transformations such as: parallel
246  // assembly, eliminating boundary conditions, applying conforming
247  // constraints for non-conforming AMR, static condensation, etc.
248  if (myid == 0) { cout << "matrix ... " << flush; }
249  if (static_cond) { a->EnableStaticCondensation(); }
250  a->Assemble();
251 
252  HypreParMatrix A;
253  Vector B, X;
254  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
255  if (myid == 0)
256  {
257  cout << "done." << endl;
258  cout << "Size of linear system: " << A.GetGlobalNumRows() << endl;
259  }
260 
261  // 14. Define and apply a parallel PCG solver for A X = B with the BoomerAMG
262  // preconditioner from hypre.
263  HypreBoomerAMG *amg = new HypreBoomerAMG(A);
264  if (amg_elast && !a->StaticCondensationIsEnabled())
265  {
266  amg->SetElasticityOptions(fespace);
267  }
268  else
269  {
270  amg->SetSystemsOptions(dim, reorder_space);
271  }
272  HyprePCG *pcg = new HyprePCG(A);
273  pcg->SetTol(1e-8);
274  pcg->SetMaxIter(500);
275  pcg->SetPrintLevel(2);
276  pcg->SetPreconditioner(*amg);
277  pcg->Mult(B, X);
278 
279  // 15. Recover the parallel grid function corresponding to X. This is the
280  // local finite element solution on each processor.
281  a->RecoverFEMSolution(X, *b, x);
282 
283  // 16. For non-NURBS meshes, make the mesh curved based on the finite element
284  // space. This means that we define the mesh elements through a fespace
285  // based transformation of the reference element. This allows us to save
286  // the displaced mesh as a curved mesh when using high-order finite
287  // element displacement field. We assume that the initial mesh (read from
288  // the file) is not higher order curved mesh compared to the chosen FE
289  // space.
290  if (!use_nodal_fespace)
291  {
292  pmesh->SetNodalFESpace(fespace);
293  }
294 
295  // 17. Save in parallel the displaced mesh and the inverted solution (which
296  // gives the backward displacements to the original grid). This output
297  // can be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
298  {
299  GridFunction *nodes = pmesh->GetNodes();
300  *nodes += x;
301  x *= -1;
302 
303  ostringstream mesh_name, sol_name;
304  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
305  sol_name << "sol." << setfill('0') << setw(6) << myid;
306 
307  ofstream mesh_ofs(mesh_name.str().c_str());
308  mesh_ofs.precision(8);
309  pmesh->Print(mesh_ofs);
310 
311  ofstream sol_ofs(sol_name.str().c_str());
312  sol_ofs.precision(8);
313  x.Save(sol_ofs);
314  }
315 
316  // 18. Send the above data by socket to a GLVis server. Use the "n" and "b"
317  // keys in GLVis to visualize the displacements.
318  if (visualization)
319  {
320  char vishost[] = "localhost";
321  int visport = 19916;
322  socketstream sol_sock(vishost, visport);
323  sol_sock << "parallel " << num_procs << " " << myid << "\n";
324  sol_sock.precision(8);
325  sol_sock << "solution\n" << *pmesh << x << flush;
326  }
327 
328  // 19. Free the used memory.
329  delete pcg;
330  delete amg;
331  delete a;
332  delete b;
333  if (fec)
334  {
335  delete fespace;
336  delete fec;
337  }
338  delete pmesh;
339 
340  MPI_Finalize();
341 
342  return 0;
343 }
void SetTol(double tol)
Definition: hypre.cpp:3569
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:790
Vector coefficient defined by an array of scalar coefficients. Coefficients that are not set will eva...
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:275
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:841
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:279
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3589
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1387
Class for parallel linear form.
Definition: plinearform.hpp:26
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
constexpr char vishost[]
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9143
constexpr int visport
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void Set(int i, Coefficient *c, bool own=true)
Sets coefficient in the vector.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4382
void SetElasticityOptions(ParFiniteElementSpace *fespace)
Definition: hypre.cpp:4569
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:4843
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:70
int Dimension() const
Definition: mesh.hpp:911
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3579
bool StaticCondensationIsEnabled() const
Check if static condensation was actually enabled by a previous call to EnableStaticCondensation().
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:204
PCG solver in hypre.
Definition: hypre.hpp:1060
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
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 &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
HYPRE_Int HYPRE_BigInt
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:206
A piecewise constant coefficient with the constants keyed off the element attribute numbers...
Definition: coefficient.hpp:94
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:573
double mu
Definition: ex25.cpp:139
void DegreeElevate(int rel_degree, int degree=16)
Definition: mesh.cpp:4665
int dim
Definition: ex24.cpp:53
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:3594
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:327
Class for parallel bilinear form.
Vector data type.
Definition: vector.hpp:60
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7343
void SetSystemsOptions(int dim, bool order_bynodes=false)
Definition: hypre.cpp:4466
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
Class for parallel grid function.
Definition: pgridfunc.hpp:32
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:121
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:277
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:3617
Class for parallel meshes.
Definition: pmesh.hpp:32
int main()
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:202
void EnableStaticCondensation()
Enable the use of static condensation. For details see the description for class StaticCondensation i...
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150