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