MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
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
39using namespace std;
40using namespace mfem;
41
42int main(int argc, char *argv[])
43{
44 // 1. Initialize MPI and HYPRE.
45 Mpi::Init(argc, argv);
46 int num_procs = Mpi::WorldSize();
47 int myid = Mpi::WorldRank();
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 const char *device_config = "cpu";
63
64 OptionsParser args(argc, argv);
65 args.AddOption(&mesh_file, "-m", "--mesh",
66 "Mesh file to use.");
67 args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
68 "Number of times to refine the mesh uniformly in serial.");
69 args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
70 "Number of times to refine the mesh uniformly in parallel.");
71 args.AddOption(&order, "-o", "--order",
72 "Finite element order (polynomial degree) or -1 for"
73 " isoparametric space.");
74 args.AddOption(&nev, "-n", "--num-eigs",
75 "Number of desired eigenmodes.");
76 args.AddOption(&seed, "-s", "--seed",
77 "Random seed used to initialize LOBPCG.");
78#ifdef MFEM_USE_SUPERLU
79 args.AddOption(&slu_solver, "-slu", "--superlu", "-no-slu",
80 "--no-superlu", "Use the SuperLU Solver.");
81#endif
82#ifdef MFEM_USE_STRUMPACK
83 args.AddOption(&sp_solver, "-sp", "--strumpack", "-no-sp",
84 "--no-strumpack", "Use the STRUMPACK Solver.");
85#endif
86 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
87 "--no-visualization",
88 "Enable or disable GLVis visualization.");
89 args.AddOption(&use_slepc, "-useslepc","--useslepc","-no-slepc",
90 "--no-slepc","Use or not SLEPc to solve the eigenvalue problem");
91 args.AddOption(&slepcrc_file, "-slepcopts", "--slepcopts",
92 "SlepcOptions file to use.");
93 args.AddOption(&device_config, "-d", "--device",
94 "Device configuration string, see Device::Configure().");
95 args.Parse();
96 if (slu_solver && sp_solver)
97 {
98 if (myid == 0)
99 cout << "WARNING: Both SuperLU and STRUMPACK have been selected,"
100 << " please choose either one." << endl
101 << " Defaulting to SuperLU." << endl;
102 sp_solver = false;
103 }
104 // The command line options are also passed to the STRUMPACK
105 // solver. So do not exit if some options are not recognized.
106 if (!sp_solver)
107 {
108 if (!args.Good())
109 {
110 if (myid == 0)
111 {
112 args.PrintUsage(cout);
113 }
114 return 1;
115 }
116 }
117 if (myid == 0)
118 {
119 args.PrintOptions(cout);
120 }
121
122 // 2b. Enable hardware devices such as GPUs, and programming models such as
123 // CUDA, OCCA, RAJA and OpenMP based on command line options.
124 Device device(device_config);
125 if (myid == 0) { device.Print(); }
126
127 // 2c. We initialize SLEPc. This internally initializes PETSc as well.
128 MFEMInitializeSlepc(NULL,NULL,slepcrc_file,NULL);
129
130 // 3. Read the (serial) mesh from the given mesh file on all processors. We
131 // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
132 // and volume meshes with the same code.
133 Mesh *mesh = new Mesh(mesh_file, 1, 1);
134 int dim = mesh->Dimension();
135
136 // 4. Refine the serial mesh on all processors to increase the resolution. In
137 // this example we do 'ref_levels' of uniform refinement (2 by default, or
138 // specified on the command line with -rs).
139 for (int lev = 0; lev < ser_ref_levels; lev++)
140 {
141 mesh->UniformRefinement();
142 }
143
144 // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
145 // this mesh further in parallel to increase the resolution (1 time by
146 // default, or specified on the command line with -rp). Once the parallel
147 // mesh is defined, the serial mesh can be deleted.
148 ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
149 delete mesh;
150 for (int lev = 0; lev < par_ref_levels; lev++)
151 {
152 pmesh->UniformRefinement();
153 }
154
155 // 6. Define a parallel finite element space on the parallel mesh. Here we
156 // use continuous Lagrange finite elements of the specified order. If
157 // order < 1, we instead use an isoparametric/isogeometric space.
159 if (order > 0)
160 {
161 fec = new H1_FECollection(order, dim);
162 }
163 else if (pmesh->GetNodes())
164 {
165 fec = pmesh->GetNodes()->OwnFEC();
166 }
167 else
168 {
169 fec = new H1_FECollection(order = 1, dim);
170 }
171 ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
172 HYPRE_BigInt size = fespace->GlobalTrueVSize();
173 if (myid == 0)
174 {
175 cout << "Number of unknowns: " << size << endl;
176 }
177
178 // 7. Set up the parallel bilinear forms a(.,.) and m(.,.) on the finite
179 // element space. The first corresponds to the Laplacian operator -Delta,
180 // while the second is a simple mass matrix needed on the right hand side
181 // of the generalized eigenvalue problem below. The boundary conditions
182 // are implemented by elimination with special values on the diagonal to
183 // shift the Dirichlet eigenvalues out of the computational range. After
184 // serial and parallel assembly we extract the corresponding parallel
185 // matrices A and M.
186 ConstantCoefficient one(1.0);
187 Array<int> ess_bdr;
188 if (pmesh->bdr_attributes.Size())
189 {
190 ess_bdr.SetSize(pmesh->bdr_attributes.Max());
191 ess_bdr = 1;
192 }
193
194 ParBilinearForm *a = new ParBilinearForm(fespace);
195 a->AddDomainIntegrator(new DiffusionIntegrator(one));
196 if (pmesh->bdr_attributes.Size() == 0)
197 {
198 // Add a mass term if the mesh has no boundary, e.g. periodic mesh or
199 // closed surface.
200 a->AddDomainIntegrator(new MassIntegrator(one));
201 }
202 a->Assemble();
203 a->EliminateEssentialBCDiag(ess_bdr, 1.0);
204 a->Finalize();
205
206 ParBilinearForm *m = new ParBilinearForm(fespace);
208 m->Assemble();
209 // shift the eigenvalue corresponding to eliminated dofs to a large value
210 m->EliminateEssentialBCDiag(ess_bdr, numeric_limits<real_t>::min());
211 m->Finalize();
212
213 PetscParMatrix *pA = NULL, *pM = NULL;
214 HypreParMatrix *A = NULL, *M = NULL;
215 Operator::Type tid =
217 OperatorHandle Ah(tid), Mh(tid);
218
219 a->ParallelAssemble(Ah);
220 if (!use_slepc) { Ah.Get(A); }
221 else { Ah.Get(pA); }
222 Ah.SetOperatorOwner(false);
223
224 m->ParallelAssemble(Mh);
225 if (!use_slepc) {Mh.Get(M); }
226 else {Mh.Get(pM); }
227 Mh.SetOperatorOwner(false);
228
229#if defined(MFEM_USE_SUPERLU) || defined(MFEM_USE_STRUMPACK)
230 Operator * Arow = NULL;
231#ifdef MFEM_USE_SUPERLU
232 if (slu_solver)
233 {
234 Arow = new SuperLURowLocMatrix(*A);
235 }
236#endif
237#ifdef MFEM_USE_STRUMPACK
238 if (sp_solver)
239 {
240 Arow = new STRUMPACKRowLocMatrix(*A);
241 }
242#endif
243#endif
244
245 delete a;
246 delete m;
247
248 // 8. Define and configure the LOBPCG eigensolver and the BoomerAMG
249 // preconditioner for A to be used within the solver. Set the matrices
250 // which define the generalized eigenproblem A x = lambda M x.
251 Solver * precond = NULL;
252 if (!use_slepc)
253 {
254 if (!slu_solver && !sp_solver)
255 {
256 HypreBoomerAMG * amg = new HypreBoomerAMG(*A);
257 amg->SetPrintLevel(0);
258 precond = amg;
259 }
260 else
261 {
262#ifdef MFEM_USE_SUPERLU
263 if (slu_solver)
264 {
265 SuperLUSolver * superlu = new SuperLUSolver(MPI_COMM_WORLD);
266 superlu->SetPrintStatistics(false);
267 superlu->SetSymmetricPattern(true);
269 superlu->SetOperator(*Arow);
270 precond = superlu;
271 }
272#endif
273#ifdef MFEM_USE_STRUMPACK
274 if (sp_solver)
275 {
276 STRUMPACKSolver * strumpack = new STRUMPACKSolver(MPI_COMM_WORLD, argc, argv);
277 strumpack->SetPrintFactorStatistics(true);
278 strumpack->SetPrintSolveStatistics(false);
279 strumpack->SetKrylovSolver(strumpack::KrylovSolver::DIRECT);
280 strumpack->SetReorderingStrategy(strumpack::ReorderingStrategy::METIS);
281 strumpack->SetMatching(strumpack::MatchingJob::NONE);
282 strumpack->SetCompression(strumpack::CompressionType::NONE);
283 strumpack->SetOperator(*Arow);
284 strumpack->SetFromCommandLine();
285 precond = strumpack;
286 }
287#endif
288 }
289 }
290
291 HypreLOBPCG * lobpcg = NULL;
292 SlepcEigenSolver * slepc = NULL;
293 if (!use_slepc)
294 {
295
296 lobpcg = new HypreLOBPCG(MPI_COMM_WORLD);
297 lobpcg->SetNumModes(nev);
298 lobpcg->SetRandomSeed(seed);
299 lobpcg->SetPreconditioner(*precond);
300 lobpcg->SetMaxIter(200);
301 lobpcg->SetTol(1e-8);
302 lobpcg->SetPrecondUsageMode(1);
303 lobpcg->SetPrintLevel(1);
304 lobpcg->SetMassMatrix(*M);
305 lobpcg->SetOperator(*A);
306 }
307 else
308 {
309 slepc = new SlepcEigenSolver(MPI_COMM_WORLD);
310 slepc->SetNumModes(nev);
312 slepc->SetTarget(0.0);
314 slepc->SetOperators(*pA,*pM);
315 }
316
317 // 9. Compute the eigenmodes and extract the array of eigenvalues. Define a
318 // parallel grid function to represent each of the eigenmodes returned by
319 // the solver.
320 Array<real_t> eigenvalues;
321 if (!use_slepc)
322 {
323 lobpcg->Solve();
324 lobpcg->GetEigenvalues(eigenvalues);
325 }
326 else
327 {
328 slepc->Solve();
329 eigenvalues.SetSize(nev);
330 for (int i=0; i<nev; i++)
331 {
332 slepc->GetEigenvalue(i,eigenvalues[i]);
333 }
334 }
335 Vector temp(fespace->GetTrueVSize());
336 ParGridFunction x(fespace);
337
338 // 10. Save the refined mesh and the modes in parallel. This output can be
339 // viewed later using GLVis: "glvis -np <np> -m mesh -g mode".
340 {
341 ostringstream mesh_name, mode_name;
342 mesh_name << "mesh." << setfill('0') << setw(6) << myid;
343
344 ofstream mesh_ofs(mesh_name.str().c_str());
345 mesh_ofs.precision(8);
346 pmesh->Print(mesh_ofs);
347
348 for (int i=0; i<nev; i++)
349 {
350 // convert eigenvector from HypreParVector to ParGridFunction
351 if (!use_slepc)
352 {
353 x = lobpcg->GetEigenvector(i);
354 }
355 else
356 {
357 slepc->GetEigenvector(i,temp);
358 x.Distribute(temp);
359 }
360
361 mode_name << "mode_" << setfill('0') << setw(2) << i << "."
362 << setfill('0') << setw(6) << myid;
363
364 ofstream mode_ofs(mode_name.str().c_str());
365 mode_ofs.precision(8);
366 x.Save(mode_ofs);
367 mode_name.str("");
368 }
369 }
370
371 // 11. Send the solution by socket to a GLVis server.
372 if (visualization)
373 {
374 char vishost[] = "localhost";
375 int visport = 19916;
376 socketstream mode_sock(vishost, visport);
377 mode_sock.precision(8);
378
379 for (int i=0; i<nev; i++)
380 {
381 if ( myid == 0 )
382 {
383 cout << "Eigenmode " << i+1 << '/' << nev
384 << ", Lambda = " << eigenvalues[i] << endl;
385 }
386
387 // convert eigenvector from HypreParVector to ParGridFunction
388 if (!use_slepc)
389 {
390 x = lobpcg->GetEigenvector(i);
391 }
392 else
393 {
394 slepc->GetEigenvector(i,temp);
395 x.Distribute(temp);
396 }
397
398 mode_sock << "parallel " << num_procs << " " << myid << "\n"
399 << "solution\n" << *pmesh << x << flush
400 << "window_title 'Eigenmode " << i+1 << '/' << nev
401 << ", Lambda = " << eigenvalues[i] << "'" << endl;
402
403 char c;
404 if (myid == 0)
405 {
406 cout << "press (q)uit or (c)ontinue --> " << flush;
407 cin >> c;
408 }
409 MPI_Bcast(&c, 1, MPI_CHAR, 0, MPI_COMM_WORLD);
410
411 if (c != 'c')
412 {
413 break;
414 }
415 }
416 mode_sock.close();
417 }
418
419 // 12. Free the used memory.
420 delete lobpcg;
421 delete slepc;
422 delete precond;
423 delete M;
424 delete A;
425 delete pA;
426 delete pM;
427#if defined(MFEM_USE_SUPERLU) || defined(MFEM_USE_STRUMPACK)
428 delete Arow;
429#endif
430 delete fespace;
431 if (order > 0)
432 {
433 delete fec;
434 }
435 delete pmesh;
436
437 // We finalize SLEPc
439
440 return 0;
441}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
int Size() const
Return the logical size of the array.
Definition array.hpp:144
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY....
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void EliminateEssentialBCDiag(const Array< int > &bdr_attr_is_ess, real_t value)
Perform elimination and set the diagonal entry to the given value.
A coefficient that is constant across space and time.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition device.hpp:123
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition device.cpp:286
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
The BoomerAMG solver in hypre.
Definition hypre.hpp:1691
void SetPrintLevel(int print_level)
Definition hypre.hpp:1771
void SetMassMatrix(Operator &M)
Definition hypre.cpp:6343
void SetPrintLevel(int logging)
Definition hypre.cpp:6272
void SetPreconditioner(Solver &precond)
Definition hypre.cpp:6287
void GetEigenvalues(Array< real_t > &eigenvalues) const
Collect the converged eigenvalues.
Definition hypre.cpp:6353
void SetTol(real_t tol)
Definition hypre.cpp:6250
void SetOperator(Operator &A)
Definition hypre.cpp:6296
void SetNumModes(int num_eigs)
Definition hypre.hpp:2103
void Solve()
Solve the eigenproblem.
Definition hypre.cpp:6410
void SetPrecondUsageMode(int pcg_mode)
Definition hypre.cpp:6281
void SetRandomSeed(int s)
Definition hypre.hpp:2105
void SetMaxIter(int max_iter)
Definition hypre.cpp:6266
const HypreParVector & GetEigenvector(unsigned int i) const
Extract a single eigenvector.
Definition hypre.cpp:6365
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.hpp:74
Mesh data type.
Definition mesh.hpp:56
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:282
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
void GetNodes(Vector &node_coord) const
Definition mesh.cpp:8973
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
Pointer to an Operator of a specified type.
Definition handle.hpp:34
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition handle.hpp:120
void Get(OpType *&A) const
Return the Operator pointer statically cast to a given OpType.
Definition handle.hpp:114
Abstract operator.
Definition operator.hpp:25
Type
Enumeration defining IDs for some classes derived from Operator.
Definition operator.hpp:284
@ Hypre_ParCSR
ID for class HypreParMatrix.
Definition operator.hpp:287
@ PETSC_MATAIJ
ID for class PetscParMatrix, MATAIJ format.
Definition operator.hpp:288
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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 'var' to receive the value. Enable/disable tags are used to set the bool...
Definition optparser.hpp:82
bool Good() const
Return true if the command line options were parsed successfully.
Class for parallel bilinear form.
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void Assemble(int skip_zeros=1)
Assemble the local matrix.
Abstract parallel finite element space.
Definition pfespace.hpp:29
HYPRE_BigInt GlobalTrueVSize() const
Definition pfespace.hpp:285
int GetTrueVSize() const override
Return the number of local vector true dofs.
Definition pfespace.hpp:289
Class for parallel grid function.
Definition pgridfunc.hpp:33
void Save(std::ostream &out) const override
void Distribute(const Vector *tv)
Class for parallel meshes.
Definition pmesh.hpp:34
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Definition pmesh.cpp:4801
Wrapper for PETSc's matrix class.
Definition petsc.hpp:319
void SetOperator(const Operator &op)
Set the operator/matrix.
void SetMatching(strumpack::MatchingJob job)
Configure static pivoting for stability.
void SetCompression(strumpack::CompressionType type)
Select compression for sparse data types.
void SetPrintFactorStatistics(bool print_stat)
Set up verbose printing during the factor step.
void SetKrylovSolver(strumpack::KrylovSolver method)
Set the Krylov solver method to use.
void SetFromCommandLine()
Set options that were captured from the command line.
void SetPrintSolveStatistics(bool print_stat)
Set up verbose printing during the solve step.
void SetReorderingStrategy(strumpack::ReorderingStrategy method)
Set matrix reordering strategy.
@ TARGET_REAL
The eigenvalues with the real component closest to the target value.
Definition slepc.hpp:146
void SetOperators(const PetscParMatrix &op, const PetscParMatrix &opB)
Set operators for generalized eigenvalue problem.
Definition slepc.cpp:93
void SetWhichEigenpairs(Which which)
Set the which eigenvalues the solver will target and the order they will be indexed in.
Definition slepc.cpp:195
void SetTarget(real_t target)
Set the target value for the eigenpairs you want when using SlepcEigenSolver::TARGET_MAGNITUDE or Sle...
Definition slepc.cpp:229
void GetEigenvalue(unsigned int i, real_t &lr) const
Get the ith eigenvalue after the system has been solved.
Definition slepc.cpp:147
void SetNumModes(int num_eigs)
Set the number of eigenmodes to compute.
Definition slepc.cpp:124
void SetSpectralTransformation(SpectralTransformation transformation)
Set the spectral transformation strategy for acceletating convergenvce. Both SlepcEigenSolver::SHIFT ...
Definition slepc.cpp:234
@ SHIFT_INVERT
Utilize the shift and invert strategy.
Definition slepc.hpp:157
void GetEigenvector(unsigned int i, Vector &vr) const
Get the ith eigenvector after the system has been solved.
Definition slepc.cpp:158
void Solve()
Solve the eigenvalue problem for the specified number of eigenvalues.
Definition slepc.cpp:130
Base class for solvers.
Definition operator.hpp:683
void SetColumnPermutation(superlu::ColPerm col_perm)
Specify how to permute the columns of the matrix.
Definition superlu.cpp:399
void SetSymmetricPattern(bool sym)
Specify whether the matrix has a symmetric pattern to avoid extra work (default false)
Definition superlu.cpp:454
void SetOperator(const Operator &op)
Set the operator/matrix.
Definition superlu.cpp:475
void SetPrintStatistics(bool print_stat)
Specify whether to print the solver statistics (default true)
Definition superlu.cpp:385
Vector data type.
Definition vector.hpp:80
int close()
Close the socketstream.
int dim
Definition ex24.cpp:53
int main()
HYPRE_Int HYPRE_BigInt
real_t a
Definition lissajous.cpp:41
@ PARMETIS
Sequential ordering on structure of using the PARMETIS package.
Definition superlu.hpp:71
const int visport
void MFEMFinalizeSlepc()
Definition slepc.cpp:53
void MFEMInitializeSlepc()
Definition slepc.cpp:29
const char vishost[]