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