Reference documentation for deal.II version 8.4.1
Enumerations | Functions
VectorTools Namespace Reference

Enumerations

enum  NormType {
  mean, L1_norm, L2_norm, Lp_norm,
  Linfty_norm, H1_seminorm, Hdiv_seminorm, H1_norm,
  W1p_seminorm, W1p_norm, W1infty_seminorm, W1infty_norm
}
 

Functions

 DeclExceptionMsg (ExcNonInterpolatingFE,"You are attempting an operation that requires the ""finite element involved to be 'interpolating', i.e., ""it needs to have support points. The finite element ""you are using here does not appear to have those.")
 
 DeclExceptionMsg (ExcPointNotAvailableHere,"The given point is inside a cell of a ""parallel::distributed::Triangulation that is not ""locally owned by this processor.")
 
Interpolation and projection
template<typename VectorType , int dim, int spacedim, template< int, int > class DoFHandlerType>
void interpolate (const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof, const Function< spacedim, double > &function, VectorType &vec)
 
template<typename VectorType , typename DoFHandlerType >
void interpolate (const DoFHandlerType &dof, const Function< DoFHandlerType::space_dimension, double > &function, VectorType &vec)
 
template<int dim, class InVector , class OutVector , int spacedim>
void interpolate (const DoFHandler< dim, spacedim > &dof_1, const DoFHandler< dim, spacedim > &dof_2, const FullMatrix< double > &transfer, const InVector &data_1, OutVector &data_2)
 
template<typename VectorType , typename DoFHandlerType >
void interpolate_based_on_material_id (const Mapping< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &mapping, const DoFHandlerType &dof_handler, const std::map< types::material_id, const Function< DoFHandlerType::space_dimension, double > * > &function_map, VectorType &dst, const ComponentMask &component_mask=ComponentMask())
 
template<int dim, int spacedim, template< int, int > class DoFHandlerType, typename VectorType >
void interpolate_to_different_mesh (const DoFHandlerType< dim, spacedim > &dof1, const VectorType &u1, const DoFHandlerType< dim, spacedim > &dof2, VectorType &u2)
 
template<int dim, int spacedim, template< int, int > class DoFHandlerType, typename VectorType >
void interpolate_to_different_mesh (const DoFHandlerType< dim, spacedim > &dof1, const VectorType &u1, const DoFHandlerType< dim, spacedim > &dof2, const ConstraintMatrix &constraints, VectorType &u2)
 
template<int dim, int spacedim, template< int, int > class DoFHandlerType, typename VectorType >
void interpolate_to_different_mesh (const InterGridMap< DoFHandlerType< dim, spacedim > > &intergridmap, const VectorType &u1, const ConstraintMatrix &constraints, VectorType &u2)
 
template<int dim, typename VectorType , int spacedim>
void project (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const ConstraintMatrix &constraints, const Quadrature< dim > &quadrature, const Function< spacedim, double > &function, VectorType &vec, const bool enforce_zero_boundary=false, const Quadrature< dim-1 > &q_boundary=(dim > 1?QGauss< dim-1 >(2):Quadrature< dim-1 >(0)), const bool project_to_boundary_first=false)
 
template<int dim, typename VectorType , int spacedim>
void project (const DoFHandler< dim, spacedim > &dof, const ConstraintMatrix &constraints, const Quadrature< dim > &quadrature, const Function< spacedim, double > &function, VectorType &vec, const bool enforce_zero_boundary=false, const Quadrature< dim-1 > &q_boundary=(dim > 1?QGauss< dim-1 >(2):Quadrature< dim-1 >(0)), const bool project_to_boundary_first=false)
 
template<int dim, typename VectorType , int spacedim>
void project (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const ConstraintMatrix &constraints, const hp::QCollection< dim > &quadrature, const Function< spacedim, double > &function, VectorType &vec, const bool enforce_zero_boundary=false, const hp::QCollection< dim-1 > &q_boundary=hp::QCollection< dim-1 >(dim > 1?QGauss< dim-1 >(2):Quadrature< dim-1 >(0)), const bool project_to_boundary_first=false)
 
template<int dim, typename VectorType , int spacedim>
void project (const hp::DoFHandler< dim, spacedim > &dof, const ConstraintMatrix &constraints, const hp::QCollection< dim > &quadrature, const Function< spacedim, double > &function, VectorType &vec, const bool enforce_zero_boundary=false, const hp::QCollection< dim-1 > &q_boundary=hp::QCollection< dim-1 >(dim > 1?QGauss< dim-1 >(2):Quadrature< dim-1 >(0)), const bool project_to_boundary_first=false)
 
template<typename DoFHandlerType >
void interpolate_boundary_values (const Mapping< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &mapping, const DoFHandlerType &dof, const typename FunctionMap< DoFHandlerType::space_dimension >::type &function_map, std::map< types::global_dof_index, double > &boundary_values, const ComponentMask &component_mask=ComponentMask())
 
template<int dim, int spacedim>
void interpolate_boundary_values (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const typename FunctionMap< spacedim >::type &function_map, std::map< types::global_dof_index, double > &boundary_values, const ComponentMask &component_mask=ComponentMask())
 
template<typename DoFHandlerType >
void interpolate_boundary_values (const Mapping< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &mapping, const DoFHandlerType &dof, const types::boundary_id boundary_component, const Function< DoFHandlerType::space_dimension, double > &boundary_function, std::map< types::global_dof_index, double > &boundary_values, const ComponentMask &component_mask=ComponentMask())
 
template<typename DoFHandlerType >
void interpolate_boundary_values (const DoFHandlerType &dof, const types::boundary_id boundary_component, const Function< DoFHandlerType::space_dimension, double > &boundary_function, std::map< types::global_dof_index, double > &boundary_values, const ComponentMask &component_mask=ComponentMask())
 
template<typename DoFHandlerType >
void interpolate_boundary_values (const DoFHandlerType &dof, const typename FunctionMap< DoFHandlerType::space_dimension >::type &function_map, std::map< types::global_dof_index, double > &boundary_values, const ComponentMask &component_mask=ComponentMask())
 
template<typename DoFHandlerType >
void interpolate_boundary_values (const Mapping< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &mapping, const DoFHandlerType &dof, const typename FunctionMap< DoFHandlerType::space_dimension >::type &function_map, ConstraintMatrix &constraints, const ComponentMask &component_mask=ComponentMask())
 
template<typename DoFHandlerType >
void interpolate_boundary_values (const Mapping< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &mapping, const DoFHandlerType &dof, const types::boundary_id boundary_component, const Function< DoFHandlerType::space_dimension, double > &boundary_function, ConstraintMatrix &constraints, const ComponentMask &component_mask=ComponentMask())
 
template<typename DoFHandlerType >
void interpolate_boundary_values (const DoFHandlerType &dof, const types::boundary_id boundary_component, const Function< DoFHandlerType::space_dimension, double > &boundary_function, ConstraintMatrix &constraints, const ComponentMask &component_mask=ComponentMask())
 
