MFEM  v3.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex1.cpp
Go to the documentation of this file.
1 // MFEM Example 1
2 //
3 // Compile with: make ex1
4 //
5 // Sample runs: ex1 -m ../data/square-disc.mesh
6 // ex1 -m ../data/star.mesh
7 // ex1 -m ../data/escher.mesh
8 // ex1 -m ../data/fichera.mesh
9 // ex1 -m ../data/square-disc-p2.vtk -o 2
10 // ex1 -m ../data/square-disc-p3.mesh -o 3
11 // ex1 -m ../data/square-disc-nurbs.mesh -o -1
12 // ex1 -m ../data/disc-nurbs.mesh -o -1
13 // ex1 -m ../data/pipe-nurbs.mesh -o -1
14 // ex1 -m ../data/star-surf.mesh
15 // ex1 -m ../data/square-disc-surf.mesh
16 //
17 // Description: This example code demonstrates the use of MFEM to define a
18 // simple finite element discretization of the Laplace problem
19 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
20 // Specifically, we discretize using a FE space of the specified
21 // order, or if order < 1 using an isoparametric/isogeometric
22 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
23 // NURBS mesh, etc.)
24 //
25 // The example highlights the use of mesh refinement, finite
26 // element grid functions, as well as linear and bilinear forms
27 // corresponding to the left-hand side and right-hand side of the
28 // discrete linear system. We also cover the explicit elimination
29 // of boundary conditions on all boundary edges, and the optional
30 // connection to the GLVis tool for visualization.
31 
32 #include "mfem.hpp"
33 #include <fstream>
34 #include <iostream>
35 
36 using namespace std;
37 using namespace mfem;
38 
39 int main(int argc, char *argv[])
40 {
41  // 1. Parse command-line options.
42  const char *mesh_file = "../data/star.mesh";
43  int order = 1;
44  bool visualization = 1;
45 
46  OptionsParser args(argc, argv);
47  args.AddOption(&mesh_file, "-m", "--mesh",
48  "Mesh file to use.");
49  args.AddOption(&order, "-o", "--order",
50  "Finite element order (polynomial degree) or -1 for"
51  " isoparametric space.");
52  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
53  "--no-visualization",
54  "Enable or disable GLVis visualization.");
55  args.Parse();
56  if (!args.Good())
57  {
58  args.PrintUsage(cout);
59  return 1;
60  }
61  args.PrintOptions(cout);
62 
63  // 2. Read the mesh from the given mesh file. We can handle triangular,
64  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
65  // the same code.
66  Mesh *mesh;
67  ifstream imesh(mesh_file);
68  if (!imesh)
69  {
70  cerr << "\nCan not open mesh file: " << mesh_file << '\n' << endl;
71  return 2;
72  }
73  mesh = new Mesh(imesh, 1, 1);
74  imesh.close();
75  int dim = mesh->Dimension();
76 
77  // 3. Refine the mesh to increase the resolution. In this example we do
78  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
79  // largest number that gives a final mesh with no more than 50,000
80  // elements.
81  {
82  int ref_levels =
83  (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
84  for (int l = 0; l < ref_levels; l++)
85  mesh->UniformRefinement();
86  }
87 
88  // 4. Define a finite element space on the mesh. Here we use continuous
89  // Lagrange finite elements of the specified order. If order < 1, we
90  // instead use an isoparametric/isogeometric space.
92  if (order > 0)
93  fec = new H1_FECollection(order, dim);
94  else if (mesh->GetNodes())
95  fec = mesh->GetNodes()->OwnFEC();
96  else
97  fec = new H1_FECollection(order = 1, dim);
98  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
99  cout << "Number of unknowns: " << fespace->GetVSize() << endl;
100 
101  // 5. Set up the linear form b(.) which corresponds to the right-hand side of
102  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
103  // the basis functions in the finite element fespace.
104  LinearForm *b = new LinearForm(fespace);
105  ConstantCoefficient one(1.0);
107  b->Assemble();
108 
109  // 6. Define the solution vector x as a finite element grid function
110  // corresponding to fespace. Initialize x with initial guess of zero,
111  // which satisfies the boundary conditions.
112  GridFunction x(fespace);
113  x = 0.0;
114 
115  // 7. Set up the bilinear form a(.,.) on the finite element space
116  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
117  // domain integrator and imposing homogeneous Dirichlet boundary
118  // conditions. The boundary conditions are implemented by marking all the
119  // boundary attributes from the mesh as essential (Dirichlet). After
120  // assembly and finalizing we extract the corresponding sparse matrix A.
121  BilinearForm *a = new BilinearForm(fespace);
123  a->Assemble();
124  Array<int> ess_bdr(mesh->bdr_attributes.Max());
125  ess_bdr = 1;
126  a->EliminateEssentialBC(ess_bdr, x, *b);
127  a->Finalize();
128  const SparseMatrix &A = a->SpMat();
129 
130 #ifndef MFEM_USE_SUITESPARSE
131  // 8. Define a simple symmetric Gauss-Seidel preconditioner and use it to
132  // solve the system Ax=b with PCG.
133  GSSmoother M(A);
134  PCG(A, M, *b, x, 1, 200, 1e-12, 0.0);
135 #else
136  // 8. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
137  UMFPackSolver umf_solver;
138  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
139  umf_solver.SetOperator(A);
140  umf_solver.Mult(*b, x);
141 #endif
142 
143  // 9. Save the refined mesh and the solution. This output can be viewed later
144  // using GLVis: "glvis -m refined.mesh -g sol.gf".
145  ofstream mesh_ofs("refined.mesh");
146  mesh_ofs.precision(8);
147  mesh->Print(mesh_ofs);
148  ofstream sol_ofs("sol.gf");
149  sol_ofs.precision(8);
150  x.Save(sol_ofs);
151 
152  // 10. Send the solution by socket to a GLVis server.
153  if (visualization)
154  {
155  char vishost[] = "localhost";
156  int visport = 19916;
157  socketstream sol_sock(vishost, visport);
158  sol_sock.precision(8);
159  sol_sock << "solution\n" << *mesh << x << flush;
160  }
161 
162  // 11. Free the used memory.
163  delete a;
164  delete b;
165  delete fespace;
166  if (order > 0)
167  delete fec;
168  delete mesh;
169 
170  return 0;
171 }
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:47
int GetVSize() const
Definition: fespace.hpp:151
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
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
Data type for Gauss-Seidel smoother of sparse matrix.
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
Data type sparse matrix.
Definition: sparsemat.hpp:38
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6225
T Max() const
Definition: array.cpp:78
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
int Dimension() const
Definition: mesh.hpp:417
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:385
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
Array< int > bdr_attributes
Definition: mesh.hpp:305
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
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:266
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.
void EliminateEssentialBC(Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs, int d=0)
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