MFEM  v3.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
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 //
19 // Description: This example code demonstrates the use of MFEM to define a
20 // simple finite element discretization of the Laplace problem
21 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
22 // Specifically, we discretize using a FE space of the specified
23 // order, or if order < 1 using an isoparametric/isogeometric
24 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
25 // NURBS mesh, etc.)
26 //
27 // The example highlights the use of mesh refinement, finite
28 // element grid functions, as well as linear and bilinear forms
29 // corresponding to the left-hand side and right-hand side of the
30 // discrete linear system. We also cover the explicit elimination
31 // of boundary conditions on all boundary edges, and the optional
32 // connection to the GLVis tool for visualization.
33 
34 #include "mfem.hpp"
35 #include <fstream>
36 #include <iostream>
37 
38 using namespace std;
39 using namespace mfem;
40 
41 int main(int argc, char *argv[])
42 {
43  // 1. Initialize MPI.
44  int num_procs, myid;
45  MPI_Init(&argc, &argv);
46  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
47  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
48 
49  // 2. Parse command-line options.
50  const char *mesh_file = "../data/star.mesh";
51  int order = 1;
52  bool visualization = 1;
53 
54  OptionsParser args(argc, argv);
55  args.AddOption(&mesh_file, "-m", "--mesh",
56  "Mesh file to use.");
57  args.AddOption(&order, "-o", "--order",
58  "Finite element order (polynomial degree) or -1 for"
59  " isoparametric space.");
60  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
61  "--no-visualization",
62  "Enable or disable GLVis visualization.");
63  args.Parse();
64  if (!args.Good())
65  {
66  if (myid == 0)
67  args.PrintUsage(cout);
68  MPI_Finalize();
69  return 1;
70  }
71  if (myid == 0)
72  args.PrintOptions(cout);
73 
74  // 3. Read the (serial) mesh from the given mesh file on all processors. We
75  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
76  // and volume meshes with the same code.
77  Mesh *mesh;
78  ifstream imesh(mesh_file);
79  if (!imesh)
80  {
81  if (myid == 0)
82  cerr << "\nCan not open mesh file: " << mesh_file << '\n' << endl;
83  MPI_Finalize();
84  return 2;
85  }
86  mesh = new Mesh(imesh, 1, 1);
87  imesh.close();
88  int dim = mesh->Dimension();
89 
90  // 4. Refine the serial mesh on all processors to increase the resolution. In
91  // this example we do 'ref_levels' of uniform refinement. We choose
92  // 'ref_levels' to be the largest number that gives a final mesh with no
93  // more than 10,000 elements.
94  {
95  int ref_levels =
96  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
97  for (int l = 0; l < ref_levels; l++)
98  mesh->UniformRefinement();
99  }
100 
101  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
102  // this mesh further in parallel to increase the resolution. Once the
103  // parallel mesh is defined, the serial mesh can be deleted.
104  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
105  delete mesh;
106  {
107  int par_ref_levels = 2;
108  for (int l = 0; l < par_ref_levels; l++)
109  pmesh->UniformRefinement();
110  }
111 
112  // 6. Define a parallel finite element space on the parallel mesh. Here we
113  // use continuous Lagrange finite elements of the specified order. If
114  // order < 1, we instead use an isoparametric/isogeometric space.
116  if (order > 0)
117  fec = new H1_FECollection(order, dim);
118  else if (pmesh->GetNodes())
119  fec = pmesh->GetNodes()->OwnFEC();
120  else
121  fec = new H1_FECollection(order = 1, dim);
122  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
123  int size = fespace->GlobalTrueVSize();
124  if (myid == 0)
125  cout << "Number of unknowns: " << size << endl;
126 
127  // 7. Set up the parallel linear form b(.) which corresponds to the
128  // right-hand side of the FEM linear system, which in this case is
129  // (1,phi_i) where phi_i are the basis functions in fespace.
130  ParLinearForm *b = new ParLinearForm(fespace);
131  ConstantCoefficient one(1.0);
133  b->Assemble();
134 
135  // 8. Define the solution vector x as a parallel finite element grid function
136  // corresponding to fespace. Initialize x with initial guess of zero,
137  // which satisfies the boundary conditions.
138  ParGridFunction x(fespace);
139  x = 0.0;
140 
141  // 9. Set up the parallel bilinear form a(.,.) on the finite element space
142  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
143  // domain integrator and imposing homogeneous Dirichlet boundary
144  // conditions. The boundary conditions are implemented by marking all the
145  // boundary attributes from the mesh as essential. After serial and
146  // parallel assembly we extract the corresponding parallel matrix A.
147  ParBilinearForm *a = new ParBilinearForm(fespace);
149  a->Assemble();
150  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
151  ess_bdr = 1;
152  a->EliminateEssentialBC(ess_bdr, x, *b);
153  a->Finalize();
154 
155  // 10. Define the parallel (hypre) matrix and vectors representing a(.,.),
156  // b(.) and the finite element approximation.
160 
161  delete a;
162  delete b;
163 
164  // 11. Define and apply a parallel PCG solver for AX=B with the BoomerAMG
165  // preconditioner from hypre.
166  HypreSolver *amg = new HypreBoomerAMG(*A);
167  HyprePCG *pcg = new HyprePCG(*A);
168  pcg->SetTol(1e-12);
169  pcg->SetMaxIter(200);
170  pcg->SetPrintLevel(2);
171  pcg->SetPreconditioner(*amg);
172  pcg->Mult(*B, *X);
173 
174  // 12. Extract the parallel grid function corresponding to the finite element
175  // approximation X. This is the local solution on each processor.
176  x = *X;
177 
178  // 13. Save the refined mesh and the solution in parallel. This output can
179  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
180  {
181  ostringstream mesh_name, sol_name;
182  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
183  sol_name << "sol." << setfill('0') << setw(6) << myid;
184 
185  ofstream mesh_ofs(mesh_name.str().c_str());
186  mesh_ofs.precision(8);
187  pmesh->Print(mesh_ofs);
188 
189  ofstream sol_ofs(sol_name.str().c_str());
190  sol_ofs.precision(8);
191  x.Save(sol_ofs);
192  }
193 
194  // 14. Send the solution by socket to a GLVis server.
195  if (visualization)
196  {
197  char vishost[] = "localhost";
198  int visport = 19916;
199  socketstream sol_sock(vishost, visport);
200  sol_sock << "parallel " << num_procs << " " << myid << "\n";
201  sol_sock.precision(8);
202  sol_sock << "solution\n" << *pmesh << x << flush;
203  }
204 
205  // 15. Free the used memory.
206  delete pcg;
207  delete amg;
208  delete X;
209  delete B;
210  delete A;
211 
212  delete fespace;
213  if (order > 0)
214  delete fec;
215  delete pmesh;
216 
217  MPI_Finalize();
218 
219  return 0;
220 }
void SetTol(double tol)
Definition: hypre.cpp:1339
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:396
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:250
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:1354
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
The BoomerAMG solver in hypre.
Definition: hypre.hpp:519
Class for parallel linear form.
Definition: plinearform.hpp:26
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6225
T Max() const
Definition: array.cpp:78
void Assemble(int skip_zeros=1)
Assemble the local matrix.
int Dimension() const
Definition: mesh.hpp:417
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:385
void SetMaxIter(int max_iter)
Definition: hypre.cpp:1344
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:38
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Array< int > bdr_attributes
Definition: mesh.hpp:305
int main(int argc, char *argv[])
Definition: ex1.cpp:39
PCG solver in hypre.
Definition: hypre.hpp:381
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 Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:1360
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:266
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:345
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:4948
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:52
Class for parallel grid function.
Definition: pgridfunc.hpp:31
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:103
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:1377
Class for parallel meshes.
Definition: pmesh.hpp:27
void ParallelAverage(Vector &tv) const
Returns the vector averaged on the true dofs.
Definition: pgridfunc.cpp:110
void EliminateEssentialBC(Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs, int d=0)
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:34
virtual void Print(std::ostream &out=std::cout) const
Definition: pmesh.cpp:2444
bool Good() const
Definition: optparser.hpp:120