template<typename DoFHandlerType >
void interpolate_boundary_values (const DoFHandlerType &dof, const typename FunctionMap< DoFHandlerType::space_dimension >::type &function_map, ConstraintMatrix &constraints, const ComponentMask &component_mask=ComponentMask())
 
template<int dim, int spacedim>
void project_boundary_values (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const typename FunctionMap< spacedim >::type &boundary_functions, const Quadrature< dim-1 > &q, std::map< types::global_dof_index, double > &boundary_values, std::vector< unsigned int > component_mapping=std::vector< unsigned int >())
 
template<int dim, int spacedim>
void project_boundary_values (const DoFHandler< dim, spacedim > &dof, const typename FunctionMap< spacedim >::type &boundary_function, const Quadrature< dim-1 > &q, std::map< types::global_dof_index, double > &boundary_values, std::vector< unsigned int > component_mapping=std::vector< unsigned int >())
 
template<int dim, int spacedim>
void project_boundary_values (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const typename FunctionMap< spacedim >::type &boundary_functions, const hp::QCollection< dim-1 > &q, std::map< types::global_dof_index, double > &boundary_values, std::vector< unsigned int > component_mapping=std::vector< unsigned int >())
 
template<int dim, int spacedim>
void project_boundary_values (const hp::DoFHandler< dim, spacedim > &dof, const typename FunctionMap< spacedim >::type &boundary_function, const hp::QCollection< dim-1 > &q, std::map< types::global_dof_index, double > &boundary_values, std::vector< unsigned int > component_mapping=std::vector< unsigned int >())
 
template<int dim, int spacedim>
void project_boundary_values (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const typename FunctionMap< spacedim >::type &boundary_functions, const Quadrature< dim-1 > &q, ConstraintMatrix &constraints, std::vector< unsigned int > component_mapping=std::vector< unsigned int >())
 
template<int dim, int spacedim>
void project_boundary_values (const DoFHandler< dim, spacedim > &dof, const typename FunctionMap< spacedim >::type &boundary_function, const Quadrature< dim-1 > &q, ConstraintMatrix &constraints, std::vector< unsigned int > component_mapping=std::vector< unsigned int >())
 
template<int dim>
void project_boundary_values_curl_conforming (const DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, double > &boundary_function, const types::boundary_id boundary_component, ConstraintMatrix &constraints, const Mapping< dim > &mapping=StaticMappingQ1< dim >::mapping)
 
template<int dim>
void project_boundary_values_curl_conforming (const hp::DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, double > &boundary_function, const types::boundary_id boundary_component, ConstraintMatrix &constraints, const hp::MappingCollection< dim, dim > &mapping_collection=hp::StaticMappingQ1< dim >::mapping_collection)
 
template<int dim>
void project_boundary_values_curl_conforming_l2 (const DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, double > &boundary_function, const types::boundary_id boundary_component, ConstraintMatrix &constraints, const Mapping< dim > &mapping=StaticMappingQ1< dim >::mapping)
 
template<int dim>
void project_boundary_values_curl_conforming_l2 (const hp::DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, double > &boundary_function, const types::boundary_id boundary_component, ConstraintMatrix &constraints, const hp::MappingCollection< dim, dim > &mapping_collection=hp::StaticMappingQ1< dim >::mapping_collection)
 
template<int dim>
void project_boundary_values_div_conforming (const DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, double > &boundary_function, const types::boundary_id boundary_component, ConstraintMatrix &constraints, const Mapping< dim > &mapping=StaticMappingQ1< dim >::mapping)
 
template<int dim>
void project_boundary_values_div_conforming (const hp::DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, double > &boundary_function, const types::boundary_id boundary_component, ConstraintMatrix &constraints, const hp::MappingCollection< dim, dim > &mapping_collection=hp::StaticMappingQ1< dim >::mapping_collection)
 
template<int dim, template< int, int > class DoFHandlerType, int spacedim>
void compute_nonzero_normal_flux_constraints (const DoFHandlerType< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, typename FunctionMap< spacedim >::type &function_map, ConstraintMatrix &constraints, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim >::mapping)
 
template<int dim, template< int, int > class DoFHandlerType, int spacedim>
void compute_no_normal_flux_constraints (const DoFHandlerType< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, ConstraintMatrix &constraints, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim >::mapping)
 
template<int dim, template< int, int > class DoFHandlerType, int spacedim>
void compute_nonzero_tangential_flux_constraints (const DoFHandlerType< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, typename FunctionMap< spacedim >::type &function_map, ConstraintMatrix &constraints, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim >::mapping)
 
template<int dim, template< int, int > class DoFHandlerType, int spacedim>
void compute_normal_flux_constraints (const DoFHandlerType< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, ConstraintMatrix &constraints, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim >::mapping)
 
