15#ifndef dealii_error_estimator_h
16#define dealii_error_estimator_h
34template <
int dim,
int spacedim>
264template <
int dim,
int spacedim = dim>
343 template <
typename Number>
364 template <
typename Number>
393 template <
typename Number>
414 template <
typename Number>
435 template <
typename Number>
460 template <
typename Number>
461 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
462 "Use the version of this function which takes a hp::MappingCollection instead.")
466 const
hp::QCollection<dim - 1> &quadrature,
483 template <
typename Number>
504 template <
typename Number>
529 template <
typename Number>
530 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
531 "Use the version of this function which takes a hp::MappingCollection instead.")
535 const
hp::QCollection<dim - 1> &quadrature,
552 template <
typename Number>
572 "You provided a ComponentMask argument that is invalid. "
573 "Component masks need to be either default constructed "
574 "(in which case they indicate that every component is "
575 "selected) or need to have a length equal to the number "
576 "of vector components of the finite element in use "
577 "by the DoFHandler object. In the latter case, at "
578 "least one component needs to be selected.");
584 "If you do specify the argument for a (possibly "
585 "spatially variable) coefficient function for this function, "
586 "then it needs to refer to a coefficient that is either "
587 "scalar (has one vector component) or has as many vector "
588 "components as there are in the finite element used by "
589 "the DoFHandler argument.");
597 <<
"You provided a function map that for boundary indicator "
598 << arg1 <<
" specifies a function with " << arg2
599 <<
" vector components. However, the finite "
600 "element in use has "
602 <<
" components, and these two numbers need to match.");
609 <<
"The number of input vectors, " << arg1
610 <<
" needs to be equal to the number of output vectors, "
612 <<
". This is not the case in your call of this function.");
617 "You need to specify at least one solution vector as "
632template <
int spacedim>
675 template <
typename Number>
686 const Function<spacedim> *coefficient =
nullptr,
698 template <
typename Number>
701 const DoFHandler<1, spacedim> &dof,
702 const Quadrature<0> &quadrature,
705 const ReadVector<Number> &solution,
706 Vector<float> &error,
707 const ComponentMask &component_mask = {},
708 const Function<spacedim> *coefficients =
nullptr,
729 template <
typename Number>
732 const Mapping<1, spacedim> &
mapping,
733 const DoFHandler<1, spacedim> &dof,
734 const Quadrature<0> &quadrature,
737 const ArrayView<
const ReadVector<Number> *> &solutions,
738 ArrayView<Vector<float> *> &errors,
739 const ComponentMask &component_mask = {},
740 const Function<spacedim> *coefficients =
nullptr,
752 template <
typename Number>
755 const DoFHandler<1, spacedim> &dof,
756 const Quadrature<0> &quadrature,
759 const ArrayView<
const ReadVector<Number> *> &solutions,
760 ArrayView<Vector<float> *> &errors,
761 const ComponentMask &component_mask = {},
762 const Function<spacedim> *coefficients =
nullptr,
774 template <
typename Number>
777 const hp::MappingCollection<1, spacedim> &
mapping,
778 const DoFHandler<1, spacedim> &dof,
779 const hp::QCollection<0> &quadrature,
782 const ReadVector<Number> &solution,
783 Vector<float> &error,
784 const ComponentMask &component_mask = {},
785 const Function<spacedim> *coefficients =
nullptr,
800 template <
typename Number>
801 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
802 "Use the version of this function which takes a hp::MappingCollection instead.")
804 const Mapping<1, spacedim> &
mapping,
805 const DoFHandler<1, spacedim> &dof,
806 const hp::QCollection<0> &quadrature,
807 const std::map<types::boundary_id, const Function<spacedim, Number> *>
809 const ReadVector<Number> &solution,
810 Vector<
float> &error,
811 const ComponentMask &component_mask = {},
812 const Function<spacedim> *coefficients =
nullptr,
824 template <
typename Number>
827 const DoFHandler<1, spacedim> &dof,
828 const hp::QCollection<0> &quadrature,
831 const ReadVector<Number> &solution,
832 Vector<float> &error,
833 const ComponentMask &component_mask = {},
834 const Function<spacedim> *coefficients =
nullptr,
846 template <
typename Number>
849 const hp::MappingCollection<1, spacedim> &
mapping,
850 const DoFHandler<1, spacedim> &dof,
851 const hp::QCollection<0> &quadrature,
854 const ArrayView<
const ReadVector<Number> *> &solutions,
855 ArrayView<Vector<float> *> &errors,
856 const ComponentMask &component_mask = {},
857 const Function<spacedim> *coefficients =
nullptr,
872 template <
typename Number>
873 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
874 "Use the version of this function which takes a hp::MappingCollection instead.")
876 const Mapping<1, spacedim> &
mapping,
877 const DoFHandler<1, spacedim> &dof,
878 const hp::QCollection<0> &quadrature,
879 const std::map<types::boundary_id, const Function<spacedim, Number> *>
881 const ArrayView<const ReadVector<Number> *> &solutions,
882 ArrayView<Vector<
float> *> &errors,
883 const ComponentMask &component_mask = {},
884 const Function<spacedim> *coefficients =
nullptr,
896 template <
typename Number>
899 const DoFHandler<1, spacedim> &dof,
900 const hp::QCollection<0> &quadrature,
903 const ArrayView<
const ReadVector<Number> *> &solutions,
904 ArrayView<Vector<float> *> &errors,
905 const ComponentMask &component_mask = {},
906 const Function<spacedim> *coefficients =
nullptr,
918 "You provided a ComponentMask argument that is invalid. "
919 "Component masks need to be either default constructed "
920 "(in which case they indicate that every component is "
921 "selected) or need to have a length equal to the number "
922 "of vector components of the finite element in use "
923 "by the DoFHandler object. In the latter case, at "
924 "least one component needs to be selected.");
931 "If you do specify the argument for a (possibly "
932 "spatially variable) coefficient function for this function, "
933 "then it needs to refer to a coefficient that is either "
934 "scalar (has one vector component) or has as many vector "
935 "components as there are in the finite element used by "
936 "the DoFHandler argument.");
945 <<
"You provided a function map that for boundary indicator "
946 << arg1 <<
" specifies a function with " << arg2
947 <<
" vector components. However, the finite "
948 "element in use has "
950 <<
" components, and these two numbers need to match.");
958 <<
"The number of input vectors, " << arg1
959 <<
" needs to be equal to the number of output vectors, "
961 <<
". This is not the case in your call of this function.");
967 "You need to specify at least one solution vector as "
@ cell_diameter
Kelly error estimator with the factor .
@ face_diameter_over_twice_max_degree
@ cell_diameter_over_24
Kelly error estimator with the factor .
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 hp::MappingCollection< 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)
@ 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)
static void estimate(const hp::MappingCollection< 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)
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 & ExcInvalidBoundaryFunction(types::boundary_id arg1, int arg2, int arg3)
static ::ExceptionBase & ExcIncompatibleNumberOfElements(int arg1, int arg2)
static ::ExceptionBase & ExcNoSolutions()
#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 & ExcInvalidCoefficient()
static ::ExceptionBase & ExcNoSolutions()
static ::ExceptionBase & ExcIncompatibleNumberOfElements(int arg1, int arg2)
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
constexpr unsigned int invalid_unsigned_int
constexpr types::material_id invalid_material_id
constexpr types::subdomain_id invalid_subdomain_id
unsigned int subdomain_id