MFEM
v4.5.2
Finite element discretization library
|
The KellyErrorEstimator class provides a fast error indication strategy for smooth scalar parallel problems. More...
#include <estimators.hpp>
Public Types | |
using | ElementCoefficientFunction = std::function< double(Mesh *, const int)> |
Function type to compute the local coefficient hₑ of an element. More... | |
using | FaceCoefficientFunction = std::function< double(Mesh *, const int, const bool)> |
Function type to compute the local coefficient hₖ of a face. The third argument is true for shared faces and false for local faces. More... | |
Public Member Functions | |
KellyErrorEstimator (BilinearFormIntegrator &di_, GridFunction &sol_, FiniteElementSpace &flux_fes_, const Array< int > &attributes_=Array< int >()) | |
Construct a new KellyErrorEstimator object for a scalar field. More... | |
KellyErrorEstimator (BilinearFormIntegrator &di_, GridFunction &sol_, FiniteElementSpace *flux_fes_, const Array< int > &attributes_=Array< int >()) | |
Construct a new KellyErrorEstimator object for a scalar field. More... | |
~KellyErrorEstimator () | |
const Vector & | GetLocalErrors () override |
Get a Vector with all element errors. More... | |
void | Reset () override |
Reset the error estimator. More... | |
virtual double | GetTotalError () const override |
Return the total error from the last error estimate. More... | |
void | SetElementCoefficientFunction (ElementCoefficientFunction compute_element_coefficient_) |
Change the method to compute hₑ on a per-element basis. More... | |
void | SetFaceCoefficientFunction (FaceCoefficientFunction compute_face_coefficient_) |
Change the method to compute hₖ on a per-element basis. More... | |
void | ResetCoefficientFunctions () |
Change the coefficients back to default as described above. More... | |
Public Member Functions inherited from mfem::ErrorEstimator | |
virtual | ~ErrorEstimator () |
Destruct the error estimator. More... | |
Public Member Functions inherited from mfem::AbstractErrorEstimator | |
virtual | ~AbstractErrorEstimator () |
The KellyErrorEstimator class provides a fast error indication strategy for smooth scalar parallel problems.
The Kelly error indicator is based on the following papers:
[1] Kelly, D. W., et al. "A posteriori error analysis and adaptive processes in the finite element method: Part I—Error analysis." International journal for numerical methods in engineering 19.11 (1983): 1593-1619.
[2] De SR Gago, J. P., et al. "A posteriori error analysis and adaptive processes in the finite element method: Part II—Adaptive mesh refinement." International journal for numerical methods in engineering 19.11 (1983): 1621-1656.
It can be roughly described by: ||∇(u-uₕ)||ₑ ≅ √( C hₑ ∑ₖ (hₖ ∫ |J[∇uₕ]|²) dS ) where "e" denotes an element, ||⋅||ₑ the corresponding local norm and k the corresponding faces. u is the analytic solution and uₕ the discretized solution. hₖ and hₑ are factors dependent on the face and element geometry. J is the jump function, i.e. the difference between the limits at each point for each side of the face. A custom method to compute hₖ can be provided. It is also possible to estimate the error only on a subspace by feeding this class an attribute array describing the subspace.
Definition at line 555 of file estimators.hpp.
using mfem::KellyErrorEstimator::ElementCoefficientFunction = std::function<double(Mesh*, const int)> |
Function type to compute the local coefficient hₑ of an element.
Definition at line 560 of file estimators.hpp.
using mfem::KellyErrorEstimator::FaceCoefficientFunction = std::function<double(Mesh*, const int, const bool)> |
Function type to compute the local coefficient hₖ of a face. The third argument is true for shared faces and false for local faces.
Definition at line 564 of file estimators.hpp.
mfem::KellyErrorEstimator::KellyErrorEstimator | ( | BilinearFormIntegrator & | di_, |
GridFunction & | sol_, | ||
FiniteElementSpace & | flux_fes_, | ||
const Array< int > & | attributes_ = Array<int>() |
||
) |
Construct a new KellyErrorEstimator object for a scalar field.
di_ | The bilinearform to compute the interface flux. |
sol_ | The solution field whose error is to be estimated. |
flux_fes_ | The finite element space for the interface flux. |
attributes_ | The attributes of the subdomain(s) for which the error should be estimated. An empty array results in estimating the error over the complete domain. |
Definition at line 64 of file estimators.cpp.
mfem::KellyErrorEstimator::KellyErrorEstimator | ( | BilinearFormIntegrator & | di_, |
GridFunction & | sol_, | ||
FiniteElementSpace * | flux_fes_, | ||
const Array< int > & | attributes_ = Array<int>() |
||
) |
Construct a new KellyErrorEstimator object for a scalar field.
di_ | The bilinearform to compute the interface flux. |
sol_ | The solution field whose error is to be estimated. |
flux_fes_ | The finite element space for the interface flux. |
attributes_ | The attributes of the subdomain(s) for which the error should be estimated. An empty array results in estimating the error over the complete domain. |
Definition at line 80 of file estimators.cpp.
mfem::KellyErrorEstimator::~KellyErrorEstimator | ( | ) |
Definition at line 96 of file estimators.cpp.
|
inlineoverridevirtual |
Get a Vector with all element errors.
Implements mfem::ErrorEstimator.
Definition at line 652 of file estimators.hpp.
|
inlineoverridevirtual |
Return the total error from the last error estimate.
Reimplemented from mfem::ErrorEstimator.
Definition at line 664 of file estimators.hpp.
|
inlineoverridevirtual |
Reset the error estimator.
Implements mfem::ErrorEstimator.
Definition at line 662 of file estimators.hpp.
void mfem::KellyErrorEstimator::ResetCoefficientFunctions | ( | ) |
Change the coefficients back to default as described above.
Definition at line 104 of file estimators.cpp.
|
inline |
Change the method to compute hₑ on a per-element basis.
compute_element_coefficient_ | A function taking a mesh and an element index to compute the local hₑ for the element. |
Definition at line 671 of file estimators.hpp.
|
inline |
Change the method to compute hₖ on a per-element basis.
compute_face_coefficient_ | A function taking a mesh and a face index to compute the local hₖ for the face. |
Definition at line 682 of file estimators.hpp.