Assembling of right hand sides
template<int dim, int spacedim>
void create_right_hand_side (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &q, const Function< spacedim, double > &rhs, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_right_hand_side (const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &q, const Function< spacedim, double > &rhs, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_right_hand_side (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim > &q, const Function< spacedim, double > &rhs, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_right_hand_side (const hp::DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim > &q, const Function< spacedim, double > &rhs, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_point_source_vector (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_point_source_vector (const DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_point_source_vector (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_point_source_vector (const hp::DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_point_source_vector (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, const Point< dim > &direction, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_point_source_vector (const DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, const Point< dim > &direction, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_point_source_vector (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, const Point< dim > &direction, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_point_source_vector (const hp::DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, const Point< dim > &direction, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_boundary_right_hand_side (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim-1 > &q, const Function< spacedim, double > &rhs, Vector< double > &rhs_vector, const std::set< types::boundary_id > &boundary_ids=std::set< types::boundary_id >())
 
template<int dim, int spacedim>
void create_boundary_right_hand_side (const DoFHandler< dim, spacedim > &dof, const Quadrature< dim-1 > &q, const Function< spacedim, double > &rhs, Vector< double > &rhs_vector, const std::set< types::boundary_id > &boundary_ids=std::set< types::boundary_id >())
 
template<int dim, int spacedim>
void create_boundary_right_hand_side (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim-1 > &q, const Function< spacedim, double > &rhs, Vector< double > &rhs_vector, const std::set< types::boundary_id > &boundary_ids=std::set< types::boundary_id >())
 
template<int dim, int spacedim>
void create_boundary_right_hand_side (const hp::DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim-1 > &q, const Function< spacedim, double > &rhs, Vector< double > &rhs_vector, const std::set< types::boundary_id > &boundary_ids=std::set< types::boundary_id >())
 
Evaluation of functions and errors
template<int dim, class InVector , class OutVector , int spacedim>
void integrate_difference (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim, double > &exact_solution, OutVector &difference, const Quadrature< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=0, const double exponent=2.)
 
template<int dim, class InVector , class OutVector , int spacedim>
void integrate_difference (const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim, double > &exact_solution, OutVector &difference, const Quadrature< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=0, const double exponent=2.)
 
template<int dim, class InVector , class OutVector , int spacedim>
void integrate_difference (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim, double > &exact_solution, OutVector &difference, const hp::QCollection< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=0, const double exponent=2.)
 
template<int dim, class InVector , class OutVector , int spacedim>
void integrate_difference (const hp::DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim, double > &exact_solution, OutVector &difference, const hp::QCollection< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=0, const double exponent=2.)
 
template<int dim, typename VectorType , int spacedim>
void point_difference (const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Function< spacedim, double > &exact_solution, Vector< double > &difference, const Point< spacedim > &point)
 
template<int dim, typename VectorType , int spacedim>
void point_difference (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Function< spacedim, double > &exact_solution, Vector< double > &difference, const Point< spacedim > &point)
 
template<int dim, typename VectorType , int spacedim>
void point_value (const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, Vector< double > &value)
 
template<int dim, typename VectorType , int spacedim>
void point_value (const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, Vector< double > &value)
 
template<int dim, typename VectorType , int spacedim>
double point_value (const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point)
 
template<int dim, typename VectorType , int spacedim>
double point_value (const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point)
 
template<int dim, typename VectorType , int spacedim>
void point_value (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, Vector< double > &value)
 
template<int dim, typename VectorType , int spacedim>
void point_value (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, Vector< double > &value)
 
template<int dim, typename VectorType , int spacedim>
double point_value (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point)
 
template<int dim, typename VectorType , int spacedim>
double point_value (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point)
 
template<int dim, typename VectorType , int spacedim>
void point_gradient (const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > &value)
 
template<int dim, typename VectorType , int spacedim>
void point_gradient (const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > &value)
 
template<int dim, typename VectorType , int spacedim>
Tensor< 1, spacedim, typename VectorType::value_type > point_gradient (const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point)
 
template<int dim, typename VectorType , int spacedim>
Tensor< 1, spacedim, typename VectorType::value_type > point_gradient (const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point)
 
template<int dim, typename VectorType , int spacedim>
void point_gradient (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > &value)
 
template<int dim, typename VectorType , int spacedim>
void point_gradient (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > &value)
 
template<int dim, typename VectorType , int spacedim>
Tensor< 1, spacedim, typename VectorType::value_type > point_gradient (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point)
 
template<int dim, typename VectorType , int spacedim>
Tensor< 1, spacedim, typename VectorType::value_type > point_gradient (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point)
 
template<typename VectorType >
void subtract_mean_value (VectorType &v, const std::vector< bool > &p_select=std::vector< bool >())
 
template<int dim, typename VectorType , int spacedim>
double compute_mean_value (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &quadrature, const VectorType &v, const unsigned int component)
 
template<int dim, typename VectorType , int spacedim>
double compute_mean_value (const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &quadrature, const VectorType &v, const unsigned int component)
 
template<typename DoFHandlerType , typename VectorType >
void get_position_vector (const DoFHandlerType &dh, VectorType &vector, const ComponentMask &mask=ComponentMask())
 

Detailed Description

Provide a namespace which offers some operations on vectors. Among these are assembling of standard vectors, integration of the difference of a finite element solution and a continuous function, interpolations and projections of continuous functions to the finite element space and other operations.

Note
There exist two versions of almost all functions, one that takes an explicit Mapping argument and one that does not. The second one generally calls the first with an implicit \(Q_1\) argument (i.e., with an argument of kind MappingQGeneric(1)). If your intend your code to use a different mapping than a (bi-/tri-)linear one, then you need to call the functions with mapping argument should be used.

Description of operations

This collection of methods offers the following operations:

All functions use the finite element given to the DoFHandler object the last time that the degrees of freedom were distributed over the triangulation. Also, if access to an object describing the exact form of the boundary is needed, the pointer stored within the triangulation object is accessed.

Note
Instantiations for this template are provided for some vector types, in particular Vector<float>, Vector<double>, BlockVector<float>, BlockVector<double>; others can be generated in application code (see the section on Template instantiations in the manual).
Author
Wolfgang Bangerth, Ralf Hartmann, Guido Kanschat, 1998, 1999, 2000, 2001

Function Documentation

template<typename VectorType , int dim, int spacedim, template< int, int > class DoFHandlerType>
void VectorTools::interpolate ( const Mapping< dim, spacedim > &  mapping,
const DoFHandlerType< dim, spacedim > &  dof,
const Function< spacedim, double > &  function,
VectorType &  vec 
)

Compute the interpolation of function at the support points to the finite element space described by the Triangulation and FiniteElement object with which the given DoFHandler argument is initialized. It is assumed that the number of components of function matches that of the finite element used by dof.

Note that you may have to call hanging_nodes.distribute(vec) with the hanging nodes from space dof afterwards, to make the result continuous again.

The template argument DoFHandlerType may either be of type DoFHandler or hp::DoFHandler.

See the general documentation of this namespace for further information.

Todo:
The mapping argument should be replaced by a hp::MappingCollection in case of a hp::DoFHandler.
template<typename VectorType , typename DoFHandlerType >
void VectorTools::interpolate ( const DoFHandlerType &  dof,
const Function< DoFHandlerType::space_dimension, double > &  function,
VectorType &  vec 
)

Calls the interpolate() function above with mapping=MappingQGeneric1<dim>@().

template<int dim, class InVector , class OutVector , int spacedim>
void VectorTools::interpolate ( const DoFHandler< dim, spacedim > &  dof_1,
const DoFHandler< dim, spacedim > &  dof_2,
const FullMatrix< double > &  transfer,
const InVector &  data_1,
OutVector &  data_2 
)

Interpolate different finite element spaces. The interpolation of vector data_1 is executed from the FE space represented by dof_1 to the vector data_2 on FE space dof_2. The interpolation on each cell is represented by the matrix transfer. Curved boundaries are neglected so far.

Note that you may have to call hanging_nodes.distribute(data_2) with the hanging nodes from space dof_2 afterwards, to make the result continuous again.

Note
Instantiations for this template are provided for some vector types (see the general documentation of the namespace), but only the same vector for InVector and OutVector. Other combinations must be instantiated by hand.
template<typename VectorType , typename DoFHandlerType >
void VectorTools::interpolate_based_on_material_id ( const Mapping< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &  mapping,
const DoFHandlerType &  dof_handler,
const std::map< types::material_id, const Function< DoFHandlerType::space_dimension, double > * > &  function_map,
VectorType &  dst,
const ComponentMask component_mask = ComponentMask() 
)

This function is a kind of generalization or modification of the very first interpolate() function in the series. It interpolations a set of functions onto the finite element space given by the DoFHandler argument where the determination which function to use is made based on the material id (see GlossMaterialId) of each cell.

Parameters
mapping- The mapping to use to determine the location of support points at which the functions are to be evaluated.
dof_handler- DoFHandler initialized with Triangulation and FiniteElement objects,
function_map- std::map reflecting the correspondence between material ids and functions,
dst- global FE vector at the support points,
component_mask- mask of components that shall be interpolated
Note
If a material id of some group of cells is missed in function_map, then dst will not be updated in the respective degrees of freedom of the output vector For example, if dst was successfully initialized to capture the degrees of freedom of the dof_handler of the problem with all zeros in it, then those zeros which correspond to the missed material ids will still remain in dst even after calling this function.
Degrees of freedom located on faces between cells of different material ids will get their value by that cell which was called last in the respective loop over cells implemented in this function. Since this process is kind of arbitrary, you cannot control it. However, if you want to have control over the order in which cells are visited, let us take a look at the following example: Let u be a variable of interest which is approximated by some CG finite element. Let 0, 1 and 2 be material ids of cells on the triangulation. Let 0: 0.0, 1: 1.0, 2: 2.0 be the whole function_map that you want to pass to this function, where key is a material id and value is a value of u. By using the whole function_map you do not really know which values will be assigned to the face DoFs. On the other hand, if you split the whole function_map into three smaller independent objects 0: 0.0 and 1: 1.0 and 2: 2.0 and make three distinct calls of this function passing each of these objects separately (the order depends on what you want to get between cells), then each subsequent call will rewrite the intercell dofs of the previous one.
Author
Valentin Zingan, 2013
template<int dim, int spacedim, template< int, int > class DoFHandlerType, typename VectorType >
void VectorTools::interpolate_to_different_mesh ( const DoFHandlerType< dim, spacedim > &  dof1,
const VectorType &  u1,
const DoFHandlerType< dim, spacedim > &  dof2,
VectorType &  u2 
)

Gives the interpolation of a dof1-function u1 to a dof2-function u2, where dof1 and dof2 represent different triangulations with a common coarse grid.

dof1 and dof2 need to have the same finite element discretization.

Note that for continuous elements on grids with hanging nodes (i.e. locally refined grids) this function does not give the expected output. Indeed, the resulting output vector does not necessarily respect continuity requirements at hanging nodes, due to local cellwise interpolation.

For this case (continuous elements on grids with hanging nodes), please use the interpolate_to_different_mesh function with an additional ConstraintMatrix argument, see below, or make the field conforming yourself by calling the ConstraintsMatrix::distribute function of your hanging node constraints object.

Note
: This function works with parallel::distributed::Triangulation, but only if the parallel partitioning is the same for both meshes (see the parallel::distributed::Triangulation<dim>::no_automatic_repartitioning flag).
template<int dim, int spacedim, template< int, int > class DoFHandlerType, typename VectorType >
void VectorTools::interpolate_to_different_mesh ( const DoFHandlerType< dim, spacedim > &  dof1,
const VectorType &  u1,
const DoFHandlerType< dim, spacedim > &  dof2,
const ConstraintMatrix constraints,
VectorType &  u2 
)

Gives the interpolation of a dof1-function u1 to a dof2-function u2, where dof1 and dof2 represent different triangulations with a common coarse grid.

dof1 and dof2 need to have the same finite element discretization.

constraints is a hanging node constraints object corresponding to dof2. This object is particularly important when interpolating onto continuous elements on grids with hanging nodes (locally refined grids): Without it - due to cellwise interpolation - the resulting output vector does not necessarily respect continuity requirements at hanging nodes.

template<int dim, int spacedim, template< int, int > class DoFHandlerType, typename VectorType >
void VectorTools::interpolate_to_different_mesh ( const InterGridMap< DoFHandlerType< dim, spacedim > > &  intergridmap,
const VectorType &  u1,
const ConstraintMatrix constraints,
VectorType &  u2 
)

The same function as above, but takes an InterGridMap object directly as a parameter. Useful for interpolating several vectors at the same time.

intergridmap has to be initialized via InterGridMap::make_mapping pointing from a source DoFHandler to a destination DoFHandler.

template<int dim, typename VectorType , int spacedim>
void VectorTools::project ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const ConstraintMatrix constraints,
const Quadrature< dim > &  quadrature,
const Function< spacedim, double > &  function,
VectorType &  vec,
const bool  enforce_zero_boundary = false,
const Quadrature< dim-1 > &  q_boundary = (dim > 1?QGauss< dim-1 >(2):Quadrature< dim-1 >(0)),
const bool  project_to_boundary_first = false 
)

Compute the projection of function to the finite element space.

By default, projection to the boundary and enforcement of zero boundary values are disabled. The ordering of arguments to this function is such that you need not give a second quadrature formula if you don't want to project to the boundary first, but that you must if you want to do so.

This function needs the mass matrix of the finite element space on the present grid. To this end, the mass matrix is assembled exactly using MatrixTools::create_mass_matrix. This function performs numerical quadrature using the given quadrature rule; you should therefore make sure that the given quadrature formula is also sufficient for the integration of the mass matrix.

See the general documentation of this namespace for further information.

In 1d, the default value of the boundary quadrature formula is an invalid object since integration on the boundary doesn't happen in 1d.

template<int dim, typename VectorType , int spacedim>
void VectorTools::project ( const DoFHandler< dim, spacedim > &  dof,
const ConstraintMatrix constraints,
const Quadrature< dim > &  quadrature,
const Function< spacedim, double > &  function,
VectorType &  vec,
const bool  enforce_zero_boundary = false,
const Quadrature< dim-1 > &  q_boundary = (dim > 1?QGauss< dim-1 >(2):Quadrature< dim-1 >(0)),
const bool  project_to_boundary_first = false 
)

Calls the project() function above, with mapping=MappingQGeneric<dim>(1).

template<int dim, typename VectorType , int spacedim>
void VectorTools::project ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const ConstraintMatrix constraints,
const hp::QCollection< dim > &  quadrature,
const Function< spacedim, double > &  function,
VectorType &  vec,
const bool  enforce_zero_boundary = false,
const hp::QCollection< dim-1 > &  q_boundary = hp::QCollection< dim-1 >(dim > 1?QGauss< dim-1 >(2):Quadrature< dim-1 >(0)),
const bool  project_to_boundary_first = false 
)

Same as above, but for arguments of type hp::DoFHandler, hp::QuadratureCollection, hp::MappingCollection

template<int dim, typename VectorType , int spacedim>
void VectorTools::project ( const hp::DoFHandler< dim, spacedim > &  dof,
const ConstraintMatrix constraints,
const hp::QCollection< dim > &  quadrature,
const Function< spacedim, double > &  function,
VectorType &  vec,
const bool  enforce_zero_boundary = false,
const hp::QCollection< dim-1 > &  q_boundary = hp::QCollection< dim-1 >(dim > 1?QGauss< dim-1 >(2):Quadrature< dim-1 >(0)),
const bool  project_to_boundary_first = false 
)

Calls the project() function above, with a collection of \(Q_1\) mapping objects, i.e., with hp::StaticMappingQ1::mapping_collection.

template<typename DoFHandlerType >
void VectorTools::interpolate_boundary_values ( const Mapping< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &  mapping,
const DoFHandlerType &  dof,
const typename FunctionMap< DoFHandlerType::space_dimension >::type &  function_map,
std::map< types::global_dof_index, double > &  boundary_values,
const ComponentMask component_mask = ComponentMask() 
)

Compute Dirichlet boundary conditions. This function makes up a map of degrees of freedom subject to Dirichlet boundary conditions and the corresponding values to be assigned to them, by interpolation around the boundary. For each degree of freedom at the boundary, if its index already exists in boundary_values then its boundary value will be overwritten, otherwise a new entry with proper index and boundary value for this degree of freedom will be inserted into boundary_values.

The parameter function_map provides a list of boundary indicators to be handled by this function and corresponding boundary value functions. The keys of this map correspond to the number boundary_id of the face. numbers::internal_face_boundary_id is an illegal value for this key since it is reserved for interior faces.

The flags in the last parameter, component_mask denote which components of the finite element space shall be interpolated. If it is left as specified by the default value (i.e. an empty array), all components are interpolated. If it is different from the default value, it is assumed that the number of entries equals the number of components in the boundary functions and the finite element, and those components in the given boundary function will be used for which the respective flag was set in the component mask. See also GlossComponentMask. As an example, assume that you are solving the Stokes equations in 2d, with variables \((u,v,p)\) and that you only want to interpolate boundary values for the velocity, then the component mask should correspond to (true,true,false).

Note
Whether a component mask has been specified or not, the number of components of the functions in function_map must match that of the finite element used by dof. In other words, for the example above, you need to provide a Function object that has 3 components (the two velocities and the pressure), even though you are only interested in the first two of them. interpolate_boundary_values() will then call this function to obtain a vector of 3 values at each interpolation point but only take the first two and discard the third. In other words, you are free to return whatever you like in the third component of the vector returned by Function::vector_value, but the Function object must state that it has 3 components.

If the finite element used has shape functions that are non-zero in more than one component (in deal.II speak: they are non-primitive), then these components can presently not be used for interpolating boundary values. Thus, the elements in the component mask corresponding to the components of these non-primitive shape functions must be false.

See the general documentation of this namespace for more information.

template<int dim, int spacedim>
void VectorTools::interpolate_boundary_values ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const typename FunctionMap< spacedim >::type &  function_map,
std::map< types::global_dof_index, double > &  boundary_values,
const ComponentMask component_mask = ComponentMask() 
)

Like the previous function, but take a mapping collection to go with the hp::DoFHandler object.

template<typename DoFHandlerType >
void VectorTools::interpolate_boundary_values ( const Mapping< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &  mapping,
const DoFHandlerType &  dof,
const types::boundary_id  boundary_component,
const Function< DoFHandlerType::space_dimension, double > &  boundary_function,
std::map< types::global_dof_index, double > &  boundary_values,
const ComponentMask component_mask = ComponentMask() 
)

Same function as above, but taking only one pair of boundary indicator and corresponding boundary function. The same comments apply as for the previous function, in particular about the use of the component mask and the requires size of the function object.

See also
Glossary entry on boundary indicators
template<typename DoFHandlerType >
void VectorTools::interpolate_boundary_values ( const DoFHandlerType &  dof,
const types::boundary_id  boundary_component,
const Function< DoFHandlerType::space_dimension, double > &  boundary_function,
std::map< types::global_dof_index, double > &  boundary_values,
const ComponentMask component_mask = ComponentMask() 
)

Calls the other interpolate_boundary_values() function, see above, with mapping=MappingQGeneric<dim,spacedim>(1). The same comments apply as for the previous function, in particular about the use of the component mask and the requires size of the function object.

See also
Glossary entry on boundary indicators
template<typename DoFHandlerType >
void VectorTools::interpolate_boundary_values ( const DoFHandlerType &  dof,
const typename FunctionMap< DoFHandlerType::space_dimension >::type &  function_map,
std::map< types::global_dof_index, double > &  boundary_values,
const ComponentMask component_mask = ComponentMask() 
)

Calls the other interpolate_boundary_values() function, see above, with mapping=MappingQGeneric<dim,spacedim>(1). The same comments apply as for the previous function, in particular about the use of the component mask and the requires size of the function object.

template<int dim, int spacedim>
void VectorTools::project_boundary_values ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const typename FunctionMap< spacedim >::type &  boundary_functions,
const Quadrature< dim-1 > &  q,
std::map< types::global_dof_index, double > &  boundary_values,
std::vector< unsigned int >  component_mapping = std::vector< unsigned int >() 
)

Project a function or a set of functions to the boundary of the domain. In other words, compute the solution of the following problem: Find \(u_h \in V_h\) (where \(V_h\) is the finite element space represented by the DoFHandler argument of this function) so that

\begin{align*} \int_{\Gamma} \varphi_i u_h = \sum_{k \in {\cal K}} \int_{\Gamma_k} \varphi_i f_k, \qquad \forall \varphi_i \in V_h \end{align*}

where \(\Gamma = \bigcup_{k \in {\cal K}} \Gamma_k\), \(\Gamma_k \subset \partial\Omega\), \(\cal K\) is the set of indices and \(f_k\) the corresponding boundary functions represented in the function map argument boundary_values to this function, and the integrals are evaluated by quadrature. This problem has a non-unique solution in the interior, but it is well defined for the degrees of freedom on the part of the boundary, \(\Gamma\), for which we do the integration. The values of \(u_h|_\Gamma\), i.e., the nodal values of the degrees of freedom of this function along the boundary, are then what is computed by this function.

Parameters
[in]mappingThe mapping that will be used in the transformations necessary to integrate along the boundary.
[in]dofThe DoFHandler that describes the finite element space and the numbering of degrees of freedom.
[in]boundary_functionsA map from boundary indicators to pointers to functions that describe the desired values on those parts of the boundary marked with this boundary indicator (see Boundary indicator). The projection happens on only those parts of the boundary whose indicators are represented in this map.
[in]qThe face quadrature used in the integration necessary to compute the mass matrix and right hand side of the projection.
[out]boundary_valuesThe result of this function. It is a map containing all indices of degrees of freedom at the boundary (as covered by the boundary parts in boundary_functions) and the computed dof value for this degree of freedom. For each degree of freedom at the boundary, if its index already exists in boundary_values then its boundary value will be overwritten, otherwise a new entry with proper index and boundary value for this degree of freedom will be inserted into boundary_values.
[in]component_mappingIt is sometimes convenient to project a vector-valued function onto only parts of a finite element space (for example, to project a function with dim components onto the velocity components of a dim+1 component DoFHandler for a Stokes problem). To allow for this, this argument allows components to be remapped. If the vector is not empty, it has to have one entry for each vector component of the finite element used in dof. This entry is the component number in boundary_functions that should be used for this component in dof. By default, no remapping is applied.
template<int dim, int spacedim>
void VectorTools::project_boundary_values ( const DoFHandler< dim, spacedim > &  dof,
const typename FunctionMap< spacedim >::type &  boundary_function,
const Quadrature< dim-1 > &  q,
std::map< types::global_dof_index, double > &  boundary_values,
std::vector< unsigned int >  component_mapping = std::vector< unsigned int >() 
)

Calls the project_boundary_values() function, see above, with mapping=MappingQGeneric<dim,spacedim>(1).

template<int dim, int spacedim>
void VectorTools::project_boundary_values ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const typename FunctionMap< spacedim >::type &  boundary_functions,
const hp::QCollection< dim-1 > &  q,
std::map< types::global_dof_index, double > &  boundary_values,
std::vector< unsigned int >  component_mapping = std::vector< unsigned int >() 
)

Same as above, but for objects of type hp::DoFHandler

template<int dim, int spacedim>
void VectorTools::project_boundary_values ( const hp::DoFHandler< dim, spacedim > &  dof,
const typename FunctionMap< spacedim >::type &  boundary_function,
const hp::QCollection< dim-1 > &  q,
std::map< types::global_dof_index, double > &  boundary_values,
std::vector< unsigned int >  component_mapping = std::vector< unsigned int >() 
)

Calls the project_boundary_values() function, see above, with mapping=MappingQGeneric<dim,spacedim>(1).

template<int dim, int spacedim>
void VectorTools::create_right_hand_side ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim > &  q,
const Function< spacedim, double > &  rhs,
Vector< double > &  rhs_vector 
)

Create a right hand side vector. Prior content of the given rhs_vector vector is deleted.

See the general documentation of this namespace for further information.

template<int dim, int spacedim>
void VectorTools::create_right_hand_side ( const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim > &  q,
const Function< spacedim, double > &  rhs,
Vector< double > &  rhs_vector 
)

Calls the create_right_hand_side() function, see above, with mapping=MappingQGeneric<dim>(1).

template<int dim, int spacedim>
void VectorTools::create_right_hand_side ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const hp::QCollection< dim > &  q,
const Function< spacedim, double > &  rhs,
Vector< double > &  rhs_vector 
)

Like the previous set of functions, but for hp objects.

template<int dim, int spacedim>
void VectorTools::create_right_hand_side ( const hp::DoFHandler< dim, spacedim > &  dof,
const hp::QCollection< dim > &  q,
const Function< spacedim, double > &  rhs,
Vector< double > &  rhs_vector 
)

Like the previous set of functions, but for hp objects.

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
Vector< double > &  rhs_vector 
)

Create a right hand side vector for a point source at point p. In other words, it creates a vector \(F\) so that \(F_i = \int_\Omega \delta(x-p) \phi_i(x) dx\). Prior content of the given rhs_vector vector is deleted.

See the general documentation of this namespace for further information.

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
Vector< double > &  rhs_vector 
)

Calls the create_point_source_vector() function, see above, with mapping=MappingQGeneric<dim>(1).

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
Vector< double > &  rhs_vector 
)

Like the previous set of functions, but for hp objects.

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const hp::DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
Vector< double > &  rhs_vector 
)

Like the previous set of functions, but for hp objects. The function uses the default Q1 mapping object. Note that if your hp::DoFHandler uses any active fe index other than zero, then you need to call the function above that provides a mapping object for each active fe index.

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
const Point< dim > &  direction,
Vector< double > &  rhs_vector 
)

