Reference documentation for deal.II version 9.2.0
|
#include <deal.II/numerics/error_estimator.h>
Public Types | |
enum | Strategy { cell_diameter_over_24 = 0, face_diameter_over_twice_max_degree, cell_diameter } |
Static Public Member Functions | |
template<typename InputVector , typename DoFHandlerType > | |
static void | estimate (const Mapping< 1, spacedim > &mapping, const DoFHandlerType &dof, const Quadrature< 0 > &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 > *coefficient=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) |
template<typename InputVector , typename DoFHandlerType > | |
static void | estimate (const DoFHandlerType &dof, const Quadrature< 0 > &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) |
template<typename InputVector , typename DoFHandlerType > | |
static void | estimate (const Mapping< 1, spacedim > &mapping, const DoFHandlerType &dof, const Quadrature< 0 > &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) |
template<typename InputVector , typename DoFHandlerType > | |
static void | estimate (const DoFHandlerType &dof, const Quadrature< 0 > &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) |
template<typename InputVector , typename DoFHandlerType > | |
static void | estimate (const Mapping< 1, spacedim > &mapping, const DoFHandlerType &dof, const hp::QCollection< 0 > &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) |
template<typename InputVector , typename DoFHandlerType > | |
static void | estimate (const DoFHandlerType &dof, const hp::QCollection< 0 > &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) |
template<typename InputVector , typename DoFHandlerType > | |
static void | estimate (const Mapping< 1, spacedim > &mapping, const DoFHandlerType &dof, const hp::QCollection< 0 > &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) |
template<typename InputVector , typename DoFHandlerType > | |
static void | estimate (const DoFHandlerType &dof, const hp::QCollection< 0 > &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 ::ExceptionBase & | ExcInvalidComponentMask () |
static ::ExceptionBase & | ExcInvalidCoefficient () |
static ::ExceptionBase & | ExcInvalidBoundaryFunction (types::boundary_id arg1, int arg2, int arg3) |
static ::ExceptionBase & | ExcIncompatibleNumberOfElements (int arg1, int arg2) |
static ::ExceptionBase & | ExcNoSolutions () |
This is a specialization of the general template for 1d. The implementation is sufficiently different for 1d to justify this specialization. The basic difference between 1d and all other space dimensions is that in 1d, there are no faces of cells, just the vertices between line segments, so we have to compute the jump terms differently. However, this class offers exactly the same public functions as the general template, so that a user will not see any difference.
Definition at line 588 of file error_estimator.h.
enum KellyErrorEstimator< 1, spacedim >::Strategy |
The enum type given to the class functions to decide on the scaling factors of the face integrals.
Definition at line 595 of file error_estimator.h.
|
static |
Implementation of the error estimator described above. You may give a coefficient, but there is a default value which denotes the constant coefficient with value one. The coefficient function may either be a scalar one, in which case it is used for all components of the finite element, or a vector-valued one with as many components as there are in the finite element; in the latter case, each component is weighted by the respective component in the coefficient.
You might give a list of components you want to evaluate, in case the finite element used by the DoFHandler object is vector-valued. You then have to set those entries to true in the bit-vector component_mask
for which the respective component is to be used in the error estimator. The default is to use all components, which is done by either providing a bit-vector with all-set entries, or an empty bit-vector. All the other parameters are as in the general case used for 2d and higher.
The parameter n_threads is no longer used and will be ignored. Multithreading is not presently implemented for 1d, but we retain the respective parameter for compatibility with the function signature in the general case.
Definition at line 62 of file error_estimator_1d.cc.
|
static |
Call the estimate
function, see above, with mapping=MappingQGeneric1<1>()
.
Definition at line 100 of file error_estimator_1d.cc.
|
static |
Same function as above, but accepts more than one solution vectors and returns one error vector for each solution vector. For the reason of existence of this function, see the general documentation of this class.
Since we do not want to force the user of this function to copy around their solution vectors, the vector of solution vectors takes pointers to the solutions, rather than being a vector of vectors. This makes it simpler to have the solution vectors somewhere in memory, rather than to have them collected somewhere special. (Note that it is not possible to construct of vector of references, so we had to use a vector of pointers.)
Definition at line 297 of file error_estimator_1d.cc.
|
static |
Call the estimate
function, see above, with mapping=MappingQGeneric1<1>()
.
Definition at line 134 of file error_estimator_1d.cc.
|
static |
Equivalent to the set of functions above, except that this one takes a quadrature collection for hp finite element dof handlers.
Definition at line 168 of file error_estimator_1d.cc.
|
static |
Equivalent to the set of functions above, except that this one takes a quadrature collection for hp finite element dof handlers.
Definition at line 205 of file error_estimator_1d.cc.
|
static |
Equivalent to the set of functions above, except that this one takes a quadrature collection for hp finite element dof handlers.
Definition at line 273 of file error_estimator_1d.cc.
|
static |
Equivalent to the set of functions above, except that this one takes a quadrature collection for hp finite element dof handlers.
Definition at line 239 of file error_estimator_1d.cc.