16#ifndef dealii_error_estimator_h
17#define dealii_error_estimator_h
33template <
int dim,
int spacedim>
260template <
int dim,
int spacedim = dim>
339 template <
typename InputVector>
348 const InputVector & solution,
361 template <
typename InputVector>
369 const InputVector & solution,
391 template <
typename InputVector>
400 const std::vector<const InputVector *> &solutions,
413 template <
typename InputVector>
421 const std::vector<const InputVector *> &solutions,
435 template <
typename InputVector>
444 const InputVector & solution,
458 template <
typename InputVector>
466 const InputVector & solution,
480 template <
typename InputVector>
489 const std::vector<const InputVector *> &solutions,
503 template <
typename InputVector>
511 const std::vector<const InputVector *> &solutions,
524 "You provided a ComponentMask argument that is invalid. "
525 "Component masks need to be either default constructed "
526 "(in which case they indicate that every component is "
527 "selected) or need to have a length equal to the number "
528 "of vector components of the finite element in use "
529 "by the DoFHandler object. In the latter case, at "
530 "least one component needs to be selected.");
536 "If you do specify the argument for a (possibly "
537 "spatially variable) coefficient function for this function, "
538 "then it needs to refer to a coefficient that is either "
539 "scalar (has one vector component) or has as many vector "
540 "components as there are in the finite element used by "
541 "the DoFHandler argument.");
549 <<
"You provided a function map that for boundary indicator "
550 << arg1 <<
" specifies a function with " << arg2
551 <<
" vector components. However, the finite "
552 "element in use has "
554 <<
" components, and these two numbers need to match.");
561 <<
"The number of input vectors, " << arg1
562 <<
" needs to be equal to the number of output vectors, "
564 <<
". This is not the case in your call of this function.");
569 "You need to specify at least one solution vector as "
584template <
int spacedim>
627 template <
typename InputVector>
636 const InputVector & solution,
651 template <
typename InputVector>
659 const InputVector & solution,
683 template <
typename InputVector>
692 const std::vector<const InputVector *> &solutions,
707 template <
typename InputVector>
715 const std::vector<const InputVector *> &solutions,
730 template <
typename InputVector>
739 const InputVector & solution,
754 template <
typename InputVector>
762 const InputVector & solution,
777 template <
typename InputVector>
786 const std::vector<const InputVector *> &solutions,
801 template <
typename InputVector>
809 const std::vector<const InputVector *> &solutions,
824 "You provided a ComponentMask argument that is invalid. "
825 "Component masks need to be either default constructed "
826 "(in which case they indicate that every component is "
827 "selected) or need to have a length equal to the number "
828 "of vector components of the finite element in use "
829 "by the DoFHandler object. In the latter case, at "
830 "least one component needs to be selected.");
837 "If you do specify the argument for a (possibly "
838 "spatially variable) coefficient function for this function, "
839 "then it needs to refer to a coefficient that is either "
840 "scalar (has one vector component) or has as many vector "
841 "components as there are in the finite element used by "
842 "the DoFHandler argument.");
851 <<
"You provided a function map that for boundary indicator "
852 << arg1 <<
" specifies a function with " << arg2
853 <<
" vector components. However, the finite "
854 "element in use has "
856 <<
" components, and these two numbers need to match.");
864 <<
"The number of input vectors, " << arg1
865 <<
" needs to be equal to the number of output vectors, "
867 <<
". This is not the case in your call of this function.");
873 "You need to specify at least one solution vector as "
@ face_diameter_over_twice_max_degree
static void estimate(const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, typename InputVector::value_type > * > &neumann_bc, const std::vector< const InputVector * > &solutions, std::vector< Vector< float > * > &errors, const ComponentMask &component_mask=ComponentMask(), 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, typename InputVector::value_type > * > &neumann_bc, const std::vector< const InputVector * > &solutions, std::vector< Vector< float > * > &errors, const ComponentMask &component_mask=ComponentMask(), 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, typename InputVector::value_type > * > &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), 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, typename InputVector::value_type > * > &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), 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, typename InputVector::value_type > * > &neumann_bc, const std::vector< const InputVector * > &solutions, std::vector< Vector< float > * > &errors, const ComponentMask &component_mask=ComponentMask(), 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, typename InputVector::value_type > * > &neumann_bc, const std::vector< const InputVector * > &solutions, std::vector< Vector< float > * > &errors, const ComponentMask &component_mask=ComponentMask(), 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, typename InputVector::value_type > * > &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), 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, typename InputVector::value_type > * > &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), 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 .
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