MFEM  v3.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex7.cpp
Go to the documentation of this file.
1 // MFEM Example 7
2 //
3 // Compile with: make ex7
4 //
5 // Sample runs: ex7 -e 0 -o 2 -r 4
6 // ex7 -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. Parse command-line options.
36  int elem_type = 1;
37  int ref_levels = 2;
38  int order = 2;
39  bool always_snap = false;
40  bool visualization = 1;
41 
42  OptionsParser args(argc, argv);
43  args.AddOption(&elem_type, "-e", "--elem",
44  "Type of elements to use: 0 - triangles, 1 - quads.");
45  args.AddOption(&order, "-o", "--order",
46  "Finite element order (polynomial degree).");
47  args.AddOption(&ref_levels, "-r", "--refine",
48  "Number of times to refine the mesh uniformly.");
49  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
50  "--no-visualization",
51  "Enable or disable GLVis visualization.");
52  args.AddOption(&always_snap, "-snap", "--always-snap", "-no-snap",
53  "--snap-at-the-end",
54  "If true, snap nodes to the sphere initially and after each refinement "
55  "otherwise, snap only after the last refinement");
56  args.Parse();
57  if (!args.Good())
58  {
59  args.PrintUsage(cout);
60  return 1;
61  }
62  args.PrintOptions(cout);
63 
64  // 2. Generate an initial high-order (surface) mesh on the unit sphere. The
65  // Mesh object represents a 2D mesh in 3 spatial dimensions. We first add
66  // the elements and the vertices of the mesh, and then make it high-order
67  // by specifying a finite element space for its nodes.
68  int Nvert = 8, Nelem = 6;
69  if (elem_type == 0)
70  {
71  Nvert = 6;
72  Nelem = 8;
73  }
74  Mesh *mesh = new Mesh(2, Nvert, Nelem, 0, 3);
75 
76  if (elem_type == 0) // inscribed octahedron
77  {
78  const double tri_v[6][3] =
79  {{ 1, 0, 0}, { 0, 1, 0}, {-1, 0, 0},
80  { 0, -1, 0}, { 0, 0, 1}, { 0, 0, -1}};
81  const int tri_e[8][3] =
82  {{0, 1, 4}, {1, 2, 4}, {2, 3, 4}, {3, 0, 4},
83  {1, 0, 5}, {2, 1, 5}, {3, 2, 5}, {0, 3, 5}};
84 
85  for (int j = 0; j < Nvert; j++)
86  {
87  mesh->AddVertex(tri_v[j]);
88  }
89  for (int j = 0; j < Nelem; j++)
90  {
91  int attribute = j + 1;
92  mesh->AddTriangle(tri_e[j], attribute);
93  }
94  mesh->FinalizeTriMesh(1, 1, true);
95  }
96  else // inscribed cube
97  {
98  const double quad_v[8][3] =
99  {{-1, -1, -1}, {+1, -1, -1}, {+1, +1, -1}, {-1, +1, -1},
100  {-1, -1, +1}, {+1, -1, +1}, {+1, +1, +1}, {-1, +1, +1}};
101  const int quad_e[6][4] =
102  {{3, 2, 1, 0}, {0, 1, 5, 4}, {1, 2, 6, 5},
103  {2, 3, 7, 6}, {3, 0, 4, 7}, {4, 5, 6, 7}};
104 
105  for (int j = 0; j < Nvert; j++)
106  {
107  mesh->AddVertex(quad_v[j]);
108  }
109  for (int j = 0; j < Nelem; j++)
110  {
111  int attribute = j + 1;
112  mesh->AddQuad(quad_e[j], attribute);
113  }
114  mesh->FinalizeQuadMesh(1, 1, true);
115  }
116 
117  // Set the space for the high-order mesh nodes.
118  H1_FECollection fec(order, mesh->Dimension());
119  FiniteElementSpace nodal_fes(mesh, &fec, mesh->SpaceDimension());
120  mesh->SetNodalFESpace(&nodal_fes);
121 
122  // 3. Refine the mesh while snapping nodes to the sphere.
123  for (int l = 0; l <= ref_levels; l++)
124  {
125  if (l > 0) // for l == 0 just perform snapping
126  mesh->UniformRefinement();
127 
128  // Snap the nodes of the refined mesh back to sphere surface.
129  if (always_snap || l == ref_levels)
130  SnapNodes(*mesh);
131  }
132 
133  // 4. Define a finite element space on the mesh. Here we use isoparametric
134  // finite elements -- the same as the mesh nodes.
135  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, &fec);
136  cout << "Number of unknowns: " << fespace->GetVSize() << endl;
137 
138  // 5. Set up the linear form b(.) which corresponds to the right-hand side of
139  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
140  // the basis functions in the finite element fespace.
141  LinearForm *b = new LinearForm(fespace);
142  ConstantCoefficient one(1.0);
145  b->AddDomainIntegrator(new DomainLFIntegrator(rhs_coef));
146  b->Assemble();
147 
148  // 6. Define the solution vector x as a finite element grid function
149  // corresponding to fespace. Initialize x with initial guess of zero,
150  // which satisfies the boundary conditions.
151  GridFunction x(fespace);
152  x = 0.0;
153 
154  // 7. Set up the bilinear form a(.,.) on the finite element space
155  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
156  // domain integrator and imposing homogeneous Dirichlet boundary
157  // conditions. The boundary conditions are implemented by marking all the
158  // boundary attributes from the mesh as essential (Dirichlet). After
159  // assembly and finalizing we extract the corresponding sparse matrix A.
160  BilinearForm *a = new BilinearForm(fespace);
162  a->AddDomainIntegrator(new MassIntegrator(one));
163  a->Assemble();
164  a->Finalize();
165  const SparseMatrix &A = a->SpMat();
166 
167 #ifndef MFEM_USE_SUITESPARSE
168  // 8. Define a simple symmetric Gauss-Seidel preconditioner and use it to
169  // solve the system Ax=b with PCG.
170  GSSmoother M(A);
171  PCG(A, M, *b, x, 1, 200, 1e-12, 0.0);
172 #else
173  // 8. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
174  UMFPackSolver umf_solver;
175  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
176  umf_solver.SetOperator(A);
177  umf_solver.Mult(*b, x);
178 #endif
179 
180  // 9. Compute and print the L^2 norm of the error.
181  cout<<"\nL2 norm of error: " << x.ComputeL2Error(sol_coef) << endl;
182 
183  // 10. Save the refined mesh and the solution. This output can be viewed
184  // later using GLVis: "glvis -m sphere_refined.mesh -g sol.gf".
185  {
186  ofstream mesh_ofs("sphere_refined.mesh");
187  mesh_ofs.precision(8);
188  mesh->Print(mesh_ofs);
189  ofstream sol_ofs("sol.gf");
190  sol_ofs.precision(8);
191  x.Save(sol_ofs);
192  }
193 
194  // 11. 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.precision(8);
201  sol_sock << "solution\n" << *mesh << x << flush;
202  }
203 
204  // 12. Free the used memory.
205  delete a;
206  delete b;
207  delete fespace;
208  delete mesh;
209 
210  return 0;
211 }
212 
214 {
215  double l2 = x(0)*x(0) + x(1)*x(1) + x(2)*x(2);
216  return x(0)*x(1)/l2;
217 }
218 
220 {
221  double l2 = x(0)*x(0) + x(1)*x(1) + x(2)*x(2);
222  return 7*x(0)*x(1)/l2;
223 }
224 
225 void SnapNodes(Mesh &mesh)
226 {
227  GridFunction &nodes = *mesh.GetNodes();
228  Vector node(mesh.SpaceDimension());
229  for (int i = 0; i < nodes.FESpace()->GetNDofs(); i++)
230  {
231  for (int d = 0; d < mesh.SpaceDimension(); d++)
232  node(d) = nodes(nodes.FESpace()->DofToVDof(i, d));
233 
234  node /= node.Norml2();
235 
236  for (int d = 0; d < mesh.SpaceDimension(); d++)
237  nodes(nodes.FESpace()->DofToVDof(i, d)) = node(d);
238  }
239 }
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:47
int GetVSize() const
Definition: fespace.hpp:151
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
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
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
Data type for Gauss-Seidel smoother of sparse matrix.
void AddVertex(const double *)
Definition: mesh.cpp:760
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:329
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:1797
double analytic_rhs(Vector &x)
Definition: ex7.cpp:219
void FinalizeTriMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Definition: mesh.cpp:892
double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:166
Data type sparse matrix.
Definition: sparsemat.hpp:38
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6225
double analytic_solution(Vector &x)
Definition: ex7.cpp:213
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
Definition: solvers.cpp:424
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
int SpaceDimension() const
Definition: mesh.hpp:418
virtual void Print(std::ostream &out=std::cout) const
Print the mesh to the given stream using the default MFEM mesh format.
Definition: mesh.cpp:7136
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
int main(int argc, char *argv[])
Definition: ex1.cpp:39
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:340
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 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 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 linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
const SparseMatrix & SpMat() const
Returns a reference to the sparse martix.
virtual void Mult(const Vector &b, Vector &x) const
Operator application.
Definition: solvers.cpp:1597
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: solvers.cpp:1505
bool Good() const
Definition: optparser.hpp:120