Create a right hand side vector for a point source at point p. This variation of the function is meant for vector-valued problems with exactly dim components (it will also work for problems with more than dim components, and in this case simply consider only the first dim components of the shape functions). It computes a right hand side that corresponds to a forcing function that is equal to a delta function times a given direction. In other words, it creates a vector \(F\) so that \(F_i = \int_\Omega [\mathbf d \delta(x-p)] \cdot \phi_i(x) dx\). Note here that \(\phi_i\) is a vector-valued function. \(\mathbf d\) is the given direction of the source term \(\mathbf d \delta(x-p)\) and corresponds to the direction argument to be passed to this function.

Prior content of the given rhs_vector vector is deleted.

See the general documentation of this namespace for further information.

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
const Point< dim > &  direction,
Vector< double > &  rhs_vector 
)

Calls the create_point_source_vector() function for vector-valued finite elements, see above, with mapping=MappingQGeneric<dim>(1).

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
const Point< dim > &  direction,
Vector< double > &  rhs_vector 
)

Like the previous set of functions, but for hp objects.

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const hp::DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
const Point< dim > &  direction,
Vector< double > &  rhs_vector 
)

Like the previous set of functions, but for hp objects. The function uses the default Q1 mapping object. Note that if your hp::DoFHandler uses any active fe index other than zero, then you need to call the function above that provides a mapping object for each active fe index.

