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