20 #ifndef LIBMESH_JUMP_ERROR_ESTIMATOR_H 21 #define LIBMESH_JUMP_ERROR_ESTIMATOR_H 24 #include "libmesh/dense_vector.h" 25 #include "libmesh/error_estimator.h" 26 #include "libmesh/fem_context.h" 89 bool estimate_parent_error =
false)
override;
174 #endif // LIBMESH_JUMP_ERROR_ESTIMATOR_H JumpErrorEstimator & operator=(const JumpErrorEstimator &)=delete
std::unique_ptr< FEMContext > fine_context
Context objects for integrating on the fine and coarse elements sharing a face.
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
This abstract base class implements utility functions for error estimators which are based on integra...
JumpErrorEstimator()
Constructor.
virtual void internal_side_integration()=0
The function, to be implemented by derived classes, which calculates an error term on an internal sid...
The libMesh namespace provides an interface to certain functionality in the library.
void reinit_sides()
A utility function to reinit the finite element data on elements sharing a side.
virtual bool boundary_side_integration()
The function, to be implemented by derived classes, which calculates an error term on a boundary side...
Manages consistently variables, degrees of freedom, and coefficient vectors.
This class holds functions that will estimate the error in a finite element solution on a given mesh...
This class provides all data required for a physics package (e.g.
std::unique_ptr< FEMContext > coarse_context
Real fine_error
The fine and coarse error values to be set by each side_integration();.
virtual ~JumpErrorEstimator()=default
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
float coarse_n_flux_faces_increment()
A utility function to correctly increase n_flux_faces for the coarse element.
bool use_unweighted_quadrature_rules
This boolean flag allows you to use "unweighted" quadrature rules (sized to exactly integrate unweigh...
unsigned int var
The variable number currently being evaluated.
bool scale_by_n_flux_faces
This boolean flag allows you to scale the error indicator result for each element by the number of "f...
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false) override
This function uses the derived class's jump error estimate formula to estimate the error on each cell...
bool integrate_boundary_sides
A boolean flag, by default false, to be set to true if integrations with boundary_side_integration() ...
virtual void init_context(FEMContext &c)
An initialization function, to give derived classes a chance to request specific data from the FE obj...