MFEM  v4.2.0
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 // PETSc Modification
3 //
4 // Compile with: make ex11p
5 //
6 // Sample runs: mpirun -np 4 ex11p -m ../../data/star.mesh
7 // mpirun -np 4 ex11p -m ../../data/star.mesh --slepcopts rc_ex11p_lobpcg
8 // mpirun -np 4 ex11p -m ../../data/star.mesh --slepcopts rc_ex11p_gd
9 //
10 // Description: This example code demonstrates the use of MFEM to solve the
11 // eigenvalue problem -Delta u = lambda u with homogeneous
12 // Dirichlet boundary conditions.
13 //
14 // We compute a number of the lowest eigenmodes by discretizing
15 // the Laplacian and Mass operators using a FE space of the
16 // specified order, or an isoparametric/isogeometric space if
17 // order < 1 (quadratic for quadratic curvilinear mesh, NURBS for
18 // NURBS mesh, etc.)
19 //
20 // The example demonstrates the use of the SLEPc eigensolver as an
21 // alternative to the LOBPCG eigenvalue solver. The shift and
22 // invert spectral transformation is used to help the convergence
23 // to the smaller eigenvalues. Alternative solver parameters can
24 // be passed in a file with "-slepcopts".
25 //
26 // Reusing a single GLVis visualization window for multiple
27 // eigenfunctions is also illustrated.
28 //
29 // We recommend viewing Example 1 before viewing this example.
30 
31 #include "mfem.hpp"
32 #include <fstream>
33 #include <iostream>
34 
35 #ifndef MFEM_USE_SLEPC
36 #error This examples requires that MFEM is build with MFEM_USE_SLEPC=YES
37 #endif
38 
39 using namespace std;
40 using namespace mfem;
41 
42 int main(int argc, char *argv[])
43 {
44  // 1. Initialize MPI.
45  int num_procs, myid;
46  MPI_Init(&argc, &argv);
47  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
48  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
49 
50  // 2. Parse command-line options.
51  const char *mesh_file = "../../data/star.mesh";
52  int ser_ref_levels = 2;
53  int par_ref_levels = 1;
54  int order = 1;
55  int nev = 5;
56  int seed = 75;
57  bool slu_solver = false;
58  bool sp_solver = false;
59  bool visualization = 1;
60  bool use_slepc = true;
61  const char *slepcrc_file = "";
62 
63  OptionsParser args(argc, argv);
64  args.AddOption(&mesh_file, "-m", "--mesh",
65  "Mesh file to use.");
66  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
67  "Number of times to refine the mesh uniformly in serial.");
68  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
69  "Number of times to refine the mesh uniformly in parallel.");
70  args.AddOption(&order, "-o", "--order",
71  "Finite element order (polynomial degree) or -1 for"
72  " isoparametric space.");
73  args.AddOption(&nev, "-n", "--num-eigs",
74  "Number of desired eigenmodes.");
75  args.AddOption(&seed, "-s", "--seed",
76  "Random seed used to initialize LOBPCG.");
77 #ifdef MFEM_USE_SUPERLU
78  args.AddOption(&slu_solver, "-slu", "--superlu", "-no-slu",
79  "--no-superlu", "Use the SuperLU Solver.");
80 #endif
81 #ifdef MFEM_USE_STRUMPACK
82  args.AddOption(&sp_solver, "-sp", "--strumpack", "-no-sp",
83  "--no-strumpack", "Use the STRUMPACK Solver.");
84 #endif
85  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
86  "--no-visualization",
87  "Enable or disable GLVis visualization.");
88  args.AddOption(&use_slepc, "-useslepc","--useslepc","-no-slepc",
89  "--no-slepc","Use or not SLEPc to solve the eigenvalue problem");
90  args.AddOption(&slepcrc_file, "-slepcopts", "--slepcopts",
91  "SlepcOptions file to use.");
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  // 2b. We initialize SLEPc. This internally initializes PETSc as well.
121  MFEMInitializeSlepc(NULL,NULL,slepcrc_file,NULL);
122 
123  // 3. Read the (serial) mesh from the given mesh file on all processors. We
124  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
125  // and volume meshes with the same code.
126  Mesh *mesh = new Mesh(mesh_file, 1, 1);
127  int dim = mesh->Dimension();
128 
129  // 4. Refine the serial mesh on all processors to increase the resolution. In
130  // this example we do 'ref_levels' of uniform refinement (2 by default, or
131  // specified on the command line with -rs).
132  for (int lev = 0; lev < ser_ref_levels; lev++)
133  {
134  mesh->UniformRefinement();
135  }
136 
137  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
138  // this mesh further in parallel to increase the resolution (1 time by
139  // default, or specified on the command line with -rp). Once the parallel
140  // mesh is defined, the serial mesh can be deleted.
141  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
142  delete mesh;
143  for (int lev = 0; lev < par_ref_levels; lev++)
144  {
145  pmesh->UniformRefinement();
146  }
147 
148  // 6. Define a parallel finite element space on the parallel mesh. Here we
149  // use continuous Lagrange finite elements of the specified order. If
150  // order < 1, we instead use an isoparametric/isogeometric space.
152  if (order > 0)
153  {
154  fec = new H1_FECollection(order, dim);
155  }
156  else if (pmesh->GetNodes())
157  {
158  fec = pmesh->GetNodes()->OwnFEC();
159  }
160  else
161  {
162  fec = new H1_FECollection(order = 1, dim);
163  }
164  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
165  HYPRE_Int size = fespace->GlobalTrueVSize();
166  if (myid == 0)
167  {
168  cout << "Number of unknowns: " << size << endl;
169  }
170 
171  // 7. Set up the parallel bilinear forms a(.,.) and m(.,.) on the finite
172  // element space. The first corresponds to the Laplacian operator -Delta,
173  // while the second is a simple mass matrix needed on the right hand side
174  // of the generalized eigenvalue problem below. The boundary conditions
175  // are implemented by elimination with special values on the diagonal to
176  // shift the Dirichlet eigenvalues out of the computational range. After
177  // serial and parallel assembly we extract the corresponding parallel
178  // matrices A and M.
179  ConstantCoefficient one(1.0);
180  Array<int> ess_bdr;
181  if (pmesh->bdr_attributes.Size())
182  {
183  ess_bdr.SetSize(pmesh->bdr_attributes.Max());
184  ess_bdr = 1;
185  }
186 
187  ParBilinearForm *a = new ParBilinearForm(fespace);
189  if (pmesh->bdr_attributes.Size() == 0)
190  {
191  // Add a mass term if the mesh has no boundary, e.g. periodic mesh or
192  // closed surface.
193  a->AddDomainIntegrator(new MassIntegrator(one));
194  }
195  a->Assemble();
196  a->EliminateEssentialBCDiag(ess_bdr, 1.0);
197  a->Finalize();
198 
199  ParBilinearForm *m = new ParBilinearForm(fespace);
200  m->AddDomainIntegrator(new MassIntegrator(one));
201  m->Assemble();
202  // shift the eigenvalue corresponding to eliminated dofs to a large value
203  m->EliminateEssentialBCDiag(ess_bdr, numeric_limits<double>::min());
204  m->Finalize();
205 
206  PetscParMatrix *pA = NULL, *pM = NULL;
207  HypreParMatrix *A = NULL, *M = NULL;
208  Operator::Type tid =
209  !use_slepc ? Operator::Hypre_ParCSR : Operator::PETSC_MATAIJ;
210  OperatorHandle Ah(tid), Mh(tid);
211 
212  a->ParallelAssemble(Ah);
213  if (!use_slepc) { Ah.Get(A); }
214  else { Ah.Get(pA); }
215  Ah.SetOperatorOwner(false);
216 
217  m->ParallelAssemble(Mh);
218  if (!use_slepc) {Mh.Get(M); }
219  else {Mh.Get(pM); }
220  Mh.SetOperatorOwner(false);
221 
222 #if defined(MFEM_USE_SUPERLU) || defined(MFEM_USE_STRUMPACK)
223  Operator * Arow = NULL;
224 #ifdef MFEM_USE_SUPERLU
225  if (slu_solver)
226  {
227  Arow = new SuperLURowLocMatrix(*A);
228  }
229 #endif
230 #ifdef MFEM_USE_STRUMPACK
231  if (sp_solver)
232  {
233  Arow = new STRUMPACKRowLocMatrix(*A);
234  }
235 #endif
236 #endif
237 
238  delete a;
239  delete m;
240 
241  // 8. Define and configure the LOBPCG eigensolver and the BoomerAMG
242  // preconditioner for A to be used within the solver. Set the matrices
243  // which define the generalized eigenproblem A x = lambda M x.
244  Solver * precond = NULL;
245  if (!use_slepc)
246  {
247  if (!slu_solver && !sp_solver)
248  {
249  HypreBoomerAMG * amg = new HypreBoomerAMG(*A);
250  amg->SetPrintLevel(0);
251  precond = amg;
252  }
253  else
254  {
255 #ifdef MFEM_USE_SUPERLU
256  if (slu_solver)
257  {
258  SuperLUSolver * superlu = new SuperLUSolver(MPI_COMM_WORLD);
259  superlu->SetPrintStatistics(false);
260  superlu->SetSymmetricPattern(true);
262  superlu->SetOperator(*Arow);
263  precond = superlu;
264  }
265 #endif
266 #ifdef MFEM_USE_STRUMPACK
267  if (sp_solver)
268  {
269  STRUMPACKSolver * strumpack = new STRUMPACKSolver(argc, argv, MPI_COMM_WORLD);
270  strumpack->SetPrintFactorStatistics(true);
271  strumpack->SetPrintSolveStatistics(false);
272  strumpack->SetKrylovSolver(strumpack::KrylovSolver::DIRECT);
273  strumpack->SetReorderingStrategy(strumpack::ReorderingStrategy::METIS);
274  strumpack->DisableMatching();
275  strumpack->SetOperator(*Arow);
276  strumpack->SetFromCommandLine();
277  precond = strumpack;
278  }
279 #endif
280  }
281  }
282 
283  HypreLOBPCG * lobpcg = NULL;
284  SlepcEigenSolver * slepc = NULL;
285  if (!use_slepc)
286  {
287 
288  lobpcg = new HypreLOBPCG(MPI_COMM_WORLD);
289  lobpcg->SetNumModes(nev);
290  lobpcg->SetRandomSeed(seed);
291  lobpcg->SetPreconditioner(*precond);
292  lobpcg->SetMaxIter(200);
293  lobpcg->SetTol(1e-8);
294  lobpcg->SetPrecondUsageMode(1);
295  lobpcg->SetPrintLevel(1);
296  lobpcg->SetMassMatrix(*M);
297  lobpcg->SetOperator(*A);
298  }
299  else
300  {
301  slepc = new SlepcEigenSolver(MPI_COMM_WORLD);
302  slepc->SetNumModes(nev);
303  slepc->SetWhichEigenpairs(SlepcEigenSolver::TARGET_REAL);
304  slepc->SetTarget(0.0);
305  slepc->SetSpectralTransformation(SlepcEigenSolver::SHIFT_INVERT);
306  slepc->SetOperators(*pA,*pM);
307  }
308 
309  // 9. Compute the eigenmodes and extract the array of eigenvalues. Define a
310  // parallel grid function to represent each of the eigenmodes returned by
311  // the solver.
312  Array<double> eigenvalues;
313  if (!use_slepc)
314  {
315  lobpcg->Solve();
316  lobpcg->GetEigenvalues(eigenvalues);
317  }
318  else
319  {
320  slepc->Solve();
321  eigenvalues.SetSize(nev);
322  for (int i=0; i<nev; i++)
323  {
324  slepc->GetEigenvalue(i,eigenvalues[i]);
325  }
326  }
327  Vector temp(fespace->GetTrueVSize());
328  ParGridFunction x(fespace);
329 
330  // 10. Save the refined mesh and the modes in parallel. This output can be
331  // viewed later using GLVis: "glvis -np <np> -m mesh -g mode".
332  {
333  ostringstream mesh_name, mode_name;
334  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
335 
336  ofstream mesh_ofs(mesh_name.str().c_str());
337  mesh_ofs.precision(8);
338  pmesh->Print(mesh_ofs);
339 
340  for (int i=0; i<nev; i++)
341  {
342  // convert eigenvector from HypreParVector to ParGridFunction
343  if (!use_slepc)
344  {
345  x = lobpcg->GetEigenvector(i);
346  }
347  else
348  {
349  slepc->GetEigenvector(i,temp);
350  x.Distribute(temp);
351 
352  }
353 
354  mode_name << "mode_" << setfill('0') << setw(2) << i << "."
355  << setfill('0') << setw(6) << myid;
356 
357  ofstream mode_ofs(mode_name.str().c_str());
358  mode_ofs.precision(8);
359  x.Save(mode_ofs);
360  mode_name.str("");
361  }
362  }
363 
364  // 11. Send the solution by socket to a GLVis server.
365  if (visualization)
366  {
367  char vishost[] = "localhost";
368  int visport = 19916;
369  socketstream mode_sock(vishost, visport);
370  mode_sock.precision(8);
371 
372  for (int i=0; i<nev; i++)
373  {
374  if ( myid == 0 )
375  {
376  cout << "Eigenmode " << i+1 << '/' << nev
377  << ", Lambda = " << eigenvalues[i] << endl;
378  }
379 
380  // convert eigenvector from HypreParVector to ParGridFunction
381  if (!use_slepc)
382  {
383  x = lobpcg->GetEigenvector(i);
384  }
385  else
386  {
387  slepc->GetEigenvector(i,temp);
388  x.Distribute(temp);
389  }
390 
391  mode_sock << "parallel " << num_procs << " " << myid << "\n"
392  << "solution\n" << *pmesh << x << flush
393  << "window_title 'Eigenmode " << i+1 << '/' << nev
394  << ", Lambda = " << eigenvalues[i] << "'" << endl;
395 
396  char c;
397  if (myid == 0)
398  {
399  cout << "press (q)uit or (c)ontinue --> " << flush;
400  cin >> c;
401  }
402  MPI_Bcast(&c, 1, MPI_CHAR, 0, MPI_COMM_WORLD);
403 
404  if (c != 'c')
405  {
406  break;
407  }
408  }
409  mode_sock.close();
410  }
411 
412  // 12. Free the used memory.
413  delete lobpcg;
414  delete slepc;
415  delete precond;
416  delete M;
417  delete A;
418  delete pA;
419  delete pM;
420 #if defined(MFEM_USE_SUPERLU) || defined(MFEM_USE_STRUMPACK)
421  delete Arow;
422 #endif
423  delete fespace;
424  if (order > 0)
425  {
426  delete fec;
427  }
428  delete pmesh;
429 
430  // We finalize SLEPc
432  MPI_Finalize();
433 
434  return 0;
435 }
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
Definition: hypre.cpp:4335
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
void MFEMInitializeSlepc()
Definition: slepc.cpp:29
Wrapper for PETSc&#39;s matrix class.
Definition: petsc.hpp:197
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
void GetEigenvalue(unsigned int i, double &lr) const
Get the corresponding eigenvalue.
Definition: slepc.cpp:141
void SetPreconditioner(Solver &precond)
Definition: hypre.cpp:4258
void SetMassMatrix(Operator &M)
Definition: hypre.cpp:4313
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetPrintLevel(int logging)
Definition: hypre.cpp:4243
int main(int argc, char *argv[])
Definition: ex1.cpp:66
void SetTarget(double target)
Definition: slepc.cpp:223
void SetSymmetricPattern(bool sym)
Definition: superlu.cpp:420
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void SetTol(double tol)
Definition: hypre.cpp:4221
int close()
Close the socketstream.
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:255
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: superlu.cpp:576
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1079
void SetSpectralTransformation(SpectralTransformation transformation)
Definition: slepc.cpp:228
void SetPrintFactorStatistics(bool print_stat)
Definition: strumpack.cpp:134
void SetMaxIter(int max_iter)
Definition: hypre.cpp:4237
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[]
void SetColumnPermutation(superlu::ColPerm col_perm)
Definition: superlu.cpp:340
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:8382
void SetPrintLevel(int print_level)
Definition: hypre.hpp:1119
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 SetNumModes(int num_eigs)
Set the number of required eigenmodes.
Definition: slepc.cpp:118
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:4169
void SetPrintSolveStatistics(bool print_stat)
Definition: strumpack.cpp:139
void SetPrintStatistics(bool print_stat)
Definition: superlu.cpp:322
void Solve()
Solve the eigenvalue problem for the specified number of eigenvalues.
Definition: slepc.cpp:124
int Dimension() const
Definition: mesh.hpp:788
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:434
void MFEMFinalizeSlepc()
Definition: slepc.cpp:46
void SetReorderingStrategy(strumpack::ReorderingStrategy method)
Definition: strumpack.cpp:149
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:237
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:201
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
Definition: hypre.cpp:4323
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
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:654
void SetRandomSeed(int s)
Definition: hypre.hpp:1337
double a
Definition: lissajous.cpp:41
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition: handle.hpp:112
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void SetWhichEigenpairs(Which which)
Definition: slepc.cpp:189
int dim
Definition: ex24.cpp:53
void SetKrylovSolver(strumpack::KrylovSolver method)
Definition: strumpack.cpp:144
void SetOperator(Operator &A)
Definition: hypre.cpp:4267
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.
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: strumpack.cpp:228
void Get(OpType *&A) const
Return the Operator pointer statically cast to a given OpType.
Definition: handle.hpp:106
Vector data type.
Definition: vector.hpp:51
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:6603
void SetPrecondUsageMode(int pcg_mode)
Definition: hypre.cpp:4252
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:159
Base class for solvers.
Definition: operator.hpp:634
Class for parallel grid function.
Definition: pgridfunc.hpp:32
void GetEigenvector(unsigned int i, Vector &vr) const
Get the corresponding eigenvector.
Definition: slepc.cpp:152
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
void SetNumModes(int num_eigs)
Definition: hypre.hpp:1335
Class for parallel meshes.
Definition: pmesh.hpp:32
void Solve()
Solve the eigenproblem.
Definition: hypre.cpp:4380
void SetOperators(const PetscParMatrix &op, const PetscParMatrix &opB)
Set operator for generalized eigenvalue problem.
Definition: slepc.cpp:87
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:145