Reference documentation for deal.II version 8.1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
Public Member Functions | Static Public Member Functions | List of all members
KellyErrorEstimator< 1, spacedim > Class Template Reference

#include <error_estimator.h>

Public Member Functions

 DeclException0 (ExcInvalidBoundaryIndicator)
 
 DeclException0 (ExcInvalidComponentMask)
 
 DeclException0 (ExcInvalidCoefficient)
 
 DeclException0 (ExcInvalidBoundaryFunction)
 
 DeclException2 (ExcIncompatibleNumberOfElements, int, int,<< "The number of elements "<< arg1<< " and "<< arg2<< " of the vectors do not match!")
 
 DeclException0 (ExcInvalidSolutionVector)
 
 DeclException0 (ExcNoSolutions)
 

Static Public Member Functions

template<typename InputVector , class DH >
static void estimate (const Mapping< 1, spacedim > &mapping, const DH &dof, const Quadrature< 0 > &quadrature, const typename FunctionMap< spacedim >::type &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=0, 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)
 
template<typename InputVector , class DH >
static void estimate (const DH &dof, const Quadrature< 0 > &quadrature, const typename FunctionMap< spacedim >::type &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=0, 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)
 
template<typename InputVector , class DH >
static void estimate (const Mapping< 1, spacedim > &mapping, const DH &dof, const Quadrature< 0 > &quadrature, const typename FunctionMap< spacedim >::type &neumann_bc, const std::vector< const InputVector * > &solutions, std::vector< Vector< float > * > &errors, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=0, 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)
 
template<typename InputVector , class DH >
static void estimate (const DH &dof, const Quadrature< 0 > &quadrature, const typename FunctionMap< spacedim >::type &neumann_bc, const std::vector< const InputVector * > &solutions, std::vector< Vector< float > * > &errors, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=0, 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)
 
template<typename InputVector , class DH >
static void estimate (const Mapping< 1, spacedim > &mapping, const DH &dof, const hp::QCollection< 0 > &quadrature, const typename FunctionMap< spacedim >::type &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=0, 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)
 
template<typename InputVector , class DH >
static void estimate (const DH &dof, const hp::QCollection< 0 > &quadrature, const typename FunctionMap< spacedim >::type &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=0, 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)
 
template<typename InputVector , class DH >
static void estimate (const Mapping< 1, spacedim > &mapping, const DH &dof, const hp::QCollection< 0 > &quadrature, const typename FunctionMap< spacedim >::type &neumann_bc, const std::vector< const InputVector * > &solutions, std::vector< Vector< float > * > &errors, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=0, 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)
 
template<typename InputVector , class DH >
static void estimate (const DH &dof, const hp::QCollection< 0 > &quadrature, const typename FunctionMap< spacedim >::type &neumann_bc, const std::vector< const InputVector * > &solutions, std::vector< Vector< float > * > &errors, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=0, 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)
 

Detailed Description

template<int spacedim>
class KellyErrorEstimator< 1, spacedim >

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.

Author
Wolfgang Bangerth, 1998, 2004.

Definition at line 481 of file error_estimator.h.

Member Function Documentation

template<int spacedim>
template<typename InputVector , class DH >
static void KellyErrorEstimator< 1, spacedim >::estimate ( const Mapping< 1, spacedim > &  mapping,
const DH &  dof,
const Quadrature< 0 > &  quadrature,
const typename FunctionMap< spacedim >::type &  neumann_bc,
const InputVector &  solution,
Vector< float > &  error,
const ComponentMask component_mask = ComponentMask(),
const Function< spacedim > *  coefficients = 0,
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 
)
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.

template<int spacedim>
template<typename InputVector , class DH >
static void KellyErrorEstimator< 1, spacedim >::estimate ( const DH &  dof,
const Quadrature< 0 > &  quadrature,
const typename FunctionMap< spacedim >::type &  neumann_bc,
const InputVector &  solution,
Vector< float > &  error,
const ComponentMask component_mask = ComponentMask(),
const Function< spacedim > *  coefficients = 0,
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 
)
static

Calls the estimate function, see above, with mapping=MappingQ1<1>().