template<int dim, int spacedim>
void VectorTools::create_boundary_right_hand_side ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim-1 > &  q,
const Function< spacedim, double > &  rhs,
Vector< double > &  rhs_vector,
const std::set< types::boundary_id > &  boundary_ids = std::set< types::boundary_id >() 
)

Create a right hand side vector from boundary forces. Prior content of the given rhs_vector vector is deleted.

See the general documentation of this namespace for further information.

See also
Glossary entry on boundary indicators
template<int dim, int spacedim>
void VectorTools::create_boundary_right_hand_side ( const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim-1 > &  q,
const Function< spacedim, double > &  rhs,
Vector< double > &  rhs_vector,
const std::set< types::boundary_id > &  boundary_ids = std::set< types::boundary_id >() 
)

Calls the create_boundary_right_hand_side() function, see above, with mapping=MappingQGeneric<dim>(1).

See also
Glossary entry on boundary indicators
template<int dim, int spacedim>
void VectorTools::create_boundary_right_hand_side ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const hp::QCollection< dim-1 > &  q,
const Function< spacedim, double > &  rhs,
Vector< double > &  rhs_vector,
const std::set< types::boundary_id > &  boundary_ids = std::set< types::boundary_id >() 
)

Same as the set of functions above, but for hp objects.

See also
Glossary entry on boundary indicators
template<int dim, int spacedim>
void VectorTools::create_boundary_right_hand_side ( const hp::DoFHandler< dim, spacedim > &  dof,
const hp::QCollection< dim-1 > &  q,
const Function< spacedim, double > &  rhs,
Vector< double > &  rhs_vector,
const std::set< types::boundary_id > &  boundary_ids = std::set< types::boundary_id >() 
)

