MFEM  v3.3
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, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
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 
81  BilinearFormIntegrator *integ; ///< Not owned.
82  GridFunction *solution; ///< Not owned.
83 
84  FiniteElementSpace *flux_space; /**< @brief Ownership based on own_flux_fes.
85  Its Update() method is called automatically by this class when needed. */
86  bool own_flux_fes; ///< Ownership flag for flux_space.
87 
88  /// Check if the mesh of the solution was modified.
90  {
91  long mesh_sequence = solution->FESpace()->GetMesh()->GetSequence();
92  MFEM_ASSERT(mesh_sequence >= current_sequence, "");
93  return (mesh_sequence > current_sequence);
94  }
95 
96  /// Compute the element error estimates.
97  void ComputeEstimates();
98 
99 public:
100  /** @brief Construct a new ZienkiewiczZhuEstimator object.
101  @param integ This BilinearFormIntegrator must implement the methods
102  ComputeElementFlux() and ComputeFluxEnergy().
103  @param sol The solution field whose error is to be estimated.
104  @param flux_fes The ZienkiewiczZhuEstimator assumes ownership of this
105  FiniteElementSpace and will call its Update() method when
106  needed. */
108  FiniteElementSpace *flux_fes)
109  : current_sequence(-1),
110  total_error(),
111  anisotropic(false),
112  integ(&integ),
113  solution(&sol),
114  flux_space(flux_fes),
115  own_flux_fes(true)
116  { }
117 
118  /** @brief Construct a new ZienkiewiczZhuEstimator object.
119  @param integ This BilinearFormIntegrator must implement the methods
120  ComputeElementFlux() and ComputeFluxEnergy().
121  @param sol The solution field whose error is to be estimated.
122  @param flux_fes The ZienkiewiczZhuEstimator does NOT assume ownership of
123  this FiniteElementSpace; will call its Update() method
124  when needed. */
126  FiniteElementSpace &flux_fes)
127  : current_sequence(-1),
128  total_error(),
129  anisotropic(false),
130  integ(&integ),
131  solution(&sol),
132  flux_space(&flux_fes),
133  own_flux_fes(false)
134  { }
135 
136  /** @brief Enable/disable anisotropic estimates. To enable this option, the
137  BilinearFormIntegrator must support the 'd_energy' parameter in its
138  ComputeFluxEnergy() method. */
139  void SetAnisotropic(bool aniso = true) { anisotropic = aniso; }
140 
141  /// Return the total error from the last error estimate.
142  double GetTotalError() const { return total_error; }
143 
144  /// Get a Vector with all element errors.
145  virtual const Vector &GetLocalErrors()
146  {
147  if (MeshIsModified()) { ComputeEstimates(); }
148  return error_estimates;
149  }
150 
151  /** @brief Get an Array<int> with anisotropic flags for all mesh elements.
152  Return an empty array when anisotropic estimates are not available or
153  enabled. */
155  {
156  if (MeshIsModified()) { ComputeEstimates(); }
157  return aniso_flags;
158  }
159 
160  /// Reset the error estimator.
161  virtual void Reset() { current_sequence = -1; }
162 
163  /** @brief Destroy a ZienkiewiczZhuEstimator object. Destroys, if owned, the
164  FiniteElementSpace, flux_space. */
166  {
167  if (own_flux_fes) { delete flux_space; }
168  }
169 };
170 
171 
172 #ifdef MFEM_USE_MPI
173 
174 /** @brief The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu
175  error estimation procedure where the flux averaging is replaced by a global
176  L2 projection (requiring a mass matrix solve).
177 
178  The required BilinearFormIntegrator must implement the methods
179  ComputeElementFlux() and ComputeFluxEnergy().
180 
181  Implemented for the parallel case only.
182  */
184 {
185 protected:
187  int local_norm_p; ///< Local L_p norm to use, default is 1.
189  double total_error;
190 
191  BilinearFormIntegrator *integ; ///< Not owned.
192  ParGridFunction *solution; ///< Not owned.
193 
194  ParFiniteElementSpace *flux_space; /**< @brief Ownership based on the flag
195  own_flux_fes. Its Update() method is called automatically by this class
196  when needed. */
197  ParFiniteElementSpace *smooth_flux_space; /**< @brief Ownership based on the
198  flag own_flux_fes. Its Update() method is called automatically by this
199  class when needed.*/
200  bool own_flux_fes; ///< Ownership flag for flux_space and smooth_flux_space.
201 
203  ParGridFunction &sol,
204  ParFiniteElementSpace *flux_fes,
205  ParFiniteElementSpace *smooth_flux_fes)
206  {
207  current_sequence = -1;
208  local_norm_p = 1;
209  total_error = 0.0;
210  this->integ = &integ;
211  solution = &sol;
212  flux_space = flux_fes;
213  smooth_flux_space = smooth_flux_fes;
214  }
215 
216  /// Check if the mesh of the solution was modified.
218  {
219  long mesh_sequence = solution->FESpace()->GetMesh()->GetSequence();
220  MFEM_ASSERT(mesh_sequence >= current_sequence, "");
221  return (mesh_sequence > current_sequence);
222  }
223 
224  /// Compute the element error estimates.
225  void ComputeEstimates();
226 
227 public:
228  /** @brief Construct a new L2ZienkiewiczZhuEstimator object.
229  @param integ This BilinearFormIntegrator must implement the methods
230  ComputeElementFlux() and ComputeFluxEnergy().
231  @param sol The solution field whose error is to be estimated.
232  @param flux_fes The L2ZienkiewiczZhuEstimator assumes ownership of this
233  FiniteElementSpace and will call its Update() method when
234  needed.
235  @param smooth_flux_fes
236  The L2ZienkiewiczZhuEstimator assumes ownership of this
237  FiniteElementSpace and will call its Update() method when
238  needed. */
240  ParGridFunction &sol,
241  ParFiniteElementSpace *flux_fes,
242  ParFiniteElementSpace *smooth_flux_fes)
243  { Init(integ, sol, flux_fes, smooth_flux_fes); own_flux_fes = true; }
244 
245  /** @brief Construct a new L2ZienkiewiczZhuEstimator object.
246  @param integ This BilinearFormIntegrator must implement the methods
247  ComputeElementFlux() and ComputeFluxEnergy().
248  @param sol The solution field whose error is to be estimated.
249  @param flux_fes The L2ZienkiewiczZhuEstimator does NOT assume ownership
250  of this FiniteElementSpace; will call its Update() method
251  when needed.
252  @param smooth_flux_fes
253  The L2ZienkiewiczZhuEstimator does NOT assume ownership
254  of this FiniteElementSpace; will call its Update() method
255  when needed. */
257  ParGridFunction &sol,
258  ParFiniteElementSpace &flux_fes,
259  ParFiniteElementSpace &smooth_flux_fes)
260  { Init(integ, sol, &flux_fes, &smooth_flux_fes); own_flux_fes = false; }
261 
262  /** @brief Set the exponent, p, of the Lp norm used for computing the local
263  element errors. Default value is 1. */
264  void SetLocalErrorNormP(int p) { local_norm_p = p; }
265 
266  /// Return the total error from the last error estimate.
267  double GetTotalError() const { return total_error; }
268 
269  /// Get a Vector with all element errors.
270  virtual const Vector &GetLocalErrors()
271  {
272  if (MeshIsModified()) { ComputeEstimates(); }
273  return error_estimates;
274  }
275 
276  /// Reset the error estimator.
277  virtual void Reset() { current_sequence = -1; }
278 
279  /** @brief Destroy a L2ZienkiewiczZhuEstimator object. Destroys, if owned,
280  the FiniteElementSpace, flux_space. */
282  {
283  if (own_flux_fes) { delete flux_space; delete smooth_flux_space; }
284  }
285 };
286 
287 #endif // MFEM_USE_MPI
288 
289 } // namespace mfem
290 
291 #endif // MFEM_ERROR_ESTIMATORS
virtual ~L2ZienkiewiczZhuEstimator()
Destroy a L2ZienkiewiczZhuEstimator object. Destroys, if owned, the FiniteElementSpace, flux_space.
Definition: estimators.hpp:281
void SetAnisotropic(bool aniso=true)
Enable/disable anisotropic estimates. To enable this option, the BilinearFormIntegrator must support ...
Definition: estimators.hpp:139
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:239
void ComputeEstimates()
Compute the element error estimates.
Definition: estimators.cpp:34
ParFiniteElementSpace * flux_space
Ownership based on the flag own_flux_fes. Its Update() method is called automatically by this class w...
Definition: estimators.hpp:194
bool MeshIsModified()
Check if the mesh of the solution was modified.
Definition: estimators.hpp:217
virtual const Vector & GetLocalErrors()
Get a Vector with all element errors.
Definition: estimators.hpp:270
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:264
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void Reset()
Reset the error estimator.
Definition: estimators.hpp:277
The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure where t...
Definition: estimators.hpp:183
virtual ~ErrorEstimator()
Definition: estimators.hpp:48
void ComputeEstimates()
Compute the element error estimates.
Definition: estimators.cpp:17
ParGridFunction * solution
Not owned.
Definition: estimators.hpp:192
void Init(BilinearFormIntegrator &integ, ParGridFunction &sol, ParFiniteElementSpace *flux_fes, ParFiniteElementSpace *smooth_flux_fes)
Definition: estimators.hpp:202
long GetSequence() const
Definition: mesh.hpp:981
bool own_flux_fes
Ownership flag for flux_space and smooth_flux_space.
Definition: estimators.hpp:200
virtual const Vector & GetLocalErrors()=0
Get a Vector with all element errors.
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:136
Base class for all error estimators.
Definition: estimators.hpp:27
GridFunction * solution
Not owned.
Definition: estimators.hpp:82
double GetTotalError() const
Return the total error from the last error estimate.
Definition: estimators.hpp:142
int local_norm_p
Local L_p norm to use, default is 1.
Definition: estimators.hpp:187
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:256
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:293
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:22
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:197
double GetTotalError() const
Return the total error from the last error estimate.
Definition: estimators.hpp:267
virtual void Reset()
Reset the error estimator.
Definition: estimators.hpp:161
FiniteElementSpace * flux_space
Ownership based on own_flux_fes. Its Update() method is called automatically by this class when neede...
Definition: estimators.hpp:84
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:154
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:191
ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, GridFunction &sol, FiniteElementSpace *flux_fes)
Construct a new ZienkiewiczZhuEstimator object.
Definition: estimators.hpp:107
The ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure...
Definition: estimators.hpp:72
Vector data type.
Definition: vector.hpp:36
virtual const Vector & GetLocalErrors()
Get a Vector with all element errors.
Definition: estimators.hpp:145
bool MeshIsModified()
Check if the mesh of the solution was modified.
Definition: estimators.hpp:89
virtual ~ZienkiewiczZhuEstimator()
Destroy a ZienkiewiczZhuEstimator object. Destroys, if owned, the FiniteElementSpace, flux_space.
Definition: estimators.hpp:165
ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, GridFunction &sol, FiniteElementSpace &flux_fes)
Construct a new ZienkiewiczZhuEstimator object.
Definition: estimators.hpp:125
Class for parallel grid function.
Definition: pgridfunc.hpp:31
BilinearFormIntegrator * integ
Not owned.
Definition: estimators.hpp:81
bool own_flux_fes
Ownership flag for flux_space.
Definition: estimators.hpp:86