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