12 #ifndef MFEM_ERROR_ESTIMATORS
13 #define MFEM_ERROR_ESTIMATORS
15 #include "../config/config.hpp"
16 #include "../linalg/vector.hpp"
46 virtual void Reset() = 0;
228 this->integ = &
integ;
305 #endif // MFEM_USE_MPI
309 #endif // MFEM_ERROR_ESTIMATORS
virtual ~L2ZienkiewiczZhuEstimator()
Destroy a L2ZienkiewiczZhuEstimator object. Destroys, if owned, the FiniteElementSpace, flux_space.
void SetAnisotropic(bool aniso=true)
Enable/disable anisotropic estimates. To enable this option, the BilinearFormIntegrator must support ...
Class for grid function - Vector with associated FE space.
L2ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, ParGridFunction &sol, ParFiniteElementSpace *flux_fes, ParFiniteElementSpace *smooth_flux_fes)
Construct a new L2ZienkiewiczZhuEstimator object.
void ComputeEstimates()
Compute the element error estimates.
ParFiniteElementSpace * flux_space
Ownership based on the flag own_flux_fes. Its Update() method is called automatically by this class w...
bool MeshIsModified()
Check if the mesh of the solution was modified.
virtual ~AbstractErrorEstimator()
virtual const Vector & GetLocalErrors()
Get a Vector with all element errors.
void SetLocalErrorNormP(int p)
Set the exponent, p, of the Lp norm used for computing the local element errors. Default value is 1...
Abstract parallel finite element space.
void SetWithCoeff(bool w_coeff=true)
Consider the coefficient in BilinearFormIntegrator to calculate the fluxes for the error estimator...
virtual void Reset()
Reset the error estimator.
The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure where t...
virtual ~ErrorEstimator()
void ComputeEstimates()
Compute the element error estimates.
ParGridFunction * solution
Not owned.
void Init(BilinearFormIntegrator &integ, ParGridFunction &sol, ParFiniteElementSpace *flux_fes, ParFiniteElementSpace *smooth_flux_fes)
bool own_flux_fes
Ownership flag for flux_space and smooth_flux_space.
virtual const Vector & GetLocalErrors()=0
Get a Vector with all element errors.
Mesh * GetMesh() const
Returns the mesh.
Base class for all error estimators.
GridFunction * solution
Not owned.
double GetTotalError() const
Return the total error from the last error estimate.
int local_norm_p
Local L_p norm to use, default is 1.
void SetFluxAveraging(int fa)
Set the way the flux is averaged (smoothed) across elements.
Base class for all element based error estimators.
L2ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, ParGridFunction &sol, ParFiniteElementSpace &flux_fes, ParFiniteElementSpace &smooth_flux_fes)
Construct a new L2ZienkiewiczZhuEstimator object.
FiniteElementSpace * FESpace()
virtual const Array< int > & GetAnisotropicFlags()=0
Get an Array<int> with anisotropic flags for all mesh elements.
virtual void Reset()=0
Force recomputation of the estimates on the next call to GetLocalErrors.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
ParFiniteElementSpace * smooth_flux_space
Ownership based on the flag own_flux_fes. Its Update() method is called automatically by this class w...
double GetTotalError() const
Return the total error from the last error estimate.
virtual void Reset()
Reset the error estimator.
FiniteElementSpace * flux_space
Ownership based on own_flux_fes. Its Update() method is called automatically by this class when neede...
virtual const Array< int > & GetAnisotropicFlags()
Get an Array<int> with anisotropic flags for all mesh elements. Return an empty array when anisotropi...
The AnisotropicErrorEstimator class is the base class for all error estimators that compute one non-n...
BilinearFormIntegrator * integ
Not owned.
ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, GridFunction &sol, FiniteElementSpace *flux_fes)
Construct a new ZienkiewiczZhuEstimator object.
The ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure...
virtual const Vector & GetLocalErrors()
Get a Vector with all element errors.
bool MeshIsModified()
Check if the mesh of the solution was modified.
virtual ~ZienkiewiczZhuEstimator()
Destroy a ZienkiewiczZhuEstimator object. Destroys, if owned, the FiniteElementSpace, flux_space.
ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, GridFunction &sol, FiniteElementSpace &flux_fes)
Construct a new ZienkiewiczZhuEstimator object.
Class for parallel grid function.
BilinearFormIntegrator * integ
Not owned.
bool own_flux_fes
Ownership flag for flux_space.