MFEM v4.8.0
Finite element discretization library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex2p.cpp
Go to the documentation of this file.
1// MFEM Example 2 - Parallel Version
2// PETSc Modification
3//
4// Compile with: make ex2p
5//
6// Sample runs:
7// mpirun -np 4 ex2p -m ../../data/beam-quad.mesh --petscopts rc_ex2p
8//
9// Description: This example code solves a simple linear elasticity problem
10// describing a multi-material cantilever beam.
11//
12// Specifically, we approximate the weak form of -div(sigma(u))=0
13// where sigma(u)=lambda*div(u)*I+mu*(grad*u+u*grad) is the stress
14// tensor corresponding to displacement field u, and lambda and mu
15// are the material Lame constants. The boundary conditions are
16// u=0 on the fixed part of the boundary with attribute 1, and
17// sigma(u).n=f on the remainder with f being a constant pull down
18// vector on boundary elements with attribute 2, and zero
19// otherwise. The geometry of the domain is assumed to be as
20// follows:
21//
22// +----------+----------+
23// boundary --->| material | material |<--- boundary
24// attribute 1 | 1 | 2 | attribute 2
25// (fixed) +----------+----------+ (pull down)
26//
27// The example demonstrates the use of high-order and NURBS vector
28// finite element spaces with the linear elasticity bilinear form,
29// meshes with curved elements, and the definition of piece-wise
30// constant and vector coefficient objects. Static condensation is
31// also illustrated. The example also shows how to form a linear
32// system using a PETSc matrix and solve with a PETSc solver.
33//
34// The example also show how to use the non-overlapping feature of
35// the ParBilinearForm class to obtain the linear operator in
36// a format suitable for the BDDC preconditioner in PETSc.
37//
38// We recommend viewing Example 1 before viewing this example.
39
40#include "mfem.hpp"
41#include <fstream>
42#include <iostream>
43
44#ifndef MFEM_USE_PETSC
45#error This example requires that MFEM is built with MFEM_USE_PETSC=YES
46#endif
47
48using namespace std;
49using namespace mfem;
50
51int main(int argc, char *argv[])
52{
53 // 1. Initialize MPI and HYPRE.
54 Mpi::Init(argc, argv);
55 int num_procs = Mpi::WorldSize();
56 int myid = Mpi::WorldRank();
58
59 // 2. Parse command-line options.
60 const char *mesh_file = "../../data/beam-tri.mesh";
61 int ser_ref_levels = -1;
62 int par_ref_levels = 1;
63 int order = 1;
64 bool static_cond = false;
65 bool visualization = 1;
66 bool amg_elast = 0;
67 bool use_petsc = true;
68 const char *petscrc_file = "";
69 bool use_nonoverlapping = false;
70 const char *device_config = "cpu";
71
72 OptionsParser args(argc, argv);
73 args.AddOption(&mesh_file, "-m", "--mesh",
74 "Mesh file to use.");
75 args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
76 "Number of times to refine the mesh uniformly in serial.");
77 args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
78 "Number of times to refine the mesh uniformly in parallel.");
79 args.AddOption(&order, "-o", "--order",
80 "Finite element order (polynomial degree).");
81 args.AddOption(&amg_elast, "-elast", "--amg-for-elasticity", "-sys",
82 "--amg-for-systems",
83 "Use the special AMG elasticity solver (GM/LN approaches), "
84 "or standard AMG for systems (unknown approach).");
85 args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
86 "--no-static-condensation", "Enable static condensation.");
87 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
88 "--no-visualization",
89 "Enable or disable GLVis visualization.");
90 args.AddOption(&use_petsc, "-usepetsc", "--usepetsc", "-no-petsc",
91 "--no-petsc",
92 "Use or not PETSc to solve the linear system.");
93 args.AddOption(&petscrc_file, "-petscopts", "--petscopts",
94 "PetscOptions file to use.");
95 args.AddOption(&use_nonoverlapping, "-nonoverlapping", "--nonoverlapping",
96 "-no-nonoverlapping", "--no-nonoverlapping",
97 "Use or not the block diagonal PETSc's matrix format "
98 "for non-overlapping domain decomposition.");
99 args.AddOption(&device_config, "-d", "--device",
100 "Device configuration string, see Device::Configure().");
101 args.Parse();
102 if (!args.Good())
103 {
104 if (myid == 0)
105 {
106 args.PrintUsage(cout);
107 }
108 return 1;
109 }
110 if (myid == 0)
111 {
112 args.PrintOptions(cout);
113 }
114
115 // 2b. Enable hardware devices such as GPUs, and programming models such as
116 // CUDA, OCCA, RAJA and OpenMP based on command line options.
117 Device device(device_config);
118 if (myid == 0) { device.Print(); }
119
120 // 2c. We initialize PETSc
121 if (use_petsc) { MFEMInitializePetsc(NULL,NULL,petscrc_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 if (mesh->attributes.Max() < 2 || mesh->bdr_attributes.Max() < 2)
130 {
131 if (myid == 0)
132 cerr << "\nInput mesh should have at least two materials and "
133 << "two boundary attributes! (See schematic in ex2.cpp)\n"
134 << endl;
135 return 3;
136 }
137
138 // 4. Select the order of the finite element discretization space. For NURBS
139 // meshes, we increase the order by degree elevation.
140 if (mesh->NURBSext)
141 {
142 mesh->DegreeElevate(order, order);
143 }
144
145 // 5. Refine the serial mesh on all processors to increase the resolution. In
146 // this example we do 'ref_levels' of uniform refinement. We choose
147 // 'ref_levels' to be the largest number that gives a final mesh with no
148 // more than 1,000 elements.
149 {
150 int ref_levels = ser_ref_levels >= 0 ? ser_ref_levels :
151 (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
152 for (int l = 0; l < ref_levels; l++)
153 {
154 mesh->UniformRefinement();
155 }
156 }
157
158 // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
159 // this mesh further in parallel to increase the resolution. Once the
160 // parallel mesh is defined, the serial mesh can be deleted.
161 ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
162 delete mesh;
163 {
164 for (int l = 0; l < par_ref_levels; l++)
165 {
166 pmesh->UniformRefinement();
167 }
168 }
169
170 // 7. Define a parallel finite element space on the parallel mesh. Here we
171 // use vector finite elements, i.e. dim copies of a scalar finite element
172 // space. We use the ordering by vector dimension (the last argument of
173 // the FiniteElementSpace constructor) which is expected in the systems
174 // version of BoomerAMG preconditioner. For NURBS meshes, we use the
175 // (degree elevated) NURBS space associated with the mesh nodes.
177 ParFiniteElementSpace *fespace;
178 const bool use_nodal_fespace = pmesh->NURBSext && !amg_elast;
179 if (use_nodal_fespace)
180 {
181 fec = NULL;
182 fespace = (ParFiniteElementSpace *)pmesh->GetNodes()->FESpace();
183 }
184 else
185 {
186 fec = new H1_FECollection(order, dim);
187 fespace = new ParFiniteElementSpace(pmesh, fec, dim, Ordering::byVDIM);
188 }
189 HYPRE_BigInt size = fespace->GlobalTrueVSize();
190 if (myid == 0)
191 {
192 cout << "Number of finite element unknowns: " << size << endl
193 << "Assembling: " << flush;
194 }
195
196 // 8. Determine the list of true (i.e. parallel conforming) essential
197 // boundary dofs. In this example, the boundary conditions are defined by
198 // marking only boundary attribute 1 from the mesh as essential and
199 // converting it to a list of true dofs.
200 Array<int> ess_tdof_list, ess_bdr(pmesh->bdr_attributes.Max());
201 ess_bdr = 0;
202 ess_bdr[0] = 1;
203 fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
204
205 // 9. Set up the parallel linear form b(.) which corresponds to the
206 // right-hand side of the FEM linear system. In this case, b_i equals the
207 // boundary integral of f*phi_i where f represents a "pull down" force on
208 // the Neumann part of the boundary and phi_i are the basis functions in
209 // the finite element fespace. The force is defined by the object f, which
210 // is a vector of Coefficient objects. The fact that f is non-zero on
211 // boundary attribute 2 is indicated by the use of piece-wise constants
212 // coefficient for its last component.
214 for (int i = 0; i < dim-1; i++)
215 {
216 f.Set(i, new ConstantCoefficient(0.0));
217 }
218 {
219 Vector pull_force(pmesh->bdr_attributes.Max());
220 pull_force = 0.0;
221 pull_force(1) = -1.0e-2;
222 f.Set(dim-1, new PWConstCoefficient(pull_force));
223 }
224
225 ParLinearForm *b = new ParLinearForm(fespace);
226 b->AddBoundaryIntegrator(new VectorBoundaryLFIntegrator(f));
227 if (myid == 0)
228 {
229 cout << "r.h.s. ... " << flush;
230 }
231 b->Assemble();
232
233 // 10. Define the solution vector x as a parallel finite element grid
234 // function corresponding to fespace. Initialize x with initial guess of
235 // zero, which satisfies the boundary conditions.
236 ParGridFunction x(fespace);
237 x = 0.0;
238
239 // 11. Set up the parallel bilinear form a(.,.) on the finite element space
240 // corresponding to the linear elasticity integrator with piece-wise
241 // constants coefficient lambda and mu.
242 Vector lambda(pmesh->attributes.Max());
243 lambda = 1.0;
244 lambda(0) = lambda(1)*50;
245 PWConstCoefficient lambda_func(lambda);
246 Vector mu(pmesh->attributes.Max());
247 mu = 1.0;
248 mu(0) = mu(1)*50;
249 PWConstCoefficient mu_func(mu);
250
251 ParBilinearForm *a = new ParBilinearForm(fespace);
252 a->AddDomainIntegrator(new ElasticityIntegrator(lambda_func, mu_func));
253
254 // 12. Assemble the parallel bilinear form and the corresponding linear
255 // system, applying any necessary transformations such as: parallel
256 // assembly, eliminating boundary conditions, applying conforming
257 // constraints for non-conforming AMR, static condensation, etc.
258 if (myid == 0) { cout << "matrix ... " << flush; }
259 if (static_cond) { a->EnableStaticCondensation(); }
260 // Here we want to try out block-size aware AMG solver in PETSc.
261 // For that to work properly, we need a fully-compliant block-size
262 // structure and we do not skip zeros when assembling.
263 a->Assemble(use_petsc ? 0 : 1);
264
265 Vector B, X;
266 if (!use_petsc)
267 {
269 a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
270 if (myid == 0)
271 {
272 cout << "done." << endl;
273 cout << "Size of linear system: " << A.GetGlobalNumRows() << endl;
274 }
275
276 // 13. Define and apply a parallel PCG solver for A X = B with the BoomerAMG
277 // preconditioner from hypre.
278 HypreBoomerAMG *amg = new HypreBoomerAMG(A);
279 if (amg_elast && !a->StaticCondensationIsEnabled())
280 {
281 amg->SetElasticityOptions(fespace);
282 }
283 else
284 {
286 }
287 HyprePCG *pcg = new HyprePCG(A);
288 pcg->SetTol(1e-8);
289 pcg->SetMaxIter(500);
290 pcg->SetPrintLevel(2);
291 pcg->SetPreconditioner(*amg);
292 pcg->Mult(B, X);
293 delete pcg;
294 delete amg;
295 }
296 else
297 {
298 // 13b. Use PETSc to solve the linear system.
299 // Assemble a PETSc matrix, so that PETSc solvers can be used natively.
301 a->SetOperatorType(use_nonoverlapping ?
303 a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
304 if (myid == 0)
305 {
306 cout << "done." << endl;
307 cout << "Size of linear system: " << A.M() << endl;
308 }
309 // Tell PETSc the matrix has a block structure
310 A.SetBlockSize(dim);
311
312 // The preconditioner for the PCG solver can be specified in the
313 // PETSc config file
314 PetscPCGSolver *pcg = new PetscPCGSolver(A);
315 PetscPreconditioner *prec = NULL;
316 if (use_nonoverlapping) // Specialized BDDC construction
317 {
318 // Compute dofs belonging to the natural boundary
319 Array<int> nat_tdof_list, nat_bdr(pmesh->bdr_attributes.Max());
320 nat_bdr = 1;
321 nat_bdr[0] = 0;
322 fespace->GetEssentialTrueDofs(nat_bdr, nat_tdof_list);
323
324 // Auxiliary class for BDDC customization
326 // Inform the solver about the finite element space
327 opts.SetSpace(fespace);
328 // Inform the solver about essential dofs
329 opts.SetEssBdrDofs(&ess_tdof_list);
330 // Inform the solver about natural dofs
331 opts.SetNatBdrDofs(&nat_tdof_list);
332 // Create a BDDC solver with parameters
333 prec = new PetscBDDCSolver(A,opts);
334 pcg->SetPreconditioner(*prec);
335 }
336
337 pcg->SetMaxIter(500);
338 pcg->SetTol(1e-8);
339 pcg->SetPrintLevel(2);
340 pcg->Mult(B, X);
341 delete pcg;
342 delete prec;
343 }
344
345 // 14. Recover the parallel grid function corresponding to X. This is the
346 // local finite element solution on each processor.
347 a->RecoverFEMSolution(X, *b, x);
348
349 // 15. For non-NURBS meshes, make the mesh curved based on the finite element
350 // space. This means that we define the mesh elements through a fespace
351 // based transformation of the reference element. This allows us to save
352 // the displaced mesh as a curved mesh when using high-order finite
353 // element displacement field. We assume that the initial mesh (read from
354 // the file) is not higher order curved mesh compared to the chosen FE
355 // space.
356 if (!use_nodal_fespace)
357 {
358 pmesh->SetNodalFESpace(fespace);
359 }
360
361 // 16. Save in parallel the displaced mesh and the inverted solution (which
362 // gives the backward displacements to the original grid). This output
363 // can be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
364 {
365 GridFunction *nodes = pmesh->GetNodes();
366 *nodes += x;
367 x *= -1;
368
369 ostringstream mesh_name, sol_name;
370 mesh_name << "mesh." << setfill('0') << setw(6) << myid;
371 sol_name << "sol." << setfill('0') << setw(6) << myid;
372
373 ofstream mesh_ofs(mesh_name.str().c_str());
374 mesh_ofs.precision(8);
375 pmesh->Print(mesh_ofs);
376
377 ofstream sol_ofs(sol_name.str().c_str());
378 sol_ofs.precision(8);
379 x.Save(sol_ofs);
380 }
381
382 // 17. Send the above data by socket to a GLVis server. Use the "n" and "b"
383 // keys in GLVis to visualize the displacements.
384 if (visualization)
385 {
386 char vishost[] = "localhost";
387 int visport = 19916;
388 socketstream sol_sock(vishost, visport);
389 sol_sock << "parallel " << num_procs << " " << myid << "\n";
390 sol_sock.precision(8);
391 sol_sock << "solution\n" << *pmesh << x << flush;
392 }
393
394 // 18. Free the used memory.
395 delete a;
396 delete b;
397 if (fec)
398 {
399 delete fespace;
400 delete fec;
401 }
402 delete pmesh;
403
404 // We finalize PETSc
405 if (use_petsc) { MFEMFinalizePetsc(); }
406
407 return 0;
408}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
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:297
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:275
The BoomerAMG solver in hypre.
Definition hypre.hpp:1717
void SetSystemsOptions(int dim, bool order_bynodes=false)
Definition hypre.cpp:5169
void SetElasticityOptions(ParFiniteElementSpace *fespace, bool interp_refine=true)
Definition hypre.cpp:5296
PCG solver in hypre.
Definition hypre.hpp:1301
void Mult(const HypreParVector &b, HypreParVector &x) const override
Solve Ax=b with hypre's PCG.
Definition hypre.cpp:4242
void SetPrintLevel(int print_lvl)
Definition hypre.cpp:4214
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition hypre.cpp:4219
void SetMaxIter(int max_iter)
Definition hypre.cpp:4204
void SetTol(real_t tol)
Definition hypre.cpp:4194
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:408
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition hypre.hpp:699
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 GetNE() const
Returns number of elements.
Definition mesh.hpp:1282
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 DegreeElevate(int rel_degree, int degree=16)
Definition mesh.cpp:6052
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:11295
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition mesh.hpp:288
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).
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0)
Form a constrained linear system using a matrix-free approach.
Definition operator.cpp:114
@ PETSC_MATIS
ID for class PetscParMatrix, MATIS format.
Definition operator.hpp:289
@ PETSC_MATAIJ
ID for class PetscParMatrix, MATAIJ format.
Definition operator.hpp:288
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear sys...
Definition operator.cpp:148
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.
A piecewise constant coefficient with the constants keyed off the element attribute numbers.
Class for parallel bilinear form.
Abstract parallel finite element space.
Definition pfespace.hpp:29
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
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 linear form.
Class for parallel meshes.
Definition pmesh.hpp:34
void SetNodalFESpace(FiniteElementSpace *nfes) override
Definition pmesh.cpp:2027
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Definition pmesh.cpp:4800
Auxiliary class for BDDC customization.
Definition petsc.hpp:831
void SetEssBdrDofs(const Array< int > *essdofs, bool loc=false)
Specify dofs on the essential boundary.
Definition petsc.hpp:850
void SetSpace(ParFiniteElementSpace *fe)
Definition petsc.hpp:845
void SetNatBdrDofs(const Array< int > *natdofs, bool loc=false)
Specify dofs on the natural boundary.
Definition petsc.hpp:858
void SetPreconditioner(Solver &precond)
Definition petsc.cpp:3129
void Mult(const Vector &b, Vector &x) const override
Application of the solver.
Definition petsc.cpp:3209
Wrapper for PETSc's matrix class.
Definition petsc.hpp:319
PetscInt M() const
Returns the global number of rows.
Definition petsc.cpp:1010
void SetBlockSize(PetscInt rbs, PetscInt cbs=-1)
Set row and column block sizes of a matrix.
Definition petsc.cpp:1031
Abstract class for PETSc's preconditioners.
Definition petsc.hpp:805
void SetTol(real_t tol)
Definition petsc.cpp:2390
void SetPrintLevel(int plev)
Definition petsc.cpp:2472
void SetMaxIter(int max_iter)
Definition petsc.cpp:2445
Vector coefficient defined by an array of scalar coefficients. Coefficients that are not set will eva...
Vector data type.
Definition vector.hpp:82
Vector & Set(const real_t a, const Vector &x)
(*this) = a * x
Definition vector.cpp:337
int dim
Definition ex24.cpp:53
real_t mu
Definition ex25.cpp:140
int main()
HYPRE_Int HYPRE_BigInt
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
Definition petsc.cpp:206
void MFEMFinalizePetsc()
Definition petsc.cpp:246
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
const char vishost[]
std::array< int, NCMesh::MaxFaceNodes > nodes