MFEM  v3.3
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
petsc/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 // We recommend viewing Example 1 before viewing this example.
35 
36 #include "mfem.hpp"
37 #include <fstream>
38 #include <iostream>
39 
40 #ifndef MFEM_USE_PETSC
41 #error This example requires that MFEM is built with MFEM_USE_PETSC=YES
42 #endif
43 
44 using namespace std;
45 using namespace mfem;
46 
47 int main(int argc, char *argv[])
48 {
49  // 1. Initialize MPI.
50  int num_procs, myid;
51  MPI_Init(&argc, &argv);
52  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
53  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
54 
55  // 2. Parse command-line options.
56  const char *mesh_file = "../../data/beam-tri.mesh";
57  int order = 1;
58  bool static_cond = false;
59  bool visualization = 1;
60  bool amg_elast = 0;
61  bool use_petsc = true;
62  const char *petscrc_file = "";
63  bool use_nonoverlapping = false;
64 
65  OptionsParser args(argc, argv);
66  args.AddOption(&mesh_file, "-m", "--mesh",
67  "Mesh file to use.");
68  args.AddOption(&order, "-o", "--order",
69  "Finite element order (polynomial degree).");
70  args.AddOption(&amg_elast, "-elast", "--amg-for-elasticity", "-sys",
71  "--amg-for-systems",
72  "Use the special AMG elasticity solver (GM/LN approaches), "
73  "or standard AMG for systems (unknown approach).");
74  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
75  "--no-static-condensation", "Enable static condensation.");
76  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
77  "--no-visualization",
78  "Enable or disable GLVis visualization.");
79  args.AddOption(&use_petsc, "-usepetsc", "--usepetsc", "-no-petsc",
80  "--no-petsc",
81  "Use or not PETSc to solve the linear system.");
82  args.AddOption(&petscrc_file, "-petscopts", "--petscopts",
83  "PetscOptions file to use.");
84  args.AddOption(&use_nonoverlapping, "-nonoverlapping", "--nonoverlapping",
85  "-no-nonoverlapping", "--no-nonoverlapping",
86  "Use or not the block diagonal PETSc's matrix format "
87  "for non-overlapping domain decomposition.");
88  args.Parse();
89  if (!args.Good())
90  {
91  if (myid == 0)
92  {
93  args.PrintUsage(cout);
94  }
95  MPI_Finalize();
96  return 1;
97  }
98  if (myid == 0)
99  {
100  args.PrintOptions(cout);
101  }
102 
103  // 2b. We initialize PETSc
104  if (use_petsc) { PetscInitialize(NULL,NULL,petscrc_file,NULL); }
105 
106  // 3. Read the (serial) mesh from the given mesh file on all processors. We
107  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
108  // and volume meshes with the same code.
109  Mesh *mesh = new Mesh(mesh_file, 1, 1);
110  int dim = mesh->Dimension();
111 
112  if (mesh->attributes.Max() < 2 || mesh->bdr_attributes.Max() < 2)
113  {
114  if (myid == 0)
115  cerr << "\nInput mesh should have at least two materials and "
116  << "two boundary attributes! (See schematic in ex2.cpp)\n"
117  << endl;
118  MPI_Finalize();
119  return 3;
120  }
121 
122  // 4. Select the order of the finite element discretization space. For NURBS
123  // meshes, we increase the order by degree elevation.
124  if (mesh->NURBSext && order > mesh->NURBSext->GetOrder())
125  {
126  mesh->DegreeElevate(order - mesh->NURBSext->GetOrder());
127  }
128 
129  // 5. Refine the serial mesh on all processors to increase the resolution. In
130  // this example we do 'ref_levels' of uniform refinement. We choose
131  // 'ref_levels' to be the largest number that gives a final mesh with no
132  // more than 1,000 elements.
133  {
134  int ref_levels =
135  (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
136  for (int l = 0; l < ref_levels; l++)
137  {
138  mesh->UniformRefinement();
139  }
140  }
141 
142  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
143  // this mesh further in parallel to increase the resolution. Once the
144  // parallel mesh is defined, the serial mesh can be deleted.
145  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
146  delete mesh;
147  {
148  int par_ref_levels = 1;
149  for (int l = 0; l < par_ref_levels; l++)
150  {
151  pmesh->UniformRefinement();
152  }
153  }
154 
155  // 7. Define a parallel finite element space on the parallel mesh. Here we
156  // use vector finite elements, i.e. dim copies of a scalar finite element
157  // space. We use the ordering by vector dimension (the last argument of
158  // the FiniteElementSpace constructor) which is expected in the systems
159  // version of BoomerAMG preconditioner. For NURBS meshes, we use the
160  // (degree elevated) NURBS space associated with the mesh nodes.
162  ParFiniteElementSpace *fespace;
163  const bool use_nodal_fespace = pmesh->NURBSext && !amg_elast;
164  if (use_nodal_fespace)
165  {
166  fec = NULL;
167  fespace = (ParFiniteElementSpace *)pmesh->GetNodes()->FESpace();
168  }
169  else
170  {
171  fec = new H1_FECollection(order, dim);
172  fespace = new ParFiniteElementSpace(pmesh, fec, dim, Ordering::byVDIM);
173  }
174  HYPRE_Int size = fespace->GlobalTrueVSize();
175  if (myid == 0)
176  {
177  cout << "Number of finite element unknowns: " << size << endl
178  << "Assembling: " << flush;
179  }
180 
181  // 8. Determine the list of true (i.e. parallel conforming) essential
182  // boundary dofs. In this example, the boundary conditions are defined by
183  // marking only boundary attribute 1 from the mesh as essential and
184  // converting it to a list of true dofs.
185  Array<int> ess_tdof_list, ess_bdr(pmesh->bdr_attributes.Max());
186  ess_bdr = 0;
187  ess_bdr[0] = 1;
188  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
189 
190  // 9. Set up the parallel linear form b(.) which corresponds to the
191  // right-hand side of the FEM linear system. In this case, b_i equals the
192  // boundary integral of f*phi_i where f represents a "pull down" force on
193  // the Neumann part of the boundary and phi_i are the basis functions in
194  // the finite element fespace. The force is defined by the object f, which
195  // is a vector of Coefficient objects. The fact that f is non-zero on
196  // boundary attribute 2 is indicated by the use of piece-wise constants
197  // coefficient for its last component.
198  VectorArrayCoefficient f(dim);
199  for (int i = 0; i < dim-1; i++)
200  {
201  f.Set(i, new ConstantCoefficient(0.0));
202  }
203  {
204  Vector pull_force(pmesh->bdr_attributes.Max());
205  pull_force = 0.0;
206  pull_force(1) = -1.0e-2;
207  f.Set(dim-1, new PWConstCoefficient(pull_force));
208  }
209 
210  ParLinearForm *b = new ParLinearForm(fespace);
212  if (myid == 0)
213  {
214  cout << "r.h.s. ... " << flush;
215  }
216  b->Assemble();
217 
218  // 10. Define the solution vector x as a parallel finite element grid
219  // function corresponding to fespace. Initialize x with initial guess of
220  // zero, which satisfies the boundary conditions.
221  ParGridFunction x(fespace);
222  x = 0.0;
223 
224  // 11. Set up the parallel bilinear form a(.,.) on the finite element space
225  // corresponding to the linear elasticity integrator with piece-wise
226  // constants coefficient lambda and mu.
227  Vector lambda(pmesh->attributes.Max());
228  lambda = 1.0;
229  lambda(0) = lambda(1)*50;
230  PWConstCoefficient lambda_func(lambda);
231  Vector mu(pmesh->attributes.Max());
232  mu = 1.0;
233  mu(0) = mu(1)*50;
234  PWConstCoefficient mu_func(mu);
235 
236  ParBilinearForm *a = new ParBilinearForm(fespace);
237  a->AddDomainIntegrator(new ElasticityIntegrator(lambda_func, mu_func));
238 
239  // 12. Assemble the parallel bilinear form and the corresponding linear
240  // system, applying any necessary transformations such as: parallel
241  // assembly, eliminating boundary conditions, applying conforming
242  // constraints for non-conforming AMR, static condensation, etc.
243  if (myid == 0) { cout << "matrix ... " << flush; }
244  if (static_cond) { a->EnableStaticCondensation(); }
245  a->Assemble();
246 
247  Vector B, X;
248  if (!use_petsc)
249  {
250  HypreParMatrix A;
251  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
252  if (myid == 0)
253  {
254  cout << "done." << endl;
255  cout << "Size of linear system: " << A.GetGlobalNumRows() << endl;
256  }
257 
258  // 13. Define and apply a parallel PCG solver for A X = B with the BoomerAMG
259  // preconditioner from hypre.
260  HypreBoomerAMG *amg = new HypreBoomerAMG(A);
261  if (amg_elast && !a->StaticCondensationIsEnabled())
262  {
263  amg->SetElasticityOptions(fespace);
264  }
265  else
266  {
267  amg->SetSystemsOptions(dim);
268  }
269  HyprePCG *pcg = new HyprePCG(A);
270  pcg->SetTol(1e-8);
271  pcg->SetMaxIter(500);
272  pcg->SetPrintLevel(2);
273  pcg->SetPreconditioner(*amg);
274  pcg->Mult(B, X);
275  delete pcg;
276  delete amg;
277  }
278  else
279  {
280  // 13b. Use PETSc to solve the linear system.
281  // Assemble a PETSc matrix, so that PETSc solvers can be used natively.
282  PetscParMatrix A;
283  a->SetOperatorType(use_nonoverlapping ?
284  Operator::PETSC_MATIS : Operator::PETSC_MATAIJ);
285  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
286  if (myid == 0)
287  {
288  cout << "done." << endl;
289  cout << "Size of linear system: " << A.M() << endl;
290  }
291  // The preconditioner for the PCG solver defined below is specified in the
292  // PETSc config file, rc_ex2p, since a Krylov solver in PETSc can also
293  // customize its preconditioner.
294  PetscPCGSolver *pcg = new PetscPCGSolver(A);
295  pcg->SetMaxIter(500);
296  pcg->SetTol(1e-8);
297  pcg->SetPrintLevel(2);
298  pcg->Mult(B, X);
299  delete pcg;
300  }
301 
302  // 14. Recover the parallel grid function corresponding to X. This is the
303  // local finite element solution on each processor.
304  a->RecoverFEMSolution(X, *b, x);
305 
306  // 15. For non-NURBS meshes, make the mesh curved based on the finite element
307  // space. This means that we define the mesh elements through a fespace
308  // based transformation of the reference element. This allows us to save
309  // the displaced mesh as a curved mesh when using high-order finite
310  // element displacement field. We assume that the initial mesh (read from
311  // the file) is not higher order curved mesh compared to the chosen FE
312  // space.
313  if (!use_nodal_fespace)
314  {
315  pmesh->SetNodalFESpace(fespace);
316  }
317 
318  // 16. Save in parallel the displaced mesh and the inverted solution (which
319  // gives the backward displacements to the original grid). This output
320  // can be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
321  {
322  GridFunction *nodes = pmesh->GetNodes();
323  *nodes += x;
324  x *= -1;
325 
326  ostringstream mesh_name, sol_name;
327  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
328  sol_name << "sol." << setfill('0') << setw(6) << myid;
329 
330  ofstream mesh_ofs(mesh_name.str().c_str());
331  mesh_ofs.precision(8);
332  pmesh->Print(mesh_ofs);
333 
334  ofstream sol_ofs(sol_name.str().c_str());
335  sol_ofs.precision(8);
336  x.Save(sol_ofs);
337  }
338 
339  // 17. Send the above data by socket to a GLVis server. Use the "n" and "b"
340  // keys in GLVis to visualize the displacements.
341  if (visualization)
342  {
343  char vishost[] = "localhost";
344  int visport = 19916;
345  socketstream sol_sock(vishost, visport);
346  sol_sock << "parallel " << num_procs << " " << myid << "\n";
347  sol_sock.precision(8);
348  sol_sock << "solution\n" << *pmesh << x << flush;
349  }
350 
351  // 18. Free the used memory.
352  delete a;
353  delete b;
354  if (fec)
355  {
356  delete fespace;
357  delete fec;
358  }
359  delete pmesh;
360 
361  // We finalize PETSc
362  if (use_petsc) { PetscFinalize(); }
363 
364  MPI_Finalize();
365 
366  return 0;
367 }
void SetTol(double tol)
Definition: hypre.cpp:1993
void SetPrintLevel(int plev)
Definition: petsc.cpp:1293
Vector coefficient defined by an array of scalar coefficients.
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
Subclass constant coefficient.
Definition: coefficient.hpp:57
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Wrapper for PETSc&#39;s matrix class.
Definition: petsc.hpp:128
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:42
void DegreeElevate(int t)
Definition: mesh.cpp:3038
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:584
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:314
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetMaxIter(int max_iter)
Definition: petsc.cpp:1266
int main(int argc, char *argv[])
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2008
PetscInt M() const
Returns the global number of rows.
Definition: petsc.cpp:338
void SetSystemsOptions(int dim)
Definition: hypre.cpp:2407
The BoomerAMG solver in hypre.
Definition: hypre.hpp:770
void SetOperatorType(Operator::Type tid)
Set the operator type id for the parallel matrix/operator.
Class for parallel linear form.
Definition: plinearform.hpp:26
HYPRE_Int GetGlobalNumRows() const
Definition: hypre.hpp:371
int dim
Definition: ex3.cpp:47
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6684
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:108
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void SetElasticityOptions(ParFiniteElementSpace *fespace)
Definition: hypre.cpp:2489
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:3261
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator.
Definition: linearform.cpp:24
int Dimension() const
Definition: mesh.hpp:611
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetMaxIter(int max_iter)
Definition: hypre.cpp:1998
bool StaticCondensationIsEnabled() const
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:142
PCG solver in hypre.
Definition: hypre.hpp:630
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list)
Definition: pfespace.cpp:521
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)
Definition: optparser.hpp:74
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:144
class for piecewise constant coefficient
Definition: coefficient.hpp:72
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2013
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
void SetTol(double tol)
Definition: petsc.cpp:1211
Class for parallel bilinear form.
Vector data type.
Definition: vector.hpp:36
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5368
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:146
Class for parallel grid function.
Definition: pgridfunc.hpp:31
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:174
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2034
virtual void Mult(const Vector &b, Vector &x) const
Application of the solver.
Definition: petsc.cpp:1624
Class for parallel meshes.
Definition: pmesh.hpp:28
void Set(int i, Coefficient *c)
Sets coefficient in the vector.
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:140
void EnableStaticCondensation()
virtual void Print(std::ostream &out=std::cout) const
Definition: pmesh.cpp:3395
bool Good() const
Definition: optparser.hpp:120