MFEM  v3.4
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 = new Mesh(mesh_file, 1, 1);
76  int dim = mesh->Dimension();
77 
78  // 3. Refine the mesh to increase the resolution. In this example we do
79  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
80  // largest number that gives a final mesh with no more than 50,000
81  // elements.
82  {
83  int ref_levels =
84  (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
85  for (int l = 0; l < ref_levels; l++)
86  {
87  mesh->UniformRefinement();
88  }
89  }
90 
91  // 4. Define a finite element space on the mesh. Here we use continuous
92  // Lagrange finite elements of the specified order. If order < 1, we
93  // instead use an isoparametric/isogeometric space.
95  if (order > 0)
96  {
97  fec = new H1_FECollection(order, dim);
98  }
99  else if (mesh->GetNodes())
100  {
101  fec = mesh->GetNodes()->OwnFEC();
102  cout << "Using isoparametric FEs: " << fec->Name() << endl;
103  }
104  else
105  {
106  fec = new H1_FECollection(order = 1, dim);
107  }
108  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
109  cout << "Number of finite element unknowns: "
110  << fespace->GetTrueVSize() << endl;
111 
112  // 5. Determine the list of true (i.e. conforming) essential boundary dofs.
113  // In this example, the boundary conditions are defined by marking all
114  // the boundary attributes from the mesh as essential (Dirichlet) and
115  // converting them to a list of true dofs.
116  Array<int> ess_tdof_list;
117  if (mesh->bdr_attributes.Size())
118  {
119  Array<int> ess_bdr(mesh->bdr_attributes.Max());
120  ess_bdr = 1;
121  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
122  }
123 
124  // 6. Set up the linear form b(.) which corresponds to the right-hand side of
125  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
126  // the basis functions in the finite element fespace.
127  LinearForm *b = new LinearForm(fespace);
128  ConstantCoefficient one(1.0);
130  b->Assemble();
131 
132  // 7. Define the solution vector x as a finite element grid function
133  // corresponding to fespace. Initialize x with initial guess of zero,
134  // which satisfies the boundary conditions.
135  GridFunction x(fespace);
136  x = 0.0;
137 
138  // 8. Set up the bilinear form a(.,.) on the finite element space
139  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
140  // domain integrator.
141  BilinearForm *a = new BilinearForm(fespace);
143 
144  // 9. Assemble the bilinear form and the corresponding linear system,
145  // applying any necessary transformations such as: eliminating boundary
146  // conditions, applying conforming constraints for non-conforming AMR,
147  // static condensation, etc.
148  if (static_cond) { a->EnableStaticCondensation(); }
149  a->Assemble();
150 
151  SparseMatrix A;
152  Vector B, X;
153  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
154 
155  cout << "Size of linear system: " << A.Height() << endl;
156 
157 #ifndef MFEM_USE_SUITESPARSE
158  // 10. Define a simple symmetric Gauss-Seidel preconditioner and use it to
159  // solve the system A X = B with PCG.
160  GSSmoother M(A);
161  PCG(A, M, B, X, 1, 200, 1e-12, 0.0);
162 #else
163  // 10. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
164  UMFPackSolver umf_solver;
165  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
166  umf_solver.SetOperator(A);
167  umf_solver.Mult(B, X);
168 #endif
169 
170  // 11. Recover the solution as a finite element grid function.
171  a->RecoverFEMSolution(X, *b, x);
172 
173  // 12. Save the refined mesh and the solution. This output can be viewed later
174  // using GLVis: "glvis -m refined.mesh -g sol.gf".
175  ofstream mesh_ofs("refined.mesh");
176  mesh_ofs.precision(8);
177  mesh->Print(mesh_ofs);
178  ofstream sol_ofs("sol.gf");
179  sol_ofs.precision(8);
180  x.Save(sol_ofs);
181 
182  // 13. Send the solution by socket to a GLVis server.
183  if (visualization)
184  {
185  char vishost[] = "localhost";
186  int visport = 19916;
187  socketstream sol_sock(vishost, visport);
188  sol_sock.precision(8);
189  sol_sock << "solution\n" << *mesh << x << flush;
190  }
191 
192  // 14. Free the used memory.
193  delete a;
194  delete b;
195  delete fespace;
196  if (order > 0) { delete fec; }
197  delete mesh;
198 
199  return 0;
200 }
int Size() const
Logical size of the array.
Definition: array.hpp:133
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1034
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:51
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Recover the solution of a linear system formed with FormLinearSystem().
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: fespace.cpp:365
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:618
int main(int argc, char *argv[])
Definition: ex1.cpp:45
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:344
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2349
int dim
Definition: ex3.cpp:47
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:36
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6741
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:109
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:457
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:256
int Dimension() const
Definition: mesh.hpp:645
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:174
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:355
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:66
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 const char * Name() const
Definition: fe_coll.hpp:50
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
Vector data type.
Definition: vector.hpp:48
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5422
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:79
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, SparseMatrix &A, Vector &X, Vector &B, int copy_interior=0)
Form a linear system, A X = B.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1785
void EnableStaticCondensation()
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:1688
bool Good() const
Definition: optparser.hpp:120