MFEM  v3.4
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex12p.cpp
Go to the documentation of this file.
1 // MFEM Example 12 - Parallel Version
2 //
3 // Compile with: make ex12p
4 //
5 // Sample runs:
6 // mpirun -np 4 ex12p -m ../data/beam-tri.mesh
7 // mpirun -np 4 ex12p -m ../data/beam-quad.mesh
8 // mpirun -np 4 ex12p -m ../data/beam-tet.mesh -n 10 -o 2 -elast
9 // mpirun -np 4 ex12p -m ../data/beam-hex.mesh -s 3876
10 // mpirun -np 4 ex12p -m ../data/beam-tri.mesh -o 2 -sys
11 // mpirun -np 4 ex12p -m ../data/beam-quad.mesh -s 4526 -n 6 -o 3 -elast
12 // mpirun -np 4 ex12p -m ../data/beam-quad-nurbs.mesh
13 // mpirun -np 4 ex12p -m ../data/beam-hex-nurbs.mesh
14 //
15 // Description: This example code solves the linear elasticity eigenvalue
16 // problem for a multi-material cantilever beam.
17 //
18 // Specifically, we compute a number of the lowest eigenmodes by
19 // approximating the weak form of -div(sigma(u)) = lambda u where
20 // 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. The geometry of the domain is
25 // assumed to be as follows:
26 //
27 // +----------+----------+
28 // boundary --->| material | material |
29 // attribute 1 | 1 | 2 |
30 // (fixed) +----------+----------+
31 //
32 // The example highlights the use of the LOBPCG eigenvalue solver
33 // together with the BoomerAMG preconditioner in HYPRE. Reusing a
34 // single GLVis visualization window for multiple eigenfunctions
35 // is also illustrated.
36 //
37 // We recommend viewing examples 2 and 11 before viewing this
38 // example.
39 
40 #include "mfem.hpp"
41 #include <fstream>
42 #include <iostream>
43 
44 using namespace std;
45 using namespace mfem;
46 
47 int main(int argc, char *argv[])
48 {
49  // 1. Initialize MPI.
50  int num_procs, myid;
51  MPI_Init(&argc, &argv);
52  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
53  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
54 
55  // 2. Parse command-line options.
56  const char *mesh_file = "../data/beam-tri.mesh";
57  int order = 1;
58  int nev = 5;
59  int seed = 75;
60  bool visualization = 1;
61  bool amg_elast = 0;
62 
63  OptionsParser args(argc, argv);
64  args.AddOption(&mesh_file, "-m", "--mesh",
65  "Mesh file to use.");
66  args.AddOption(&order, "-o", "--order",
67  "Finite element order (polynomial degree).");
68  args.AddOption(&nev, "-n", "--num-eigs",
69  "Number of desired eigenmodes.");
70  args.AddOption(&seed, "-s", "--seed",
71  "Random seed used to initialize LOBPCG.");
72  args.AddOption(&amg_elast, "-elast", "--amg-for-elasticity", "-sys",
73  "--amg-for-systems",
74  "Use the special AMG elasticity solver (GM/LN approaches), "
75  "or standard AMG for systems (unknown approach).");
76  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
77  "--no-visualization",
78  "Enable or disable GLVis visualization.");
79  args.Parse();
80  if (!args.Good())
81  {
82  if (myid == 0)
83  {
84  args.PrintUsage(cout);
85  }
86  MPI_Finalize();
87  return 1;
88  }
89  if (myid == 0)
90  {
91  args.PrintOptions(cout);
92  }
93 
94  // 3. Read the (serial) mesh from the given mesh file on all processors. We
95  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
96  // and volume meshes with the same code.
97  Mesh *mesh = new Mesh(mesh_file, 1, 1);
98  int dim = mesh->Dimension();
99 
100  if (mesh->attributes.Max() < 2)
101  {
102  if (myid == 0)
103  cerr << "\nInput mesh should have at least two materials!"
104  << " (See schematic in ex12p.cpp)\n"
105  << endl;
106  MPI_Finalize();
107  return 3;
108  }
109 
110  // 4. Select the order of the finite element discretization space. For NURBS
111  // meshes, we increase the order by degree elevation.
112  if (mesh->NURBSext)
113  {
114  mesh->DegreeElevate(order, order);
115  }
116 
117  // 5. Refine the serial mesh on all processors to increase the resolution. In
118  // this example we do 'ref_levels' of uniform refinement. We choose
119  // 'ref_levels' to be the largest number that gives a final mesh with no
120  // more than 1,000 elements.
121  {
122  int ref_levels =
123  (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
124  for (int l = 0; l < ref_levels; l++)
125  {
126  mesh->UniformRefinement();
127  }
128  }
129 
130  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
131  // this mesh further in parallel to increase the resolution. Once the
132  // parallel mesh is defined, the serial mesh can be deleted.
133  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
134  delete mesh;
135  {
136  int par_ref_levels = 1;
137  for (int l = 0; l < par_ref_levels; l++)
138  {
139  pmesh->UniformRefinement();
140  }
141  }
142 
143  // 7. Define a parallel finite element space on the parallel mesh. Here we
144  // use vector finite elements, i.e. dim copies of a scalar finite element
145  // space. We use the ordering by vector dimension (the last argument of
146  // the FiniteElementSpace constructor) which is expected in the systems
147  // version of BoomerAMG preconditioner. For NURBS meshes, we use the
148  // (degree elevated) NURBS space associated with the mesh nodes.
150  ParFiniteElementSpace *fespace;
151  const bool use_nodal_fespace = pmesh->NURBSext && !amg_elast;
152  if (use_nodal_fespace)
153  {
154  fec = NULL;
155  fespace = (ParFiniteElementSpace *)pmesh->GetNodes()->FESpace();
156  }
157  else
158  {
159  fec = new H1_FECollection(order, dim);
160  fespace = new ParFiniteElementSpace(pmesh, fec, dim, Ordering::byVDIM);
161  }
162  HYPRE_Int size = fespace->GlobalTrueVSize();
163  if (myid == 0)
164  {
165  cout << "Number of unknowns: " << size << endl
166  << "Assembling: " << flush;
167  }
168 
169  // 8. Set up the parallel bilinear forms a(.,.) and m(.,.) on the finite
170  // element space corresponding to the linear elasticity integrator with
171  // piece-wise constants coefficient lambda and mu, a simple mass matrix
172  // needed on the right hand side of the generalized eigenvalue problem
173  // below. The boundary conditions are implemented by marking only boundary
174  // attribute 1 as essential. We use special values on the diagonal to
175  // shift the Dirichlet eigenvalues out of the computational range. After
176  // serial/parallel assembly we extract the corresponding parallel matrices
177  // A and M.
178  Vector lambda(pmesh->attributes.Max());
179  lambda = 1.0;
180  lambda(0) = lambda(1)*50;
181  PWConstCoefficient lambda_func(lambda);
182  Vector mu(pmesh->attributes.Max());
183  mu = 1.0;
184  mu(0) = mu(1)*50;
185  PWConstCoefficient mu_func(mu);
186 
187  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
188  ess_bdr = 0;
189  ess_bdr[0] = 1;
190 
191  ParBilinearForm *a = new ParBilinearForm(fespace);
192  a->AddDomainIntegrator(new ElasticityIntegrator(lambda_func, mu_func));
193  if (myid == 0)
194  {
195  cout << "matrix ... " << flush;
196  }
197  a->Assemble();
198  a->EliminateEssentialBCDiag(ess_bdr, 1.0);
199  a->Finalize();
200 
201  ParBilinearForm *m = new ParBilinearForm(fespace);
203  m->Assemble();
204  // shift the eigenvalue corresponding to eliminated dofs to a large value
205  m->EliminateEssentialBCDiag(ess_bdr, numeric_limits<double>::min());
206  m->Finalize();
207  if (myid == 0)
208  {
209  cout << "done." << endl;
210  }
211 
214 
215  delete a;
216  delete m;
217 
218  // 9. Define and configure the LOBPCG eigensolver and the BoomerAMG
219  // preconditioner for A to be used within the solver. Set the matrices
220  // which define the generalized eigenproblem A x = lambda M x.
221  HypreBoomerAMG * amg = new HypreBoomerAMG(*A);
222  amg->SetPrintLevel(0);
223  if (amg_elast)
224  {
225  amg->SetElasticityOptions(fespace);
226  }
227  else
228  {
229  amg->SetSystemsOptions(dim);
230  }
231 
232  HypreLOBPCG * lobpcg = new HypreLOBPCG(MPI_COMM_WORLD);
233  lobpcg->SetNumModes(nev);
234  lobpcg->SetRandomSeed(seed);
235  lobpcg->SetPreconditioner(*amg);
236  lobpcg->SetMaxIter(100);
237  lobpcg->SetTol(1e-8);
238  lobpcg->SetPrecondUsageMode(1);
239  lobpcg->SetPrintLevel(1);
240  lobpcg->SetMassMatrix(*M);
241  lobpcg->SetOperator(*A);
242 
243  // 10. Compute the eigenmodes and extract the array of eigenvalues. Define a
244  // parallel grid function to represent each of the eigenmodes returned by
245  // the solver.
246  Array<double> eigenvalues;
247  lobpcg->Solve();
248  lobpcg->GetEigenvalues(eigenvalues);
249  ParGridFunction x(fespace);
250 
251  // 11. For non-NURBS meshes, make the mesh curved based on the finite element
252  // space. This means that we define the mesh elements through a fespace
253  // based transformation of the reference element. This allows us to save
254  // the displaced mesh as a curved mesh when using high-order finite
255  // element displacement field. We assume that the initial mesh (read from
256  // the file) is not higher order curved mesh compared to the chosen FE
257  // space.
258  if (!use_nodal_fespace)
259  {
260  pmesh->SetNodalFESpace(fespace);
261  }
262 
263  // 12. Save the refined mesh and the modes in parallel. This output can be
264  // viewed later using GLVis: "glvis -np <np> -m mesh -g mode".
265  {
266  ostringstream mesh_name, mode_name;
267  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
268 
269  ofstream mesh_ofs(mesh_name.str().c_str());
270  mesh_ofs.precision(8);
271  pmesh->Print(mesh_ofs);
272 
273  for (int i=0; i<nev; i++)
274  {
275  // convert eigenvector from HypreParVector to ParGridFunction
276  x = lobpcg->GetEigenvector(i);
277 
278  mode_name << "mode_" << setfill('0') << setw(2) << i << "."
279  << setfill('0') << setw(6) << myid;
280 
281  ofstream mode_ofs(mode_name.str().c_str());
282  mode_ofs.precision(8);
283  x.Save(mode_ofs);
284  mode_name.str("");
285  }
286  }
287 
288  // 13. Send the above data by socket to a GLVis server. Use the "n" and "b"
289  // keys in GLVis to visualize the displacements.
290  if (visualization)
291  {
292  char vishost[] = "localhost";
293  int visport = 19916;
294  socketstream mode_sock(vishost, visport);
295 
296  for (int i=0; i<nev; i++)
297  {
298  if ( myid == 0 )
299  {
300  cout << "Eigenmode " << i+1 << '/' << nev
301  << ", Lambda = " << eigenvalues[i] << endl;
302  }
303 
304  // convert eigenvector from HypreParVector to ParGridFunction
305  x = lobpcg->GetEigenvector(i);
306 
307  mode_sock << "parallel " << num_procs << " " << myid << "\n"
308  << "solution\n" << *pmesh << x << flush
309  << "window_title 'Eigenmode " << i+1 << '/' << nev
310  << ", Lambda = " << eigenvalues[i] << "'" << endl;
311 
312  char c;
313  if (myid == 0)
314  {
315  cout << "press (q)uit or (c)ontinue --> " << flush;
316  cin >> c;
317  }
318  MPI_Bcast(&c, 1, MPI_CHAR, 0, MPI_COMM_WORLD);
319 
320  if (c != 'c')
321  {
322  break;
323  }
324  }
325  mode_sock.close();
326  }
327 
328  // 14. Free the used memory.
329  delete lobpcg;
330  delete amg;
331  delete M;
332  delete A;
333 
334  if (fec)
335  {
336  delete fespace;
337  delete fec;
338  }
339  delete pmesh;
340 
341  MPI_Finalize();
342 
343  return 0;
344 }
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
Definition: hypre.cpp:3371
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:618
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:398
void SetPreconditioner(Solver &precond)
Definition: hypre.cpp:3294
void SetMassMatrix(Operator &M)
Definition: hypre.cpp:3349
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetPrintLevel(int logging)
Definition: hypre.cpp:3279
int main(int argc, char *argv[])
Definition: ex1.cpp:45
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void SetTol(double tol)
Definition: hypre.cpp:3257
void SetSystemsOptions(int dim)
Definition: hypre.cpp:2584
The BoomerAMG solver in hypre.
Definition: hypre.hpp:796
int dim
Definition: ex3.cpp:47
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3273
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6741
void SetPrintLevel(int print_level)
Definition: hypre.hpp:837
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:109
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:247
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void EliminateEssentialBCDiag(const Array< int > &bdr_attr_is_ess, double value)
Perform elimination and set the diagonal entry to the given value.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:3428
void SetElasticityOptions(ParFiniteElementSpace *fespace)
Definition: hypre.cpp:2666
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:3338
int Dimension() const
Definition: mesh.hpp:645
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:174
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
Definition: hypre.cpp:3359
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
void SetRandomSeed(int s)
Definition: hypre.hpp:1041
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:176
class for piecewise constant coefficient
Definition: coefficient.hpp:72
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void DegreeElevate(int rel_degree, int degree=16)
Definition: mesh.cpp:3120
void SetOperator(Operator &A)
Definition: hypre.cpp:3303
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
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:5422
void SetPrecondUsageMode(int pcg_mode)
Definition: hypre.cpp:3288
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:79
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:175
void SetNumModes(int num_eigs)
Definition: hypre.hpp:1039
Class for parallel meshes.
Definition: pmesh.hpp:32
void Solve()
Solve the eigenproblem.
Definition: hypre.cpp:3416
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:172
bool Good() const
Definition: optparser.hpp:120