Calls the create_boundary_right_hand_side() function, see above, with a single Q1 mapping as collection. This function therefore will only work if the only active fe index in use is zero.

See also
Glossary entry on boundary indicators
template<int dim, class InVector , class OutVector , int spacedim>
void VectorTools::integrate_difference ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Function< spacedim, double > &  exact_solution,
OutVector &  difference,
const Quadrature< dim > &  q,
const NormType norm,
const Function< spacedim, double > *  weight = 0,
const double  exponent = 2. 
)

Compute the error of the finite element solution. Integrate the difference between a reference function which is given as a continuous function object, and a finite element function. The result of this function is the vector difference that contains one value per active cell \(K\) of the triangulation. Each of the values of this vector \(d\) equals

\begin{align*} d_K = \| u-u_h \|_X \end{align*}

where \(X\) denotes the norm chosen and \(u\) represents the exact solution.

It is assumed that the number of components of the function exact_solution matches that of the finite element used by dof.

Parameters
[in]mappingThe mapping that is used when integrating the difference \(u-u_h\).
[in]dofThe DoFHandler object that describes the finite element space in which the solution vector lives.
[in]fe_functionA vector with nodal values representing the numerical approximation \(u_h\). This vector needs to correspond to the finite element space represented by dof.
[in]exact_solutionThe exact solution that is used to compute the error.
[out]differenceThe vector of values \(d_K\) computed as above.
[in]qThe quadrature formula used to approximate the integral shown above. Note that some quadrature formulas are more useful than other in integrating \(u-u_h\). For example, it is known that the \(Q_1\) approximation \(u_h\) to the exact solution \(u\) of a Laplace equation is particularly accurate (in fact, superconvergent, i.e. accurate to higher order) at the 4 Gauss points of a cell in 2d (or 8 points in 3d) that correspond to a QGauss(2) object. Consequently, because a QGauss(2) formula only evaluates the two solutions at these particular points, choosing this quadrature formula may indicate an error far smaller than it actually is.
[in]normThe norm \(X\) shown above that should be computed. If the norm is NormType::Hdiv_seminorm, then the finite element on which this function is called needs to have at least dim vector components, and the divergence will be computed on the first div components. This works, for example, on the finite elements used for the mixed Laplace (step-20) and the Stokes equations (step-22).
[in]weightThe additional argument weight allows to evaluate weighted norms. The weight function may be scalar, establishing a spatially variable weight in the domain for all components equally. This may be used, for instance, to only integrate over parts of the domain. The weight function may also be vector-valued, with as many components as the finite element: Then, different components get different weights. A typical application is when the error with respect to only one or a subset of the solution variables is to be computed, in which case the other components would have weight values equal to zero. The ComponentSelectFunction class is particularly useful for this purpose as it provides such a "mask" weight. The weight function is expected to be positive, but negative values are not filtered. The default value of this function, a null pointer, is interpreted as "no weighting function", i.e., weight=1 in the whole domain for all vector components uniformly.
[in]exponentThis value denotes the \(p\) used in computing \(L^p\)-norms and \(W^{1,p}\)-norms. The value is ignores if a norm other than NormType::Lp_norm or NormType::W1p_norm is chosen.

