MFEM  v3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
examples/ex1p.cpp
Go to the documentation of this file.
1 // MFEM Example 1 - Parallel Version
2 //
3 // Compile with: make ex1p
4 //
5 // Sample runs: mpirun -np 4 ex1p -m ../data/square-disc.mesh
6 // mpirun -np 4 ex1p -m ../data/star.mesh
7 // mpirun -np 4 ex1p -m ../data/escher.mesh
8 // mpirun -np 4 ex1p -m ../data/fichera.mesh
9 // mpirun -np 4 ex1p -m ../data/square-disc-p2.vtk -o 2
10 // mpirun -np 4 ex1p -m ../data/square-disc-p3.mesh -o 3
11 // mpirun -np 4 ex1p -m ../data/square-disc-nurbs.mesh -o -1
12 // mpirun -np 4 ex1p -m ../data/disc-nurbs.mesh -o -1
13 // mpirun -np 4 ex1p -m ../data/pipe-nurbs.mesh -o -1
14 // mpirun -np 4 ex1p -m ../data/ball-nurbs.mesh -o 2
15 // mpirun -np 4 ex1p -m ../data/star-surf.mesh
16 // mpirun -np 4 ex1p -m ../data/square-disc-surf.mesh
17 // mpirun -np 4 ex1p -m ../data/inline-segment.mesh
18 // mpirun -np 4 ex1p -m ../data/amr-quad.mesh
19 // mpirun -np 4 ex1p -m ../data/amr-hex.mesh
20 // mpirun -np 4 ex1p -m ../data/mobius-strip.mesh
21 // mpirun -np 4 ex1p -m ../data/mobius-strip.mesh -o -1 -sc
22 //
23 // Description: This example code demonstrates the use of MFEM to define a
24 // simple finite element discretization of the Laplace problem
25 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
26 // Specifically, we discretize using a FE space of the specified
27 // order, or if order < 1 using an isoparametric/isogeometric
28 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
29 // NURBS mesh, etc.)
30 //
31 // The example highlights the use of mesh refinement, finite
32 // element grid functions, as well as linear and bilinear forms
33 // corresponding to the left-hand side and right-hand side of the
34 // discrete linear system. We also cover the explicit elimination
35 // of essential boundary conditions, static condensation, and the
36 // optional connection to the GLVis tool for visualization.
37 
38 #include "mfem.hpp"
39 #include <fstream>
40 #include <iostream>
41 
42 using namespace std;
43 using namespace mfem;
44 
45 int main(int argc, char *argv[])
46 {
47  // 1. Initialize MPI.
48  int num_procs, myid;
49  MPI_Init(&argc, &argv);
50  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
51  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
52 
53  // 2. Parse command-line options.
54  const char *mesh_file = "../data/star.mesh";
55  int order = 1;
56  bool static_cond = false;
57  bool visualization = 1;
58 
59  OptionsParser args(argc, argv);
60  args.AddOption(&mesh_file, "-m", "--mesh",
61  "Mesh file to use.");
62  args.AddOption(&order, "-o", "--order",
63  "Finite element order (polynomial degree) or -1 for"
64  " isoparametric space.");
65  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
66  "--no-static-condensation", "Enable static condensation.");
67  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
68  "--no-visualization",
69  "Enable or disable GLVis visualization.");
70  args.Parse();
71  if (!args.Good())
72  {
73  if (myid == 0)
74  {
75  args.PrintUsage(cout);
76  }
77  MPI_Finalize();
78  return 1;
79  }
80  if (myid == 0)
81  {
82  args.PrintOptions(cout);
83  }
84 
85  // 3. Read the (serial) mesh from the given mesh file on all processors. We
86  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
87  // and volume meshes with the same code.
88  Mesh *mesh = new Mesh(mesh_file, 1, 1);
89  int dim = mesh->Dimension();
90 
91  // 4. Refine the serial mesh on all processors to increase the resolution. In
92  // this example we do 'ref_levels' of uniform refinement. We choose
93  // 'ref_levels' to be the largest number that gives a final mesh with no
94  // more than 10,000 elements.
95  {
96  int ref_levels =
97  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
98  for (int l = 0; l < ref_levels; l++)
99  {
100  mesh->UniformRefinement();
101  }
102  }
103 
104  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
105  // this mesh further in parallel to increase the resolution. Once the
106  // parallel mesh is defined, the serial mesh can be deleted.
107  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
108  delete mesh;
109  {
110  int par_ref_levels = 2;
111  for (int l = 0; l < par_ref_levels; l++)
112  {
113  pmesh->UniformRefinement();
114  }
115  }
116 
117  // 6. Define a parallel finite element space on the parallel mesh. Here we
118  // use continuous Lagrange finite elements of the specified order. If
119  // order < 1, we instead use an isoparametric/isogeometric space.
121  if (order > 0)
122  {
123  fec = new H1_FECollection(order, dim);
124  }
125  else if (pmesh->GetNodes())
126  {
127  fec = pmesh->GetNodes()->OwnFEC();
128  if (myid == 0)
129  {
130  cout << "Using isoparametric FEs: " << fec->Name() << endl;
131  }
132  }
133  else
134  {
135  fec = new H1_FECollection(order = 1, dim);
136  }
137  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
138  HYPRE_Int size = fespace->GlobalTrueVSize();
139  if (myid == 0)
140  {
141  cout << "Number of finite element unknowns: " << size << endl;
142  }
143 
144  // 7. Determine the list of true (i.e. parallel conforming) essential
145  // boundary dofs. In this example, the boundary conditions are defined
146  // by marking all the boundary attributes from the mesh as essential
147  // (Dirichlet) and converting them to a list of true dofs.
148  Array<int> ess_tdof_list;
149  if (pmesh->bdr_attributes.Size())
150  {
151  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
152  ess_bdr = 1;
153  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
154  }
155 
156  // 8. Set up the parallel linear form b(.) which corresponds to the
157  // right-hand side of the FEM linear system, which in this case is
158  // (1,phi_i) where phi_i are the basis functions in fespace.
159  ParLinearForm *b = new ParLinearForm(fespace);
160  ConstantCoefficient one(1.0);
162  b->Assemble();
163 
164  // 9. Define the solution vector x as a parallel finite element grid function
165  // corresponding to fespace. Initialize x with initial guess of zero,
166  // which satisfies the boundary conditions.
167  ParGridFunction x(fespace);
168  x = 0.0;
169 
170  // 10. Set up the parallel bilinear form a(.,.) on the finite element space
171  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
172  // domain integrator.
173  ParBilinearForm *a = new ParBilinearForm(fespace);
175 
176  // 11. Assemble the parallel bilinear form and the corresponding linear
177  // system, applying any necessary transformations such as: parallel
178  // assembly, eliminating boundary conditions, applying conforming
179  // constraints for non-conforming AMR, static condensation, etc.
180  if (static_cond) { a->EnableStaticCondensation(); }
181  a->Assemble();
182 
183  HypreParMatrix A;
184  Vector B, X;
185  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
186 
187  if (myid == 0)
188  {
189  cout << "Size of linear system: " << A.GetGlobalNumRows() << endl;
190  }
191 
192  // 12. Define and apply a parallel PCG solver for AX=B with the BoomerAMG
193  // preconditioner from hypre.
194  HypreSolver *amg = new HypreBoomerAMG(A);
195  HyprePCG *pcg = new HyprePCG(A);
196  pcg->SetTol(1e-12);
197  pcg->SetMaxIter(200);
198  pcg->SetPrintLevel(2);
199  pcg->SetPreconditioner(*amg);
200  pcg->Mult(B, X);
201 
202  // 13. Recover the parallel grid function corresponding to X. This is the
203  // local finite element solution on each processor.
204  a->RecoverFEMSolution(X, *b, x);
205 
206  // 14. Save the refined mesh and the solution in parallel. This output can
207  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
208  {
209  ostringstream mesh_name, sol_name;
210  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
211  sol_name << "sol." << setfill('0') << setw(6) << myid;
212 
213  ofstream mesh_ofs(mesh_name.str().c_str());
214  mesh_ofs.precision(8);
215  pmesh->Print(mesh_ofs);
216 
217  ofstream sol_ofs(sol_name.str().c_str());
218  sol_ofs.precision(8);
219  x.Save(sol_ofs);
220  }
221 
222  // 15. Send the solution by socket to a GLVis server.
223  if (visualization)
224  {
225  char vishost[] = "localhost";
226  int visport = 19916;
227  socketstream sol_sock(vishost, visport);
228  sol_sock << "parallel " << num_procs << " " << myid << "\n";
229  sol_sock.precision(8);
230  sol_sock << "solution\n" << *pmesh << x << flush;
231  }
232 
233  // 16. Free the used memory.
234  delete pcg;
235  delete amg;
236  delete a;
237  delete b;
238  delete fespace;
239  if (order > 0) { delete fec; }
240  delete pmesh;
241 
242  MPI_Finalize();
243 
244  return 0;
245 }
void SetTol(double tol)
Definition: hypre.cpp:1964
int Size() const
Logical size of the array.
Definition: array.hpp:109
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:47
Subclass constant coefficient.
Definition: coefficient.hpp:57
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:34
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:496
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:312
Abstract parallel finite element space.
Definition: pfespace.hpp:28
int main(int argc, char *argv[])
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:1979
The BoomerAMG solver in hypre.
Definition: hypre.hpp:705
Class for parallel linear form.
Definition: plinearform.hpp:26
HYPRE_Int GetGlobalNumRows() const
Definition: hypre.hpp:340
int dim
Definition: ex3.cpp:47
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6402
T Max() const
Definition: array.cpp:90
void Assemble(int skip_zeros=1)
Assemble the local matrix.
int Dimension() const
Definition: mesh.hpp:523
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetMaxIter(int max_iter)
Definition: hypre.cpp:1969
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Array< int > bdr_attributes
Definition: mesh.hpp:141
PCG solver in hypre.
Definition: hypre.hpp:565
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list)
Definition: pfespace.cpp:479
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
void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
virtual const char * Name() const
Definition: fe_coll.hpp:50
void FormLinearSystem(Array< int > &ess_tdof_list, Vector &x, Vector &b, HypreParMatrix &A, Vector &X, Vector &B, int copy_interior=0)
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:1984
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:529
Vector data type.
Definition: vector.hpp:33
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5114
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:77
Class for parallel grid function.
Definition: pgridfunc.hpp:31
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:150
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2005
Class for parallel meshes.
Definition: pmesh.hpp:28
void EnableStaticCondensation()
virtual void Print(std::ostream &out=std::cout) const
Definition: pmesh.cpp:2927
bool Good() const
Definition: optparser.hpp:120