MFEM  v3.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex7p.cpp
Go to the documentation of this file.
1 // MFEM Example 7 - Parallel Version
2 //
3 // Compile with: make ex7p
4 //
5 // Sample runs: mpirun -np 4 ex7p -e 0 -o 2 -r 4
6 // mpirun -np 4 ex7p -e 1 -o 2 -r 4 -snap
7 //
8 // Description: This example code demonstrates the use of MFEM to define a
9 // triangulation of a unit sphere and a simple isoparametric
10 // finite element discretization of the Laplace problem with mass
11 // term, -Delta u + u = f.
12 //
13 // The example highlights mesh generation, the use of mesh
14 // refinement, high-order meshes and finite elements, as well as
15 // surface-based linear and bilinear forms corresponding to the
16 // left-hand side and right-hand side of the discrete linear
17 // system.
18 //
19 // We recommend viewing Example 1 before viewing this example.
20 
21 #include "mfem.hpp"
22 #include <fstream>
23 #include <iostream>
24 
25 using namespace std;
26 using namespace mfem;
27 
28 // Exact solution and r.h.s., see below for implementation.
29 double analytic_solution(Vector &x);
30 double analytic_rhs(Vector &x);
31 void SnapNodes(Mesh &mesh);
32 
33 int main(int argc, char *argv[])
34 {
35  // 1. Initialize MPI.
36  int num_procs, myid;
37  MPI_Init(&argc, &argv);
38  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
39  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
40 
41  // 2. Parse command-line options.
42  int elem_type = 1;
43  int ref_levels = 2;
44  int order = 2;
45  bool always_snap = false;
46  bool visualization = 1;
47 
48  OptionsParser args(argc, argv);
49  args.AddOption(&elem_type, "-e", "--elem",
50  "Type of elements to use: 0 - triangles, 1 - quads.");
51  args.AddOption(&order, "-o", "--order",
52  "Finite element order (polynomial degree).");
53  args.AddOption(&ref_levels, "-r", "--refine",
54  "Number of times to refine the mesh uniformly.");
55  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
56  "--no-visualization",
57  "Enable or disable GLVis visualization.");
58  args.AddOption(&always_snap, "-snap", "--always-snap", "-no-snap",
59  "--snap-at-the-end",
60  "If true, snap nodes to the sphere initially and after each refinement "
61  "otherwise, snap only after the last refinement");
62  args.Parse();
63  if (!args.Good())
64  {
65  if (myid == 0)
66  args.PrintUsage(cout);
67  MPI_Finalize();
68  return 1;
69  }
70  if (myid == 0)
71  args.PrintOptions(cout);
72 
73  // 3. Generate an initial high-order (surface) mesh on the unit sphere. The
74  // Mesh object represents a 2D mesh in 3 spatial dimensions. We first add
75  // the elements and the vertices of the mesh, and then make it high-order
76  // by specifying a finite element space for its nodes.
77  int Nvert = 8, Nelem = 6;
78  if (elem_type == 0)
79  {
80  Nvert = 6;
81  Nelem = 8;
82  }
83  Mesh *mesh = new Mesh(2, Nvert, Nelem, 0, 3);
84 
85  if (elem_type == 0) // inscribed octahedron
86  {
87  const double tri_v[6][3] =
88  {{ 1, 0, 0}, { 0, 1, 0}, {-1, 0, 0},
89  { 0, -1, 0}, { 0, 0, 1}, { 0, 0, -1}};
90  const int tri_e[8][3] =
91  {{0, 1, 4}, {1, 2, 4}, {2, 3, 4}, {3, 0, 4},
92  {1, 0, 5}, {2, 1, 5}, {3, 2, 5}, {0, 3, 5}};
93 
94  for (int j = 0; j < Nvert; j++)
95  {
96  mesh->AddVertex(tri_v[j]);
97  }
98  for (int j = 0; j < Nelem; j++)
99  {
100  int attribute = j + 1;
101  mesh->AddTriangle(tri_e[j], attribute);
102  }
103  mesh->FinalizeTriMesh(1, 1, true);
104  }
105  else // inscribed cube
106  {
107  const double quad_v[8][3] =
108  {{-1, -1, -1}, {+1, -1, -1}, {+1, +1, -1}, {-1, +1, -1},
109  {-1, -1, +1}, {+1, -1, +1}, {+1, +1, +1}, {-1, +1, +1}};
110  const int quad_e[6][4] =
111  {{3, 2, 1, 0}, {0, 1, 5, 4}, {1, 2, 6, 5},
112  {2, 3, 7, 6}, {3, 0, 4, 7}, {4, 5, 6, 7}};
113 
114  for (int j = 0; j < Nvert; j++)
115  {
116  mesh->AddVertex(quad_v[j]);
117  }
118  for (int j = 0; j < Nelem; j++)
119  {
120  int attribute = j + 1;
121  mesh->AddQuad(quad_e[j], attribute);
122  }
123  mesh->FinalizeQuadMesh(1, 1, true);
124  }
125 
126  // Set the space for the high-order mesh nodes.
127  H1_FECollection fec(order, mesh->Dimension());
128  FiniteElementSpace nodal_fes(mesh, &fec, mesh->SpaceDimension());
129  mesh->SetNodalFESpace(&nodal_fes);
130 
131  // 4. Refine the mesh while snapping nodes to the sphere. Number of parallel
132  // refinements is fixed to 2.
133  for (int l = 0; l <= ref_levels; l++)
134  {
135  if (l > 0) // for l == 0 just perform snapping
136  mesh->UniformRefinement();
137 
138  // Snap the nodes of the refined mesh back to sphere surface.
139  if (always_snap)
140  SnapNodes(*mesh);
141  }
142 
143  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
144  delete mesh;
145  {
146  int par_ref_levels = 2;
147  for (int l = 0; l < par_ref_levels; l++)
148  {
149  pmesh->UniformRefinement();
150 
151  // Snap the nodes of the refined mesh back to sphere surface.
152  if (always_snap)
153  SnapNodes(*pmesh);
154  }
155  if (!always_snap || par_ref_levels < 1)
156  SnapNodes(*pmesh);
157  }
158 
159  // 5. Define a finite element space on the mesh. Here we use isoparametric
160  // finite elements -- the same as the mesh nodes.
161  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, &fec);
162  int size = fespace->GlobalTrueVSize();
163  if (myid == 0)
164  cout << "Number of unknowns: " << size << endl;
165 
166  // 6. Set up the linear form b(.) which corresponds to the right-hand side of
167  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
168  // the basis functions in the finite element fespace.
169  ParLinearForm *b = new ParLinearForm(fespace);
170  ConstantCoefficient one(1.0);
173  b->AddDomainIntegrator(new DomainLFIntegrator(rhs_coef));
174  b->Assemble();
175 
176  // 7. Define the solution vector x as a finite element grid function
177  // corresponding to fespace. Initialize x with initial guess of zero,
178  // which satisfies the boundary conditions.
179  ParGridFunction x(fespace);
180  x = 0.0;
181 
182  // 8. Set up the bilinear form a(.,.) on the finite element space
183  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
184  // domain integrator and imposing homogeneous Dirichlet boundary
185  // conditions. The boundary conditions are implemented by marking all the
186  // boundary attributes from the mesh as essential (Dirichlet).
187  ParBilinearForm *a = new ParBilinearForm(fespace);
189  a->AddDomainIntegrator(new MassIntegrator(one));
190  a->Assemble();
191  a->Finalize();
192 
193  // 9. Define the parallel (hypre) matrix and vectors representing a(.,.),
194  // b(.) and the finite element approximation.
195  HypreParMatrix * A = a->ParallelAssemble();
196  HypreParVector * B = b->ParallelAssemble();
198 
199  delete a;
200  delete b;
201 
202  // 10. Define and apply a parallel PCG solver for AX=B with the BoomerAMG
203  // preconditioner from hypre. Extract the parallel grid function x
204  // corresponding to the finite element approximation X. This is the local
205  // solution on each processor.
206  HypreSolver *amg = new HypreBoomerAMG(*A);
207  HyprePCG *pcg = new HyprePCG(*A);
208  pcg->SetTol(1e-12);
209  pcg->SetMaxIter(200);
210  pcg->SetPrintLevel(2);
211  pcg->SetPreconditioner(*amg);
212  pcg->Mult(*B, *X);
213  x = *X;
214 
215  // 11. Compute and print the L^2 norm of the error.
216  double err = x.ComputeL2Error(sol_coef);
217  if (myid == 0)
218  cout << "\nL2 norm of error: " << err << endl;
219 
220  // 12. Save the refined mesh and the solution. This output can be viewed
221  // later using GLVis: "glvis -np <np> -m sphere_refined -g sol".
222  {
223  ostringstream mesh_name, sol_name;
224  mesh_name << "sphere_refined." << setfill('0') << setw(6) << myid;
225  sol_name << "sol." << setfill('0') << setw(6) << myid;
226 
227  ofstream mesh_ofs(mesh_name.str().c_str());
228  mesh_ofs.precision(8);
229  pmesh->Print(mesh_ofs);
230 
231  ofstream sol_ofs(sol_name.str().c_str());
232  sol_ofs.precision(8);
233  x.Save(sol_ofs);
234  }
235 
236  // 13. Send the solution by socket to a GLVis server.
237  if (visualization)
238  {
239  char vishost[] = "localhost";
240  int visport = 19916;
241  socketstream sol_sock(vishost, visport);
242  sol_sock << "parallel " << num_procs << " " << myid << "\n";
243  sol_sock.precision(8);
244  sol_sock << "solution\n" << *pmesh << x << flush;
245  }
246 
247  // 14. Free the used memory.
248  delete pcg;
249  delete amg;
250  delete X;
251  delete B;
252  delete A;
253  delete fespace;
254  delete pmesh;
255 
256  MPI_Finalize();
257 
258  return 0;
259 }
260 
261 double analytic_solution(Vector &x)
262 {
263  double l2 = x(0)*x(0) + x(1)*x(1) + x(2)*x(2);
264  return x(0)*x(1)/l2;
265 }
266 
267 double analytic_rhs(Vector &x)
268 {
269  double l2 = x(0)*x(0) + x(1)*x(1) + x(2)*x(2);
270  return 7*x(0)*x(1)/l2;
271 }
272 
273 void SnapNodes(Mesh &mesh)
274 {
275  GridFunction &nodes = *mesh.GetNodes();
276  Vector node(mesh.SpaceDimension());
277  for (int i = 0; i < nodes.FESpace()->GetNDofs(); i++)
278  {
279  for (int d = 0; d < mesh.SpaceDimension(); d++)
280  node(d) = nodes(nodes.FESpace()->DofToVDof(i, d));
281 
282  node /= node.Norml2();
283 
284  for (int d = 0; d < mesh.SpaceDimension(); d++)
285  nodes(nodes.FESpace()->DofToVDof(i, d)) = node(d);
286  }
287 }
void SetTol(double tol)
Definition: hypre.cpp:1339
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:47
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:149
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:26
double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:128
Subclass constant coefficient.
Definition: coefficient.hpp:57
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:532
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:34
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.
void AddVertex(const double *)
Definition: mesh.cpp:760
double analytic_rhs(Vector &x)
Definition: ex7.cpp:219
The BoomerAMG solver in hypre.
Definition: hypre.hpp:519
Class for parallel linear form.
Definition: plinearform.hpp:26
void FinalizeTriMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Definition: mesh.cpp:892
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6225
void Assemble(int skip_zeros=1)
Assemble the local matrix.
double analytic_solution(Vector &x)
Definition: ex7.cpp:213
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:3066
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:246
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
int SpaceDimension() const
Definition: mesh.hpp:418
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:38
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
int main(int argc, char *argv[])
Definition: ex1.cpp:39
PCG solver in hypre.
Definition: hypre.hpp:381
Abstract finite element space.
Definition: fespace.hpp:61
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 AddQuad(const int *vi, int attr=1)
Definition: mesh.cpp:779
void AddTriangle(const int *vi, int attr=1)
Definition: mesh.cpp:774
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Definition: mesh.cpp:916
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
void SnapNodes(Mesh &mesh)
Definition: ex7.cpp:225
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:345
class for C-function coefficient
int DofToVDof(int dof, int vd) const
Definition: fespace.cpp:88
Vector data type.
Definition: vector.hpp:29
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 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