See the general documentation of this namespace for more information.

Note
If the integration here happens over the cells of a parallel::distribute::Triangulation object, then this function computes the vector elements \(d_K\) for an output vector with as many cells as there are active cells of the triangulation object of the current processor. However, not all active cells are in fact locally owned: some may be ghost or artificial cells (see here and here). The vector computed will, in the case of a distributed triangulation, contain zeros for cells that are not locally owned. As a consequence, in order to compute the global \(L_2\) error (for example), the errors from different processors need to be combined, but this is simple because every processor only computes contributions for those cells of the global triangulation it locally owns (and these sets are, by definition, mutually disjoint). Consequently, the following piece of code computes the global \(L_2\) error across multiple processors sharing a parallel::distribute::Triangulation:
Vector<double> local_errors (tria.n_active_cells());
solution, exact_solution,
local_errors,
QGauss<dim>(fe.degree+2),
const double total_local_error = local_errors.l2_norm();
const double total_global_error
= std::sqrt (Utilities::MPI::sum (total_local_error * total_local_error, MPI_COMM_WORLD));
The squaring and taking the square root is necessary in order to compute the sum of squares of norms over all all cells in the definition of the \(L_2\) norm:

\begin{align*} \textrm{error} = \sqrt{\sum_K \|u-u_h\|_{L_2(K)}^2} \end{align*}

Obviously, if you are interested in computing the \(L_1\) norm of the error, the correct form of the last two lines would have been
const double total_local_error = local_errors.l1_norm();
const double total_global_error
= Utilities::MPI::sum (total_local_error, MPI_COMM_WORLD);
instead, and similar considerations hold when computing the \(L_\infty\) norm of the error.

Instantiations for this template are provided for some vector types (see the general documentation of the namespace), but only for InVectors as in the documentation of the namespace, OutVector only Vector<double> and Vector<float>.

template<int dim, class InVector , class OutVector , int spacedim>
void VectorTools::integrate_difference ( const DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Function< spacedim, double > &  exact_solution,
OutVector &  difference,
const Quadrature< dim > &  q,
const NormType norm,
const Function< spacedim, double > *  weight = 0,
const double  exponent = 2. 
)

Calls the integrate_difference() function, see above, with mapping=MappingQGeneric<dim>(1).

template<int dim, class InVector , class OutVector , int spacedim>
void VectorTools::integrate_difference ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Function< spacedim, double > &  exact_solution,
OutVector &  difference,
const hp::QCollection< dim > &  q,
const NormType norm,
const Function< spacedim, double > *  weight = 0,
const double  exponent = 2. 
)

Same as above for hp.

template<int dim, class InVector , class OutVector , int spacedim>
void VectorTools::integrate_difference ( const hp::DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Function< spacedim, double > &  exact_solution,
OutVector &  difference,
const hp::QCollection< dim > &  q,
const NormType norm,
const Function< spacedim, double > *  weight = 0,
const double  exponent = 2. 
)

Calls the integrate_difference() function, see above, with mapping=MappingQGeneric<dim>(1).

template<int dim, typename VectorType , int spacedim>
void VectorTools::point_difference ( const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Function< spacedim, double > &  exact_solution,
Vector< double > &  difference,
const Point< spacedim > &  point 
)

Point error evaluation. Find the first cell containing the given point and compute the difference of a (possibly vector-valued) finite element function and a continuous function (with as many vector components as the finite element) at this point.

This is a wrapper function using a Q1-mapping for cell boundaries to call the other point_difference() function.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
void VectorTools::point_difference ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Function< spacedim, double > &  exact_solution,
Vector< double > &  difference,
const Point< spacedim > &  point 
)

Point error evaluation. Find the first cell containing the given point and compute the difference of a (possibly vector-valued) finite element function and a continuous function (with as many vector components as the finite element) at this point.

Compared with the other function of the same name, this function uses an arbitrary mapping to evaluate the difference.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
void VectorTools::point_value ( const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point,
Vector< double > &  value 
)

Evaluate a possibly vector-valued finite element function defined by the given DoFHandler and nodal vector at the given point, and return the (vector) value of this function through the last argument.

This is a wrapper function using a Q1-mapping for cell boundaries to call the other point_difference() function.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
void VectorTools::point_value ( const hp::DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point,
Vector< double > &  value 
)

Same as above for hp.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
double VectorTools::point_value ( const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point 
)

Evaluate a scalar finite element function defined by the given DoFHandler and nodal vector at the given point, and return the value of this function.

Compared with the other function of the same name, this is a wrapper function using a Q1-mapping for cells.

This function is used in the "Possibilities for extensions" part of the results section of step-3.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
double VectorTools::point_value ( const hp::DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point 
)

Same as above for hp.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
void VectorTools::point_value ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point,
Vector< double > &  value 
)

Evaluate a possibly vector-valued finite element function defined by the given DoFHandler and nodal vector at the given point, and return the (vector) value of this function through the last argument.

Compared with the other function of the same name, this function uses an arbitrary mapping to evaluate the difference.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
void VectorTools::point_value ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point,
Vector< double > &  value 
)

Same as above for hp.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
double VectorTools::point_value ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point 
)

