MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
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. More...
 
 ZienkiewiczZhuEstimator (BilinearFormIntegrator &integ, GridFunction &sol, FiniteElementSpace &flux_fes)
 Construct a new ZienkiewiczZhuEstimator object. More...
 
void SetWithCoeff (bool w_coeff=true)
 Consider the coefficient in BilinearFormIntegrator to calculate the fluxes for the error estimator. More...
 
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. More...
 
void SetFluxAveraging (int fa)
 Set the way the flux is averaged (smoothed) across elements. More...
 
double GetTotalError () const
 Return the total error from the last error estimate. More...
 
virtual const VectorGetLocalErrors ()
 Get a Vector with all element errors. More...
 
virtual const Array< int > & GetAnisotropicFlags ()
 Get an Array<int> with anisotropic flags for all mesh elements. Return an empty array when anisotropic estimates are not available or enabled. More...
 
virtual void Reset ()
 Reset the error estimator. More...
 
virtual ~ZienkiewiczZhuEstimator ()
 Destroy a ZienkiewiczZhuEstimator object. Destroys, if owned, the FiniteElementSpace, flux_space. More...
 
- Public Member Functions inherited from mfem::ErrorEstimator
virtual ~ErrorEstimator ()
 Destruct the error estimator. More...
 
- Public Member Functions inherited from mfem::AbstractErrorEstimator
virtual ~AbstractErrorEstimator ()
 

Protected Member Functions

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

Protected Attributes

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

Detailed Description

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

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).

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 81 of file estimators.hpp.

Constructor & Destructor Documentation

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 118 of file estimators.hpp.

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 138 of file estimators.hpp.

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

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

Definition at line 192 of file estimators.hpp.

Member Function Documentation

void mfem::ZienkiewiczZhuEstimator::ComputeEstimates ( )
protected

Compute the element error estimates.

Definition at line 17 of file estimators.cpp.

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

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 181 of file estimators.hpp.

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

Get a Vector with all element errors.

Implements mfem::ErrorEstimator.

Definition at line 172 of file estimators.hpp.

double mfem::ZienkiewiczZhuEstimator::GetTotalError ( ) const
inline

Return the total error from the last error estimate.

Definition at line 169 of file estimators.hpp.

bool mfem::ZienkiewiczZhuEstimator::MeshIsModified ( )
inlineprotected

Check if the mesh of the solution was modified.

Definition at line 100 of file estimators.hpp.

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

Reset the error estimator.

Implements mfem::ErrorEstimator.

Definition at line 188 of file estimators.hpp.

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 158 of file estimators.hpp.

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 globally. When fa is non-zero, the flux averaging is performed locally for each mesh attribute, i.e. the flux is not averaged across interfaces between different mesh attributes.

Definition at line 166 of file estimators.hpp.

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 153 of file estimators.hpp.

Member Data Documentation

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

Definition at line 88 of file estimators.hpp.

bool mfem::ZienkiewiczZhuEstimator::anisotropic
protected

Definition at line 87 of file estimators.hpp.

long mfem::ZienkiewiczZhuEstimator::current_sequence
protected

Definition at line 84 of file estimators.hpp.

Vector mfem::ZienkiewiczZhuEstimator::error_estimates
protected

Definition at line 85 of file estimators.hpp.

int mfem::ZienkiewiczZhuEstimator::flux_averaging
protected

Definition at line 89 of file estimators.hpp.

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 94 of file estimators.hpp.

BilinearFormIntegrator* mfem::ZienkiewiczZhuEstimator::integ
protected

Not owned.

Definition at line 91 of file estimators.hpp.

bool mfem::ZienkiewiczZhuEstimator::own_flux_fes
protected

Ownership flag for flux_space.

Definition at line 97 of file estimators.hpp.

GridFunction* mfem::ZienkiewiczZhuEstimator::solution
protected

Not owned.

Definition at line 92 of file estimators.hpp.

double mfem::ZienkiewiczZhuEstimator::total_error
protected

Definition at line 86 of file estimators.hpp.

bool mfem::ZienkiewiczZhuEstimator::with_coeff
protected

Definition at line 96 of file estimators.hpp.


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