MFEM  v4.1.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 -pa -d ceed-cpu
35 // ex1 -pa -d ceed-cuda
36 // ex1 -m ../data/beam-hex.mesh -pa -d cuda
37 // ex1 -m ../data/beam-tet.mesh -pa -d ceed-cpu
38 // ex1 -m ../data/beam-tet.mesh -pa -d ceed-cuda:/gpu/cuda/ref
39 //
40 // Description: This example code demonstrates the use of MFEM to define a
41 // simple finite element discretization of the Laplace problem
42 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
43 // Specifically, we discretize using a FE space of the specified
44 // order, or if order < 1 using an isoparametric/isogeometric
45 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
46 // NURBS mesh, etc.)
47 //
48 // The example highlights the use of mesh refinement, finite
49 // element grid functions, as well as linear and bilinear forms
50 // corresponding to the left-hand side and right-hand side of the
51 // discrete linear system. We also cover the explicit elimination
52 // of essential boundary conditions, static condensation, and the
53 // optional connection to the GLVis tool for visualization.
54 
55 #include "mfem.hpp"
56 #include <fstream>
57 #include <iostream>
58 
59 using namespace std;
60 using namespace mfem;
61 
62 int main(int argc, char *argv[])
63 {
64  // 1. Parse command-line options.
65  const char *mesh_file = "../data/star.mesh";
66  int order = 1;
67  bool static_cond = false;
68  bool pa = false;
69  const char *device_config = "cpu";
70  bool visualization = true;
71 
72  OptionsParser args(argc, argv);
73  args.AddOption(&mesh_file, "-m", "--mesh",
74  "Mesh file to use.");
75  args.AddOption(&order, "-o", "--order",
76  "Finite element order (polynomial degree) or -1 for"
77  " isoparametric space.");
78  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
79  "--no-static-condensation", "Enable static condensation.");
80  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
81  "--no-partial-assembly", "Enable Partial Assembly.");
82  args.AddOption(&device_config, "-d", "--device",
83  "Device configuration string, see Device::Configure().");
84  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
85  "--no-visualization",
86  "Enable or disable GLVis visualization.");
87  args.Parse();
88  if (!args.Good())
89  {
90  args.PrintUsage(cout);
91  return 1;
92  }
93  args.PrintOptions(cout);
94 
95  // 2. Enable hardware devices such as GPUs, and programming models such as
96  // CUDA, OCCA, RAJA and OpenMP based on command line options.
97  Device device(device_config);
98  device.Print();
99 
100  // 3. Read the mesh from the given mesh file. We can handle triangular,
101  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
102  // the same code.
103  Mesh *mesh = new Mesh(mesh_file, 1, 1);
104  int dim = mesh->Dimension();
105 
106  // 4. Refine the mesh to increase the resolution. In this example we do
107  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
108  // largest number that gives a final mesh with no more than 50,000
109  // elements.
110  {
111  int ref_levels =
112  (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
113  for (int l = 0; l < ref_levels; l++)
114  {
115  mesh->UniformRefinement();
116  }
117  }
118 
119  // 5. Define a finite element space on the mesh. Here we use continuous
120  // Lagrange finite elements of the specified order. If order < 1, we
121  // instead use an isoparametric/isogeometric space.
123  if (order > 0)
124  {
125  fec = new H1_FECollection(order, dim);
126  }
127  else if (mesh->GetNodes())
128  {
129  fec = mesh->GetNodes()->OwnFEC();
130  cout << "Using isoparametric FEs: " << fec->Name() << endl;
131  }
132  else
133  {
134  fec = new H1_FECollection(order = 1, dim);
135  }
136  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
137  cout << "Number of finite element unknowns: "
138  << fespace->GetTrueVSize() << endl;
139 
140  // 6. Determine the list of true (i.e. conforming) essential boundary dofs.
141  // In this example, the boundary conditions are defined by marking all
142  // the boundary attributes from the mesh as essential (Dirichlet) and
143  // converting them to a list of true dofs.
144  Array<int> ess_tdof_list;
145  if (mesh->bdr_attributes.Size())
146  {
147  Array<int> ess_bdr(mesh->bdr_attributes.Max());
148  ess_bdr = 1;
149  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
150  }
151 
152  // 7. Set up the linear form b(.) which corresponds to the right-hand side of
153  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
154  // the basis functions in the finite element fespace.
155  LinearForm *b = new LinearForm(fespace);
156  ConstantCoefficient one(1.0);
158  b->Assemble();
159 
160  // 8. Define the solution vector x as a finite element grid function
161  // corresponding to fespace. Initialize x with initial guess of zero,
162  // which satisfies the boundary conditions.
163  GridFunction x(fespace);
164  x = 0.0;
165 
166  // 9. Set up the bilinear form a(.,.) on the finite element space
167  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
168  // domain integrator.
169  BilinearForm *a = new BilinearForm(fespace);
170  if (pa) { a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
172 
173  // 10. Assemble the bilinear form and the corresponding linear system,
174  // applying any necessary transformations such as: eliminating boundary
175  // conditions, applying conforming constraints for non-conforming AMR,
176  // static condensation, etc.
177  if (static_cond) { a->EnableStaticCondensation(); }
178  a->Assemble();
179 
180  OperatorPtr A;
181  Vector B, X;
182  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
183 
184  cout << "Size of linear system: " << A->Height() << endl;
185 
186  // 11. Solve the linear system A X = B.
187  if (!pa)
188  {
189 #ifndef MFEM_USE_SUITESPARSE
190  // Use a simple symmetric Gauss-Seidel preconditioner with PCG.
191  GSSmoother M((SparseMatrix&)(*A));
192  PCG(*A, M, B, X, 1, 200, 1e-12, 0.0);
193 #else
194  // If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
195  UMFPackSolver umf_solver;
196  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
197  umf_solver.SetOperator(*A);
198  umf_solver.Mult(B, X);
199 #endif
200  }
201  else // Jacobi preconditioning in partial assembly mode
202  {
203  if (UsesTensorBasis(*fespace))
204  {
205  OperatorJacobiSmoother M(*a, ess_tdof_list);
206  PCG(*A, M, B, X, 1, 400, 1e-12, 0.0);
207  }
208  else
209  {
210  CG(*A, B, X, 1, 400, 1e-12, 0.0);
211  }
212  }
213 
214  // 12. Recover the solution as a finite element grid function.
215  a->RecoverFEMSolution(X, *b, x);
216 
217  // 13. Save the refined mesh and the solution. This output can be viewed later
218  // using GLVis: "glvis -m refined.mesh -g sol.gf".
219  ofstream mesh_ofs("refined.mesh");
220  mesh_ofs.precision(8);
221  mesh->Print(mesh_ofs);
222  ofstream sol_ofs("sol.gf");
223  sol_ofs.precision(8);
224  x.Save(sol_ofs);
225 
226  // 14. Send the solution by socket to a GLVis server.
227  if (visualization)
228  {
229  char vishost[] = "localhost";
230  int visport = 19916;
231  socketstream sol_sock(vishost, visport);
232  sol_sock.precision(8);
233  sol_sock << "solution\n" << *mesh << x << flush;
234  }
235 
236  // 15. Free the used memory.
237  delete a;
238  delete b;
239  delete fespace;
240  if (order > 0) { delete fec; }
241  delete mesh;
242 
243  return 0;
244 }
int Size() const
Logical size of the array.
Definition: array.hpp:124
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1188
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:693
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:236
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::FULL.
bool UsesTensorBasis(const FiniteElementSpace &fes)
Definition: fespace.hpp:919
int main(int argc, char *argv[])
Definition: ex1.cpp:62
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:580
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2731
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:57
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:84
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7982
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:516
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:529
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:384
int Dimension() const
Definition: mesh.hpp:744
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:189
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:591
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
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
double a
Definition: lissajous.cpp:41
virtual const char * Name() const
Definition: fe_coll.hpp:53
int dim
Definition: ex24.cpp:43
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:6203
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:114
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:2382
void EnableStaticCondensation()
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:2285
bool Good() const
Definition: optparser.hpp:125