MFEM  v4.4.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 
48 using namespace std;
49 using namespace mfem;
50 
51 int 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();
57  Hypre::Init();
58 
59  // 2. Parse command-line options.
60  const char *mesh_file = "../../data/beam-tri.mesh";
61  int order = 1;
62  bool static_cond = false;
63  bool visualization = 1;
64  bool amg_elast = 0;
65  bool use_petsc = true;
66  const char *petscrc_file = "";
67  bool use_nonoverlapping = false;
68  int ser_ref_levels = -1, par_ref_levels = 1;
69 
70  OptionsParser args(argc, argv);
71  args.AddOption(&mesh_file, "-m", "--mesh",
72  "Mesh file to use.");
73  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
74  "Number of times to refine the mesh uniformly in serial.");
75  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
76  "Number of times to refine the mesh uniformly in parallel.");
77  args.AddOption(&order, "-o", "--order",
78  "Finite element order (polynomial degree).");
79  args.AddOption(&amg_elast, "-elast", "--amg-for-elasticity", "-sys",
80  "--amg-for-systems",
81  "Use the special AMG elasticity solver (GM/LN approaches), "
82  "or standard AMG for systems (unknown approach).");
83  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
84  "--no-static-condensation", "Enable static condensation.");
85  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
86  "--no-visualization",
87  "Enable or disable GLVis visualization.");
88  args.AddOption(&use_petsc, "-usepetsc", "--usepetsc", "-no-petsc",
89  "--no-petsc",
90  "Use or not PETSc to solve the linear system.");
91  args.AddOption(&petscrc_file, "-petscopts", "--petscopts",
92  "PetscOptions file to use.");
93  args.AddOption(&use_nonoverlapping, "-nonoverlapping", "--nonoverlapping",
94  "-no-nonoverlapping", "--no-nonoverlapping",
95  "Use or not the block diagonal PETSc's matrix format "
96  "for non-overlapping domain decomposition.");
97  args.Parse();
98  if (!args.Good())
99  {
100  if (myid == 0)
101  {
102  args.PrintUsage(cout);
103  }
104  return 1;
105  }
106  if (myid == 0)
107  {
108  args.PrintOptions(cout);
109  }
110 
111  // 2b. We initialize PETSc
112  if (use_petsc) { MFEMInitializePetsc(NULL,NULL,petscrc_file,NULL); }
113 
114  // 3. Read the (serial) mesh from the given mesh file on all processors. We
115  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
116  // and volume meshes with the same code.
117  Mesh *mesh = new Mesh(mesh_file, 1, 1);
118  int dim = mesh->Dimension();
119 
120  if (mesh->attributes.Max() < 2 || mesh->bdr_attributes.Max() < 2)
121  {
122  if (myid == 0)
123  cerr << "\nInput mesh should have at least two materials and "
124  << "two boundary attributes! (See schematic in ex2.cpp)\n"
125  << endl;
126  return 3;
127  }
128 
129  // 4. Select the order of the finite element discretization space. For NURBS
130  // meshes, we increase the order by degree elevation.
131  if (mesh->NURBSext)
132  {
133  mesh->DegreeElevate(order, order);
134  }
135 
136  // 5. Refine the serial mesh on all processors to increase the resolution. In
137  // this example we do 'ref_levels' of uniform refinement. We choose
138  // 'ref_levels' to be the largest number that gives a final mesh with no
139  // more than 1,000 elements.
140  {
141  int ref_levels = ser_ref_levels >= 0 ? ser_ref_levels :
142  (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
143  for (int l = 0; l < ref_levels; l++)
144  {
145  mesh->UniformRefinement();
146  }
147  }
148 
149  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
150  // this mesh further in parallel to increase the resolution. Once the
151  // parallel mesh is defined, the serial mesh can be deleted.
152  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
153  delete mesh;
154  {
155  for (int l = 0; l < par_ref_levels; l++)
156  {
157  pmesh->UniformRefinement();
158  }
159  }
160 
161  // 7. Define a parallel finite element space on the parallel mesh. Here we
162  // use vector finite elements, i.e. dim copies of a scalar finite element
163  // space. We use the ordering by vector dimension (the last argument of
164  // the FiniteElementSpace constructor) which is expected in the systems
165  // version of BoomerAMG preconditioner. For NURBS meshes, we use the
166  // (degree elevated) NURBS space associated with the mesh nodes.
168  ParFiniteElementSpace *fespace;
169  const bool use_nodal_fespace = pmesh->NURBSext && !amg_elast;
170  if (use_nodal_fespace)
171  {
172  fec = NULL;
173  fespace = (ParFiniteElementSpace *)pmesh->GetNodes()->FESpace();
174  }
175  else
176  {
177  fec = new H1_FECollection(order, dim);
178  fespace = new ParFiniteElementSpace(pmesh, fec, dim, Ordering::byVDIM);
179  }
180  HYPRE_BigInt size = fespace->GlobalTrueVSize();
181  if (myid == 0)
182  {
183  cout << "Number of finite element unknowns: " << size << endl
184  << "Assembling: " << flush;
185  }
186 
187  // 8. Determine the list of true (i.e. parallel conforming) essential
188  // boundary dofs. In this example, the boundary conditions are defined by
189  // marking only boundary attribute 1 from the mesh as essential and
190  // converting it to a list of true dofs.
191  Array<int> ess_tdof_list, ess_bdr(pmesh->bdr_attributes.Max());
192  ess_bdr = 0;
193  ess_bdr[0] = 1;
194  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
195 
196  // 9. Set up the parallel linear form b(.) which corresponds to the
197  // right-hand side of the FEM linear system. In this case, b_i equals the
198  // boundary integral of f*phi_i where f represents a "pull down" force on
199  // the Neumann part of the boundary and phi_i are the basis functions in
200  // the finite element fespace. The force is defined by the object f, which
201  // is a vector of Coefficient objects. The fact that f is non-zero on
202  // boundary attribute 2 is indicated by the use of piece-wise constants
203  // coefficient for its last component.
205  for (int i = 0; i < dim-1; i++)
206  {
207  f.Set(i, new ConstantCoefficient(0.0));
208  }
209  {
210  Vector pull_force(pmesh->bdr_attributes.Max());
211  pull_force = 0.0;
212  pull_force(1) = -1.0e-2;
213  f.Set(dim-1, new PWConstCoefficient(pull_force));
214  }
215 
216  ParLinearForm *b = new ParLinearForm(fespace);
218  if (myid == 0)
219  {
220  cout << "r.h.s. ... " << flush;
221  }
222  b->Assemble();
223 
224  // 10. Define the solution vector x as a parallel finite element grid
225  // function corresponding to fespace. Initialize x with initial guess of
226  // zero, which satisfies the boundary conditions.
227  ParGridFunction x(fespace);
228  x = 0.0;
229 
230  // 11. Set up the parallel bilinear form a(.,.) on the finite element space
231  // corresponding to the linear elasticity integrator with piece-wise
232  // constants coefficient lambda and mu.
233  Vector lambda(pmesh->attributes.Max());
234  lambda = 1.0;
235  lambda(0) = lambda(1)*50;
236  PWConstCoefficient lambda_func(lambda);
237  Vector mu(pmesh->attributes.Max());
238  mu = 1.0;
239  mu(0) = mu(1)*50;
240  PWConstCoefficient mu_func(mu);
241 
242  ParBilinearForm *a = new ParBilinearForm(fespace);
243  a->AddDomainIntegrator(new ElasticityIntegrator(lambda_func, mu_func));
244 
245  // 12. Assemble the parallel bilinear form and the corresponding linear
246  // system, applying any necessary transformations such as: parallel
247  // assembly, eliminating boundary conditions, applying conforming
248  // constraints for non-conforming AMR, static condensation, etc.
249  if (myid == 0) { cout << "matrix ... " << flush; }
250  if (static_cond) { a->EnableStaticCondensation(); }
251  a->Assemble();
252 
253  Vector B, X;
254  if (!use_petsc)
255  {
256  HypreParMatrix A;
257  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
258  if (myid == 0)
259  {
260  cout << "done." << endl;
261  cout << "Size of linear system: " << A.GetGlobalNumRows() << endl;
262  }
263 
264  // 13. Define and apply a parallel PCG solver for A X = B with the BoomerAMG
265  // preconditioner from hypre.
266  HypreBoomerAMG *amg = new HypreBoomerAMG(A);
267  if (amg_elast && !a->StaticCondensationIsEnabled())
268  {
269  amg->SetElasticityOptions(fespace);
270  }
271  else
272  {
273  amg->SetSystemsOptions(dim);
274  }
275  HyprePCG *pcg = new HyprePCG(A);
276  pcg->SetTol(1e-8);
277  pcg->SetMaxIter(500);
278  pcg->SetPrintLevel(2);
279  pcg->SetPreconditioner(*amg);
280  pcg->Mult(B, X);
281  delete pcg;
282  delete amg;
283  }
284  else
285  {
286  // 13b. Use PETSc to solve the linear system.
287  // Assemble a PETSc matrix, so that PETSc solvers can be used natively.
288  PetscParMatrix A;
289  a->SetOperatorType(use_nonoverlapping ?
290  Operator::PETSC_MATIS : Operator::PETSC_MATAIJ);
291  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
292  if (myid == 0)
293  {
294  cout << "done." << endl;
295  cout << "Size of linear system: " << A.M() << endl;
296  }
297  PetscPCGSolver *pcg = new PetscPCGSolver(A);
298 
299  // The preconditioner for the PCG solver defined below is specified in the
300  // PETSc config file, rc_ex2p, since a Krylov solver in PETSc can also
301  // customize its preconditioner.
302  PetscPreconditioner *prec = NULL;
303  if (use_nonoverlapping)
304  {
305  // Compute dofs belonging to the natural boundary
306  Array<int> nat_tdof_list, nat_bdr(pmesh->bdr_attributes.Max());
307  nat_bdr = 1;
308  nat_bdr[0] = 0;
309  fespace->GetEssentialTrueDofs(nat_bdr, nat_tdof_list);
310 
311  // Auxiliary class for BDDC customization
313  // Inform the solver about the finite element space
314  opts.SetSpace(fespace);
315  // Inform the solver about essential dofs
316  opts.SetEssBdrDofs(&ess_tdof_list);
317  // Inform the solver about natural dofs
318  opts.SetNatBdrDofs(&nat_tdof_list);
319  // Create a BDDC solver with parameters
320  prec = new PetscBDDCSolver(A,opts);
321  pcg->SetPreconditioner(*prec);
322  }
323 
324  pcg->SetMaxIter(500);
325  pcg->SetTol(1e-8);
326  pcg->SetPrintLevel(2);
327  pcg->Mult(B, X);
328  delete pcg;
329  delete prec;
330  }
331 
332  // 14. Recover the parallel grid function corresponding to X. This is the
333  // local finite element solution on each processor.
334  a->RecoverFEMSolution(X, *b, x);
335 
336  // 15. For non-NURBS meshes, make the mesh curved based on the finite element
337  // space. This means that we define the mesh elements through a fespace
338  // based transformation of the reference element. This allows us to save
339  // the displaced mesh as a curved mesh when using high-order finite
340  // element displacement field. We assume that the initial mesh (read from
341  // the file) is not higher order curved mesh compared to the chosen FE
342  // space.
343  if (!use_nodal_fespace)
344  {
345  pmesh->SetNodalFESpace(fespace);
346  }
347 
348  // 16. Save in parallel the displaced mesh and the inverted solution (which
349  // gives the backward displacements to the original grid). This output
350  // can be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
351  {
352  GridFunction *nodes = pmesh->GetNodes();
353  *nodes += x;
354  x *= -1;
355 
356  ostringstream mesh_name, sol_name;
357  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
358  sol_name << "sol." << setfill('0') << setw(6) << myid;
359 
360  ofstream mesh_ofs(mesh_name.str().c_str());
361  mesh_ofs.precision(8);
362  pmesh->Print(mesh_ofs);
363 
364  ofstream sol_ofs(sol_name.str().c_str());
365  sol_ofs.precision(8);
366  x.Save(sol_ofs);
367  }
368 
369  // 17. Send the above data by socket to a GLVis server. Use the "n" and "b"
370  // keys in GLVis to visualize the displacements.
371  if (visualization)
372  {
373  char vishost[] = "localhost";
374  int visport = 19916;
375  socketstream sol_sock(vishost, visport);
376  sol_sock << "parallel " << num_procs << " " << myid << "\n";
377  sol_sock.precision(8);
378  sol_sock << "solution\n" << *pmesh << x << flush;
379  }
380 
381  // 18. Free the used memory.
382  delete a;
383  delete b;
384  if (fec)
385  {
386  delete fespace;
387  delete fec;
388  }
389  delete pmesh;
390 
391  // We finalize PETSc
392  if (use_petsc) { MFEMFinalizePetsc(); }
393 
394  return 0;
395 }
void SetTol(double tol)
Definition: hypre.cpp:3775
void SetPreconditioner(Solver &precond)
Definition: petsc.cpp:3021
void SetPrintLevel(int plev)
Definition: petsc.cpp:2370
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:1032
Vector coefficient defined by an array of scalar coefficients. Coefficients that are not set will eva...
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
Abstract class for PETSc&#39;s preconditioners.
Definition: petsc.hpp:775
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
Wrapper for PETSc&#39;s matrix class.
Definition: petsc.hpp:307
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:285
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:928
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:873
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetMaxIter(int max_iter)
Definition: petsc.cpp:2343
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3795
PetscInt M() const
Returns the global number of rows.
Definition: petsc.cpp:924
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1443
void SetOperatorType(Operator::Type tid)
Set the operator type id for the parallel matrix/operator when using AssemblyLevel::LEGACY.
Class for parallel linear form.
Definition: plinearform.hpp:26
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
constexpr char vishost[]
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9737
constexpr int visport
Auxiliary class for BDDC customization.
Definition: petsc.hpp:801
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void Set(int i, Coefficient *c, bool own=true)
Sets coefficient in the vector.
void SetElasticityOptions(ParFiniteElementSpace *fespace)
Definition: hypre.cpp:4793
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:5273
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:70
int Dimension() const
Definition: mesh.hpp:999
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3785
bool StaticCondensationIsEnabled() const
Check if static condensation was actually enabled by a previous call to EnableStaticCondensation().
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:270
PCG solver in hypre.
Definition: hypre.hpp:1116
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
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 &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
HYPRE_Int HYPRE_BigInt
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
FDualNumber< tbase > log(const FDualNumber< tbase > &f)
log([dual number])
Definition: fdual.hpp:515
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:272
void MFEMFinalizePetsc()
Definition: petsc.cpp:192
A piecewise constant coefficient with the constants keyed off the element attribute numbers...
Definition: coefficient.hpp:94
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:629
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
Definition: petsc.cpp:168
double mu
Definition: ex25.cpp:139
void SetSpace(ParFiniteElementSpace *fe)
Definition: petsc.hpp:816
void DegreeElevate(int rel_degree, int degree=16)
Definition: mesh.cpp:5095
int dim
Definition: ex24.cpp:53
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:3800
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:324
void SetTol(double tol)
Definition: petsc.cpp:2288
Class for parallel bilinear form.
Vector data type.
Definition: vector.hpp:60
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7841
void SetSystemsOptions(int dim, bool order_bynodes=false)
Definition: hypre.cpp:4672
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
void Print(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4684
Class for parallel grid function.
Definition: pgridfunc.hpp:32
void SetNatBdrDofs(const Array< int > *natdofs, bool loc=false)
Specify dofs on the natural boundary.
Definition: petsc.hpp:829
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:337
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:3823
virtual void Mult(const Vector &b, Vector &x) const
Application of the solver.
Definition: petsc.cpp:3101
Class for parallel meshes.
Definition: pmesh.hpp:32
void SetEssBdrDofs(const Array< int > *essdofs, bool loc=false)
Specify dofs on the essential boundary.
Definition: petsc.hpp:821
int main()
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:268
void EnableStaticCondensation()
Enable the use of static condensation. For details see the description for class StaticCondensation i...
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150