MFEM  v4.0
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/star-mixed.mesh
8 // ex1 -m ../data/escher.mesh
9 // ex1 -m ../data/fichera.mesh
10 // ex1 -m ../data/fichera-mixed.mesh
11 // ex1 -m ../data/toroid-wedge.mesh
12 // ex1 -m ../data/square-disc-p2.vtk -o 2
13 // ex1 -m ../data/square-disc-p3.mesh -o 3
14 // ex1 -m ../data/square-disc-nurbs.mesh -o -1
15 // ex1 -m ../data/star-mixed-p2.mesh -o 2
16 // ex1 -m ../data/disc-nurbs.mesh -o -1
17 // ex1 -m ../data/pipe-nurbs.mesh -o -1
18 // ex1 -m ../data/fichera-mixed-p2.mesh -o 2
19 // ex1 -m ../data/star-surf.mesh
20 // ex1 -m ../data/square-disc-surf.mesh
21 // ex1 -m ../data/inline-segment.mesh
22 // ex1 -m ../data/amr-quad.mesh
23 // ex1 -m ../data/amr-hex.mesh
24 // ex1 -m ../data/fichera-amr.mesh
25 // ex1 -m ../data/mobius-strip.mesh
26 // ex1 -m ../data/mobius-strip.mesh -o -1 -sc
27 //
28 // Device sample runs:
29 // ex1 -pa -d cuda
30 // ex1 -pa -d raja-cuda
31 // ex1 -pa -d occa-cuda
32 // ex1 -pa -d raja-omp
33 // ex1 -pa -d occa-omp
34 // ex1 -m ../data/beam-hex.mesh -pa -d cuda
35 //
36 // Description: This example code demonstrates the use of MFEM to define a
37 // simple finite element discretization of the Laplace problem
38 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
39 // Specifically, we discretize using a FE space of the specified
40 // order, or if order < 1 using an isoparametric/isogeometric
41 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
42 // NURBS mesh, etc.)
43 //
44 // The example highlights the use of mesh refinement, finite
45 // element grid functions, as well as linear and bilinear forms
46 // corresponding to the left-hand side and right-hand side of the
47 // discrete linear system. We also cover the explicit elimination
48 // of essential boundary conditions, static condensation, and the
49 // optional connection to the GLVis tool for visualization.
50 
51 #include "mfem.hpp"
52 #include <fstream>
53 #include <iostream>
54 
55 using namespace std;
56 using namespace mfem;
57 
58 int main(int argc, char *argv[])
59 {
60  // 1. Parse command-line options.
61  const char *mesh_file = "../data/star.mesh";
62  int order = 1;
63  bool static_cond = false;
64  bool pa = false;
65  const char *device_config = "cpu";
66  bool visualization = true;
67 
68  OptionsParser args(argc, argv);
69  args.AddOption(&mesh_file, "-m", "--mesh",
70  "Mesh file to use.");
71  args.AddOption(&order, "-o", "--order",
72  "Finite element order (polynomial degree) or -1 for"
73  " isoparametric space.");
74  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
75  "--no-static-condensation", "Enable static condensation.");
76  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
77  "--no-partial-assembly", "Enable Partial Assembly.");
78  args.AddOption(&device_config, "-d", "--device",
79  "Device configuration string, see Device::Configure().");
80  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
81  "--no-visualization",
82  "Enable or disable GLVis visualization.");
83  args.Parse();
84  if (!args.Good())
85  {
86  args.PrintUsage(cout);
87  return 1;
88  }
89  args.PrintOptions(cout);
90 
91  // 2. Enable hardware devices such as GPUs, and programming models such as
92  // CUDA, OCCA, RAJA and OpenMP based on command line options.
93  Device device(device_config);
94  device.Print();
95 
96  // 3. Read the mesh from the given mesh file. We can handle triangular,
97  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
98  // the same code.
99  Mesh *mesh = new Mesh(mesh_file, 1, 1);
100  int dim = mesh->Dimension();
101 
102  // 4. Refine the mesh to increase the resolution. In this example we do
103  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
104  // largest number that gives a final mesh with no more than 50,000
105  // elements.
106  {
107  int ref_levels =
108  (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
109  for (int l = 0; l < ref_levels; l++)
110  {
111  mesh->UniformRefinement();
112  }
113  }
114 
115  // 5. Define a finite element space on the mesh. Here we use continuous
116  // Lagrange finite elements of the specified order. If order < 1, we
117  // instead use an isoparametric/isogeometric space.
119  if (order > 0)
120  {
121  fec = new H1_FECollection(order, dim);
122  }
123  else if (mesh->GetNodes())
124  {
125  fec = mesh->GetNodes()->OwnFEC();
126  cout << "Using isoparametric FEs: " << fec->Name() << endl;
127  }
128  else
129  {
130  fec = new H1_FECollection(order = 1, dim);
131  }
132  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
133  cout << "Number of finite element unknowns: "
134  << fespace->GetTrueVSize() << endl;
135 
136  // 6. Determine the list of true (i.e. conforming) essential boundary dofs.
137  // In this example, the boundary conditions are defined by marking all
138  // the boundary attributes from the mesh as essential (Dirichlet) and
139  // converting them to a list of true dofs.
140  Array<int> ess_tdof_list;
141  if (mesh->bdr_attributes.Size())
142  {
143  Array<int> ess_bdr(mesh->bdr_attributes.Max());
144  ess_bdr = 1;
145  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
146  }
147 
148  // 7. Set up the linear form b(.) which corresponds to the right-hand side of
149  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
150  // the basis functions in the finite element fespace.
151  LinearForm *b = new LinearForm(fespace);
152  ConstantCoefficient one(1.0);
154  b->Assemble();
155 
156  // 8. Define the solution vector x as a finite element grid function
157  // corresponding to fespace. Initialize x with initial guess of zero,
158  // which satisfies the boundary conditions.
159  GridFunction x(fespace);
160  x = 0.0;
161 
162  // 9. Set up the bilinear form a(.,.) on the finite element space
163  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
164  // domain integrator.
165  BilinearForm *a = new BilinearForm(fespace);
166  if (pa) { a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
168 
169  // 10. Assemble the bilinear form and the corresponding linear system,
170  // applying any necessary transformations such as: eliminating boundary
171  // conditions, applying conforming constraints for non-conforming AMR,
172  // static condensation, etc.
173  if (static_cond) { a->EnableStaticCondensation(); }
174  a->Assemble();
175 
176  OperatorPtr A;
177  Vector B, X;
178  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
179 
180  cout << "Size of linear system: " << A->Height() << endl;
181 
182  // 11. Solve the linear system A X = B.
183  if (!pa)
184  {
185 #ifndef MFEM_USE_SUITESPARSE
186  // Use a simple symmetric Gauss-Seidel preconditioner with PCG.
187  GSSmoother M((SparseMatrix&)(*A));
188  PCG(*A, M, B, X, 1, 200, 1e-12, 0.0);
189 #else
190  // If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
191  UMFPackSolver umf_solver;
192  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
193  umf_solver.SetOperator(*A);
194  umf_solver.Mult(B, X);
195 #endif
196  }
197  else // No preconditioning for now in partial assembly mode.
198  {
199  CG(*A, B, X, 1, 2000, 1e-12, 0.0);
200  }
201 
202  // 12. Recover the solution as a finite element grid function.
203  a->RecoverFEMSolution(X, *b, x);
204 
205  // 13. Save the refined mesh and the solution. This output can be viewed later
206  // using GLVis: "glvis -m refined.mesh -g sol.gf".
207  ofstream mesh_ofs("refined.mesh");
208  mesh_ofs.precision(8);
209  mesh->Print(mesh_ofs);
210  ofstream sol_ofs("sol.gf");
211  sol_ofs.precision(8);
212  x.Save(sol_ofs);
213 
214  // 14. Send the solution by socket to a GLVis server.
215  if (visualization)
216  {
217  char vishost[] = "localhost";
218  int visport = 19916;
219  socketstream sol_sock(vishost, visport);
220  sol_sock.precision(8);
221  sol_sock << "solution\n" << *mesh << x << flush;
222  }
223 
224  // 15. Free the used memory.
225  delete a;
226  delete b;
227  delete fespace;
228  if (order > 0) { delete fec; }
229  delete mesh;
230 
231  return 0;
232 }
int Size() const
Logical size of the array.
Definition: array.hpp:118
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1156
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:67
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
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:366
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:676
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:98
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::FULL.
int main(int argc, char *argv[])
Definition: ex1.cpp:58
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:2627
int dim
Definition: ex3.cpp:48
Data type sparse matrix.
Definition: sparsemat.hpp:40
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:7591
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
void CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Conjugate gradient method. (tolerances are squared)
Definition: solvers.cpp:448
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:461
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:350
int Dimension() const
Definition: mesh.hpp:713
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:179
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:85
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:76
virtual const char * Name() const
Definition: fe_coll.hpp:53
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
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:5837
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:83
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:96
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1794
void EnableStaticCondensation()
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:1697
bool Good() const
Definition: optparser.hpp:122