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

The LSZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure [1,2] using face-based patches [3]. More...

#include <estimators.hpp>

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

Public Member Functions

 LSZienkiewiczZhuEstimator (BilinearFormIntegrator &integ, GridFunction &sol)
 Construct a new LSZienkiewiczZhuEstimator object.
 
void SetWithCoeff (bool w_coeff=true)
 Consider the coefficient in BilinearFormIntegrator to calculate the fluxes for the error estimator.
 
void DisableReconstructionAcrossSubdomains ()
 Disable reconstructing the flux in patches spanning different subdomains.
 
void SetTichonovRegularization (real_t tcoeff=1.0e-8)
 Solve a Tichonov-regularized least-squares problem for the reconstructed fluxes. This is especially helpful for when not using tensor product elements, which typically require fewer integration points and, therefore, may lead to an ill-conditioned linear system.
 
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 void Reset () override
 Reset the error estimator.
 
virtual ~LSZienkiewiczZhuEstimator ()
 
- 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 subdomain_reconstruction = true
 
real_t tichonov_coeff
 
BilinearFormIntegratorinteg
 
GridFunctionsolution
 
bool with_coeff
 

Detailed Description

The LSZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure [1,2] using face-based patches [3].

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

[3] Bartels, S. and Carstensen, C., Each averaging technique yields reliable a posteriori error control in FEM on unstructured grids. Part II: Higher order FEM. Math. Comp. 71(239), 971-994 (2002)

The required BilinearFormIntegrator must implement the method ComputeElementFlux().

Note
  • The present implementation ignores all single-element patches corresponding to boundary faces. This is appropriate for Dirichlet boundaries, but suboptimal for Neumann boundaries. Reference 3 shows that a constrained least-squares problem, where the reconstructed flux is constrained by the Neumann boundary data, is appropriate to handle this case. NOTE THAT THIS CONSTRAINED LS PROBLEM IS NOT YET IMPLEMENTED, so it is possible that the local error estimates for elements on a Neumann boundary are suboptimal.
  • The global polynomial basis used for the flux reconstruction is, by default, aligned with the physical Cartesian axis. For patches with 2D elements, this has been improved on so that the basis is aligned with the physical patch orientation. Reorientation of the flux reconstruction basis is helpful to maintain symmetry in the refinement pattern and could be extended to 3D.
  • This estimator is ONLY implemented IN SERIAL.
  • Anisotropic refinement is NOT YET SUPPORTED.

Definition at line 241 of file estimators.hpp.

Constructor & Destructor Documentation

◆ LSZienkiewiczZhuEstimator()

mfem::LSZienkiewiczZhuEstimator::LSZienkiewiczZhuEstimator ( BilinearFormIntegrator & integ,
GridFunction & sol )
inline

Construct a new LSZienkiewiczZhuEstimator object.

Parameters
integThis BilinearFormIntegrator must implement only the method ComputeElementFlux().
solThe solution field whose error is to be estimated.

Definition at line 271 of file estimators.hpp.

◆ ~LSZienkiewiczZhuEstimator()

virtual mfem::LSZienkiewiczZhuEstimator::~LSZienkiewiczZhuEstimator ( )
inlinevirtual

Definition at line 313 of file estimators.hpp.

Member Function Documentation

◆ ComputeEstimates()

void mfem::LSZienkiewiczZhuEstimator::ComputeEstimates ( )
protected

Compute the element error estimates.

Definition at line 33 of file estimators.cpp.

◆ DisableReconstructionAcrossSubdomains()

void mfem::LSZienkiewiczZhuEstimator::DisableReconstructionAcrossSubdomains ( )
inline

Disable reconstructing the flux in patches spanning different subdomains.

Definition at line 287 of file estimators.hpp.

◆ GetLocalErrors()

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

Get a Vector with all element errors.

Implements mfem::ErrorEstimator.

Definition at line 304 of file estimators.hpp.

◆ GetTotalError()

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

Return the total error from the last error estimate.

Reimplemented from mfem::ErrorEstimator.

Definition at line 301 of file estimators.hpp.

◆ MeshIsModified()

bool mfem::LSZienkiewiczZhuEstimator::MeshIsModified ( )
inlineprotected

Check if the mesh of the solution was modified.

Definition at line 255 of file estimators.hpp.

◆ Reset()

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

Reset the error estimator.

Implements mfem::ErrorEstimator.

Definition at line 311 of file estimators.hpp.

◆ SetTichonovRegularization()

void mfem::LSZienkiewiczZhuEstimator::SetTichonovRegularization ( real_t tcoeff = 1.0e-8)
inline

Solve a Tichonov-regularized least-squares problem for the reconstructed fluxes. This is especially helpful for when not using tensor product elements, which typically require fewer integration points and, therefore, may lead to an ill-conditioned linear system.

Definition at line 294 of file estimators.hpp.

◆ SetWithCoeff()

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

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

Definition at line 283 of file estimators.hpp.

Member Data Documentation

◆ current_sequence

long mfem::LSZienkiewiczZhuEstimator::current_sequence
protected

Definition at line 244 of file estimators.hpp.

◆ error_estimates

Vector mfem::LSZienkiewiczZhuEstimator::error_estimates
protected

Definition at line 245 of file estimators.hpp.

◆ integ

BilinearFormIntegrator& mfem::LSZienkiewiczZhuEstimator::integ
protected

Definition at line 250 of file estimators.hpp.

◆ solution

GridFunction& mfem::LSZienkiewiczZhuEstimator::solution
protected

Definition at line 251 of file estimators.hpp.

◆ subdomain_reconstruction

bool mfem::LSZienkiewiczZhuEstimator::subdomain_reconstruction = true
protected

Definition at line 247 of file estimators.hpp.

◆ tichonov_coeff

real_t mfem::LSZienkiewiczZhuEstimator::tichonov_coeff
protected

Definition at line 248 of file estimators.hpp.

◆ total_error

real_t mfem::LSZienkiewiczZhuEstimator::total_error
protected

Definition at line 246 of file estimators.hpp.

◆ with_coeff

bool mfem::LSZienkiewiczZhuEstimator::with_coeff
protected

Definition at line 252 of file estimators.hpp.


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