MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
nurbs_ex11p.cpp
Go to the documentation of this file.
1// MFEM Example 11 - Parallel NURBS Version
2//
3// Compile with: make nurbs_ex11p
4//
5// Sample runs: mpirun -np 4 nurbs_ex11p -m ../../data/square-disc.mesh
6// mpirun -np 4 nurbs_ex11p -m ../../data/star.mesh
7// mpirun -np 4 nurbs_ex11p -m ../../data/escher.mesh
8// mpirun -np 4 nurbs_ex11p -m ../../data/fichera.mesh
9// mpirun -np 4 nurbs_ex11p -m ../../data/square-disc-p2.vtk -o 2
10// mpirun -np 4 nurbs_ex11p -m ../../data/square-disc-p3.mesh -o 3
11// mpirun -np 4 nurbs_ex11p -m ../../data/square-disc-nurbs.mesh -o -1
12// mpirun -np 4 nurbs_ex11p -m ../../data/disc-nurbs.mesh -o -1 -n 20
13// mpirun -np 4 nurbs_ex11p -m ../../data/pipe-nurbs.mesh -o -1
14// mpirun -np 4 nurbs_ex11p -m ../../data/ball-nurbs.mesh -o 2
15// mpirun -np 4 nurbs_ex11p -m ../../data/star-surf.mesh
16// mpirun -np 4 nurbs_ex11p -m ../../data/square-disc-surf.mesh
17// mpirun -np 4 nurbs_ex11p -m ../../data/inline-segment.mesh
18// mpirun -np 4 nurbs_ex11p -m ../../data/amr-quad.mesh
19// mpirun -np 4 nurbs_ex11p -m ../../data/amr-hex.mesh
20// mpirun -np 4 nurbs_ex11p -m ../../data/mobius-strip.mesh -n 8
21// mpirun -np 4 nurbs_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
45using namespace std;
46using namespace mfem;
47
48int main(int argc, char *argv[])
49{
50 // 1. Initialize MPI and HYPRE.
51 Mpi::Init(argc, argv);
52 int num_procs = Mpi::WorldSize();
53 int myid = Mpi::WorldRank();
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 Array<int> order(1);
61 order[0] = 0;
62 int nev = 5;
63 int seed = 75;
64 bool slu_solver = false;
65 bool sp_solver = false;
66 int visport = 19916;
67 bool visualization = 1;
68
69 OptionsParser args(argc, argv);
70 args.AddOption(&mesh_file, "-m", "--mesh",
71 "Mesh file to use.");
72 args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
73 "Number of times to refine the mesh uniformly in serial.");
74 args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
75 "Number of times to refine the mesh uniformly in parallel.");
76 args.AddOption(&order, "-o", "--order",
77 "Finite element order (polynomial degree) or -1 for"
78 " isoparametric space.");
79 args.AddOption(&nev, "-n", "--num-eigs",
80 "Number of desired eigenmodes.");
81 args.AddOption(&seed, "-s", "--seed",
82 "Random seed used to initialize LOBPCG.");
83#ifdef MFEM_USE_SUPERLU
84 args.AddOption(&slu_solver, "-slu", "--superlu", "-no-slu",
85 "--no-superlu", "Use the SuperLU Solver.");
86#endif
87#ifdef MFEM_USE_STRUMPACK
88 args.AddOption(&sp_solver, "-sp", "--strumpack", "-no-sp",
89 "--no-strumpack", "Use the STRUMPACK Solver.");
90#endif
91 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
92 "--no-visualization",
93 "Enable or disable GLVis visualization.");
94 args.AddOption(&visport, "-p", "--send-port", "Socket for GLVis.");
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 // 3. Read the (serial) mesh from the given mesh file on all processors. We
123 // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
124 // and volume meshes with the same code.
125 Mesh *mesh = new Mesh(mesh_file, 1, 1);
126 int dim = mesh->Dimension();
127
128 // 4. Refine the serial mesh on all processors to increase the resolution. In
129 // this example we do 'ref_levels' of uniform refinement (2 by default, or
130 // specified on the command line with -rs).
131 for (int lev = 0; lev < ser_ref_levels; lev++)
132 {
133 mesh->UniformRefinement();
134 }
135
136 // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
137 // this mesh further in parallel to increase the resolution (1 time by
138 // default, or specified on the command line with -rp). Once the parallel
139 // mesh is defined, the serial mesh can be deleted.
140 ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
141 delete mesh;
142 for (int lev = 0; lev < par_ref_levels; lev++)
143 {
144 pmesh->UniformRefinement();
145 }
146
147 // 6. Define a parallel finite element space on the parallel mesh. Here we
148 // use continuous Lagrange finite elements of the specified order. If
149 // order < 1, we instead use an isoparametric/isogeometric space.
151 NURBSExtension *NURBSext = NULL;
152 int own_fec = 0;
153
154 if (order[0] == 0) // Isoparametric
155 {
156 if (pmesh->GetNodes())
157 {
158 fec = pmesh->GetNodes()->OwnFEC();
159 own_fec = 0;
160 if (myid == 0)
161 {
162 cout << "Using isoparametric FEs: " << fec->Name() << endl;
163 }
164 }
165 else
166 {
167 if (myid == 0)
168 {
169 cout <<"Mesh does not have FEs --> Assume order 1.\n";
170 }
171 fec = new H1_FECollection(1, dim);
172 own_fec = 1;
173 }
174 }
175 else if (pmesh->NURBSext && (order[0] > 0) ) // Subparametric NURBS
176 {
177 fec = new NURBSFECollection(order[0]);
178 own_fec = 1;
179 int nkv = pmesh->NURBSext->GetNKV();
180
181 if (order.Size() == 1)
182 {
183 int tmp = order[0];
184 order.SetSize(nkv);
185 order = tmp;
186 }
187 if (order.Size() != nkv ) { mfem_error("Wrong number of orders set."); }
188 NURBSext = new NURBSExtension(pmesh->NURBSext, order);
189 }
190 else
191 {
192 if (order.Size() > 1) { cout <<"Wrong number of orders set, needs one.\n"; }
193 fec = new H1_FECollection(abs(order[0]), dim);
194 own_fec = 1;
195 }
196 ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh,NURBSext,fec);
197 HYPRE_BigInt size = fespace->GlobalTrueVSize();
198 if (myid == 0)
199 {
200 cout << "Number of unknowns: " << size << endl;
201 }
202
203 // 7. Set up the parallel bilinear forms a(.,.) and m(.,.) on the finite
204 // element space. The first corresponds to the Laplacian operator -Delta,
205 // while the second is a simple mass matrix needed on the right hand side
206 // of the generalized eigenvalue problem below. The boundary conditions
207 // are implemented by elimination with special values on the diagonal to
208 // shift the Dirichlet eigenvalues out of the computational range. After
209 // serial and parallel assembly we extract the corresponding parallel
210 // matrices A and M.
211 ConstantCoefficient one(1.0);
212 Array<int> ess_bdr;
213 if (pmesh->bdr_attributes.Size())
214 {
215 ess_bdr.SetSize(pmesh->bdr_attributes.Max());
216 ess_bdr = 1;
217 }
218
219 ParBilinearForm *a = new ParBilinearForm(fespace);
220 a->AddDomainIntegrator(new DiffusionIntegrator(one));
221 if (pmesh->bdr_attributes.Size() == 0)
222 {
223 // Add a mass term if the mesh has no boundary, e.g. periodic mesh or
224 // closed surface.
225 a->AddDomainIntegrator(new MassIntegrator(one));
226 }
227 a->Assemble();
228 a->EliminateEssentialBCDiag(ess_bdr, 1.0);
229 a->Finalize();
230
231 ParBilinearForm *m = new ParBilinearForm(fespace);
233 m->Assemble();
234 // shift the eigenvalue corresponding to eliminated dofs to a large value
235 m->EliminateEssentialBCDiag(ess_bdr, numeric_limits<real_t>::min());
236 m->Finalize();
237
238 HypreParMatrix *A = a->ParallelAssemble();
240
241#if defined(MFEM_USE_SUPERLU) || defined(MFEM_USE_STRUMPACK)
242 Operator * Arow = NULL;
243#ifdef MFEM_USE_SUPERLU
244 if (slu_solver)
245 {
246 Arow = new SuperLURowLocMatrix(*A);
247 }
248#endif
249#ifdef MFEM_USE_STRUMPACK
250 if (sp_solver)
251 {
252 Arow = new STRUMPACKRowLocMatrix(*A);
253 }
254#endif
255#endif
256
257 delete a;
258 delete m;
259
260 // 8. Define and configure the LOBPCG eigensolver and the BoomerAMG
261 // preconditioner for A to be used within the solver. Set the matrices
262 // which define the generalized eigenproblem A x = lambda M x.
263 Solver * precond = NULL;
264 if (!slu_solver && !sp_solver)
265 {
266 HypreBoomerAMG * amg = new HypreBoomerAMG(*A);
267 amg->SetPrintLevel(0);
268 precond = amg;
269 }
270 else
271 {
272#ifdef MFEM_USE_SUPERLU
273 if (slu_solver)
274 {
275 SuperLUSolver * superlu = new SuperLUSolver(MPI_COMM_WORLD);
276 superlu->SetPrintStatistics(false);
277 superlu->SetSymmetricPattern(true);
279 superlu->SetOperator(*Arow);
280 precond = superlu;
281 }
282#endif
283#ifdef MFEM_USE_STRUMPACK
284 if (sp_solver)
285 {
286 STRUMPACKSolver * strumpack = new STRUMPACKSolver(MPI_COMM_WORLD, argc, argv);
287 strumpack->SetPrintFactorStatistics(true);
288 strumpack->SetPrintSolveStatistics(false);
289 strumpack->SetKrylovSolver(strumpack::KrylovSolver::DIRECT);
290 strumpack->SetReorderingStrategy(strumpack::ReorderingStrategy::METIS);
291 strumpack->SetMatching(strumpack::MatchingJob::NONE);
292 strumpack->SetCompression(strumpack::CompressionType::NONE);
293 strumpack->SetOperator(*Arow);
294 strumpack->SetFromCommandLine();
295 precond = strumpack;
296 }
297#endif
298 }
299
300
301 HypreLOBPCG * lobpcg = new HypreLOBPCG(MPI_COMM_WORLD);
302 lobpcg->SetNumModes(nev);
303 lobpcg->SetRandomSeed(seed);
304 lobpcg->SetPreconditioner(*precond);
305 lobpcg->SetMaxIter(200);
306 lobpcg->SetTol(1e-8);
307 lobpcg->SetPrecondUsageMode(1);
308 lobpcg->SetPrintLevel(1);
309 lobpcg->SetMassMatrix(*M);
310 lobpcg->SetOperator(*A);
311
312 // 9. Compute the eigenmodes and extract the array of eigenvalues. Define a
313 // parallel grid function to represent each of the eigenmodes returned by
314 // the solver.
315 Array<real_t> eigenvalues;
316 lobpcg->Solve();
317 lobpcg->GetEigenvalues(eigenvalues);
318 ParGridFunction x(fespace);
319
320 // 10. Save the refined mesh and the modes in parallel. This output can be
321 // viewed later using GLVis: "glvis -np <np> -m mesh -g mode".
322 {
323 ostringstream mesh_name, mode_name;
324 mesh_name << "mesh." << setfill('0') << setw(6) << myid;
325
326 ofstream mesh_ofs(mesh_name.str().c_str());
327 mesh_ofs.precision(8);
328 pmesh->Print(mesh_ofs);
329
330 for (int i=0; i<nev; i++)
331 {
332 // convert eigenvector from HypreParVector to ParGridFunction
333 x = lobpcg->GetEigenvector(i);
334
335 mode_name << "mode_" << setfill('0') << setw(2) << i << "."
336 << setfill('0') << setw(6) << myid;
337
338 ofstream mode_ofs(mode_name.str().c_str());
339 mode_ofs.precision(8);
340 x.Save(mode_ofs);
341 mode_name.str("");
342 }
343 }
344
345 // 11. Send the solution by socket to a GLVis server.
346 if (visualization)
347 {
348 char vishost[] = "localhost";
349 socketstream mode_sock(vishost, visport);
350 mode_sock.precision(8);
351
352 for (int i=0; i<nev; i++)
353 {
354 if ( myid == 0 )
355 {
356 cout << "Eigenmode " << i+1 << '/' << nev
357 << ", Lambda = " << eigenvalues[i] << endl;
358 }
359
360 // convert eigenvector from HypreParVector to ParGridFunction
361 x = lobpcg->GetEigenvector(i);
362
363 mode_sock << "parallel " << num_procs << " " << myid << "\n"
364 << "solution\n" << *pmesh << x << flush
365 << "window_title 'Eigenmode " << i+1 << '/' << nev
366 << ", Lambda = " << eigenvalues[i] << "'" << endl;
367
368 char c;
369 if (myid == 0)
370 {
371 cout << "press (q)uit or (c)ontinue --> " << flush;
372 cin >> c;
373 }
374 MPI_Bcast(&c, 1, MPI_CHAR, 0, MPI_COMM_WORLD);
375
376 if (c != 'c')
377 {
378 break;
379 }
380 }
381 mode_sock.close();
382 }
383
384 // 12. Free the used memory.
385 delete lobpcg;
386 delete precond;
387 delete M;
388 delete A;
389#if defined(MFEM_USE_SUPERLU) || defined(MFEM_USE_STRUMPACK)
390 delete Arow;
391#endif
392
393 delete fespace;
394 if (own_fec)
395 {
396 delete fec;
397 }
398 delete pmesh;
399
400 return 0;
401}
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.
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
virtual const char * Name() const
Definition fe_coll.hpp:79
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
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
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition mesh.hpp:298
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).
NURBSExtension generally contains multiple NURBSPatch objects spanning an entire Mesh....
Definition nurbs.hpp:449
int GetNKV() const
Return the number of KnotVectors.
Definition nurbs.hpp:758
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition fe_coll.hpp:699
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 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
void mfem_error(const char *msg)
Definition error.cpp:154
const char vishost[]