MFEM  v3.1
Finite element discretization library
 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 // ex1 -m ../data/inline-segment.mesh
17 // ex1 -m ../data/amr-quad.mesh
18 // ex1 -m ../data/amr-hex.mesh
19 // ex1 -m ../data/fichera-amr.mesh
20 // ex1 -m ../data/mobius-strip.mesh
21 // ex1 -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. Parse command-line options.
48  const char *mesh_file = "../data/star.mesh";
49  int order = 1;
50  bool static_cond = false;
51  bool visualization = 1;
52 
53  OptionsParser args(argc, argv);
54  args.AddOption(&mesh_file, "-m", "--mesh",
55  "Mesh file to use.");
56  args.AddOption(&order, "-o", "--order",
57  "Finite element order (polynomial degree) or -1 for"
58  " isoparametric space.");
59  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
60  "--no-static-condensation", "Enable static condensation.");
61  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
62  "--no-visualization",
63  "Enable or disable GLVis visualization.");
64  args.Parse();
65  if (!args.Good())
66  {
67  args.PrintUsage(cout);
68  return 1;
69  }
70  args.PrintOptions(cout);
71 
72  // 2. Read the mesh from the given mesh file. We can handle triangular,
73  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
74  // the same code.
75  Mesh *mesh;
76  ifstream imesh(mesh_file);
77  if (!imesh)
78  {
79  cerr << "\nCan not open mesh file: " << mesh_file << '\n' << endl;
80  return 2;
81  }
82  mesh = new Mesh(imesh, 1, 1);
83  imesh.close();
84  int dim = mesh->Dimension();
85 
86  // 3. Refine the mesh to increase the resolution. In this example we do
87  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
88  // largest number that gives a final mesh with no more than 50,000
89  // elements.
90  {
91  int ref_levels =
92  (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
93  for (int l = 0; l < ref_levels; l++)
94  {
95  mesh->UniformRefinement();
96  }
97  }
98 
99  // 4. Define a finite element space on the mesh. Here we use continuous
100  // Lagrange finite elements of the specified order. If order < 1, we
101  // instead use an isoparametric/isogeometric space.
103  if (order > 0)
104  {
105  fec = new H1_FECollection(order, dim);
106  }
107  else if (mesh->GetNodes())
108  {
109  fec = mesh->GetNodes()->OwnFEC();
110  cout << "Using isoparametric FEs: " << fec->Name() << endl;
111  }
112  else
113  {
114  fec = new H1_FECollection(order = 1, dim);
115  }
116  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
117  cout << "Number of finite element unknowns: "
118  << fespace->GetTrueVSize() << endl;
119 
120  // 5. Determine the list of true (i.e. conforming) essential boundary dofs.
121  // In this example, the boundary conditions are defined by marking all
122  // the boundary attributes from the mesh as essential (Dirichlet) and
123  // converting them to a list of true dofs.
124  Array<int> ess_tdof_list;
125  if (mesh->bdr_attributes.Size())
126  {
127  Array<int> ess_bdr(mesh->bdr_attributes.Max());
128  ess_bdr = 1;
129  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
130  }
131 
132  // 6. Set up the linear form b(.) which corresponds to the right-hand side of
133  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
134  // the basis functions in the finite element fespace.
135  LinearForm *b = new LinearForm(fespace);
136  ConstantCoefficient one(1.0);
138  b->Assemble();
139 
140  // 7. Define the solution vector x as a finite element grid function
141  // corresponding to fespace. Initialize x with initial guess of zero,
142  // which satisfies the boundary conditions.
143  GridFunction x(fespace);
144  x = 0.0;
145 
146  // 8. Set up the bilinear form a(.,.) on the finite element space
147  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
148  // domain integrator.
149  BilinearForm *a = new BilinearForm(fespace);
151 
152  // 9. Assemble the bilinear form and the corresponding linear system,
153  // applying any necessary transformations such as: eliminating boundary
154  // conditions, applying conforming constraints for non-conforming AMR,
155  // static condensation, etc.
156  if (static_cond) { a->EnableStaticCondensation(); }
157  a->Assemble();
158 
159  SparseMatrix A;
160  Vector B, X;
161  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
162 
163  cout << "Size of linear system: " << A.Height() << endl;
164 
165 #ifndef MFEM_USE_SUITESPARSE
166  // 10. Define a simple symmetric Gauss-Seidel preconditioner and use it to
167  // solve the system A X = B with PCG.
168  GSSmoother M(A);
169  PCG(A, M, B, X, 1, 200, 1e-12, 0.0);
170 #else
171  // 10. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
172  UMFPackSolver umf_solver;
173  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
174  umf_solver.SetOperator(A);
175  umf_solver.Mult(B, X);
176 #endif
177 
178  // 11. Recover the solution as a finite element grid function.
179  a->RecoverFEMSolution(X, *b, x);
180 
181  // 12. Save the refined mesh and the solution. This output can be viewed later
182  // using GLVis: "glvis -m refined.mesh -g sol.gf".
183  ofstream mesh_ofs("refined.mesh");
184  mesh_ofs.precision(8);
185  mesh->Print(mesh_ofs);
186  ofstream sol_ofs("sol.gf");
187  sol_ofs.precision(8);
188  x.Save(sol_ofs);
189 
190  // 13. Send the solution by socket to a GLVis server.
191  if (visualization)
192  {
193  char vishost[] = "localhost";
194  int visport = 19916;
195  socketstream sol_sock(vishost, visport);
196  sol_sock.precision(8);
197  sol_sock << "solution\n" << *mesh << x << flush;
198  }
199 
200  // 14. Free the used memory.
201  delete a;
202  delete b;
203  delete fespace;
204  if (order > 0) { delete fec; }
205  delete mesh;
206 
207  return 0;
208 }
int Size() const
Logical size of the array.
Definition: array.hpp:109
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:47
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
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
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list)
Definition: fespace.cpp:469
void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:454
void FormLinearSystem(Array< int > &ess_tdof_list, Vector &x, Vector &b, SparseMatrix &A, Vector &X, Vector &B, int copy_interior=0)
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:332
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2138
int dim
Definition: ex3.cpp:48
Data type sparse matrix.
Definition: sparsemat.hpp:38
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows.
Definition: operator.hpp:35
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7316
T Max() const
Definition: array.cpp:90
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:440
int Dimension() const
Definition: mesh.hpp:475
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
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:8332
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Array< int > bdr_attributes
Definition: mesh.hpp:140
int main(int argc, char *argv[])
Definition: ex1.cpp:45
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:343
Abstract finite element space.
Definition: fespace.hpp:62
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 int GetTrueVSize()
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:167
virtual const char * Name() const
Definition: fe_coll.hpp:36
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
Vector data type.
Definition: vector.hpp:33
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5844
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:54
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
virtual void Mult(const Vector &b, Vector &x) const
Operator application.
Definition: solvers.cpp:1725
void EnableStaticCondensation()
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: solvers.cpp:1625
bool Good() const
Definition: optparser.hpp:120