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