MFEM v4.8.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//
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/star-mixed.mesh
8// mpirun -np 4 ex11p -m ../data/escher.mesh
9// mpirun -np 4 ex11p -m ../data/fichera.mesh
10// mpirun -np 4 ex11p -m ../data/fichera-mixed.mesh
11// mpirun -np 4 ex11p -m ../data/periodic-annulus-sector.msh
12// mpirun -np 4 ex11p -m ../data/periodic-torus-sector.msh -rs 1
13// mpirun -np 4 ex11p -m ../data/toroid-wedge.mesh -o 2
14// mpirun -np 4 ex11p -m ../data/square-disc-p2.vtk -o 2
15// mpirun -np 4 ex11p -m ../data/square-disc-p3.mesh -o 3
16// mpirun -np 4 ex11p -m ../data/square-disc-nurbs.mesh -o -1
17// mpirun -np 4 ex11p -m ../data/disc-nurbs.mesh -o -1 -n 20
18// mpirun -np 4 ex11p -m ../data/pipe-nurbs.mesh -o -1
19// mpirun -np 4 ex11p -m ../data/ball-nurbs.mesh -o 2
20// mpirun -np 4 ex11p -m ../data/star-surf.mesh
21// mpirun -np 4 ex11p -m ../data/square-disc-surf.mesh
22// mpirun -np 4 ex11p -m ../data/inline-segment.mesh
23// mpirun -np 4 ex11p -m ../data/inline-quad.mesh
24// mpirun -np 4 ex11p -m ../data/inline-tri.mesh
25// mpirun -np 4 ex11p -m ../data/inline-hex.mesh
26// mpirun -np 4 ex11p -m ../data/inline-tet.mesh
27// mpirun -np 4 ex11p -m ../data/inline-wedge.mesh -s 83
28// mpirun -np 4 ex11p -m ../data/amr-quad.mesh
29// mpirun -np 4 ex11p -m ../data/amr-hex.mesh
30// mpirun -np 4 ex11p -m ../data/mobius-strip.mesh -n 8
31// mpirun -np 4 ex11p -m ../data/klein-bottle.mesh -n 10
32//
33// Description: This example code demonstrates the use of MFEM to solve the
34// eigenvalue problem -Delta u = lambda u with homogeneous
35// Dirichlet boundary conditions.
36//
37// We compute a number of the lowest eigenmodes by discretizing
38// the Laplacian and Mass operators using a FE space of the
39// specified order, or an isoparametric/isogeometric space if
40// order < 1 (quadratic for quadratic curvilinear mesh, NURBS for
41// NURBS mesh, etc.)
42//
43// The example highlights the use of the LOBPCG eigenvalue solver
44// together with the BoomerAMG preconditioner in HYPRE, as well as
45// optionally the SuperLU or STRUMPACK parallel direct solvers.
46// Reusing a single GLVis visualization window for multiple
47// eigenfunctions is also illustrated.
48//
49// We recommend viewing Example 1 before viewing this example.
50
51#include "mfem.hpp"
52#include <fstream>
53#include <iostream>
54
55using namespace std;
56using namespace mfem;
57
58int main(int argc, char *argv[])
59{
60 // 1. Initialize MPI and HYPRE.
61 Mpi::Init(argc, argv);
62 int num_procs = Mpi::WorldSize();
63 int myid = Mpi::WorldRank();
65
66 // 2. Parse command-line options.
67 const char *mesh_file = "../data/star.mesh";
68 int ser_ref_levels = 2;
69 int par_ref_levels = 1;
70 int order = 1;
71 int nev = 5;
72 int seed = 75;
73 bool slu_solver = false;
74 bool sp_solver = false;
75 bool cpardiso_solver = false;
76 bool visualization = 1;
77
78 OptionsParser args(argc, argv);
79 args.AddOption(&mesh_file, "-m", "--mesh",
80 "Mesh file to use.");
81 args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
82 "Number of times to refine the mesh uniformly in serial.");
83 args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
84 "Number of times to refine the mesh uniformly in parallel.");
85 args.AddOption(&order, "-o", "--order",
86 "Finite element order (polynomial degree) or -1 for"
87 " isoparametric space.");
88 args.AddOption(&nev, "-n", "--num-eigs",
89 "Number of desired eigenmodes.");
90 args.AddOption(&seed, "-s", "--seed",
91 "Random seed used to initialize LOBPCG.");
92#ifdef MFEM_USE_SUPERLU
93 args.AddOption(&slu_solver, "-slu", "--superlu", "-no-slu",
94 "--no-superlu", "Use the SuperLU Solver.");
95#endif
96#ifdef MFEM_USE_STRUMPACK
97 args.AddOption(&sp_solver, "-sp", "--strumpack", "-no-sp",
98 "--no-strumpack", "Use the STRUMPACK Solver.");
99#endif
100#ifdef MFEM_USE_MKL_CPARDISO
101 args.AddOption(&cpardiso_solver, "-cpardiso", "--cpardiso", "-no-cpardiso",
102 "--no-cpardiso", "Use the MKL CPardiso Solver.");
103#endif
104 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
105 "--no-visualization",
106 "Enable or disable GLVis visualization.");
107 args.Parse();
108 if (slu_solver && sp_solver)
109 {
110 if (myid == 0)
111 cout << "WARNING: Both SuperLU and STRUMPACK have been selected,"
112 << " please choose either one." << endl
113 << " Defaulting to SuperLU." << endl;
114 sp_solver = false;
115 }
116 // The command line options are also passed to the STRUMPACK
117 // solver. So do not exit if some options are not recognized.
118 if (!sp_solver)
119 {
120 if (!args.Good())
121 {
122 if (myid == 0)
123 {
124 args.PrintUsage(cout);
125 }
126 return 1;
127 }
128 }
129 if (myid == 0)
130 {
131 args.PrintOptions(cout);
132 }
133
134 // 3. Read the (serial) mesh from the given mesh file on all processors. We
135 // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
136 // and volume meshes with the same code.
137 Mesh *mesh = new Mesh(mesh_file, 1, 1);
138 int dim = mesh->Dimension();
139
140 // 4. Refine the serial mesh on all processors to increase the resolution. In
141 // this example we do 'ref_levels' of uniform refinement (2 by default, or
142 // specified on the command line with -rs).
143 for (int lev = 0; lev < ser_ref_levels; lev++)
144 {
145 mesh->UniformRefinement();
146 }
147
148 // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
149 // this mesh further in parallel to increase the resolution (1 time by
150 // default, or specified on the command line with -rp). Once the parallel
151 // mesh is defined, the serial mesh can be deleted.
152 ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
153 delete mesh;
154 for (int lev = 0; lev < par_ref_levels; lev++)
155 {
156 pmesh->UniformRefinement();
157 }
158
159 // 6. Define a parallel finite element space on the parallel mesh. Here we
160 // use continuous Lagrange finite elements of the specified order. If
161 // order < 1, we instead use an isoparametric/isogeometric space.
163 if (order > 0)
164 {
165 fec = new H1_FECollection(order, dim);
166 }
167 else if (pmesh->GetNodes())
168 {
169 fec = pmesh->GetNodes()->OwnFEC();
170 }
171 else
172 {
173 fec = new H1_FECollection(order = 1, dim);
174 }
175 ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
176 HYPRE_BigInt size = fespace->GlobalTrueVSize();
177 if (myid == 0)
178 {
179 cout << "Number of unknowns: " << size << endl;
180 }
181
182 // 7. Set up the parallel bilinear forms a(.,.) and m(.,.) on the finite
183 // element space. The first corresponds to the Laplacian operator -Delta,
184 // while the second is a simple mass matrix needed on the right hand side
185 // of the generalized eigenvalue problem below. The boundary conditions
186 // are implemented by elimination with special values on the diagonal to
187 // shift the Dirichlet eigenvalues out of the computational range. After
188 // serial and parallel assembly we extract the corresponding parallel
189 // matrices A and M.
190 ConstantCoefficient one(1.0);
191 Array<int> ess_bdr;
192 if (pmesh->bdr_attributes.Size())
193 {
194 ess_bdr.SetSize(pmesh->bdr_attributes.Max());
195 ess_bdr = 0;
196 // Apply boundary conditions on all external boundaries:
197 pmesh->MarkExternalBoundaries(ess_bdr);
198 // Boundary conditions can also be applied based on named attributes:
199 // pmesh->MarkNamedBoundaries(set_name, ess_bdr)
200 }
201
202 ParBilinearForm *a = new ParBilinearForm(fespace);
203 a->AddDomainIntegrator(new DiffusionIntegrator(one));
204 if (pmesh->bdr_attributes.Size() == 0)
205 {
206 // Add a mass term if the mesh has no boundary, e.g. periodic mesh or
207 // closed surface.
208 a->AddDomainIntegrator(new MassIntegrator(one));
209 }
210 a->Assemble();
211 a->EliminateEssentialBCDiag(ess_bdr, 1.0);
212 a->Finalize();
213
214 ParBilinearForm *m = new ParBilinearForm(fespace);
216 m->Assemble();
217 // shift the eigenvalue corresponding to eliminated dofs to a large value
218 m->EliminateEssentialBCDiag(ess_bdr, numeric_limits<real_t>::min());
219 m->Finalize();
220
221 HypreParMatrix *A = a->ParallelAssemble();
223
224#if defined(MFEM_USE_SUPERLU) || defined(MFEM_USE_STRUMPACK)
225 Operator * Arow = NULL;
226#ifdef MFEM_USE_SUPERLU
227 if (slu_solver)
228 {
229 Arow = new SuperLURowLocMatrix(*A);
230 }
231#endif
232#ifdef MFEM_USE_STRUMPACK
233 if (sp_solver)
234 {
235 Arow = new STRUMPACKRowLocMatrix(*A);
236 }
237#endif
238#endif
239
240 delete a;
241 delete m;
242
243 // 8. Define and configure the LOBPCG eigensolver and the BoomerAMG
244 // preconditioner for A to be used within the solver. Set the matrices
245 // which define the generalized eigenproblem A x = lambda M x.
246 Solver * precond = NULL;
247 if (!slu_solver && !sp_solver && !cpardiso_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(MPI_COMM_WORLD, argc, argv);
270 strumpack->SetPrintFactorStatistics(true);
271 strumpack->SetPrintSolveStatistics(false);
272 strumpack->SetKrylovSolver(strumpack::KrylovSolver::DIRECT);
273 strumpack->SetReorderingStrategy(strumpack::ReorderingStrategy::METIS);
274 strumpack->SetMatching(strumpack::MatchingJob::NONE);
275 strumpack->SetCompression(strumpack::CompressionType::NONE);
276 strumpack->SetOperator(*Arow);
277 strumpack->SetFromCommandLine();
278 precond = strumpack;
279 }
280#endif
281#ifdef MFEM_USE_MKL_CPARDISO
282 if (cpardiso_solver)
283 {
284 auto cpardiso = new CPardisoSolver(A->GetComm());
285 cpardiso->SetMatrixType(CPardisoSolver::MatType::REAL_STRUCTURE_SYMMETRIC);
286 cpardiso->SetPrintLevel(1);
287 cpardiso->SetOperator(*A);
288 precond = cpardiso;
289 }
290#endif
291 }
292
293 HypreLOBPCG * lobpcg = new HypreLOBPCG(MPI_COMM_WORLD);
294 lobpcg->SetNumModes(nev);
295 lobpcg->SetRandomSeed(seed);
296 lobpcg->SetPreconditioner(*precond);
297 lobpcg->SetMaxIter(200);
298 lobpcg->SetTol(1e-8);
299 lobpcg->SetPrecondUsageMode(1);
300 lobpcg->SetPrintLevel(1);
301 lobpcg->SetMassMatrix(*M);
302 lobpcg->SetOperator(*A);
303
304 // 9. Compute the eigenmodes and extract the array of eigenvalues. Define a
305 // parallel grid function to represent each of the eigenmodes returned by
306 // the solver.
307 Array<real_t> eigenvalues;
308 lobpcg->Solve();
309 lobpcg->GetEigenvalues(eigenvalues);
310 ParGridFunction x(fespace);
311
312 // 10. Save the refined mesh and the modes in parallel. This output can be
313 // viewed later using GLVis: "glvis -np <np> -m mesh -g mode".
314 {
315 ostringstream mesh_name, mode_name;
316 mesh_name << "mesh." << setfill('0') << setw(6) << myid;
317
318 ofstream mesh_ofs(mesh_name.str().c_str());
319 mesh_ofs.precision(8);
320 pmesh->Print(mesh_ofs);
321
322 for (int i=0; i<nev; i++)
323 {
324 // convert eigenvector from HypreParVector to ParGridFunction
325 x = lobpcg->GetEigenvector(i);
326
327 mode_name << "mode_" << setfill('0') << setw(2) << i << "."
328 << setfill('0') << setw(6) << myid;
329
330 ofstream mode_ofs(mode_name.str().c_str());
331 mode_ofs.precision(8);
332 x.Save(mode_ofs);
333 mode_name.str("");
334 }
335 }
336
337 // 11. Send the solution by socket to a GLVis server.
338 if (visualization)
339 {
340 char vishost[] = "localhost";
341 int visport = 19916;
342 socketstream mode_sock(vishost, visport);
343 mode_sock.precision(8);
344
345 for (int i=0; i<nev; i++)
346 {
347 if ( myid == 0 )
348 {
349 cout << "Eigenmode " << i+1 << '/' << nev
350 << ", Lambda = " << eigenvalues[i] << endl;
351 }
352
353 // convert eigenvector from HypreParVector to ParGridFunction
354 x = lobpcg->GetEigenvector(i);
355
356 mode_sock << "parallel " << num_procs << " " << myid << "\n"
357 << "solution\n" << *pmesh << x << flush
358 << "window_title 'Eigenmode " << i+1 << '/' << nev
359 << ", Lambda = " << eigenvalues[i] << "'" << endl;
360
361 char c;
362 if (myid == 0)
363 {
364 cout << "press (q)uit or (c)ontinue --> " << flush;
365 cin >> c;
366 }
367 MPI_Bcast(&c, 1, MPI_CHAR, 0, MPI_COMM_WORLD);
368
369 if (c != 'c')
370 {
371 break;
372 }
373 }
374 mode_sock.close();
375 }
376
377 // 12. Free the used memory.
378 delete lobpcg;
379 delete precond;
380 delete M;
381 delete A;
382#if defined(MFEM_USE_SUPERLU) || defined(MFEM_USE_STRUMPACK)
383 delete Arow;
384#endif
385
386 delete fespace;
387 if (order > 0)
388 {
389 delete fec;
390 }
391 delete pmesh;
392
393 return 0;
394}
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:758
int Size() const
Return the logical size of the array.
Definition array.hpp:147
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void Finalize(int skip_zeros=1) override
Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY....
void EliminateEssentialBCDiag(const Array< int > &bdr_attr_is_ess, real_t value)
Perform elimination and set the diagonal entry to the given value.
MKL Parallel Direct Sparse Solver for Clusters.
Definition cpardiso.hpp:31
A coefficient that is constant across space and time.
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:275
The BoomerAMG solver in hypre.
Definition hypre.hpp:1717
void SetPrintLevel(int print_level)
Definition hypre.hpp:1797
void SetMassMatrix(Operator &M)
Definition hypre.cpp:6401
void SetPrintLevel(int logging)
Definition hypre.cpp:6330
void SetPreconditioner(Solver &precond)
Definition hypre.cpp:6345
void GetEigenvalues(Array< real_t > &eigenvalues) const
Collect the converged eigenvalues.
Definition hypre.cpp:6411
void SetTol(real_t tol)
Definition hypre.cpp:6308
void SetOperator(Operator &A)
Definition hypre.cpp:6354
void SetNumModes(int num_eigs)
Definition hypre.hpp:2129
void Solve()
Solve the eigenproblem.
Definition hypre.cpp:6468
void SetPrecondUsageMode(int pcg_mode)
Definition hypre.cpp:6339
void SetRandomSeed(int s)
Definition hypre.hpp:2131
void SetMaxIter(int max_iter)
Definition hypre.cpp:6324
const HypreParVector & GetEigenvector(unsigned int i) const
Extract a single eigenvector.
Definition hypre.cpp:6423
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:408
MPI_Comm GetComm() const
MPI communicator.
Definition hypre.hpp:598
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.cpp:33
Mesh data type.
Definition mesh.hpp:64
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:290
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1216
void GetNodes(Vector &node_coord) const
Definition mesh.cpp:9294
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:11295
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).
Abstract operator.
Definition operator.hpp:25
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:346
Class for parallel grid function.
Definition pgridfunc.hpp:50
void Save(std::ostream &out) const override
Class for parallel meshes.
Definition pmesh.hpp:34
void MarkExternalBoundaries(Array< int > &bdr_marker, bool excl=true) const override
Mark boundary attributes of external boundaries.
Definition pmesh.cpp:6744
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Definition pmesh.cpp:4800
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 SetOperator(const Operator &op) override
Set the operator/matrix.
void SetReorderingStrategy(strumpack::ReorderingStrategy method)
Set matrix reordering strategy.
Base class for solvers.
Definition operator.hpp:780
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
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 char vishost[]