MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::ZienkiewiczZhuEstimator Class Reference

The ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure. More...

#include <estimators.hpp>

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

Public Member Functions

 ZienkiewiczZhuEstimator (BilinearFormIntegrator &integ, GridFunction &sol, FiniteElementSpace *flux_fes)
 Construct a new ZienkiewiczZhuEstimator object.
 
 ZienkiewiczZhuEstimator (BilinearFormIntegrator &integ, GridFunction &sol, FiniteElementSpace &flux_fes)
 Construct a new ZienkiewiczZhuEstimator object.
 
void SetWithCoeff (bool w_coeff=true)
 Consider the coefficient in BilinearFormIntegrator to calculate the fluxes for the error estimator.
 
void SetAnisotropic (bool aniso=true)
 Enable/disable anisotropic estimates. To enable this option, the BilinearFormIntegrator must support the 'd_energy' parameter in its ComputeFluxEnergy() method.
 
void SetFluxAveraging (int fa)
 Set the way the flux is averaged (smoothed) across elements.
 
virtual real_t GetTotalError () const override
 Return the total error from the last error estimate.
 
virtual const VectorGetLocalErrors () override
 Get a Vector with all element errors.
 
virtual const Array< int > & GetAnisotropicFlags () override
 Get an Array<int> with anisotropic flags for all mesh elements. Return an empty array when anisotropic estimates are not available or enabled.
 
virtual void Reset () override
 Reset the error estimator.
 
virtual ~ZienkiewiczZhuEstimator ()
 Destroy a ZienkiewiczZhuEstimator object. Destroys, if owned, the FiniteElementSpace, flux_space.
 
- Public Member Functions inherited from mfem::ErrorEstimator
virtual ~ErrorEstimator ()
 Destruct the error estimator.
 
- Public Member Functions inherited from mfem::AbstractErrorEstimator
virtual ~AbstractErrorEstimator ()
 

Protected Member Functions

bool MeshIsModified ()
 Check if the mesh of the solution was modified.
 
void ComputeEstimates ()
 Compute the element error estimates.
 

Protected Attributes

long current_sequence
 
Vector error_estimates
 
real_t total_error
 
bool anisotropic
 
Array< int > aniso_flags
 
int flux_averaging
 
BilinearFormIntegratorinteg
 
GridFunctionsolution
 
FiniteElementSpaceflux_space
 Ownership based on own_flux_fes. Its Update() method is called automatically by this class when needed.
 
bool with_coeff
 
bool own_flux_fes
 Ownership flag for flux_space.
 

Detailed Description

The ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure.

[1] Zienkiewicz, O.C. and Zhu, J.Z., The superconvergent patch recovery and a posteriori error estimates. Part 1: The recovery technique. Int. J. Num. Meth. Engng. 33, 1331-1364 (1992).

[2] Zienkiewicz, O.C. and Zhu, J.Z., The superconvergent patch recovery and a posteriori error estimates. Part 2: Error estimates and adaptivity. Int. J. Num. Meth. Engng. 33, 1365-1382 (1992).

The required BilinearFormIntegrator must implement the methods ComputeElementFlux() and ComputeFluxEnergy().

Definition at line 88 of file estimators.hpp.

Constructor & Destructor Documentation

◆ ZienkiewiczZhuEstimator() [1/2]

mfem::ZienkiewiczZhuEstimator::ZienkiewiczZhuEstimator ( BilinearFormIntegrator & integ,
GridFunction & sol,
FiniteElementSpace * flux_fes )
inline

Construct a new ZienkiewiczZhuEstimator object.

Parameters
integThis BilinearFormIntegrator must implement the methods ComputeElementFlux() and ComputeFluxEnergy().
solThe solution field whose error is to be estimated.
flux_fesThe ZienkiewiczZhuEstimator assumes ownership of this FiniteElementSpace and will call its Update() method when needed.

Definition at line 125 of file estimators.hpp.

◆ ZienkiewiczZhuEstimator() [2/2]

mfem::ZienkiewiczZhuEstimator::ZienkiewiczZhuEstimator ( BilinearFormIntegrator & integ,
GridFunction & sol,
FiniteElementSpace & flux_fes )
inline

Construct a new ZienkiewiczZhuEstimator object.

Parameters
integThis BilinearFormIntegrator must implement the methods ComputeElementFlux() and ComputeFluxEnergy().
solThe solution field whose error is to be estimated.
flux_fesThe ZienkiewiczZhuEstimator does NOT assume ownership of this FiniteElementSpace; will call its Update() method when needed.

Definition at line 145 of file estimators.hpp.

◆ ~ZienkiewiczZhuEstimator()

virtual mfem::ZienkiewiczZhuEstimator::~ZienkiewiczZhuEstimator ( )
inlinevirtual

