MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Types | Public Member Functions | List of all members
mfem::KellyErrorEstimator Class Referencefinal

The KellyErrorEstimator class provides a fast error indication strategy for smooth scalar parallel problems. More...

#include <estimators.hpp>

Inheritance diagram for mfem::KellyErrorEstimator:
[legend]
Collaboration diagram for mfem::KellyErrorEstimator:
[legend]

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 VectorGetLocalErrors () 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 ()
 

Detailed Description

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.

Note
This algorithm is appropriate only for problems with scalar diffusion coefficients (e.g. Poisson problems), because it measures only the flux of the gradient of the discrete solution. The current implementation also does not include the volume term present in Equation 75 of Kelly et al [1]. Instead, it includes only the flux term recommended in Equation 82. The current implementation also does not include the constant factor "C". It furthermore assumes that the approximation error at the boundary is small enough, as the implementation ignores boundary faces.

Definition at line 443 of file estimators.hpp.

Member Typedef Documentation

using mfem::KellyErrorEstimator::ElementCoefficientFunction = std::function<double(Mesh*, const int)>

Function type to compute the local coefficient hₑ of an element.

Definition at line 448 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 452 of file estimators.hpp.

Constructor & Destructor Documentation

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.

Parameters
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 53 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.

Parameters
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 69 of file estimators.cpp.

mfem::KellyErrorEstimator::~KellyErrorEstimator ( )

Definition at line 85 of file estimators.cpp.

Member Function Documentation

const Vector& mfem::KellyErrorEstimator::GetLocalErrors ( )
inlineoverridevirtual

Get a Vector with all element errors.

Implements mfem::ErrorEstimator.

Definition at line 540 of file estimators.hpp.

virtual double mfem::KellyErrorEstimator::GetTotalError ( ) const
inlineoverridevirtual

Return the total error from the last error estimate.

Note
This method is optional for derived classes to override and the base class implementation simply returns 0.

Reimplemented from mfem::ErrorEstimator.

Definition at line 552 of file estimators.hpp.

void mfem::KellyErrorEstimator::Reset ( )
inlineoverridevirtual

Reset the error estimator.

Implements mfem::ErrorEstimator.

Definition at line 550 of file estimators.hpp.

void mfem::KellyErrorEstimator::ResetCoefficientFunctions ( )

Change the coefficients back to default as described above.

Definition at line 93 of file estimators.cpp.

void mfem::KellyErrorEstimator::SetElementCoefficientFunction ( ElementCoefficientFunction  compute_element_coefficient_)
inline

Change the method to compute hₑ on a per-element basis.

Parameters
compute_element_coefficient_A function taking a mesh and an element index to compute the local hₑ for the element.

Definition at line 559 of file estimators.hpp.

void mfem::KellyErrorEstimator::SetFaceCoefficientFunction ( FaceCoefficientFunction  compute_face_coefficient_)
inline

Change the method to compute hₖ on a per-element basis.

Parameters
compute_face_coefficient_A function taking a mesh and a face index to compute the local hₖ for the face.

Definition at line 570 of file estimators.hpp.


The documentation for this class was generated from the following files: