MFEM  v4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex15.cpp
Go to the documentation of this file.
1 // MFEM Example 15
2 //
3 // Compile with: make ex15
4 //
5 // Sample runs: ex15
6 // ex15 -o 1 -y 0.4
7 // ex15 -o 4 -y 0.1
8 // ex15 -n 5
9 // ex15 -p 1 -n 3
10 //
11 // Other meshes:
12 //
13 // ex15 -m ../data/square-disc-nurbs.mesh
14 // ex15 -m ../data/disc-nurbs.mesh
15 // ex15 -m ../data/fichera.mesh -tf 0.3
16 // ex15 -m ../data/ball-nurbs.mesh -tf 0.3
17 // ex15 -m ../data/mobius-strip.mesh
18 // ex15 -m ../data/amr-quad.mesh
19 //
20 // Conforming meshes (no derefinement):
21 //
22 // ex15 -m ../data/square-disc.mesh
23 // ex15 -m ../data/escher.mesh -r 2 -tf 0.3
24 //
25 // Description: Building on Example 6, this example demonstrates dynamic AMR.
26 // The mesh is adapted to a time-dependent solution by refinement
27 // as well as by derefinement. For simplicity, the solution is
28 // prescribed and no time integration is done. However, the error
29 // estimation and refinement/derefinement decisions are realistic.
30 //
31 // At each outer iteration the right hand side function is changed
32 // to mimic a time dependent problem. Within each inner iteration
33 // the problem is solved on a sequence of meshes which are locally
34 // refined according to a simple ZZ error estimator. At the end
35 // of the inner iteration the error estimates are also used to
36 // identify any elements which may be over-refined and a single
37 // derefinement step is performed.
38 //
39 // The example demonstrates MFEM's capability to refine and
40 // derefine nonconforming meshes, in 2D and 3D, and on linear,
41 // curved and surface meshes. Interpolation of functions between
42 // coarse and fine meshes, persistent GLVis visualization, and
43 // saving of time-dependent fields for external visualization with
44 // VisIt (visit.llnl.gov) are also illustrated.
45 //
46 // We recommend viewing Examples 1, 6 and 9 before viewing this
47 // example.
48 
49 #include "mfem.hpp"
50 #include <fstream>
51 #include <iostream>
52 
53 using namespace std;
54 using namespace mfem;
55 
56 // Choices for the problem setup. Affect bdr_func and rhs_func.
57 int problem;
59 
60 // Prescribed time-dependent boundary and right-hand side functions.
61 double bdr_func(const Vector &pt, double t);
62 double rhs_func(const Vector &pt, double t);
63 
64 // Update the finite element space, interpolate the solution and perform
65 // parallel load balancing.
66 void UpdateProblem(Mesh &mesh, FiniteElementSpace &fespace,
68 
69 
70 int main(int argc, char *argv[])
71 {
72  // 1. Parse command-line options.
73  problem = 0;
74  nfeatures = 1;
75  const char *mesh_file = "../data/star-hilbert.mesh";
76  int order = 2;
77  double t_final = 1.0;
78  double max_elem_error = 5.0e-3;
79  double hysteresis = 0.15; // derefinement safety coefficient
80  int ref_levels = 0;
81  int nc_limit = 3; // maximum level of hanging nodes
82  bool visualization = true;
83  bool visit = false;
84 
85  OptionsParser args(argc, argv);
86  args.AddOption(&mesh_file, "-m", "--mesh",
87  "Mesh file to use.");
88  args.AddOption(&problem, "-p", "--problem",
89  "Problem setup to use: 0 = spherical front, 1 = ball.");
90  args.AddOption(&nfeatures, "-n", "--nfeatures",
91  "Number of solution features (fronts/balls).");
92  args.AddOption(&order, "-o", "--order",
93  "Finite element order (polynomial degree).");
94  args.AddOption(&max_elem_error, "-e", "--max-err",
95  "Maximum element error");
96  args.AddOption(&hysteresis, "-y", "--hysteresis",
97  "Derefinement safety coefficient.");
98  args.AddOption(&ref_levels, "-r", "--ref-levels",
99  "Number of initial uniform refinement levels.");
100  args.AddOption(&nc_limit, "-l", "--nc-limit",
101  "Maximum level of hanging nodes.");
102  args.AddOption(&t_final, "-tf", "--t-final",
103  "Final time; start time is 0.");
104  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
105  "--no-visualization",
106  "Enable or disable GLVis visualization.");
107  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
108  "--no-visit-datafiles",
109  "Save data files for VisIt (visit.llnl.gov) visualization.");
110  args.Parse();
111  if (!args.Good())
112  {
113  args.PrintUsage(cout);
114  return 1;
115  }
116  args.PrintOptions(cout);
117 
118  // 2. Read the mesh from the given mesh file on all processors. We can handle
119  // triangular, quadrilateral, tetrahedral, hexahedral, surface and volume
120  // meshes with the same code.
121  Mesh mesh(mesh_file, 1, 1);
122  int dim = mesh.Dimension();
123  int sdim = mesh.SpaceDimension();
124 
125  // 3. Project a NURBS mesh to a piecewise-quadratic curved mesh. Make sure
126  // that the mesh is non-conforming if it has quads or hexes and refine it.
127  if (mesh.NURBSext)
128  {
129  mesh.UniformRefinement();
130  if (ref_levels > 0) { ref_levels--; }
131  mesh.SetCurvature(2);
132  }
133  mesh.EnsureNCMesh();
134  for (int l = 0; l < ref_levels; l++)
135  {
136  mesh.UniformRefinement();
137  }
138  // Make sure tet-only meshes are marked for local refinement.
139  mesh.Finalize(true);
140 
141  // 4. All boundary attributes will be used for essential (Dirichlet) BC.
142  MFEM_VERIFY(mesh.bdr_attributes.Size() > 0,
143  "Boundary attributes required in the mesh.");
144  Array<int> ess_bdr(mesh.bdr_attributes.Max());
145  ess_bdr = 1;
146 
147  // 5. Define a finite element space on the mesh. The polynomial order is one
148  // (linear) by default, but this can be changed on the command line.
149  H1_FECollection fec(order, dim);
150  FiniteElementSpace fespace(&mesh, &fec);
151 
152  // 6. As in Example 1p, we set up bilinear and linear forms corresponding to
153  // the Laplace problem -\Delta u = 1. We don't assemble the discrete
154  // problem yet, this will be done in the inner loop.
155  BilinearForm a(&fespace);
156  LinearForm b(&fespace);
157 
158  ConstantCoefficient one(1.0);
161 
163  a.AddDomainIntegrator(integ);
165 
166  // 7. The solution vector x and the associated finite element grid function
167  // will be maintained over the AMR iterations.
168  GridFunction x(&fespace);
169 
170  // 8. Connect to GLVis. Prepare for VisIt output.
171  char vishost[] = "localhost";
172  int visport = 19916;
173 
174  socketstream sout;
175  if (visualization)
176  {
177  sout.open(vishost, visport);
178  if (!sout)
179  {
180  cout << "Unable to connect to GLVis server at "
181  << vishost << ':' << visport << endl;
182  cout << "GLVis visualization disabled.\n";
183  visualization = false;
184  }
185  sout.precision(8);
186  }
187 
188  VisItDataCollection visit_dc("Example15", &mesh);
189  visit_dc.RegisterField("solution", &x);
190  int vis_cycle = 0;
191 
192  // 9. As in Example 6, we set up a Zienkiewicz-Zhu estimator that will be
193  // used to obtain element error indicators. The integrator needs to
194  // provide the method ComputeElementFlux. The smoothed flux space is a
195  // vector valued H1 space here.
196  FiniteElementSpace flux_fespace(&mesh, &fec, sdim);
197  ZienkiewiczZhuEstimator estimator(*integ, x, flux_fespace);
198 
199  // 10. As in Example 6, we also need a refiner. This time the refinement
200  // strategy is based on a fixed threshold that is applied locally to each
201  // element. The global threshold is turned off by setting the total error
202  // fraction to zero. We also enforce a maximum refinement ratio between
203  // adjacent elements.
204  ThresholdRefiner refiner(estimator);
205  refiner.SetTotalErrorFraction(0.0); // use purely local threshold
206  refiner.SetLocalErrorGoal(max_elem_error);
207  refiner.PreferConformingRefinement();
208  refiner.SetNCLimit(nc_limit);
209 
210  // 11. A derefiner selects groups of elements that can be coarsened to form
211  // a larger element. A conservative enough threshold needs to be set to
212  // prevent derefining elements that would immediately be refined again.
213  ThresholdDerefiner derefiner(estimator);
214  derefiner.SetThreshold(hysteresis * max_elem_error);
215  derefiner.SetNCLimit(nc_limit);
216 
217  // 12. The outer time loop. In each iteration we update the right hand side,
218  // solve the problem on the current mesh, visualize the solution and
219  // refine the mesh as many times as necessary. Then we derefine any
220  // elements which have very small errors.
221  x = 0.0;
222  for (double time = 0.0; time < t_final + 1e-10; time += 0.01)
223  {
224  cout << "\nTime " << time << "\n\nRefinement:" << endl;
225 
226  // Set the current time in the coefficients.
227  bdr.SetTime(time);
228  rhs.SetTime(time);
229 
230  // Make sure errors will be recomputed in the following.
231  refiner.Reset();
232  derefiner.Reset();
233 
234  // 13. The inner refinement loop. At the end we want to have the current
235  // time step resolved to the prescribed tolerance in each element.
236  for (int ref_it = 1; ; ref_it++)
237  {
238  cout << "Iteration: " << ref_it << ", number of unknowns: "
239  << fespace.GetVSize() << endl;
240 
241  // 14. Recompute the field on the current mesh: assemble the stiffness
242  // matrix and the right-hand side.
243  a.Assemble();
244  b.Assemble();
245 
246  // 15. Project the exact solution to the essential boundary DOFs.
247  x.ProjectBdrCoefficient(bdr, ess_bdr);
248 
249  // 16. Create and solve the linear system.
250  Array<int> ess_tdof_list;
251  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
252 
253  SparseMatrix A;
254  Vector B, X;
255  a.FormLinearSystem(ess_tdof_list, x, b, A, X, B);
256 
257 #ifndef MFEM_USE_SUITESPARSE
258  GSSmoother M(A);
259  PCG(A, M, B, X, 0, 500, 1e-12, 0.0);
260 #else
261  UMFPackSolver umf_solver;
262  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
263  umf_solver.SetOperator(A);
264  umf_solver.Mult(B, X);
265 #endif
266 
267  // 17. Extract the local solution on each processor.
268  a.RecoverFEMSolution(X, b, x);
269 
270  // 18. Send the solution by socket to a GLVis server and optionally
271  // save it in VisIt format.
272  if (visualization)
273  {
274  sout.precision(8);
275  sout << "solution\n" << mesh << x << flush;
276  }
277  if (visit)
278  {
279  visit_dc.SetCycle(vis_cycle++);
280  visit_dc.SetTime(time);
281  visit_dc.Save();
282  }
283 
284  // 19. Apply the refiner on the mesh. The refiner calls the error
285  // estimator to obtain element errors, then it selects elements to
286  // be refined and finally it modifies the mesh. The Stop() method
287  // determines if all elements satisfy the local threshold.
288  refiner.Apply(mesh);
289  if (refiner.Stop())
290  {
291  break;
292  }
293 
294  // 20. Update the space and interpolate the solution.
295  UpdateProblem(mesh, fespace, x, a, b);
296  }
297 
298  // 21. Use error estimates from the last inner iteration to check for
299  // possible derefinements. The derefiner works similarly as the
300  // refiner. The errors are not recomputed because the mesh did not
301  // change (and also the estimator was not Reset() at this time).
302  if (derefiner.Apply(mesh))
303  {
304  cout << "\nDerefined elements." << endl;
305 
306  // 22. Update the space and interpolate the solution.
307  UpdateProblem(mesh, fespace, x, a, b);
308  }
309 
310  a.Update();
311  b.Update();
312  }
313 
314  return 0;
315 }
316 
317 
318 void UpdateProblem(Mesh &mesh, FiniteElementSpace &fespace,
320 {
321  // Update the space: recalculate the number of DOFs and construct a matrix
322  // that will adjust any GridFunctions to the new mesh state.
323  fespace.Update();
324 
325  // Interpolate the solution on the new mesh by applying the transformation
326  // matrix computed in the finite element space. Multiple GridFunctions could
327  // be updated here.
328  x.Update();
329 
330  // Free any transformation matrices to save memory.
331  fespace.UpdatesFinished();
332 
333  // Inform the linear and bilinear forms that the space has changed.
334  a.Update();
335  b.Update();
336 }
337 
338 
339 const double alpha = 0.02;
340 
341 // Spherical front with a Gaussian cross section and radius t
342 double front(double x, double y, double z, double t, int)
343 {
344  double r = sqrt(x*x + y*y + z*z);
345  return exp(-0.5*pow((r - t)/alpha, 2));
346 }
347 
348 double front_laplace(double x, double y, double z, double t, int dim)
349 {
350  double x2 = x*x, y2 = y*y, z2 = z*z, t2 = t*t;
351  double r = sqrt(x2 + y2 + z2);
352  double a2 = alpha*alpha, a4 = a2*a2;
353  return -exp(-0.5*pow((r - t)/alpha, 2)) / a4 *
354  (-2*t*(x2 + y2 + z2 - (dim-1)*a2/2)/r + x2 + y2 + z2 + t2 - dim*a2);
355 }
356 
357 // Smooth spherical step function with radius t
358 double ball(double x, double y, double z, double t, int)
359 {
360  double r = sqrt(x*x + y*y + z*z);
361  return -atan(2*(r - t)/alpha);
362 }
363 
364 double ball_laplace(double x, double y, double z, double t, int dim)
365 {
366  double x2 = x*x, y2 = y*y, z2 = z*z, t2 = 4*t*t;
367  double r = sqrt(x2 + y2 + z2);
368  double a2 = alpha*alpha;
369  double den = pow(-a2 - 4*(x2 + y2 + z2 - 2*r*t) - t2, 2.0);
370  return (dim == 2) ? 2*alpha*(a2 + t2 - 4*x2 - 4*y2)/r/den
371  /* */ : 4*alpha*(a2 + t2 - 4*r*t)/r/den;
372 }
373 
374 // Composes several features into one function
375 template<typename F0, typename F1>
376 double composite_func(const Vector &pt, double t, F0 f0, F1 f1)
377 {
378  int dim = pt.Size();
379  double x = pt(0), y = pt(1), z = 0.0;
380  if (dim == 3) { z = pt(2); }
381 
382  if (problem == 0)
383  {
384  if (nfeatures <= 1)
385  {
386  return f0(x, y, z, t, dim);
387  }
388  else
389  {
390  double sum = 0.0;
391  for (int i = 0; i < nfeatures; i++)
392  {
393  double x0 = 0.5*cos(2*M_PI * i / nfeatures);
394  double y0 = 0.5*sin(2*M_PI * i / nfeatures);
395  sum += f0(x - x0, y - y0, z, t, dim);
396  }
397  return sum;
398  }
399  }
400  else
401  {
402  double sum = 0.0;
403  for (int i = 0; i < nfeatures; i++)
404  {
405  double x0 = 0.5*cos(2*M_PI * i / nfeatures + M_PI*t);
406  double y0 = 0.5*sin(2*M_PI * i / nfeatures + M_PI*t);
407  sum += f1(x - x0, y - y0, z, 0.25, dim);
408  }
409  return sum;
410  }
411 }
412 
413 // Exact solution, used for the Dirichlet BC.
414 double bdr_func(const Vector &pt, double t)
415 {
416  return composite_func(pt, t, front, ball);
417 }
418 
419 // Laplace of the exact solution, used for the right hand side.
420 double rhs_func(const Vector &pt, double t)
421 {
422  return composite_func(pt, t, front_laplace, ball_laplace);
423 }
424 
double front_laplace(double x, double y, double z, double t, int dim)
Definition: ex15.cpp:348
int Size() const
Logical size of the array.
Definition: array.hpp:118
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:347
double bdr_func(const Vector &pt, double t)
Definition: ex15.cpp:414
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 SetCycle(int c)
Set time cycle (for time-dependent simulations)
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
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().
void SetNCLimit(int nc_limit)
Set the maximum ratio of refinement levels of adjacent elements (0 = unlimited).
int Size() const
Returns the size of the vector.
Definition: vector.hpp:150
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: fespace.cpp:366
virtual void Save()
Save the collection and a VisIt root file.
bool Stop() const
Check if STOP action is requested, e.g. stopping criterion is satisfied.
int main(int argc, char *argv[])
Definition: ex1.cpp:58
double rhs_func(const Vector &pt, double t)
Definition: ex15.cpp:420
void UpdateProblem(Mesh &mesh, FiniteElementSpace &fespace, GridFunction &x, BilinearForm &a, LinearForm &b)
Definition: ex15.cpp:318
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
virtual void Reset()
Reset the associated estimator.
Data type sparse matrix.
Definition: sparsemat.hpp:40
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7591
virtual void Reset()
Reset the associated estimator.
Mesh refinement operator using an error threshold.
Data collection with VisIt I/O routines.
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
int nfeatures
Definition: ex15.cpp:58
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
int Dimension() const
Definition: mesh.hpp:713
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetTime(double t)
Set physical time (for time-dependent simulations)
void SetTime(double t)
Definition: coefficient.hpp:39
void SetNCLimit(int nc_limit)
Set the maximum ratio of refinement levels of adjacent elements (0 = unlimited).
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 ball_laplace(double x, double y, double z, double t, int dim)
Definition: ex15.cpp:364
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
int problem
Definition: ex15.cpp:57
void PreferConformingRefinement()
Use conforming refinement, if possible (triangles, tetrahedra) – this is the default.
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
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:181
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: mesh.cpp:2098
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
void EnsureNCMesh(bool triangles_nonconforming=false)
Definition: mesh.cpp:7253
double ball(double x, double y, double z, double t, int)
Definition: ex15.cpp:358
int open(const char hostname[], int port)
void SetThreshold(double thresh)
Set the de-refinement threshold. The default value is zero.
The ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure...
Definition: estimators.hpp:72
void SetLocalErrorGoal(double err_goal)
Set the local stopping criterion: stop when local_err_i &lt;= local_err_goal. The default value is zero...
const double alpha
Definition: ex15.cpp:339
class for C-function coefficient
Vector data type.
Definition: vector.hpp:48
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
Definition: fespace.hpp:611
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...
double front(double x, double y, double z, double t, int)
Definition: ex15.cpp:342
double composite_func(const Vector &pt, double t, F0 f0, F1 f1)
Definition: ex15.cpp:376
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
De-refinement operator using an error threshold.
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