MFEM  v3.4
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex11p.cpp
Go to the documentation of this file.
1 // MFEM Example 11 - Parallel Version
2 //
3 // Compile with: make ex11p
4 //
5 // Sample runs: mpirun -np 4 ex11p -m ../data/square-disc.mesh
6 // mpirun -np 4 ex11p -m ../data/star.mesh
7 // mpirun -np 4 ex11p -m ../data/escher.mesh
8 // mpirun -np 4 ex11p -m ../data/fichera.mesh
9 // mpirun -np 4 ex11p -m ../data/square-disc-p2.vtk -o 2
10 // mpirun -np 4 ex11p -m ../data/square-disc-p3.mesh -o 3
11 // mpirun -np 4 ex11p -m ../data/square-disc-nurbs.mesh -o -1
12 // mpirun -np 4 ex11p -m ../data/disc-nurbs.mesh -o -1 -n 20
13 // mpirun -np 4 ex11p -m ../data/pipe-nurbs.mesh -o -1
14 // mpirun -np 4 ex11p -m ../data/ball-nurbs.mesh -o 2
15 // mpirun -np 4 ex11p -m ../data/star-surf.mesh
16 // mpirun -np 4 ex11p -m ../data/square-disc-surf.mesh
17 // mpirun -np 4 ex11p -m ../data/inline-segment.mesh
18 // mpirun -np 4 ex11p -m ../data/amr-quad.mesh
19 // mpirun -np 4 ex11p -m ../data/amr-hex.mesh
20 // mpirun -np 4 ex11p -m ../data/mobius-strip.mesh -n 8
21 // mpirun -np 4 ex11p -m ../data/klein-bottle.mesh -n 10
22 //
23 // Description: This example code demonstrates the use of MFEM to solve the
24 // eigenvalue problem -Delta u = lambda u with homogeneous
25 // Dirichlet boundary conditions.
26 //
27 // We compute a number of the lowest eigenmodes by discretizing
28 // the Laplacian and Mass operators using a FE space of the
29 // specified order, or an isoparametric/isogeometric space if
30 // order < 1 (quadratic for quadratic curvilinear mesh, NURBS for
31 // NURBS mesh, etc.)
32 //
33 // The example highlights the use of the LOBPCG eigenvalue solver
34 // together with the BoomerAMG preconditioner in HYPRE, as well as
35 // optionally the SuperLU or STRUMPACK parallel direct solvers.
36 // Reusing a single GLVis visualization window for multiple
37 // eigenfunctions is also illustrated.
38 //
39 // We recommend viewing Example 1 before viewing this example.
40 
41 #include "mfem.hpp"
42 #include <fstream>
43 #include <iostream>
44 
45 using namespace std;
46 using namespace mfem;
47 
48 int main(int argc, char *argv[])
49 {
50  // 1. Initialize MPI.
51  int num_procs, myid;
52  MPI_Init(&argc, &argv);
53  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
54  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
55 
56  // 2. Parse command-line options.
57  const char *mesh_file = "../data/star.mesh";
58  int ser_ref_levels = 2;
59  int par_ref_levels = 1;
60  int order = 1;
61  int nev = 5;
62  int seed = 75;
63  bool slu_solver = false;
64  bool sp_solver = false;
65  bool visualization = 1;
66 
67  OptionsParser args(argc, argv);
68  args.AddOption(&mesh_file, "-m", "--mesh",
69  "Mesh file to use.");
70  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
71  "Number of times to refine the mesh uniformly in serial.");
72  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
73  "Number of times to refine the mesh uniformly in parallel.");
74  args.AddOption(&order, "-o", "--order",
75  "Finite element order (polynomial degree) or -1 for"
76  " isoparametric space.");
77  args.AddOption(&nev, "-n", "--num-eigs",
78  "Number of desired eigenmodes.");
79  args.AddOption(&seed, "-s", "--seed",
80  "Random seed used to initialize LOBPCG.");
81 #ifdef MFEM_USE_SUPERLU
82  args.AddOption(&slu_solver, "-slu", "--superlu", "-no-slu",
83  "--no-superlu", "Use the SuperLU Solver.");
84 #endif
85 #ifdef MFEM_USE_STRUMPACK
86  args.AddOption(&sp_solver, "-sp", "--strumpack", "-no-sp",
87  "--no-strumpack", "Use the STRUMPACK Solver.");
88 #endif
89  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
90  "--no-visualization",
91  "Enable or disable GLVis visualization.");
92  args.Parse();
93  if (slu_solver && sp_solver)
94  {
95  if (myid == 0)
96  cout << "WARNING: Both SuperLU and STRUMPACK have been selected,"
97  << " please choose either one." << endl
98  << " Defaulting to SuperLU." << endl;
99  sp_solver = false;
100  }
101  // The command line options are also passed to the STRUMPACK
102  // solver. So do not exit if some options are not recognized.
103  if (!sp_solver)
104  {
105  if (!args.Good())
106  {
107  if (myid == 0)
108  {
109  args.PrintUsage(cout);
110  }
111  MPI_Finalize();
112  return 1;
113  }
114  }
115  if (myid == 0)
116  {
117  args.PrintOptions(cout);
118  }
119 
120  // 3. Read the (serial) mesh from the given mesh file on all processors. We
121  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
122  // and volume meshes with the same code.
123  Mesh *mesh = new Mesh(mesh_file, 1, 1);
124  int dim = mesh->Dimension();
125 
126  // 4. Refine the serial mesh on all processors to increase the resolution. In
127  // this example we do 'ref_levels' of uniform refinement (2 by default, or
128  // specified on the command line with -rs).
129  for (int lev = 0; lev < ser_ref_levels; lev++)
130  {
131  mesh->UniformRefinement();
132  }
133 
134  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
135  // this mesh further in parallel to increase the resolution (1 time by
136  // default, or specified on the command line with -rp). Once the parallel
137  // mesh is defined, the serial mesh can be deleted.
138  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
139  delete mesh;
140  for (int lev = 0; lev < par_ref_levels; lev++)
141  {
142  pmesh->UniformRefinement();
143  }
144 
145  // 6. Define a parallel finite element space on the parallel mesh. Here we
146  // use continuous Lagrange finite elements of the specified order. If
147  // order < 1, we instead use an isoparametric/isogeometric space.
149  if (order > 0)
150  {
151  fec = new H1_FECollection(order, dim);
152  }
153  else if (pmesh->GetNodes())
154  {
155  fec = pmesh->GetNodes()->OwnFEC();
156  }
157  else
158  {
159  fec = new H1_FECollection(order = 1, dim);
160  }
161  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
162  HYPRE_Int size = fespace->GlobalTrueVSize();
163  if (myid == 0)
164  {
165  cout << "Number of unknowns: " << size << endl;
166  }
167 
168  // 7. Set up the parallel bilinear forms a(.,.) and m(.,.) on the finite
169  // element space. The first corresponds to the Laplacian operator -Delta,
170  // while the second is a simple mass matrix needed on the right hand side
171  // of the generalized eigenvalue problem below. The boundary conditions
172  // are implemented by elimination with special values on the diagonal to
173  // shift the Dirichlet eigenvalues out of the computational range. After
174  // serial and parallel assembly we extract the corresponding parallel
175  // matrices A and M.
176  ConstantCoefficient one(1.0);
177  Array<int> ess_bdr;
178  if (pmesh->bdr_attributes.Size())
179  {
180  ess_bdr.SetSize(pmesh->bdr_attributes.Max());
181  ess_bdr = 1;
182  }
183 
184  ParBilinearForm *a = new ParBilinearForm(fespace);
186  if (pmesh->bdr_attributes.Size() == 0)
187  {
188  // Add a mass term if the mesh has no boundary, e.g. periodic mesh or
189  // closed surface.
190  a->AddDomainIntegrator(new MassIntegrator(one));
191  }
192  a->Assemble();
193  a->EliminateEssentialBCDiag(ess_bdr, 1.0);
194  a->Finalize();
195 
196  ParBilinearForm *m = new ParBilinearForm(fespace);
197  m->AddDomainIntegrator(new MassIntegrator(one));
198  m->Assemble();
199  // shift the eigenvalue corresponding to eliminated dofs to a large value
200  m->EliminateEssentialBCDiag(ess_bdr, numeric_limits<double>::min());
201  m->Finalize();
202 
205 
206 #if defined(MFEM_USE_SUPERLU) || defined(MFEM_USE_STRUMPACK)
207  Operator * Arow = NULL;
208 #ifdef MFEM_USE_SUPERLU
209  if (slu_solver)
210  {
211  Arow = new SuperLURowLocMatrix(*A);
212  }
213 #endif
214 #ifdef MFEM_USE_STRUMPACK
215  if (sp_solver)
216  {
217  Arow = new STRUMPACKRowLocMatrix(*A);
218  }
219 #endif
220 #endif
221 
222  delete a;
223  delete m;
224 
225  // 8. Define and configure the LOBPCG eigensolver and the BoomerAMG
226  // preconditioner for A to be used within the solver. Set the matrices
227  // which define the generalized eigenproblem A x = lambda M x.
228  Solver * precond = NULL;
229  if (!slu_solver && !sp_solver)
230  {
231  HypreBoomerAMG * amg = new HypreBoomerAMG(*A);
232  amg->SetPrintLevel(0);
233  precond = amg;
234  }
235  else
236  {
237 #ifdef MFEM_USE_SUPERLU
238  if (slu_solver)
239  {
240  SuperLUSolver * superlu = new SuperLUSolver(MPI_COMM_WORLD);
241  superlu->SetPrintStatistics(false);
242  superlu->SetSymmetricPattern(true);
244  superlu->SetOperator(*Arow);
245  precond = superlu;
246  }
247 #endif
248 #ifdef MFEM_USE_STRUMPACK
249  if (sp_solver)
250  {
251  STRUMPACKSolver * strumpack = new STRUMPACKSolver(argc, argv, MPI_COMM_WORLD);
252  strumpack->SetPrintFactorStatistics(true);
253  strumpack->SetPrintSolveStatistics(false);
254  strumpack->SetKrylovSolver(strumpack::KrylovSolver::DIRECT);
255  strumpack->SetReorderingStrategy(strumpack::ReorderingStrategy::METIS);
256  strumpack->SetMC64Job(strumpack::MC64Job::NONE);
257  // strumpack->SetSymmetricPattern(true);
258  strumpack->SetOperator(*Arow);
259  strumpack->SetFromCommandLine();
260  precond = strumpack;
261  }
262 #endif
263  }
264 
265 
266  HypreLOBPCG * lobpcg = new HypreLOBPCG(MPI_COMM_WORLD);
267  lobpcg->SetNumModes(nev);
268  lobpcg->SetRandomSeed(seed);
269  lobpcg->SetPreconditioner(*precond);
270  lobpcg->SetMaxIter(200);
271  lobpcg->SetTol(1e-8);
272  lobpcg->SetPrecondUsageMode(1);
273  lobpcg->SetPrintLevel(1);
274  lobpcg->SetMassMatrix(*M);
275  lobpcg->SetOperator(*A);
276 
277  // 9. Compute the eigenmodes and extract the array of eigenvalues. Define a
278  // parallel grid function to represent each of the eigenmodes returned by
279  // the solver.
280  Array<double> eigenvalues;
281  lobpcg->Solve();
282  lobpcg->GetEigenvalues(eigenvalues);
283  ParGridFunction x(fespace);
284 
285  // 10. Save the refined mesh and the modes in parallel. This output can be
286  // viewed later using GLVis: "glvis -np <np> -m mesh -g mode".
287  {
288  ostringstream mesh_name, mode_name;
289  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
290 
291  ofstream mesh_ofs(mesh_name.str().c_str());
292  mesh_ofs.precision(8);
293  pmesh->Print(mesh_ofs);
294 
295  for (int i=0; i<nev; i++)
296  {
297  // convert eigenvector from HypreParVector to ParGridFunction
298  x = lobpcg->GetEigenvector(i);
299 
300  mode_name << "mode_" << setfill('0') << setw(2) << i << "."
301  << setfill('0') << setw(6) << myid;
302 
303  ofstream mode_ofs(mode_name.str().c_str());
304  mode_ofs.precision(8);
305  x.Save(mode_ofs);
306  mode_name.str("");
307  }
308  }
309 
310  // 11. Send the solution by socket to a GLVis server.
311  if (visualization)
312  {
313  char vishost[] = "localhost";
314  int visport = 19916;
315  socketstream mode_sock(vishost, visport);
316  mode_sock.precision(8);
317 
318  for (int i=0; i<nev; i++)
319  {
320  if ( myid == 0 )
321  {
322  cout << "Eigenmode " << i+1 << '/' << nev
323  << ", Lambda = " << eigenvalues[i] << endl;
324  }
325 
326  // convert eigenvector from HypreParVector to ParGridFunction
327  x = lobpcg->GetEigenvector(i);
328 
329  mode_sock << "parallel " << num_procs << " " << myid << "\n"
330  << "solution\n" << *pmesh << x << flush
331  << "window_title 'Eigenmode " << i+1 << '/' << nev
332  << ", Lambda = " << eigenvalues[i] << "'" << endl;
333 
334  char c;
335  if (myid == 0)
336  {
337  cout << "press (q)uit or (c)ontinue --> " << flush;
338  cin >> c;
339  }
340  MPI_Bcast(&c, 1, MPI_CHAR, 0, MPI_COMM_WORLD);
341 
342  if (c != 'c')
343  {
344  break;
345  }
346  }
347  mode_sock.close();
348  }
349 
350  // 12. Free the used memory.
351  delete lobpcg;
352  delete precond;
353  delete M;
354  delete A;
355 #if defined(MFEM_USE_SUPERLU) || defined(MFEM_USE_STRUMPACK)
356  delete Arow;
357 #endif
358 
359  delete fespace;
360  if (order > 0)
361  {
362  delete fec;
363  }
364  delete pmesh;
365 
366  MPI_Finalize();
367 
368  return 0;
369 }
int Size() const
Logical size of the array.
Definition: array.hpp:133
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
Definition: hypre.cpp:3371
void SetMC64Job(strumpack::MC64Job job)
Definition: strumpack.cpp:148
Subclass constant coefficient.
Definition: coefficient.hpp:57
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
void SetSymmetricPattern(bool sym)
Definition: superlu.cpp:395
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 SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: superlu.cpp:529
The BoomerAMG solver in hypre.
Definition: hypre.hpp:796
int dim
Definition: ex3.cpp:47
void SetPrintFactorStatistics(bool print_stat)
Definition: strumpack.cpp:127
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3273
void SetColumnPermutation(superlu::ColPerm col_perm)
Definition: superlu.cpp:315
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 SetPrintSolveStatistics(bool print_stat)
Definition: strumpack.cpp:132
void SetPrintStatistics(bool print_stat)
Definition: superlu.cpp:297
int Dimension() const
Definition: mesh.hpp:645
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetReorderingStrategy(strumpack::ReorderingStrategy method)
Definition: strumpack.cpp:142
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 SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:569
void SetRandomSeed(int s)
Definition: hypre.hpp:1041
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void SetKrylovSolver(strumpack::KrylovSolver method)
Definition: strumpack.cpp:137
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.
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: strumpack.cpp:198
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
Base class for solvers.
Definition: operator.hpp:268
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Abstract operator.
Definition: operator.hpp:21
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
bool Good() const
Definition: optparser.hpp:120