MFEM v4.7.0
Finite element discretization library
|
The LSZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure [1,2] using face-based patches [3]. More...
#include <estimators.hpp>
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 Vector & | GetLocalErrors () 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 |
BilinearFormIntegrator & | integ |
GridFunction & | solution |
bool | with_coeff |
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().
Definition at line 241 of file estimators.hpp.
|
inline |
Construct a new LSZienkiewiczZhuEstimator object.
integ | This BilinearFormIntegrator must implement only the method ComputeElementFlux(). |
sol | The solution field whose error is to be estimated. |
Definition at line 271 of file estimators.hpp.
|
inlinevirtual |
Definition at line 313 of file estimators.hpp.
|
protected |
Compute the element error estimates.
Definition at line 33 of file estimators.cpp.
|
inline |
Disable reconstructing the flux in patches spanning different subdomains.
Definition at line 287 of file estimators.hpp.
|
inlineoverridevirtual |
Get a Vector with all element errors.
Implements mfem::ErrorEstimator.
Definition at line 304 of file estimators.hpp.
|
inlineoverridevirtual |
Return the total error from the last error estimate.
Reimplemented from mfem::ErrorEstimator.
Definition at line 301 of file estimators.hpp.
|
inlineprotected |
Check if the mesh of the solution was modified.
Definition at line 255 of file estimators.hpp.
|
inlineoverridevirtual |
Reset the error estimator.
Implements mfem::ErrorEstimator.
Definition at line 311 of file estimators.hpp.
|
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.
|
inline |
Consider the coefficient in BilinearFormIntegrator to calculate the fluxes for the error estimator.
Definition at line 283 of file estimators.hpp.
|
protected |
Definition at line 244 of file estimators.hpp.
|
protected |
Definition at line 245 of file estimators.hpp.
|
protected |
Definition at line 250 of file estimators.hpp.
|
protected |
Definition at line 251 of file estimators.hpp.
|
protected |
Definition at line 247 of file estimators.hpp.
|
protected |
Definition at line 248 of file estimators.hpp.
|
protected |
Definition at line 246 of file estimators.hpp.
|
protected |
Definition at line 252 of file estimators.hpp.