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