Destroy a ZienkiewiczZhuEstimator object. Destroys, if owned, the FiniteElementSpace, flux_space.

Definition at line 198 of file estimators.hpp.

Member Function Documentation

◆ ComputeEstimates()

void mfem::ZienkiewiczZhuEstimator::ComputeEstimates ( )
protected

Compute the element error estimates.

Definition at line 17 of file estimators.cpp.

◆ GetAnisotropicFlags()

virtual const Array< int > & mfem::ZienkiewiczZhuEstimator::GetAnisotropicFlags ( )
inlineoverridevirtual

Get an Array<int> with anisotropic flags for all mesh elements. Return an empty array when anisotropic estimates are not available or enabled.

Implements mfem::AnisotropicErrorEstimator.

Definition at line 187 of file estimators.hpp.

◆ GetLocalErrors()

virtual const Vector & mfem::ZienkiewiczZhuEstimator::GetLocalErrors ( )
inlineoverridevirtual

Get a Vector with all element errors.

Implements mfem::ErrorEstimator.

Definition at line 178 of file estimators.hpp.

◆ GetTotalError()

virtual real_t mfem::ZienkiewiczZhuEstimator::GetTotalError ( ) const
inlineoverridevirtual

Return the total error from the last error estimate.

Reimplemented from mfem::ErrorEstimator.

Definition at line 175 of file estimators.hpp.

◆ MeshIsModified()

bool mfem::ZienkiewiczZhuEstimator::MeshIsModified ( )
inlineprotected

Check if the mesh of the solution was modified.

Definition at line 107 of file estimators.hpp.

◆ Reset()

virtual void mfem::ZienkiewiczZhuEstimator::Reset ( )
inlineoverridevirtual

Reset the error estimator.

Implements mfem::ErrorEstimator.

Definition at line 194 of file estimators.hpp.

◆ SetAnisotropic()

void mfem::ZienkiewiczZhuEstimator::SetAnisotropic ( bool aniso = true)
inline

Enable/disable anisotropic estimates. To enable this option, the BilinearFormIntegrator must support the 'd_energy' parameter in its ComputeFluxEnergy() method.

Definition at line 165 of file estimators.hpp.

◆ SetFluxAveraging()

void mfem::ZienkiewiczZhuEstimator::SetFluxAveraging ( int fa)
inline

Set the way the flux is averaged (smoothed) across elements.

When fa is zero (default), averaging is performed across interfaces between different mesh attributes. When fa is non-zero, the flux is not averaged across interfaces between different mesh attributes.

Definition at line 172 of file estimators.hpp.

◆ SetWithCoeff()

void mfem::ZienkiewiczZhuEstimator::SetWithCoeff ( bool w_coeff = true)
inline

Consider the coefficient in BilinearFormIntegrator to calculate the fluxes for the error estimator.

Definition at line 160 of file estimators.hpp.

Member Data Documentation

◆ aniso_flags

Array<int> mfem::ZienkiewiczZhuEstimator::aniso_flags
protected

Definition at line 95 of file estimators.hpp.

◆ anisotropic

bool mfem::ZienkiewiczZhuEstimator::anisotropic
protected

Definition at line 94 of file estimators.hpp.

◆ current_sequence

long mfem::ZienkiewiczZhuEstimator::current_sequence
protected

Definition at line 91 of file estimators.hpp.

◆ error_estimates

Vector mfem::ZienkiewiczZhuEstimator::error_estimates
protected

Definition at line 92 of file estimators.hpp.

◆ flux_averaging

int mfem::ZienkiewiczZhuEstimator::flux_averaging
protected

Definition at line 96 of file estimators.hpp.

◆ flux_space

FiniteElementSpace* mfem::ZienkiewiczZhuEstimator::flux_space
protected

Ownership based on own_flux_fes. Its Update() method is called automatically by this class when needed.

Definition at line 101 of file estimators.hpp.

◆ integ

BilinearFormIntegrator& mfem::ZienkiewiczZhuEstimator::integ
protected

Definition at line 98 of file estimators.hpp.

◆ own_flux_fes

bool mfem::ZienkiewiczZhuEstimator::own_flux_fes
protected

Ownership flag for flux_space.

Definition at line 104 of file estimators.hpp.

◆ solution

GridFunction& mfem::ZienkiewiczZhuEstimator::solution
protected

Definition at line 99 of file estimators.hpp.

◆ total_error

real_t mfem::ZienkiewiczZhuEstimator::total_error
protected

Definition at line 93 of file estimators.hpp.

◆ with_coeff

bool mfem::ZienkiewiczZhuEstimator::with_coeff
protected

Definition at line 103 of file estimators.hpp.


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