Go to the documentation of this file.
26 #include "libmesh/libmesh_common.h"
27 #include "libmesh/kelly_error_estimator.h"
28 #include "libmesh/error_vector.h"
29 #include "libmesh/fe_base.h"
30 #include "libmesh/libmesh_logging.h"
31 #include "libmesh/elem.h"
32 #include "libmesh/system.h"
33 #include "libmesh/dense_vector.h"
34 #include "libmesh/tensor_tools.h"
35 #include "libmesh/enum_error_estimator_type.h"
36 #include "libmesh/enum_norm_type.h"
62 for (
unsigned int v=0; v<
n_vars; v++)
68 FEBase * side_fe =
nullptr;
70 const std::set<unsigned char> & elem_dims =
73 for (
const auto &
dim : elem_dims)
95 FEBase * fe_fine =
nullptr;
98 FEBase * fe_coarse =
nullptr;
104 std::vector<std::vector<RealGradient>> dphi_coarse = fe_coarse->
get_dphi();
105 std::vector<std::vector<RealGradient>> dphi_fine = fe_fine->
get_dphi();
106 std::vector<Point> face_normals = fe_fine->
get_normals();
107 std::vector<Real> JxW_face = fe_fine->
get_JxW();
109 for (
unsigned int qp=0; qp != n_qp; ++qp)
119 const Number jump = (grad_fine - grad_coarse)*face_normals[qp];
123 error += JxW_face[qp] * jump2;
139 FEBase * fe_fine =
nullptr;
142 const std::string & var_name =
145 std::vector<std::vector<RealGradient>> dphi_fine = fe_fine->
get_dphi();
146 std::vector<Point> face_normals = fe_fine->
get_normals();
147 std::vector<Real> JxW_face = fe_fine->
get_JxW();
148 std::vector<Point> qface_point = fe_fine->
get_xyz();
157 qface_point[0], var_name).first)
168 for (
unsigned int qp=0; qp<n_qp; qp++)
171 const std::pair<bool,Real> flux_bc =
173 qface_point[qp], var_name);
177 libmesh_assert_equal_to (flux_bc.first,
true);
183 const Number jump = flux_bc.second - grad_fine*face_normals[qp];
193 error += JxW_face[qp]*jump2;
209 const std::string & var_name))
Manages consistently variables, degrees of freedom, and coefficient vectors.
std::unique_ptr< FEMContext > coarse_context
bool integrate_boundary_sides
A boolean flag, by default false, to be set to true if integrations with boundary_side_integration() ...
const std::vector< Point > & get_xyz() const
virtual bool boundary_side_integration() override
The function which calculates a normal derivative jump based error term on a boundary side.
SystemNorm error_norm
When estimating the error in a single system, the error_norm is used to control the scaling and norm ...
ErrorEstimatorType
Defines an enum for the different types of error estimators which are available.
const std::vector< Point > & get_normals() const
The libMesh namespace provides an interface to certain functionality in the library.
virtual unsigned short dim() const =0
unsigned int n_vars() const
Number of variables in solution.
unsigned int var
The variable number currently being evaluated.
This class forms the foundation from which generic finite elements may be derived.
Real fine_error
The fine and coarse error values to be set by each side_integration();.
const std::vector< Real > & get_JxW() const
Real weight(unsigned int var) const
const std::vector< std::vector< OutputGradient > > & get_dphi() const
virtual void internal_side_integration() override
The function which calculates a normal derivative jump based error term on an internal side.
std::unique_ptr< FEMContext > fine_context
Context objects for integrating on the fine and coarse elements sharing a face.
A Point defines a location in LIBMESH_DIM dimensional Real space.
virtual unsigned int n_quadrature_points() const =0
This abstract base class implements utility functions for error estimators which are based on integra...
virtual ErrorEstimatorType type() const override
virtual void init_context(FEMContext &c) override
An initialization function, for requesting specific data from the FE objects.
virtual Real hmax() const
This is the base class from which all geometric element types are derived.
const std::set< unsigned char > & elem_dimensions() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void attach_flux_bc_function(std::pair< bool, Real > fptr(const System &system, const Point &p, const std::string &var_name))
Register a user function to use in computing the flux BCs.
Number fptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
std::pair< bool, Real >(* _bc_function)(const System &system, const Point &p, const std::string &var_name)
Pointer to function that provides BC information.
KellyErrorEstimator()
Constructor.
This class provides all data required for a physics package (e.g.