MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
estimators.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_ERROR_ESTIMATORS
13 #define MFEM_ERROR_ESTIMATORS
14 
15 #include "../config/config.hpp"
16 #include "../linalg/vector.hpp"
17 #include "bilinearform.hpp"
18 #ifdef MFEM_USE_MPI
19 #include "pgridfunc.hpp"
20 #endif
21 
22 namespace mfem
23 {
24 
25 /** @brief Base class for all error estimators.
26  */
28 {
29 public:
31 };
32 
33 
34 /** @brief Base class for all element based error estimators.
35 
36  At a minimum, an ErrorEstimator must be able compute one non-negative real
37  (double) number for each element in the Mesh.
38  */
40 {
41 public:
42  /// Get a Vector with all element errors.
43  virtual const Vector &GetLocalErrors() = 0;
44 
45  /// Force recomputation of the estimates on the next call to GetLocalErrors.
46  virtual void Reset() = 0;
47 
48  /// Destruct the error estimator
49  virtual ~ErrorEstimator() { }
50 };
51 
52 
53 /** @brief The AnisotropicErrorEstimator class is the base class for all error
54  estimators that compute one non-negative real (double) number and an
55  anisotropic flag for every element in the Mesh.
56  */
58 {
59 public:
60  /** @brief Get an Array<int> with anisotropic flags for all mesh elements.
61  @return An empty array when anisotropic estimates are not available or
62  enabled. */
63  virtual const Array<int> &GetAnisotropicFlags() = 0;
64 };
65 
66 
67 /** @brief The ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu
68  error estimation procedure.
69 
70  Zienkiewicz, O.C. and Zhu, J.Z., The superconvergent patch recovery
71  and a posteriori error estimates. Part 1: The recovery technique.
72  Int. J. Num. Meth. Engng. 33, 1331-1364 (1992).
73 
74  Zienkiewicz, O.C. and Zhu, J.Z., The superconvergent patch recovery
75  and a posteriori error estimates. Part 2: Error estimates and adaptivity.
76  Int. J. Num. Meth. Engng. 33, 1365-1382 (1992).
77 
78  The required BilinearFormIntegrator must implement the methods
79  ComputeElementFlux() and ComputeFluxEnergy().
80  */
82 {
83 protected:
86  double total_error;
89  int flux_averaging; // see SetFluxAveraging()
90 
91  BilinearFormIntegrator *integ; ///< Not owned.
92  GridFunction *solution; ///< Not owned.
93 
94  FiniteElementSpace *flux_space; /**< @brief Ownership based on own_flux_fes.
95  Its Update() method is called automatically by this class when needed. */
96  bool with_coeff;
97  bool own_flux_fes; ///< Ownership flag for flux_space.
98 
99  /// Check if the mesh of the solution was modified.
101  {
102  long mesh_sequence = solution->FESpace()->GetMesh()->GetSequence();
103  MFEM_ASSERT(mesh_sequence >= current_sequence, "");
104  return (mesh_sequence > current_sequence);
105  }
106 
107  /// Compute the element error estimates.
108  void ComputeEstimates();
109 
110 public:
111  /** @brief Construct a new ZienkiewiczZhuEstimator object.
112  @param integ This BilinearFormIntegrator must implement the methods
113  ComputeElementFlux() and ComputeFluxEnergy().
114  @param sol The solution field whose error is to be estimated.
115  @param flux_fes The ZienkiewiczZhuEstimator assumes ownership of this
116  FiniteElementSpace and will call its Update() method when
117  needed.*/
119  FiniteElementSpace *flux_fes)
120  : current_sequence(-1),
121  total_error(),
122  anisotropic(false),
123  flux_averaging(0),
124  integ(&integ),
125  solution(&sol),
126  flux_space(flux_fes),
127  with_coeff(false),
128  own_flux_fes(true)
129  { }
130 
131  /** @brief Construct a new ZienkiewiczZhuEstimator object.
132  @param integ This BilinearFormIntegrator must implement the methods
133  ComputeElementFlux() and ComputeFluxEnergy().
134  @param sol The solution field whose error is to be estimated.
135  @param flux_fes The ZienkiewiczZhuEstimator does NOT assume ownership of
136  this FiniteElementSpace; will call its Update() method
137  when needed. */
139  FiniteElementSpace &flux_fes)
140  : current_sequence(-1),
141  total_error(),
142  anisotropic(false),
143  flux_averaging(0),
144  integ(&integ),
145  solution(&sol),
146  flux_space(&flux_fes),
147  with_coeff(false),
148  own_flux_fes(false)
149  { }
150 
151  /** @brief Consider the coefficient in BilinearFormIntegrator to calculate the
152  fluxes for the error estimator.*/
153  void SetWithCoeff(bool w_coeff = true) { with_coeff = w_coeff; }
154 
155  /** @brief Enable/disable anisotropic estimates. To enable this option, the
156  BilinearFormIntegrator must support the 'd_energy' parameter in its
157  ComputeFluxEnergy() method. */
158  void SetAnisotropic(bool aniso = true) { anisotropic = aniso; }
159 
160  /** @brief Set the way the flux is averaged (smoothed) across elements.
161 
162  When @a fa is zero (default), averaging is performed globally. When @a fa
163  is non-zero, the flux averaging is performed locally for each mesh
164  attribute, i.e. the flux is not averaged across interfaces between
165  different mesh attributes. */
166  void SetFluxAveraging(int fa) { flux_averaging = fa; }
167 
168  /// Return the total error from the last error estimate.
169  double GetTotalError() const { return total_error; }
170 
171  /// Get a Vector with all element errors.
172  virtual const Vector &GetLocalErrors()
173  {
174  if (MeshIsModified()) { ComputeEstimates(); }
175  return error_estimates;
176  }
177 
178  /** @brief Get an Array<int> with anisotropic flags for all mesh elements.
179  Return an empty array when anisotropic estimates are not available or
180  enabled. */
182  {
183  if (MeshIsModified()) { ComputeEstimates(); }
184  return aniso_flags;
185  }
186 
187  /// Reset the error estimator.
188  virtual void Reset() { current_sequence = -1; }
189 
190  /** @brief Destroy a ZienkiewiczZhuEstimator object. Destroys, if owned, the
191  FiniteElementSpace, flux_space. */
193  {
194  if (own_flux_fes) { delete flux_space; }
195  }
196 };
197 
198 
199 #ifdef MFEM_USE_MPI
200 
201 /** @brief The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu
202  error estimation procedure where the flux averaging is replaced by a global
203  L2 projection (requiring a mass matrix solve).
204 
205  The required BilinearFormIntegrator must implement the methods
206  ComputeElementFlux() and ComputeFluxEnergy().
207 
208  Implemented for the parallel case only.
209  */
211 {
212 protected:
214  int local_norm_p; ///< Local L_p norm to use, default is 1.
216  double total_error;
217 
218  BilinearFormIntegrator *integ; ///< Not owned.
219  ParGridFunction *solution; ///< Not owned.
220 
221  ParFiniteElementSpace *flux_space; /**< @brief Ownership based on the flag
222  own_flux_fes. Its Update() method is called automatically by this class
223  when needed. */
224  ParFiniteElementSpace *smooth_flux_space; /**< @brief Ownership based on the
225  flag own_flux_fes. Its Update() method is called automatically by this
226  class when needed.*/
227  bool own_flux_fes; ///< Ownership flag for flux_space and smooth_flux_space.
228 
229  /// Initialize with the integrator, solution, and flux finite element spaces.
231  ParGridFunction &sol,
232  ParFiniteElementSpace *flux_fes,
233  ParFiniteElementSpace *smooth_flux_fes)
234  {
235  current_sequence = -1;
236  local_norm_p = 1;
237  total_error = 0.0;
238  this->integ = &integ;
239  solution = &sol;
240  flux_space = flux_fes;
241  smooth_flux_space = smooth_flux_fes;
242  }
243 
244  /// Check if the mesh of the solution was modified.
246  {
247  long mesh_sequence = solution->FESpace()->GetMesh()->GetSequence();
248  MFEM_ASSERT(mesh_sequence >= current_sequence, "");
249  return (mesh_sequence > current_sequence);
250  }
251 
252  /// Compute the element error estimates.
253  void ComputeEstimates();
254 
255 public:
256  /** @brief Construct a new L2ZienkiewiczZhuEstimator object.
257  @param integ This BilinearFormIntegrator must implement the methods
258  ComputeElementFlux() and ComputeFluxEnergy().
259  @param sol The solution field whose error is to be estimated.
260  @param flux_fes The L2ZienkiewiczZhuEstimator assumes ownership of this
261  FiniteElementSpace and will call its Update() method when
262  needed.
263  @param smooth_flux_fes
264  The L2ZienkiewiczZhuEstimator assumes ownership of this
265  FiniteElementSpace and will call its Update() method when
266  needed. */
268  ParGridFunction &sol,
269  ParFiniteElementSpace *flux_fes,
270  ParFiniteElementSpace *smooth_flux_fes)
271  { Init(integ, sol, flux_fes, smooth_flux_fes); own_flux_fes = true; }
272 
273  /** @brief Construct a new L2ZienkiewiczZhuEstimator object.
274  @param integ This BilinearFormIntegrator must implement the methods
275  ComputeElementFlux() and ComputeFluxEnergy().
276  @param sol The solution field whose error is to be estimated.
277  @param flux_fes The L2ZienkiewiczZhuEstimator does NOT assume ownership
278  of this FiniteElementSpace; will call its Update() method
279  when needed.
280  @param smooth_flux_fes
281  The L2ZienkiewiczZhuEstimator does NOT assume ownership
282  of this FiniteElementSpace; will call its Update() method
283  when needed. */
285  ParGridFunction &sol,
286  ParFiniteElementSpace &flux_fes,
287  ParFiniteElementSpace &smooth_flux_fes)
288  { Init(integ, sol, &flux_fes, &smooth_flux_fes); own_flux_fes = false; }
289 
290  /** @brief Set the exponent, p, of the Lp norm used for computing the local
291  element errors. Default value is 1. */
293 
294  /// Return the total error from the last error estimate.
295  double GetTotalError() const { return total_error; }
296 
297  /// Get a Vector with all element errors.
298  virtual const Vector &GetLocalErrors()
299  {
300  if (MeshIsModified()) { ComputeEstimates(); }
301  return error_estimates;
302  }
303 
304  /// Reset the error estimator.
305  virtual void Reset() { current_sequence = -1; }
306 
307  /** @brief Destroy a L2ZienkiewiczZhuEstimator object. Destroys, if owned,
308  the FiniteElementSpace, flux_space. */
310  {
311  if (own_flux_fes) { delete flux_space; delete smooth_flux_space; }
312  }
313 };
314 
315 #endif // MFEM_USE_MPI
316 
317 /** @brief The LpErrorEstimator class compares the solution to a known
318  coefficient.
319 
320  This class can be used, for example, to adapt a mesh to a non-trivial
321  initial condition in a time-dependent simulation. It can also be used to
322  force refinement in the neighborhood of small features before switching to a
323  more traditional error estimator.
324 
325  The LpErrorEstimator supports either scalar or vector coefficients and works
326  both in serial and in parallel.
327 */
329 {
330 protected:
334 
338 
339  /// Check if the mesh of the solution was modified.
341  {
342  long mesh_sequence = sol->FESpace()->GetMesh()->GetSequence();
343  MFEM_ASSERT(mesh_sequence >= current_sequence, "");
344  return (mesh_sequence > current_sequence);
345  }
346 
347  /// Compute the element error estimates.
348  void ComputeEstimates();
349 
350 public:
351  /** @brief Construct a new LpErrorEstimator object for a scalar field.
352  @param p Integer which selects which Lp norm to use.
353  @param sol The GridFunction representation of the scalar field.
354  Note: the coefficient must be set before use with the SetCoef method.
355  */
357  : current_sequence(-1), local_norm_p(p),
358  error_estimates(0), coef(NULL), vcoef(NULL), sol(&sol) { }
359 
360  /** @brief Construct a new LpErrorEstimator object for a scalar field.
361  @param p Integer which selects which Lp norm to use.
362  @param coef The scalar Coefficient to compare to the solution.
363  @param sol The GridFunction representation of the scalar field.
364  */
366  : current_sequence(-1), local_norm_p(p),
367  error_estimates(0), coef(&coef), vcoef(NULL), sol(&sol) { }
368 
369  /** @brief Construct a new LpErrorEstimator object for a vector field.
370  @param p Integer which selects which Lp norm to use.
371  @param coef The vector VectorCoefficient to compare to the solution.
372  @param sol The GridFunction representation of the vector field.
373  */
375  : current_sequence(-1), local_norm_p(p),
376  error_estimates(0), coef(NULL), vcoef(&coef), sol(&sol) { }
377 
378  /** @brief Set the exponent, p, of the Lp norm used for computing the local
379  element errors. */
381 
382  void SetCoef(Coefficient &A) { coef = &A; }
383  void SetCoef(VectorCoefficient &A) { vcoef = &A; }
384 
385  /// Reset the error estimator.
386  virtual void Reset() { current_sequence = -1; }
387 
388  /// Get a Vector with all element errors.
389  virtual const Vector &GetLocalErrors()
390  {
391  if (MeshIsModified()) { ComputeEstimates(); }
392  return error_estimates;
393  }
394 
395  /// Destructor
396  virtual ~LpErrorEstimator() {}
397 };
398 
399 } // namespace mfem
400 
401 #endif // MFEM_ERROR_ESTIMATORS
virtual const Vector & GetLocalErrors()
Get a Vector with all element errors.
Definition: estimators.hpp:389
virtual ~L2ZienkiewiczZhuEstimator()
Destroy a L2ZienkiewiczZhuEstimator object. Destroys, if owned, the FiniteElementSpace, flux_space.
Definition: estimators.hpp:309
void SetAnisotropic(bool aniso=true)
Enable/disable anisotropic estimates. To enable this option, the BilinearFormIntegrator must support ...
Definition: estimators.hpp:158
The LpErrorEstimator class compares the solution to a known coefficient.
Definition: estimators.hpp:328
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
L2ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, ParGridFunction &sol, ParFiniteElementSpace *flux_fes, ParFiniteElementSpace *smooth_flux_fes)
Construct a new L2ZienkiewiczZhuEstimator object.
Definition: estimators.hpp:267
Base class for vector Coefficients that optionally depend on time and space.
void ComputeEstimates()
Compute the element error estimates.
Definition: estimators.cpp:36
ParFiniteElementSpace * flux_space
Ownership based on the flag own_flux_fes. Its Update() method is called automatically by this class w...
Definition: estimators.hpp:221
void SetLocalErrorNormP(int p)
Set the exponent, p, of the Lp norm used for computing the local element errors.
Definition: estimators.hpp:380
VectorCoefficient * vcoef
Definition: estimators.hpp:336
LpErrorEstimator(int p, Coefficient &coef, GridFunction &sol)
Construct a new LpErrorEstimator object for a scalar field.
Definition: estimators.hpp:365
bool MeshIsModified()
Check if the mesh of the solution was modified.
Definition: estimators.hpp:245
virtual const Vector & GetLocalErrors()
Get a Vector with all element errors.
Definition: estimators.hpp:298
void SetLocalErrorNormP(int p)
Set the exponent, p, of the Lp norm used for computing the local element errors. Default value is 1...
Definition: estimators.hpp:292
Abstract parallel finite element space.
Definition: pfespace.hpp:28
GridFunction * sol
Definition: estimators.hpp:337
void SetWithCoeff(bool w_coeff=true)
Consider the coefficient in BilinearFormIntegrator to calculate the fluxes for the error estimator...
Definition: estimators.hpp:153
bool MeshIsModified()
Check if the mesh of the solution was modified.
Definition: estimators.hpp:340
virtual void Reset()
Reset the error estimator.
Definition: estimators.hpp:305
virtual void Reset()
Reset the error estimator.
Definition: estimators.hpp:386
The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure where t...
Definition: estimators.hpp:210
virtual ~ErrorEstimator()
Destruct the error estimator.
Definition: estimators.hpp:49
void ComputeEstimates()
Compute the element error estimates.
Definition: estimators.cpp:17
ParGridFunction * solution
Not owned.
Definition: estimators.hpp:219
void Init(BilinearFormIntegrator &integ, ParGridFunction &sol, ParFiniteElementSpace *flux_fes, ParFiniteElementSpace *smooth_flux_fes)
Initialize with the integrator, solution, and flux finite element spaces.
Definition: estimators.hpp:230
long GetSequence() const
Definition: mesh.hpp:1227
bool own_flux_fes
Ownership flag for flux_space and smooth_flux_space.
Definition: estimators.hpp:227
LpErrorEstimator(int p, GridFunction &sol)
Construct a new LpErrorEstimator object for a scalar field.
Definition: estimators.hpp:356
virtual const Vector & GetLocalErrors()=0
Get a Vector with all element errors.
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:314
Base class for all error estimators.
Definition: estimators.hpp:27
GridFunction * solution
Not owned.
Definition: estimators.hpp:92
double GetTotalError() const
Return the total error from the last error estimate.
Definition: estimators.hpp:169
int local_norm_p
Local L_p norm to use, default is 1.
Definition: estimators.hpp:214
void SetFluxAveraging(int fa)
Set the way the flux is averaged (smoothed) across elements.
Definition: estimators.hpp:166
LpErrorEstimator(int p, VectorCoefficient &coef, GridFunction &sol)
Construct a new LpErrorEstimator object for a vector field.
Definition: estimators.hpp:374
Base class for all element based error estimators.
Definition: estimators.hpp:39
void SetCoef(Coefficient &A)
Definition: estimators.hpp:382
L2ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, ParGridFunction &sol, ParFiniteElementSpace &flux_fes, ParFiniteElementSpace &smooth_flux_fes)
Construct a new L2ZienkiewiczZhuEstimator object.
Definition: estimators.hpp:284
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:595
virtual const Array< int > & GetAnisotropicFlags()=0
Get an Array&lt;int&gt; with anisotropic flags for all mesh elements.
virtual void Reset()=0
Force recomputation of the estimates on the next call to GetLocalErrors.
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:31
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
ParFiniteElementSpace * smooth_flux_space
Ownership based on the flag own_flux_fes. Its Update() method is called automatically by this class w...
Definition: estimators.hpp:224
double GetTotalError() const
Return the total error from the last error estimate.
Definition: estimators.hpp:295
virtual void Reset()
Reset the error estimator.
Definition: estimators.hpp:188
FiniteElementSpace * flux_space
Ownership based on own_flux_fes. Its Update() method is called automatically by this class when neede...
Definition: estimators.hpp:94
virtual const Array< int > & GetAnisotropicFlags()
Get an Array&lt;int&gt; with anisotropic flags for all mesh elements. Return an empty array when anisotropi...
Definition: estimators.hpp:181
The AnisotropicErrorEstimator class is the base class for all error estimators that compute one non-n...
Definition: estimators.hpp:57
BilinearFormIntegrator * integ
Not owned.
Definition: estimators.hpp:218
void ComputeEstimates()
Compute the element error estimates.
Definition: estimators.cpp:53
virtual ~LpErrorEstimator()
Destructor.
Definition: estimators.hpp:396
ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, GridFunction &sol, FiniteElementSpace *flux_fes)
Construct a new ZienkiewiczZhuEstimator object.
Definition: estimators.hpp:118
The ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure...
Definition: estimators.hpp:81
Vector data type.
Definition: vector.hpp:51
virtual const Vector & GetLocalErrors()
Get a Vector with all element errors.
Definition: estimators.hpp:172
bool MeshIsModified()
Check if the mesh of the solution was modified.
Definition: estimators.hpp:100
virtual ~ZienkiewiczZhuEstimator()
Destroy a ZienkiewiczZhuEstimator object. Destroys, if owned, the FiniteElementSpace, flux_space.
Definition: estimators.hpp:192
ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, GridFunction &sol, FiniteElementSpace &flux_fes)
Construct a new ZienkiewiczZhuEstimator object.
Definition: estimators.hpp:138
Class for parallel grid function.
Definition: pgridfunc.hpp:32
BilinearFormIntegrator * integ
Not owned.
Definition: estimators.hpp:91
void SetCoef(VectorCoefficient &A)
Definition: estimators.hpp:383
bool own_flux_fes
Ownership flag for flux_space.
Definition: estimators.hpp:97