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