MFEM  v4.2.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 
65  OptionsParser args(argc, argv);
66  args.AddOption(&mesh_file, "-m", "--mesh",
67  "Mesh file to use.");
68  args.AddOption(&order, "-o", "--order",
69  "Finite element order (polynomial degree).");
70  args.AddOption(&amg_elast, "-elast", "--amg-for-elasticity", "-sys",
71  "--amg-for-systems",
72  "Use the special AMG elasticity solver (GM/LN approaches), "
73  "or standard AMG for systems (unknown approach).");
74  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
75  "--no-static-condensation", "Enable static condensation.");
76  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
77  "--no-visualization",
78  "Enable or disable GLVis visualization.");
79  args.AddOption(&reorder_space, "-nodes", "--by-nodes", "-vdim", "--by-vdim",
80  "Use byNODES ordering of vector space instead of byVDIM");
81  args.Parse();
82  if (!args.Good())
83  {
84  if (myid == 0)
85  {
86  args.PrintUsage(cout);
87  }
88  MPI_Finalize();
89  return 1;
90  }
91  if (myid == 0)
92  {
93  args.PrintOptions(cout);
94  }
95 
96  // 3. Read the (serial) mesh from the given mesh file on all processors. We
97  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
98  // and volume meshes with the same code.
99  Mesh *mesh = new Mesh(mesh_file, 1, 1);
100  int dim = mesh->Dimension();
101 
102  if (mesh->attributes.Max() < 2 || mesh->bdr_attributes.Max() < 2)
103  {
104  if (myid == 0)
105  cerr << "\nInput mesh should have at least two materials and "
106  << "two boundary attributes! (See schematic in ex2.cpp)\n"
107  << endl;
108  MPI_Finalize();
109  return 3;
110  }
111 
112  // 4. Select the order of the finite element discretization space. For NURBS
113  // meshes, we increase the order by degree elevation.
114  if (mesh->NURBSext)
115  {
116  mesh->DegreeElevate(order, order);
117  }
118 
119  // 5. Refine the serial mesh on all processors to increase the resolution. In
120  // this example we do 'ref_levels' of uniform refinement. We choose
121  // 'ref_levels' to be the largest number that gives a final mesh with no
122  // more than 1,000 elements.
123  {
124  int ref_levels =
125  (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
126  for (int l = 0; l < ref_levels; l++)
127  {
128  mesh->UniformRefinement();
129  }
130  }
131 
132  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
133  // this mesh further in parallel to increase the resolution. Once the
134  // parallel mesh is defined, the serial mesh can be deleted.
135  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
136  delete mesh;
137  {
138  int par_ref_levels = 1;
139  for (int l = 0; l < par_ref_levels; l++)
140  {
141  pmesh->UniformRefinement();
142  }
143  }
144 
145  // 7. Define a parallel finite element space on the parallel mesh. Here we
146  // use vector finite elements, i.e. dim copies of a scalar finite element
147  // space. We use the ordering by vector dimension (the last argument of
148  // the FiniteElementSpace constructor) which is expected in the systems
149  // version of BoomerAMG preconditioner. For NURBS meshes, we use the
150  // (degree elevated) NURBS space associated with the mesh nodes.
152  ParFiniteElementSpace *fespace;
153  const bool use_nodal_fespace = pmesh->NURBSext && !amg_elast;
154  if (use_nodal_fespace)
155  {
156  fec = NULL;
157  fespace = (ParFiniteElementSpace *)pmesh->GetNodes()->FESpace();
158  }
159  else
160  {
161  fec = new H1_FECollection(order, dim);
162  if (reorder_space)
163  {
164  fespace = new ParFiniteElementSpace(pmesh, fec, dim, Ordering::byNODES);
165  }
166  else
167  {
168  fespace = new ParFiniteElementSpace(pmesh, fec, dim, Ordering::byVDIM);
169  }
170  }
171  HYPRE_Int size = fespace->GlobalTrueVSize();
172  if (myid == 0)
173  {
174  cout << "Number of finite element unknowns: " << size << endl
175  << "Assembling: " << flush;
176  }
177 
178  // 8. Determine the list of true (i.e. parallel conforming) essential
179  // boundary dofs. In this example, the boundary conditions are defined by
180  // marking only boundary attribute 1 from the mesh as essential and
181  // converting it to a list of true dofs.
182  Array<int> ess_tdof_list, ess_bdr(pmesh->bdr_attributes.Max());
183  ess_bdr = 0;
184  ess_bdr[0] = 1;
185  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
186 
187  // 9. Set up the parallel linear form b(.) which corresponds to the
188  // right-hand side of the FEM linear system. In this case, b_i equals the
189  // boundary integral of f*phi_i where f represents a "pull down" force on
190  // the Neumann part of the boundary and phi_i are the basis functions in
191  // the finite element fespace. The force is defined by the object f, which
192  // is a vector of Coefficient objects. The fact that f is non-zero on
193  // boundary attribute 2 is indicated by the use of piece-wise constants
194  // coefficient for its last component.
196  for (int i = 0; i < dim-1; i++)
197  {
198  f.Set(i, new ConstantCoefficient(0.0));
199  }
200  {
201  Vector pull_force(pmesh->bdr_attributes.Max());
202  pull_force = 0.0;
203  pull_force(1) = -1.0e-2;
204  f.Set(dim-1, new PWConstCoefficient(pull_force));
205  }
206 
207  ParLinearForm *b = new ParLinearForm(fespace);
209  if (myid == 0)
210  {
211  cout << "r.h.s. ... " << flush;
212  }
213  b->Assemble();
214 
215  // 10. Define the solution vector x as a parallel finite element grid
216  // function corresponding to fespace. Initialize x with initial guess of
217  // zero, which satisfies the boundary conditions.
218  ParGridFunction x(fespace);
219  x = 0.0;
220 
221  // 11. Set up the parallel bilinear form a(.,.) on the finite element space
222  // corresponding to the linear elasticity integrator with piece-wise
223  // constants coefficient lambda and mu.
224  Vector lambda(pmesh->attributes.Max());
225  lambda = 1.0;
226  lambda(0) = lambda(1)*50;
227  PWConstCoefficient lambda_func(lambda);
228  Vector mu(pmesh->attributes.Max());
229  mu = 1.0;
230  mu(0) = mu(1)*50;
231  PWConstCoefficient mu_func(mu);
232 
233  ParBilinearForm *a = new ParBilinearForm(fespace);
234  a->AddDomainIntegrator(new ElasticityIntegrator(lambda_func, mu_func));
235 
236  // 12. Assemble the parallel bilinear form and the corresponding linear
237  // system, applying any necessary transformations such as: parallel
238  // assembly, eliminating boundary conditions, applying conforming
239  // constraints for non-conforming AMR, static condensation, etc.
240  if (myid == 0) { cout << "matrix ... " << flush; }
241  if (static_cond) { a->EnableStaticCondensation(); }
242  a->Assemble();
243 
244  HypreParMatrix A;
245  Vector B, X;
246  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
247  if (myid == 0)
248  {
249  cout << "done." << endl;
250  cout << "Size of linear system: " << A.GetGlobalNumRows() << endl;
251  }
252 
253  // 13. Define and apply a parallel PCG solver for A X = B with the BoomerAMG
254  // preconditioner from hypre.
255  HypreBoomerAMG *amg = new HypreBoomerAMG(A);
256  if (amg_elast && !a->StaticCondensationIsEnabled())
257  {
258  amg->SetElasticityOptions(fespace);
259  }
260  else
261  {
262  amg->SetSystemsOptions(dim, reorder_space);
263  }
264  HyprePCG *pcg = new HyprePCG(A);
265  pcg->SetTol(1e-8);
266  pcg->SetMaxIter(500);
267  pcg->SetPrintLevel(2);
268  pcg->SetPreconditioner(*amg);
269  pcg->Mult(B, X);
270 
271  // 14. Recover the parallel grid function corresponding to X. This is the
272  // local finite element solution on each processor.
273  a->RecoverFEMSolution(X, *b, x);
274 
275  // 15. For non-NURBS meshes, make the mesh curved based on the finite element
276  // space. This means that we define the mesh elements through a fespace
277  // based transformation of the reference element. This allows us to save
278  // the displaced mesh as a curved mesh when using high-order finite
279  // element displacement field. We assume that the initial mesh (read from
280  // the file) is not higher order curved mesh compared to the chosen FE
281  // space.
282  if (!use_nodal_fespace)
283  {
284  pmesh->SetNodalFESpace(fespace);
285  }
286 
287  // 16. Save in parallel the displaced mesh and the inverted solution (which
288  // gives the backward displacements to the original grid). This output
289  // can be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
290  {
291  GridFunction *nodes = pmesh->GetNodes();
292  *nodes += x;
293  x *= -1;
294 
295  ostringstream mesh_name, sol_name;
296  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
297  sol_name << "sol." << setfill('0') << setw(6) << myid;
298 
299  ofstream mesh_ofs(mesh_name.str().c_str());
300  mesh_ofs.precision(8);
301  pmesh->Print(mesh_ofs);
302 
303  ofstream sol_ofs(sol_name.str().c_str());
304  sol_ofs.precision(8);
305  x.Save(sol_ofs);
306  }
307 
308  // 17. Send the above data by socket to a GLVis server. Use the "n" and "b"
309  // keys in GLVis to visualize the displacements.
310  if (visualization)
311  {
312  char vishost[] = "localhost";
313  int visport = 19916;
314  socketstream sol_sock(vishost, visport);
315  sol_sock << "parallel " << num_procs << " " << myid << "\n";
316  sol_sock.precision(8);
317  sol_sock << "solution\n" << *pmesh << x << flush;
318  }
319 
320  // 18. Free the used memory.
321  delete pcg;
322  delete amg;
323  delete a;
324  delete b;
325  if (fec)
326  {
327  delete fespace;
328  delete fec;
329  }
330  delete pmesh;
331 
332  MPI_Finalize();
333 
334  return 0;
335 }
void SetTol(double tol)
Definition: hypre.cpp:2619
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:775
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(...
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:819
Abstract parallel finite element space.
Definition: pfespace.hpp:28
int main(int argc, char *argv[])
Definition: ex1.cpp:66
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2634
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1079
Class for parallel linear form.
Definition: plinearform.hpp:26
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:415
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:8382
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
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:251
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:4169
void SetElasticityOptions(ParFiniteElementSpace *fespace)
Definition: hypre.cpp:3561
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:4126
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:53
int Dimension() const
Definition: mesh.hpp:788
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:434
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2624
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:201
PCG solver in hypre.
Definition: hypre.hpp:759
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
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:203
A piecewise constant coefficient with the constants keyed off the element attribute numbers...
Definition: coefficient.hpp:94
double mu
Definition: ex25.cpp:135
void DegreeElevate(int rel_degree, int degree=16)
Definition: mesh.cpp:3948
int dim
Definition: ex24.cpp:53
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2639
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:304
Class for parallel bilinear form.
Vector data type.
Definition: vector.hpp:51
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:6603
void SetSystemsOptions(int dim, bool order_bynodes=false)
Definition: hypre.cpp:3458
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:159
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2662
Class for parallel meshes.
Definition: pmesh.hpp:32
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:199
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:145