16#ifndef dealii_error_estimator_h
17#define dealii_error_estimator_h
259template <
int dim,
int spacedim = dim>
338 template <
typename InputVector>
347 const InputVector & solution,
360 template <
typename InputVector>
368 const InputVector & solution,
390 template <
typename InputVector>
399 const std::vector<const InputVector *> &solutions,
412 template <
typename InputVector>
420 const std::vector<const InputVector *> &solutions,
434 template <
typename InputVector>
443 const InputVector & solution,
457 template <
typename InputVector>
465 const InputVector & solution,
479 template <
typename InputVector>
488 const std::vector<const InputVector *> &solutions,
502 template <
typename InputVector>
510 const std::vector<const InputVector *> &solutions,
523 "You provided a ComponentMask argument that is invalid. "
524 "Component masks need to be either default constructed "
525 "(in which case they indicate that every component is "
526 "selected) or need to have a length equal to the number "
527 "of vector components of the finite element in use "
528 "by the DoFHandler object. In the latter case, at "
529 "least one component needs to be selected.");
535 "If you do specify the argument for a (possibly "
536 "spatially variable) coefficient function for this function, "
537 "then it needs to refer to a coefficient that is either "
538 "scalar (has one vector component) or has as many vector "
539 "components as there are in the finite element used by "
540 "the DoFHandler argument.");
548 <<
"You provided a function map that for boundary indicator "
549 << arg1 <<
" specifies a function with " << arg2
550 <<
" vector components. However, the finite "
551 "element in use has "
553 <<
" components, and these two numbers need to match.");
560 <<
"The number of input vectors, " << arg1
561 <<
" needs to be equal to the number of output vectors, "
563 <<
". This is not the case in your call of this function.");
568 "You need to specify at least one solution vector as "
583template <
int spacedim>
624 template <
typename InputVector>
633 const InputVector & solution,
646 template <
typename InputVector>
654 const InputVector & solution,
676 template <
typename InputVector>
685 const std::vector<const InputVector *> &solutions,
698 template <
typename InputVector>
706 const std::vector<const InputVector *> &solutions,
720 template <
typename InputVector>
729 const InputVector & solution,
743 template <
typename InputVector>
751 const InputVector & solution,
765 template <
typename InputVector>
774 const std::vector<const InputVector *> &solutions,
788 template <
typename InputVector>
796 const std::vector<const InputVector *> &solutions,
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.");
821 "If you do specify the argument for a (possibly "
822 "spatially variable) coefficient function for this function, "
823 "then it needs to refer to a coefficient that is either "
824 "scalar (has one vector component) or has as many vector "
825 "components as there are in the finite element used by "
826 "the DoFHandler argument.");
834 <<
"You provided a function map that for boundary indicator "
835 << arg1 <<
" specifies a function with " << arg2
836 <<
" vector components. However, the finite "
837 "element in use has "
839 <<
" components, and these two numbers need to match.");
846 <<
"The number of input vectors, " << arg1
847 <<
" needs to be equal to the number of output vectors, "
849 <<
". This is not the case in your call of this function.");
854 "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_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