template<int spacedim>
template<typename InputVector , class DH >
static void KellyErrorEstimator< 1, spacedim >::estimate ( const Mapping< 1, spacedim > &  mapping,
const DH &  dof,
const Quadrature< 0 > &  quadrature,
const typename FunctionMap< spacedim >::type &  neumann_bc,
const std::vector< const InputVector * > &  solutions,
std::vector< Vector< float > * > &  errors,
const ComponentMask component_mask = ComponentMask(),
const Function< spacedim > *  coefficients = 0,
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 
)
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.)

template<int spacedim>
template<typename InputVector , class DH >
static void KellyErrorEstimator< 1, spacedim >::estimate ( const DH &  dof,
const Quadrature< 0 > &  quadrature,
const typename FunctionMap< spacedim >::type &  neumann_bc,
const std::vector< const InputVector * > &  solutions,
std::vector< Vector< float > * > &  errors,
const ComponentMask component_mask = ComponentMask(),
const Function< spacedim > *  coefficients = 0,
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 
)
static

Calls the estimate function, see above, with mapping=MappingQ1<1>().

template<int spacedim>
template<typename InputVector , class DH >
static void KellyErrorEstimator< 1, spacedim >::estimate ( const Mapping< 1, spacedim > &  mapping,
const DH &  dof,
const hp::QCollection< 0 > &  quadrature,
const typename FunctionMap< spacedim >::type &  neumann_bc,
const InputVector &  solution,
Vector< float > &  error,
const ComponentMask component_mask = ComponentMask(),
const Function< spacedim > *  coefficients = 0,
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 
)
static

Equivalent to the set of functions above, except that this one takes a quadrature collection for hp finite element dof handlers.

template<int spacedim>
template<typename InputVector , class DH >
static void KellyErrorEstimator< 1, spacedim >::estimate ( const DH &  dof,
const hp::QCollection< 0 > &  quadrature,
const typename FunctionMap< spacedim >::type &  neumann_bc,
const InputVector &  solution,
Vector< float > &  error,
const ComponentMask component_mask = ComponentMask(),
const Function< spacedim > *  coefficients = 0,
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 
)
static

Equivalent to the set of functions above, except that this one takes a quadrature collection for hp finite element dof handlers.

template<int spacedim>
template<typename InputVector , class DH >
static void KellyErrorEstimator< 1, spacedim >::estimate ( const Mapping< 1, spacedim > &  mapping,
const DH &  dof,
const hp::QCollection< 0 > &  quadrature,
const typename FunctionMap< spacedim >::type &  neumann_bc,
const std::vector< const InputVector * > &  solutions,
std::vector< Vector< float > * > &  errors,
const ComponentMask component_mask = ComponentMask(),
const Function< spacedim > *  coefficients = 0,
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 
)
static

Equivalent to the set of functions above, except that this one takes a quadrature collection for hp finite element dof handlers.

template<int spacedim>
template<typename InputVector , class DH >
static void KellyErrorEstimator< 1, spacedim >::estimate ( const DH &  dof,
const hp::QCollection< 0 > &  quadrature,
const typename FunctionMap< spacedim >::type &  neumann_bc,
const std::vector< const InputVector * > &  solutions,
std::vector< Vector< float > * > &  errors,
const ComponentMask component_mask = ComponentMask(),
const Function< spacedim > *  coefficients = 0,
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 
)
static

Equivalent to the set of functions above, except that this one takes a quadrature collection for hp finite element dof handlers.

template<int spacedim>
KellyErrorEstimator< 1, spacedim >::DeclException0 ( ExcInvalidBoundaryIndicator  )

Exception

template<int spacedim>
KellyErrorEstimator< 1, spacedim >::DeclException0 ( ExcInvalidComponentMask  )

Exception

template<int spacedim>
KellyErrorEstimator< 1, spacedim >::DeclException0 ( ExcInvalidCoefficient  )

Exception

template<int spacedim>
KellyErrorEstimator< 1, spacedim >::DeclException0 ( ExcInvalidBoundaryFunction  )

Exception

template<int spacedim>
KellyErrorEstimator< 1, spacedim >::DeclException2 ( ExcIncompatibleNumberOfElements  ,
int  ,
int  ,
<< "The number of elements "<< arg1<< " and "<< arg2<< " of the vectors do not match!"   
)

Exception

template<int spacedim>
KellyErrorEstimator< 1, spacedim >::DeclException0 ( ExcInvalidSolutionVector  )

Exception

template<int spacedim>
KellyErrorEstimator< 1, spacedim >::DeclException0 ( ExcNoSolutions  )

Exception


The documentation for this class was generated from the following file: