15#ifndef dealii_error_estimator_h
16#define dealii_error_estimator_h
34template <
int dim,
int spacedim>
261template <
int dim,
int spacedim = dim>
340 template <
typename Number>
361 template <
typename Number>
390 template <
typename Number>
411 template <
typename Number>
432 template <
typename Number>
454 template <
typename Number>
475 template <
typename Number>
497 template <
typename Number>
517 "You provided a ComponentMask argument that is invalid. "
518 "Component masks need to be either default constructed "
519 "(in which case they indicate that every component is "
520 "selected) or need to have a length equal to the number "
521 "of vector components of the finite element in use "
522 "by the DoFHandler object. In the latter case, at "
523 "least one component needs to be selected.");
529 "If you do specify the argument for a (possibly "
530 "spatially variable) coefficient function for this function, "
531 "then it needs to refer to a coefficient that is either "
532 "scalar (has one vector component) or has as many vector "
533 "components as there are in the finite element used by "
534 "the DoFHandler argument.");
542 <<
"You provided a function map that for boundary indicator "
543 << arg1 <<
" specifies a function with " << arg2
544 <<
" vector components. However, the finite "
545 "element in use has "
547 <<
" components, and these two numbers need to match.");
554 <<
"The number of input vectors, " << arg1
555 <<
" needs to be equal to the number of output vectors, "
557 <<
". This is not the case in your call of this function.");
562 "You need to specify at least one solution vector as "
577template <
int spacedim>
620 template <
typename Number>
643 template <
typename Number>
674 template <
typename Number>
697 template <
typename Number>
719 template <
typename Number>
742 template <
typename Number>
764 template <
typename Number>
787 template <
typename Number>
809 "You provided a ComponentMask argument that is invalid. "
810 "Component masks need to be either default constructed "
811 "(in which case they indicate that every component is "
812 "selected) or need to have a length equal to the number "
813 "of vector components of the finite element in use "
814 "by the DoFHandler object. In the latter case, at "
815 "least one component needs to be selected.");
822 "If you do specify the argument for a (possibly "
823 "spatially variable) coefficient function for this function, "
824 "then it needs to refer to a coefficient that is either "
825 "scalar (has one vector component) or has as many vector "
826 "components as there are in the finite element used by "
827 "the DoFHandler argument.");
836 <<
"You provided a function map that for boundary indicator "
837 << arg1 <<
" specifies a function with " << arg2
838 <<
" vector components. However, the finite "
839 "element in use has "
841 <<
" components, and these two numbers need to match.");
849 <<
"The number of input vectors, " << arg1
850 <<
" needs to be equal to the number of output vectors, "
852 <<
". This is not the case in your call of this function.");
858 "You need to specify at least one solution vector as "
@ face_diameter_over_twice_max_degree
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ReadVector< Number > &solution, Vector< float > &error, const ComponentMask &component_mask={}, const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ReadVector< Number > &solution, Vector< float > &error, const ComponentMask &component_mask={}, const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
static void estimate(const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ArrayView< const ReadVector< Number > * > &solutions, ArrayView< Vector< float > * > &errors, const ComponentMask &component_mask={}, const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ArrayView< const ReadVector< Number > * > &solutions, ArrayView< Vector< float > * > &errors, const ComponentMask &component_mask={}, const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
static void estimate(const DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ArrayView< const ReadVector< Number > * > &solutions, ArrayView< Vector< float > * > &errors, const ComponentMask &component_mask={}, const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ArrayView< const ReadVector< Number > * > &solutions, ArrayView< Vector< float > * > &errors, const ComponentMask &component_mask={}, const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
@ face_diameter_over_twice_max_degree
@ cell_diameter_over_24
Kelly error estimator with the factor .
@ cell_diameter
Kelly error estimator with the factor .
static void estimate(const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ReadVector< Number > &solution, Vector< float > &error, const ComponentMask &component_mask={}, const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
static void estimate(const DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ReadVector< Number > &solution, Vector< float > &error, const ComponentMask &component_mask={}, const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
Abstract base class for mapping classes.
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_CXX20_REQUIRES(condition)
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcInvalidBoundaryFunction(types::boundary_id arg1, int arg2, int arg3)
static ::ExceptionBase & ExcIncompatibleNumberOfElements(int arg1, int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInvalidComponentMask()
#define DeclException3(Exception3, type1, type2, type3, outsequence)
static ::ExceptionBase & ExcInvalidCoefficient()
static ::ExceptionBase & ExcNoSolutions()
const types::material_id invalid_material_id
const types::subdomain_id invalid_subdomain_id
static const unsigned int invalid_unsigned_int
unsigned int subdomain_id