MFEM  v4.1.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  virtual ~ErrorEstimator() { }
49 };
50 
51 
52 /** @brief The AnisotropicErrorEstimator class is the base class for all error
53  estimators that compute one non-negative real (double) number and an
54  anisotropic flag for every element in the Mesh.
55  */
57 {
58 public:
59  /** @brief Get an Array<int> with anisotropic flags for all mesh elements.
60  @return An empty array when anisotropic estimates are not available or
61  enabled. */
62  virtual const Array<int> &GetAnisotropicFlags() = 0;
63 };
64 
65 
66 /** @brief The ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu
67  error estimation procedure.
68 
69  The required BilinearFormIntegrator must implement the methods
70  ComputeElementFlux() and ComputeFluxEnergy().
71  */
73 {
74 protected:
77  double total_error;
80  int flux_averaging; // see SetFluxAveraging()
81 
82  BilinearFormIntegrator *integ; ///< Not owned.
83  GridFunction *solution; ///< Not owned.
84 
85  FiniteElementSpace *flux_space; /**< @brief Ownership based on own_flux_fes.
86  Its Update() method is called automatically by this class when needed. */
87  bool with_coeff;
88  bool own_flux_fes; ///< Ownership flag for flux_space.
89 
90  /// Check if the mesh of the solution was modified.
92  {
93  long mesh_sequence = solution->FESpace()->GetMesh()->GetSequence();
94  MFEM_ASSERT(mesh_sequence >= current_sequence, "");
95  return (mesh_sequence > current_sequence);
96  }
97 
98  /// Compute the element error estimates.
99  void ComputeEstimates();
100 
101 public:
102  /** @brief Construct a new ZienkiewiczZhuEstimator object.
103  @param integ This BilinearFormIntegrator must implement the methods
104  ComputeElementFlux() and ComputeFluxEnergy().
105  @param sol The solution field whose error is to be estimated.
106  @param flux_fes The ZienkiewiczZhuEstimator assumes ownership of this
107  FiniteElementSpace and will call its Update() method when
108  needed.*/
110  FiniteElementSpace *flux_fes)
111  : current_sequence(-1),
112  total_error(),
113  anisotropic(false),
114  flux_averaging(0),
115  integ(&integ),
116  solution(&sol),
117  flux_space(flux_fes),
118  with_coeff(false),
119  own_flux_fes(true)
120  { }
121 
122  /** @brief Construct a new ZienkiewiczZhuEstimator object.
123  @param integ This BilinearFormIntegrator must implement the methods
124  ComputeElementFlux() and ComputeFluxEnergy().
125  @param sol The solution field whose error is to be estimated.
126  @param flux_fes The ZienkiewiczZhuEstimator does NOT assume ownership of
127  this FiniteElementSpace; will call its Update() method
128  when needed. */
130  FiniteElementSpace &flux_fes)
131  : current_sequence(-1),
132  total_error(),
133  anisotropic(false),
134  flux_averaging(0),
135  integ(&integ),
136  solution(&sol),
137  flux_space(&flux_fes),
138  with_coeff(false),
139  own_flux_fes(false)
140  { }
141 
142  /** @brief Consider the coefficient in BilinearFormIntegrator to calculate the
143  fluxes for the error estimator.*/
144  void SetWithCoeff(bool w_coeff = true) { with_coeff = w_coeff; }
145 
146  /** @brief Enable/disable anisotropic estimates. To enable this option, the
147  BilinearFormIntegrator must support the 'd_energy' parameter in its
148  ComputeFluxEnergy() method. */
149  void SetAnisotropic(bool aniso = true) { anisotropic = aniso; }
150 
151  /** @brief Set the way the flux is averaged (smoothed) across elements.
152 
153  When @a fa is zero (default), averaging is performed globally. When @a fa
154  is non-zero, the flux averaging is performed locally for each mesh
155  attribute, i.e. the flux is not averaged across interfaces between
156  different mesh attributes. */
157  void SetFluxAveraging(int fa) { flux_averaging = fa; }
158 
159  /// Return the total error from the last error estimate.
160  double GetTotalError() const { return total_error; }
161 
162  /// Get a Vector with all element errors.
163  virtual const Vector &GetLocalErrors()
164  {
165  if (MeshIsModified()) { ComputeEstimates(); }
166  return error_estimates;
167  }
168 
169  /** @brief Get an Array<int> with anisotropic flags for all mesh elements.
170  Return an empty array when anisotropic estimates are not available or
171  enabled. */
173  {
174  if (MeshIsModified()) { ComputeEstimates(); }
175  return aniso_flags;
176  }
177 
178  /// Reset the error estimator.
179  virtual void Reset() { current_sequence = -1; }
180 
181  /** @brief Destroy a ZienkiewiczZhuEstimator object. Destroys, if owned, the
182  FiniteElementSpace, flux_space. */
184  {
185  if (own_flux_fes) { delete flux_space; }
186  }
187 };
188 
189 
190 #ifdef MFEM_USE_MPI
191 
192 /** @brief The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu
193  error estimation procedure where the flux averaging is replaced by a global
194  L2 projection (requiring a mass matrix solve).
195 
196  The required BilinearFormIntegrator must implement the methods
197  ComputeElementFlux() and ComputeFluxEnergy().
198 
199  Implemented for the parallel case only.
200  */
202 {
203 protected:
205  int local_norm_p; ///< Local L_p norm to use, default is 1.
207  double total_error;
208 
209  BilinearFormIntegrator *integ; ///< Not owned.
210  ParGridFunction *solution; ///< Not owned.
211 
212  ParFiniteElementSpace *flux_space; /**< @brief Ownership based on the flag
213  own_flux_fes. Its Update() method is called automatically by this class
214  when needed. */
215  ParFiniteElementSpace *smooth_flux_space; /**< @brief Ownership based on the
216  flag own_flux_fes. Its Update() method is called automatically by this
217  class when needed.*/
218  bool own_flux_fes; ///< Ownership flag for flux_space and smooth_flux_space.
219 
221  ParGridFunction &sol,
222  ParFiniteElementSpace *flux_fes,
223  ParFiniteElementSpace *smooth_flux_fes)
224  {
225  current_sequence = -1;
226  local_norm_p = 1;
227  total_error = 0.0;
228  this->integ = &integ;
229  solution = &sol;
230  flux_space = flux_fes;
231  smooth_flux_space = smooth_flux_fes;
232  }
233 
234  /// Check if the mesh of the solution was modified.
236  {
237  long mesh_sequence = solution->FESpace()->GetMesh()->GetSequence();
238  MFEM_ASSERT(mesh_sequence >= current_sequence, "");
239  return (mesh_sequence > current_sequence);
240  }
241 
242  /// Compute the element error estimates.
243  void ComputeEstimates();
244 
245 public:
246  /** @brief Construct a new L2ZienkiewiczZhuEstimator object.
247  @param integ This BilinearFormIntegrator must implement the methods
248  ComputeElementFlux() and ComputeFluxEnergy().
249  @param sol The solution field whose error is to be estimated.
250  @param flux_fes The L2ZienkiewiczZhuEstimator assumes ownership of this
251  FiniteElementSpace and will call its Update() method when
252  needed.
253  @param smooth_flux_fes
254  The L2ZienkiewiczZhuEstimator assumes ownership of this
255  FiniteElementSpace and will call its Update() method when
256  needed. */
258  ParGridFunction &sol,
259  ParFiniteElementSpace *flux_fes,
260  ParFiniteElementSpace *smooth_flux_fes)
261  { Init(integ, sol, flux_fes, smooth_flux_fes); own_flux_fes = true; }
262 
263  /** @brief Construct a new L2ZienkiewiczZhuEstimator object.
264  @param integ This BilinearFormIntegrator must implement the methods
265  ComputeElementFlux() and ComputeFluxEnergy().
266  @param sol The solution field whose error is to be estimated.
267  @param flux_fes The L2ZienkiewiczZhuEstimator does NOT assume ownership
268  of this FiniteElementSpace; will call its Update() method
269  when needed.
270  @param smooth_flux_fes
271  The L2ZienkiewiczZhuEstimator does NOT assume ownership
272  of this FiniteElementSpace; will call its Update() method
273  when needed. */
275  ParGridFunction &sol,
276  ParFiniteElementSpace &flux_fes,
277  ParFiniteElementSpace &smooth_flux_fes)
278  { Init(integ, sol, &flux_fes, &smooth_flux_fes); own_flux_fes = false; }
279 
280  /** @brief Set the exponent, p, of the Lp norm used for computing the local
281  element errors. Default value is 1. */
282  void SetLocalErrorNormP(int p) { local_norm_p = p; }
283 
284  /// Return the total error from the last error estimate.
285  double GetTotalError() const { return total_error; }
286 
287  /// Get a Vector with all element errors.
288  virtual const Vector &GetLocalErrors()
289  {
290  if (MeshIsModified()) { ComputeEstimates(); }
291  return error_estimates;
292  }
293 
294  /// Reset the error estimator.
295  virtual void Reset() { current_sequence = -1; }
296 
297  /** @brief Destroy a L2ZienkiewiczZhuEstimator object. Destroys, if owned,
298  the FiniteElementSpace, flux_space. */
300  {
301  if (own_flux_fes) { delete flux_space; delete smooth_flux_space; }
302  }
303 };
304 
305 #endif // MFEM_USE_MPI
306 
307 } // namespace mfem
308 
309 #endif // MFEM_ERROR_ESTIMATORS
virtual ~L2ZienkiewiczZhuEstimator()
Destroy a L2ZienkiewiczZhuEstimator object. Destroys, if owned, the FiniteElementSpace, flux_space.
Definition: estimators.hpp:299
void SetAnisotropic(bool aniso=true)
Enable/disable anisotropic estimates. To enable this option, the BilinearFormIntegrator must support ...
Definition: estimators.hpp:149
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
L2ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, ParGridFunction &sol, ParFiniteElementSpace *flux_fes, ParFiniteElementSpace *smooth_flux_fes)
Construct a new L2ZienkiewiczZhuEstimator object.
Definition: estimators.hpp:257
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:212
bool MeshIsModified()
Check if the mesh of the solution was modified.
Definition: estimators.hpp:235
virtual const Vector & GetLocalErrors()
Get a Vector with all element errors.
Definition: estimators.hpp:288
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:282
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetWithCoeff(bool w_coeff=true)
Consider the coefficient in BilinearFormIntegrator to calculate the fluxes for the error estimator...
Definition: estimators.hpp:144
virtual void Reset()
Reset the error estimator.
Definition: estimators.hpp:295
The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure where t...
Definition: estimators.hpp:201
virtual ~ErrorEstimator()
Definition: estimators.hpp:48
void ComputeEstimates()
Compute the element error estimates.
Definition: estimators.cpp:17
ParGridFunction * solution
Not owned.
Definition: estimators.hpp:210
void Init(BilinearFormIntegrator &integ, ParGridFunction &sol, ParFiniteElementSpace *flux_fes, ParFiniteElementSpace *smooth_flux_fes)
Definition: estimators.hpp:220
long GetSequence() const
Definition: mesh.hpp:1181
bool own_flux_fes
Ownership flag for flux_space and smooth_flux_space.
Definition: estimators.hpp:218
virtual const Vector & GetLocalErrors()=0
Get a Vector with all element errors.
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:295
Base class for all error estimators.
Definition: estimators.hpp:27
GridFunction * solution
Not owned.
Definition: estimators.hpp:83
double GetTotalError() const
Return the total error from the last error estimate.
Definition: estimators.hpp:160
int local_norm_p
Local L_p norm to use, default is 1.
Definition: estimators.hpp:205
void SetFluxAveraging(int fa)
Set the way the flux is averaged (smoothed) across elements.
Definition: estimators.hpp:157
Base class for all element based error estimators.
Definition: estimators.hpp:39
L2ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, ParGridFunction &sol, ParFiniteElementSpace &flux_fes, ParFiniteElementSpace &smooth_flux_fes)
Construct a new L2ZienkiewiczZhuEstimator object.
Definition: estimators.hpp:274
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:442
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:24
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
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:215
double GetTotalError() const
Return the total error from the last error estimate.
Definition: estimators.hpp:285
virtual void Reset()
Reset the error estimator.
Definition: estimators.hpp:179
FiniteElementSpace * flux_space
Ownership based on own_flux_fes. Its Update() method is called automatically by this class when neede...
Definition: estimators.hpp:85
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:172
The AnisotropicErrorEstimator class is the base class for all error estimators that compute one non-n...
Definition: estimators.hpp:56
BilinearFormIntegrator * integ
Not owned.
Definition: estimators.hpp:209
ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, GridFunction &sol, FiniteElementSpace *flux_fes)
Construct a new ZienkiewiczZhuEstimator object.
Definition: estimators.hpp:109
The ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure...
Definition: estimators.hpp:72
Vector data type.
Definition: vector.hpp:48
virtual const Vector & GetLocalErrors()
Get a Vector with all element errors.
Definition: estimators.hpp:163
bool MeshIsModified()
Check if the mesh of the solution was modified.
Definition: estimators.hpp:91
virtual ~ZienkiewiczZhuEstimator()
Destroy a ZienkiewiczZhuEstimator object. Destroys, if owned, the FiniteElementSpace, flux_space.
Definition: estimators.hpp:183
ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, GridFunction &sol, FiniteElementSpace &flux_fes)
Construct a new ZienkiewiczZhuEstimator object.
Definition: estimators.hpp:129
Class for parallel grid function.
Definition: pgridfunc.hpp:32
BilinearFormIntegrator * integ
Not owned.
Definition: estimators.hpp:82
bool own_flux_fes
Ownership flag for flux_space.
Definition: estimators.hpp:88