Evaluate a scalar finite element function defined by the given DoFHandler and nodal vector at the given point, and return the value of this function.

Compared with the other function of the same name, this function uses an arbitrary mapping to evaluate the difference.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
double VectorTools::point_value ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point 
)

Same as above for hp.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
void VectorTools::point_gradient ( const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point,
std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > &  value 
)

Evaluate a possibly vector-valued finite element function defined by the given DoFHandler and nodal vector at the given point, and return the (vector) gradient of this function through the last argument.

This is a wrapper function using a Q1-mapping for cell boundaries to call the other point_gradient() function.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
void VectorTools::point_gradient ( const hp::DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point,
std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > &  value 
)

Same as above for hp.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
Tensor<1, spacedim, typename VectorType::value_type> VectorTools::point_gradient ( const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point 
)

Evaluate a scalar finite element function defined by the given DoFHandler and nodal vector at the given point, and return the gradient of this function.

Compared with the other function of the same name, this is a wrapper function using a Q1-mapping for cells.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
Tensor<1, spacedim, typename VectorType::value_type> VectorTools::point_gradient ( const hp::DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point 
)

Same as above for hp.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
void VectorTools::point_gradient ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point,
std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > &  value 
)

Evaluate a possibly vector-valued finite element function defined by the given DoFHandler and nodal vector at the given point, and return the gradients of this function through the last argument.

Compared with the other function of the same name, this function uses an arbitrary mapping for evaluation.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
void VectorTools::point_gradient ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point,
std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > &  value 
)

Same as above for hp.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
Tensor<1, spacedim, typename VectorType::value_type> VectorTools::point_gradient ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point 
)

Evaluate a scalar finite element function defined by the given DoFHandler and nodal vector at the given point, and return the gradient of this function.

Compared with the other function of the same name, this function uses an arbitrary mapping for evaluation.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<int dim, typename VectorType , int spacedim>
Tensor<1, spacedim, typename VectorType::value_type> VectorTools::point_gradient ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim > &  point 
)

Same as above for hp.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
template<typename VectorType >
void VectorTools::subtract_mean_value ( VectorType &  v,
const std::vector< bool > &  p_select = std::vector< bool >() 
)

Mean value operations Subtract the (algebraic) mean value from a vector.

This function is most frequently used as a mean-value filter for Stokes: The pressure in Stokes' equations with only Dirichlet boundaries for the velocities is only determined up to a constant. This function allows to subtract the mean value of the pressure. It is usually called in a preconditioner and generates updates with mean value zero. The mean value is computed as the mean value of the degrees of freedom values as given by the input vector; they are not weighted by the area of cells, i.e. the mean is computed as \(\sum_i v_i\), rather than as \(\int_\Omega v(x) = \int_\Omega \sum_i v_i \phi_i(x)\). The latter can be obtained from the VectorTools::compute_mean_function, however.

Apart from the vector v to operate on, this function takes a boolean mask p_select that has a true entry for every element of the vector for which the mean value shall be computed and later subtracted. The argument is used to denote which components of the solution vector correspond to the pressure, and avoid touching all other components of the vector, such as the velocity components. (Note, however, that the mask is not a GlossComponentMask operating on the vector components of the finite element the solution vector v may be associated with; rather, it is a mask on the entire vector, without reference to what the vector elements mean.)

The boolean mask p_select has an empty vector as default value, which will be interpreted as selecting all vector elements, hence, subtracting the algebraic mean value on the whole vector. This allows to call this function without a boolean mask if the whole vector should be processed.

Note
In the context of using this function to filter out the kernel of an operator (such as the null space of the Stokes operator that consists of the constant pressures), this function only makes sense for finite elements for which the null space indeed consists of the vector \((1,1,\ldots,1)^T\). This is the case for example for the usual Lagrange elements where the sum of all shape functions equals the function that is constant one. However, it is not true for some other functions: for example, for the FE_DGP element (another valid choice for the pressure in Stokes discretizations), the first shape function on each cell is constant while further elements are \(L_2\) orthogonal to it (on the reference cell); consequently, the sum of all shape functions is not equal to one, and the vector that is associated with the constant mode is not equal to \((1,1,\ldots,1)^T\). For such elements, a different procedure has to be used when subtracting the mean value.
template<int dim, typename VectorType , int spacedim>
double VectorTools::compute_mean_value ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim > &  quadrature,
const VectorType &  v,
const unsigned int  component 
)

Compute the mean value of one component of the solution.

This function integrates the chosen component over the whole domain and returns the result, i.e. it computes \(\frac{1}{|\Omega|}\int_\Omega [u_h(x)]_c \; dx\) where \(c\) is the vector component and \(u_h\) is the function representation of the nodal vector given as fourth argument. The integral is evaluated numerically using the quadrature formula given as third argument.

This function is used in the "Possibilities for extensions" part of the results section of step-3.

Note
The function is most often used when solving a problem whose solution is only defined up to a constant, for example a pure Neumann problem or the pressure in a Stokes or Navier-Stokes problem. In both cases, subtracting the mean value as computed by the current function, from the nodal vector does not generally yield the desired result of a finite element function with mean value zero. In fact, it only works for Lagrangian elements. For all other elements, you will need to compute the mean value and subtract it right inside the evaluation routine.
template<int dim, typename VectorType , int spacedim>
double VectorTools::compute_mean_value ( const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim > &  quadrature,
const VectorType &  v,
const unsigned int  component 
)

Calls the other compute_mean_value() function, see above, with mapping=MappingQGeneric<dim>(1).

template<typename DoFHandlerType , typename VectorType >
void VectorTools::get_position_vector ( const DoFHandlerType &  dh,
VectorType &  vector,
const ComponentMask mask = ComponentMask() 
)

Geometrical interpolation Given a DoFHandler containing at least a spacedim vector field, this function interpolates the Triangulation at the support points of a FE_Q() finite element of the same degree as the degree of the required components.

Curved manifold are respected, and the resulting VectorType will be geometrically consistent. The resulting map is guaranteed to be interpolatory at the support points of a FE_Q() finite element of the same degree as the degree of the required components.

If the underlying finite element is an FE_Q(1)^spacedim, then the resulting VectorType is a finite element field representation of the vertices of the Triangulation.

The optional ComponentMask argument can be used to specify what components of the FiniteElement to use to describe the geometry. If no mask is specified at construction time, then a default one is used, i.e., the first spacedim components of the FiniteElement are assumed to represent the geometry of the problem.

This function is only implemented for FiniteElements where the specified components are primitive.

Author
Luca Heltai, 2015
VectorTools::DeclExceptionMsg ( ExcNonInterpolatingFE  ,
"You are attempting an operation that requires the ""finite element involved to be 'interpolating'  ,
i.  e.,
""it needs to have support points.The finite element""you are using here does not appear to have those."   
)

Exception

VectorTools::DeclExceptionMsg ( ExcPointNotAvailableHere  ,
"The given point is inside a cell of a ""parallel::distributed::Triangulation that is not ""locally owned by this processor."   
)

Exception