MFEM  v4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex6.cpp
Go to the documentation of this file.
1 // MFEM Example 6
2 //
3 // Compile with: make ex6
4 //
5 // Sample runs: ex6 -m ../data/square-disc.mesh -o 1
6 // ex6 -m ../data/square-disc.mesh -o 2
7 // ex6 -m ../data/square-disc-nurbs.mesh -o 2
8 // ex6 -m ../data/star.mesh -o 3
9 // ex6 -m ../data/escher.mesh -o 2
10 // ex6 -m ../data/fichera.mesh -o 2
11 // ex6 -m ../data/disc-nurbs.mesh -o 2
12 // ex6 -m ../data/ball-nurbs.mesh
13 // ex6 -m ../data/pipe-nurbs.mesh
14 // ex6 -m ../data/star-surf.mesh -o 2
15 // ex6 -m ../data/square-disc-surf.mesh -o 2
16 // ex6 -m ../data/amr-quad.mesh
17 //
18 // Device sample runs:
19 // ex6 -pa -d cuda
20 // ex6 -pa -d occa-cuda
21 // ex6 -pa -d raja-omp
22 //
23 // Description: This is a version of Example 1 with a simple adaptive mesh
24 // refinement loop. The problem being solved is again the Laplace
25 // equation -Delta u = 1 with homogeneous Dirichlet boundary
26 // conditions. The problem is solved on a sequence of meshes which
27 // are locally refined in a conforming (triangles, tetrahedrons)
28 // or non-conforming (quadrilaterals, hexahedra) manner according
29 // to a simple ZZ error estimator.
30 //
31 // The example demonstrates MFEM's capability to work with both
32 // conforming and nonconforming refinements, in 2D and 3D, on
33 // linear, curved and surface meshes. Interpolation of functions
34 // from coarse to fine meshes, as well as persistent GLVis
35 // visualization are also illustrated.
36 //
37 // We recommend viewing Example 1 before viewing this example.
38 
39 #include "mfem.hpp"
40 #include <fstream>
41 #include <iostream>
42 
43 using namespace std;
44 using namespace mfem;
45 
46 int main(int argc, char *argv[])
47 {
48  // 1. Parse command-line options.
49  const char *mesh_file = "../data/star.mesh";
50  int order = 1;
51  bool pa = false;
52  const char *device_config = "cpu";
53  bool visualization = true;
54 
55  OptionsParser args(argc, argv);
56  args.AddOption(&mesh_file, "-m", "--mesh",
57  "Mesh file to use.");
58  args.AddOption(&order, "-o", "--order",
59  "Finite element order (polynomial degree).");
60  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
61  "--no-partial-assembly", "Enable Partial Assembly.");
62  args.AddOption(&device_config, "-d", "--device",
63  "Device configuration string, see Device::Configure().");
64  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
65  "--no-visualization",
66  "Enable or disable GLVis visualization.");
67  args.Parse();
68  if (!args.Good())
69  {
70  args.PrintUsage(cout);
71  return 1;
72  }
73  args.PrintOptions(cout);
74 
75  // 2. Enable hardware devices such as GPUs, and programming models such as
76  // CUDA, OCCA, RAJA and OpenMP based on command line options.
77  Device device(device_config);
78  device.Print();
79 
80  // 3. Read the mesh from the given mesh file. We can handle triangular,
81  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
82  // the same code.
83  Mesh mesh(mesh_file, 1, 1);
84  int dim = mesh.Dimension();
85  int sdim = mesh.SpaceDimension();
86 
87  // 4. Since a NURBS mesh can currently only be refined uniformly, we need to
88  // convert it to a piecewise-polynomial curved mesh. First we refine the
89  // NURBS mesh a bit more and then project the curvature to quadratic Nodes.
90  if (mesh.NURBSext)
91  {
92  for (int i = 0; i < 2; i++)
93  {
94  mesh.UniformRefinement();
95  }
96  mesh.SetCurvature(2);
97  }
98 
99  // 5. Define a finite element space on the mesh. The polynomial order is
100  // one (linear) by default, but this can be changed on the command line.
101  H1_FECollection fec(order, dim);
102  FiniteElementSpace fespace(&mesh, &fec);
103 
104  // 6. As in Example 1, we set up bilinear and linear forms corresponding to
105  // the Laplace problem -\Delta u = 1. We don't assemble the discrete
106  // problem yet, this will be done in the main loop.
107  BilinearForm a(&fespace);
108  if (pa) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
109  LinearForm b(&fespace);
110 
111  ConstantCoefficient one(1.0);
112  ConstantCoefficient zero(0.0);
113 
115  a.AddDomainIntegrator(integ);
117 
118  // 7. The solution vector x and the associated finite element grid function
119  // will be maintained over the AMR iterations. We initialize it to zero.
120  GridFunction x(&fespace);
121  x = 0.0;
122 
123  // 8. All boundary attributes will be used for essential (Dirichlet) BC.
124  MFEM_VERIFY(mesh.bdr_attributes.Size() > 0,
125  "Boundary attributes required in the mesh.");
126  Array<int> ess_bdr(mesh.bdr_attributes.Max());
127  ess_bdr = 1;
128 
129  // 9. Connect to GLVis.
130  char vishost[] = "localhost";
131  int visport = 19916;
132  socketstream sol_sock;
133  if (visualization)
134  {
135  sol_sock.open(vishost, visport);
136  }
137 
138  // 10. Set up an error estimator. Here we use the Zienkiewicz-Zhu estimator
139  // that uses the ComputeElementFlux method of the DiffusionIntegrator to
140  // recover a smoothed flux (gradient) that is subtracted from the element
141  // flux to get an error indicator. We need to supply the space for the
142  // smoothed flux: an (H1)^sdim (i.e., vector-valued) space is used here.
143  FiniteElementSpace flux_fespace(&mesh, &fec, sdim);
144  ZienkiewiczZhuEstimator estimator(*integ, x, flux_fespace);
145  estimator.SetAnisotropic();
146 
147  // 11. A refiner selects and refines elements based on a refinement strategy.
148  // The strategy here is to refine elements with errors larger than a
149  // fraction of the maximum element error. Other strategies are possible.
150  // The refiner will call the given error estimator.
151  ThresholdRefiner refiner(estimator);
152  refiner.SetTotalErrorFraction(0.7);
153 
154  // 12. The main AMR loop. In each iteration we solve the problem on the
155  // current mesh, visualize the solution, and refine the mesh.
156  const int max_dofs = 50000;
157  for (int it = 0; ; it++)
158  {
159  int cdofs = fespace.GetTrueVSize();
160  cout << "\nAMR iteration " << it << endl;
161  cout << "Number of unknowns: " << cdofs << endl;
162 
163  // 13. Assemble the right-hand side.
164  b.Assemble();
165 
166  // 14. Set Dirichlet boundary values in the GridFunction x.
167  // Determine the list of Dirichlet true DOFs in the linear system.
168  Array<int> ess_tdof_list;
169  x.ProjectBdrCoefficient(zero, ess_bdr);
170  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
171 
172  // 15. Assemble the stiffness matrix.
173  a.Assemble();
174 
175  // 16. Create the linear system: eliminate boundary conditions, constrain
176  // hanging nodes and possibly apply other transformations. The system
177  // will be solved for true (unconstrained) DOFs only.
178  OperatorPtr A;
179  Vector B, X;
180 
181  const int copy_interior = 1;
182  a.FormLinearSystem(ess_tdof_list, x, b, A, X, B, copy_interior);
183 
184  // 17. Solve the linear system A X = B.
185  if (!pa)
186  {
187 #ifndef MFEM_USE_SUITESPARSE
188  // Use a simple symmetric Gauss-Seidel preconditioner with PCG.
189  GSSmoother M((SparseMatrix&)(*A));
190  PCG(*A, M, B, X, 3, 200, 1e-12, 0.0);
191 #else
192  // If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
193  UMFPackSolver umf_solver;
194  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
195  umf_solver.SetOperator(*A);
196  umf_solver.Mult(B, X);
197 #endif
198  }
199  else // No preconditioning for now in partial assembly mode.
200  {
201  CG(*A, B, X, 3, 2000, 1e-12, 0.0);
202  }
203 
204  // 18. After solving the linear system, reconstruct the solution as a
205  // finite element GridFunction. Constrained nodes are interpolated
206  // from true DOFs (it may therefore happen that x.Size() >= X.Size()).
207  a.RecoverFEMSolution(X, b, x);
208 
209  // 19. Send solution by socket to the GLVis server.
210  if (visualization && sol_sock.good())
211  {
212  sol_sock.precision(8);
213  sol_sock << "solution\n" << mesh << x << flush;
214  }
215 
216  if (cdofs > max_dofs)
217  {
218  cout << "Reached the maximum number of dofs. Stop." << endl;
219  break;
220  }
221 
222  // 20. Call the refiner to modify the mesh. The refiner calls the error
223  // estimator to obtain element errors, then it selects elements to be
224  // refined and finally it modifies the mesh. The Stop() method can be
225  // used to determine if a stopping criterion was met.
226  refiner.Apply(mesh);
227  if (refiner.Stop())
228  {
229  cout << "Stopping criterion satisfied. Stop." << endl;
230  break;
231  }
232 
233  // 21. Update the space to reflect the new state of the mesh. Also,
234  // interpolate the solution x so that it lies in the new space but
235  // represents the same function. This saves solver iterations later
236  // since we'll have a good initial guess of x in the next step.
237  // Internally, FiniteElementSpace::Update() calculates an
238  // interpolation matrix which is then used by GridFunction::Update().
239  fespace.Update();
240  x.Update();
241 
242  // 22. Inform also the bilinear and linear forms that the space has
243  // changed.
244  a.Update();
245  b.Update();
246  }
247 
248  return 0;
249 }
int Size() const
Logical size of the array.
Definition: array.hpp:118
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
void SetAnisotropic(bool aniso=true)
Enable/disable anisotropic estimates. To enable this option, the BilinearFormIntegrator must support ...
Definition: estimators.hpp:142
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
virtual void Update(bool want_transform=true)
Definition: fespace.cpp:1942
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
bool Stop() const
Check if STOP action is requested, e.g. stopping criterion is satisfied.
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:98
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::FULL.
int main(int argc, char *argv[])
Definition: ex1.cpp:58
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:344
bool Apply(Mesh &mesh)
Perform the mesh operation.
int dim
Definition: ex3.cpp:48
Data type sparse matrix.
Definition: sparsemat.hpp:40
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7591
Mesh refinement operator using an error threshold.
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3644
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:448
virtual void Update(FiniteElementSpace *nfes=NULL)
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:461
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:350
int Dimension() const
Definition: mesh.hpp:713
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
int SpaceDimension() const
Definition: mesh.hpp:714
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:23
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:179
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:355
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:85
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
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: gridfunc.cpp:152
void Update()
Update the object according to the associated FE space fes.
Definition: linearform.hpp:158
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:181
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
int open(const char hostname[], int port)
The ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure...
Definition: estimators.hpp:72
Vector data type.
Definition: vector.hpp:48
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
void SetTotalErrorFraction(double fraction)
Set the total fraction used in the computation of the threshold. The default value is 1/2...
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:96
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1794
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
Definition: gridfunc.hpp:275
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:1697
bool Good() const
Definition: optparser.hpp:122