Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Namespaces | Enumerations | Functions
VectorTools Namespace Reference

Namespaces

namespace  EvaluationFlags
 

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

static ::ExceptionBaseExcPointNotAvailableHere ()
 
template<>
void create_boundary_right_hand_side (const Mapping< 1, 1 > &, const DoFHandler< 1, 1 > &, const Quadrature< 0 > &, const Function< 1 > &, Vector< double > &, const std::set< types::boundary_id > &)
 
template<>
void create_boundary_right_hand_side (const Mapping< 1, 2 > &, const DoFHandler< 1, 2 > &, const Quadrature< 0 > &, const Function< 2 > &, Vector< double > &, const std::set< types::boundary_id > &)
 
Evaluation of functions and errors
template<int dim, typename Number , class OutVector , int spacedim>
void integrate_difference (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const ReadVector< Number > &fe_function, const Function< spacedim, Number > &exact_solution, OutVector &difference, const Quadrature< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=nullptr, const double exponent=2.)
 
template<int dim, typename Number , class OutVector , int spacedim>
void integrate_difference (const DoFHandler< dim, spacedim > &dof, const ReadVector< Number > &fe_function, const Function< spacedim, Number > &exact_solution, OutVector &difference, const Quadrature< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=nullptr, const double exponent=2.)
 
template<int dim, typename Number , class OutVector , int spacedim>
void integrate_difference (const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const ReadVector< Number > &fe_function, const Function< spacedim, Number > &exact_solution, OutVector &difference, const hp::QCollection< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=nullptr, const double exponent=2.)
 
template<int dim, typename Number , class OutVector , int spacedim>
void integrate_difference (const DoFHandler< dim, spacedim > &dof, const ReadVector< Number > &fe_function, const Function< spacedim, Number > &exact_solution, OutVector &difference, const hp::QCollection< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=nullptr, const double exponent=2.)
 
template<int dim, int spacedim, class InVector >
double compute_global_error (const Triangulation< dim, spacedim > &tria, const InVector &cellwise_error, const NormType &norm, const double exponent=2.)
 
template<int dim, typename VectorType , int spacedim>
void point_gradient (const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim, double > &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, double > &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, double > &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 DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim, double > &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, double > &point)
 
template<int dim, typename VectorType , int spacedim>
Tensor< 1, spacedim, typename VectorType::value_type > point_gradient (const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim, double > &point)
 
template<int dim, typename VectorType , int spacedim>
void point_difference (const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Function< spacedim, typename VectorType::value_type > &exact_solution, Vector< typename VectorType::value_type > &difference, const Point< spacedim, double > &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, typename VectorType::value_type > &exact_solution, Vector< typename VectorType::value_type > &difference, const Point< spacedim, double > &point)
 
template<int dim, typename VectorType , int spacedim>
void point_value (const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim, double > &point, Vector< typename VectorType::value_type > &value)
 
template<int dim, typename VectorType , int spacedim>
VectorType::value_type point_value (const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim, double > &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, double > &point, Vector< typename VectorType::value_type > &value)
 
template<int dim, typename VectorType , int spacedim>
void point_value (const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim, double > &point, Vector< typename VectorType::value_type > &value)
 
template<int dim, typename VectorType , int spacedim>
VectorType::value_type point_value (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim, double > &point)
 
template<int dim, typename VectorType , int spacedim>
VectorType::value_type point_value (const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim, double > &point)
 
template<int dim, int spacedim, typename VectorType >
void get_position_vector (const DoFHandler< dim, spacedim > &dh, VectorType &vector, const ComponentMask &mask={})
 
template<int dim, int spacedim, typename VectorType >
void get_position_vector (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dh, VectorType &vector, const ComponentMask &mask={})
 
template<typename VectorType >
void subtract_mean_value (VectorType &v, const std::vector< bool > &p_select={})
 
template<typename VectorType , int dim, int spacedim = dim>
void add_constant (VectorType &solution, const DoFHandler< dim, spacedim > &dof_handler, const unsigned int component, const typename VectorType::value_type constant_adjustment)
 
template<int dim, typename Number , int spacedim>
Number compute_mean_value (const hp::MappingCollection< dim, spacedim > &mapping_collection, const DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim > &q_collection, const ReadVector< Number > &v, const unsigned int component)
 
template<int dim, typename Number , int spacedim>
Number compute_mean_value (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &quadrature, const ReadVector< Number > &v, const unsigned int component)
 
template<int dim, typename Number , int spacedim>
Number compute_mean_value (const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &quadrature, const ReadVector< Number > &v, const unsigned int component)
 
Assembling of right hand sides
template<int dim, int spacedim>
void create_point_source_vector (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim, double > &p, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_point_source_vector (const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim, double > &p, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_point_source_vector (const DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim, double > &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_handler, const Point< spacedim, double > &p, const Point< dim, double > &direction, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_point_source_vector (const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim, double > &p, const Point< dim, double > &direction, Vector< double > &rhs_vector)
 
template<int dim, int spacedim>
void create_point_source_vector (const DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim, double > &p, const Point< dim, double > &direction, Vector< double > &rhs_vector)
 
template<int dim, int spacedim, typename VectorType >
void create_right_hand_side (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &q, const Function< spacedim, typename VectorType::value_type > &rhs, VectorType &rhs_vector, const AffineConstraints< typename VectorType::value_type > &constraints=AffineConstraints< typename VectorType::value_type >())
 
template<int dim, int spacedim, typename VectorType >
void create_right_hand_side (const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &q, const Function< spacedim, typename VectorType::value_type > &rhs, VectorType &rhs_vector, const AffineConstraints< typename VectorType::value_type > &constraints=AffineConstraints< typename VectorType::value_type >())
 
template<int dim, int spacedim, typename VectorType >
void create_right_hand_side (const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim > &q, const Function< spacedim, typename VectorType::value_type > &rhs, VectorType &rhs_vector, const AffineConstraints< typename VectorType::value_type > &constraints=AffineConstraints< typename VectorType::value_type >())
 
template<int dim, int spacedim, typename VectorType >
void create_right_hand_side (const DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim > &q, const Function< spacedim, typename VectorType::value_type > &rhs, VectorType &rhs_vector, const AffineConstraints< typename VectorType::value_type > &constraints=AffineConstraints< typename VectorType::value_type >())
 
template<int dim, int spacedim, typename VectorType >
void create_boundary_right_hand_side (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &q, const Function< spacedim, typename VectorType::value_type > &rhs, VectorType &rhs_vector, const std::set< types::boundary_id > &boundary_ids=std::set< types::boundary_id >())
 
template<int dim, int spacedim, typename VectorType >
void create_boundary_right_hand_side (const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &q, const Function< spacedim, typename VectorType::value_type > &rhs, VectorType &rhs_vector, const std::set< types::boundary_id > &boundary_ids=std::set< types::boundary_id >())
 
template<int dim, int spacedim, typename VectorType >
void create_boundary_right_hand_side (const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim - 1 > &q, const Function< spacedim, typename VectorType::value_type > &rhs, VectorType &rhs_vector, const std::set< types::boundary_id > &boundary_ids=std::set< types::boundary_id >())
 
template<int dim, int spacedim, typename VectorType >
void create_boundary_right_hand_side (const DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim - 1 > &q, const Function< spacedim, typename VectorType::value_type > &rhs, VectorType &rhs_vector, const std::set< types::boundary_id > &boundary_ids=std::set< types::boundary_id >())
 

Interpolation and projection

spacedim & dof_1
 
spacedim const DoFHandler< dim, spacedim > & dof_2
 
spacedim const DoFHandler< dim, spacedim > const FullMatrix< double > & transfer
 
spacedim const DoFHandler< dim, spacedim > const FullMatrix< double > const InVector & data_1
 
spacedim const DoFHandler< dim, spacedim > const FullMatrix< double > const InVector OutVector & data_2
 
template<int dim, int spacedim, typename number >
void interpolate_boundary_values (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &function_map, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask={})
 
template<int dim, int spacedim, typename number >
void interpolate_boundary_values (const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &function_map, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask={})
 
template<int dim, int spacedim, typename number >
void interpolate_boundary_values (const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &function_map, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask={})
 
template<int dim, int spacedim, typename number >
void interpolate_boundary_values (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const types::boundary_id boundary_indicator, const Function< spacedim, number > &boundary_function, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask={})
 
template<int dim, int spacedim, typename number >
void interpolate_boundary_values (const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const types::boundary_id boundary_indicator, const Function< spacedim, number > &boundary_function, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask={})
 
template<int dim, int spacedim, typename number >
void interpolate_boundary_values (const DoFHandler< dim, spacedim > &dof, const types::boundary_id boundary_indicator, const Function< spacedim, number > &boundary_function, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask={})
 
template<int dim, int spacedim, typename number >
void interpolate_boundary_values (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &function_map, AffineConstraints< number > &constraints, const ComponentMask &component_mask={})
 
template<int dim, int spacedim, typename number >
void interpolate_boundary_values (const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &function_map, AffineConstraints< number > &constraints, const ComponentMask &component_mask={})
 
template<int dim, int spacedim, typename number >
void interpolate_boundary_values (const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &function_map, AffineConstraints< number > &constraints, const ComponentMask &component_mask={})
 
template<int dim, int spacedim, typename number >
void interpolate_boundary_values (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const types::boundary_id boundary_indicator, const Function< spacedim, number > &boundary_function, AffineConstraints< number > &constraints, const ComponentMask &component_mask={})
 
template<int dim, int spacedim, typename number >
void interpolate_boundary_values (const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const types::boundary_id boundary_indicator, const Function< spacedim, number > &boundary_function, AffineConstraints< number > &constraints, const ComponentMask &component_mask={})
 
template<int dim, int spacedim, typename number >
void interpolate_boundary_values (const DoFHandler< dim, spacedim > &dof, const types::boundary_id boundary_indicator, const Function< spacedim, number > &boundary_function, AffineConstraints< number > &constraints, const ComponentMask &component_mask={})
 
template<int dim, int spacedim, typename number >
void project_boundary_values (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &boundary_functions, const Quadrature< dim - 1 > &q, std::map< types::global_dof_index, number > &boundary_values, std::vector< unsigned int > component_mapping={})
 
template<int dim, int spacedim, typename number >
void project_boundary_values (const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &boundary_function, const Quadrature< dim - 1 > &q, std::map< types::global_dof_index, number > &boundary_values, std::vector< unsigned int > component_mapping={})
 
template<int dim, int spacedim, typename number >
void project_boundary_values (const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &boundary_functions, const hp::QCollection< dim - 1 > &q, std::map< types::global_dof_index, number > &boundary_values, std::vector< unsigned int > component_mapping={})
 
template<int dim, int spacedim, typename number >
void project_boundary_values (const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &boundary_function, const hp::QCollection< dim - 1 > &q, std::map< types::global_dof_index, number > &boundary_values, std::vector< unsigned int > component_mapping={})
 
template<int dim, int spacedim, typename number >
void project_boundary_values (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &boundary_functions, const Quadrature< dim - 1 > &q, AffineConstraints< number > &constraints, std::vector< unsigned int > component_mapping={})
 
template<int dim, int spacedim, typename number >
void project_boundary_values (const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &boundary_function, const Quadrature< dim - 1 > &q, AffineConstraints< number > &constraints, std::vector< unsigned int > component_mapping={})
 
template<int dim, typename number >
void project_boundary_values_curl_conforming_l2 (const DoFHandler< dim, dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, number > &boundary_function, const types::boundary_id boundary_component, AffineConstraints< number > &constraints, const Mapping< dim > &mapping)
 
template<int dim, typename number >
void project_boundary_values_curl_conforming_l2 (const DoFHandler< dim, dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, number > &boundary_function, const types::boundary_id boundary_component, AffineConstraints< number > &constraints, const hp::MappingCollection< dim, dim > &mapping_collection=hp::StaticMappingQ1< dim >::mapping_collection)
 
template<int dim, typename number , typename number2 = number>
void project_boundary_values_div_conforming (const DoFHandler< dim, dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, number2 > &boundary_function, const types::boundary_id boundary_component, AffineConstraints< number > &constraints, const Mapping< dim > &mapping)
 
template<int dim, typename number , typename number2 = number>
void project_boundary_values_div_conforming (const DoFHandler< dim, dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, number2 > &boundary_function, const types::boundary_id boundary_component, AffineConstraints< number > &constraints, const hp::MappingCollection< dim, dim > &mapping_collection=hp::StaticMappingQ1< dim >::mapping_collection)
 
template<int dim, int spacedim>
void compute_nonzero_normal_flux_constraints (const DoFHandler< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, const std::map< types::boundary_id, const Function< spacedim, double > * > &function_map, AffineConstraints< double > &constraints, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
 
template<int dim, int spacedim>
void compute_nonzero_normal_flux_constraints_on_level (const DoFHandler< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, const std::map< types::boundary_id, const Function< spacedim, double > * > &function_map, AffineConstraints< double > &constraints, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()), const IndexSet &refinement_edge_indices=IndexSet(), const unsigned int level=numbers::invalid_unsigned_int)
 
template<int dim, int spacedim>
void compute_no_normal_flux_constraints (const DoFHandler< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, AffineConstraints< double > &constraints, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
 
template<int dim, int spacedim>
void compute_no_normal_flux_constraints_on_level (const DoFHandler< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, AffineConstraints< double > &constraints, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()), const IndexSet &refinement_edge_indices=IndexSet(), const unsigned int level=numbers::invalid_unsigned_int)
 
template<int dim, int spacedim>
void compute_nonzero_tangential_flux_constraints (const DoFHandler< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, const std::map< types::boundary_id, const Function< spacedim, double > * > &function_map, AffineConstraints< double > &constraints, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
 
template<int dim, int spacedim>
void compute_normal_flux_constraints (const DoFHandler< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, AffineConstraints< double > &constraints, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
 
template<int dim, int spacedim, typename VectorType >
void interpolate (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Function< spacedim, typename VectorType::value_type > &function, VectorType &vec, const ComponentMask &component_mask={})
 
template<int dim, int spacedim, typename VectorType >
void interpolate (const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Function< spacedim, typename VectorType::value_type > &function, VectorType &vec, const ComponentMask &component_mask={})
 
template<int dim, int spacedim, typename VectorType >
void interpolate (const DoFHandler< dim, spacedim > &dof, const Function< spacedim, typename VectorType::value_type > &function, VectorType &vec, const ComponentMask &component_mask={})
 
template<int dim, class InVector , class OutVector , int spacedim>
 DEAL_II_CXX20_REQUIRES (concepts::is_dealii_vector_type< InVector > &&concepts::is_writable_dealii_vector_type< OutVector >) void interpolate(const DoFHandler< dim
 
template<int dim, int spacedim, typename VectorType >
void interpolate_based_on_material_id (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const std::map< types::material_id, const Function< spacedim, typename VectorType::value_type > * > &function_map, VectorType &dst, const ComponentMask &component_mask={})
 
template<int dim, int spacedim, typename VectorType >
void interpolate_to_different_mesh (const DoFHandler< dim, spacedim > &dof_handler_1, const VectorType &u1, const DoFHandler< dim, spacedim > &dof_handler_2, VectorType &u2)
 
template<int dim, int spacedim, typename VectorType >
void interpolate_to_different_mesh (const DoFHandler< dim, spacedim > &dof_handler_1, const VectorType &u1, const DoFHandler< dim, spacedim > &dof_handler_2, const AffineConstraints< typename VectorType::value_type > &constraints, VectorType &u2)
 
template<int dim, int spacedim, typename VectorType >
void interpolate_to_different_mesh (const InterGridMap< DoFHandler< dim, spacedim > > &intergridmap, const VectorType &u1, const AffineConstraints< typename VectorType::value_type > &constraints, VectorType &u2)
 
template<int dim, int spacedim, typename VectorType >
void interpolate_to_coarser_mesh (const DoFHandler< dim, spacedim > &dof_handler_fine, const VectorType &u_fine, const DoFHandler< dim, spacedim > &dof_handler_coarse, const AffineConstraints< typename VectorType::value_type > &constraints_coarse, VectorType &u_coarse)
 
template<int dim, int spacedim, typename VectorType >
void interpolate_to_finer_mesh (const DoFHandler< dim, spacedim > &dof_handler_coarse, const VectorType &u_coarse, const DoFHandler< dim, spacedim > &dof_handler_fine, const AffineConstraints< typename VectorType::value_type > &constraints_fine, VectorType &u_fine)
 
template<int dim, typename VectorType , int spacedim>
void project (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const AffineConstraints< typename VectorType::value_type > &constraints, const Quadrature< dim > &quadrature, const Function< spacedim, typename VectorType::value_type > &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 AffineConstraints< typename VectorType::value_type > &constraints, const Quadrature< dim > &quadrature, const Function< spacedim, typename VectorType::value_type > &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 DoFHandler< dim, spacedim > &dof, const AffineConstraints< typename VectorType::value_type > &constraints, const hp::QCollection< dim > &quadrature, const Function< spacedim, typename VectorType::value_type > &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 DoFHandler< dim, spacedim > &dof, const AffineConstraints< typename VectorType::value_type > &constraints, const hp::QCollection< dim > &quadrature, const Function< spacedim, typename VectorType::value_type > &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 Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const AffineConstraints< typename VectorType::value_type > &constraints, const Quadrature< dim > &quadrature, const std::function< typename VectorType::value_type(const typename DoFHandler< dim, spacedim >::active_cell_iterator &, const unsigned int)> &func, VectorType &vec_result)
 
template<int dim, typename VectorType >
void project (std::shared_ptr< const MatrixFree< dim, typename VectorType::value_type, VectorizedArray< typename VectorType::value_type > > > data, const AffineConstraints< typename VectorType::value_type > &constraints, const unsigned int n_q_points_1d, const std::function< VectorizedArray< typename VectorType::value_type >(const unsigned int, const unsigned int)> &func, VectorType &vec_result, const unsigned int fe_component=0)
 
template<int dim, typename VectorType >
void project (std::shared_ptr< const MatrixFree< dim, typename VectorType::value_type, VectorizedArray< typename VectorType::value_type > > > data, const AffineConstraints< typename VectorType::value_type > &constraints, const std::function< VectorizedArray< typename VectorType::value_type >(const unsigned int, const unsigned int)> &func, VectorType &vec_result, const unsigned int fe_component=0)
 

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 MappingQ(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).

Enumeration Type Documentation

◆ NormType

Denote which norm/integral is to be computed by the integrate_difference() function on each cell and compute_global_error() for the whole domain. Let \(f:\Omega \rightarrow \mathbb{R}^c\) be a finite element function with \(c\) components where component \(c\) is denoted by \(f_c\) and \(\hat{f}\) be the reference function (the fe_function and exact_solution arguments to integrate_difference()). Let \(e_c = \hat{f}_c - f_c\) be the difference or error between the two. Further, let \(w:\Omega \rightarrow \mathbb{R}^c\) be the weight function of integrate_difference(), which is assumed to be equal to one if not supplied. Finally, let \(p\) be the exponent argument (for \(L_p\)-norms).

In the following,we denote by \(E_K\) the local error computed by integrate_difference() on cell \(K\), whereas \(E\) is the global error computed by compute_global_error(). Note that integrals are approximated by quadrature in the usual way:

\[ \int_A f(x) dx \approx \sum_q f(x_q) \omega_q. \]

Similarly for suprema over a cell \(T\):

\[ \sup_{x\in T} |f(x)| dx \approx \max_q |f(x_q)|. \]

Enumerator
mean 

The function or difference of functions is integrated on each cell \(K\):

\[ E_K = \int_K \sum_c (\hat{f}_c - f_c) \, w_c = \int_K \sum_c e_c \, w_c \]

and summed up to get

\[ E = \sum_K E_K = \int_\Omega \sum_c (\hat{f}_c - f_c) \, w_c \]

or, for \(w \equiv 1\):

\[ E = \int_\Omega (\hat{f} - f) = \int_\Omega e. \]

Note: This differs from what is typically known as the mean of a function by a factor of \(\frac{1}{|\Omega|}\). To compute the mean you can also use compute_mean_value(). Finally, pay attention to the sign: if \(\hat{f}=0\), this will compute the negative of the mean of \(f\).

L1_norm 

The absolute value of the function is integrated:

\[ E_K = \int_K \sum_c |e_c| \, w_c \]

and

\[ E = \sum_K E_K = \int_\Omega \sum_c |e_c| w_c, \]

or, for \(w \equiv 1\):

\[ E = \| e \|_{L^1}. \]

L2_norm 

The square of the function is integrated and the square root of the result is computed on each cell:

\[ E_K = \sqrt{ \int_K \sum_c e_c^2 \, w_c } \]

and

\[ E = \sqrt{\sum_K E_K^2} = \sqrt{ \int_\Omega \sum_c e_c^2 \, w_c } \]

or, for \(w \equiv 1\):

\[ E = \sqrt{ \int_\Omega e^2 } = \| e \|_{L^2} \]

Lp_norm 

The absolute value to the \(p\)-th power is integrated and the \(p\)-th root is computed on each cell. The exponent \(p\) is the exponent argument of integrate_difference() and compute_mean_value():

\[ E_K = \left( \int_K \sum_c |e_c|^p \, w_c \right)^{1/p} \]

and

\[ E = \left( \sum_K E_K^p \right)^{1/p} \]

or, for \(w \equiv 1\):

\[ E = \| e \|_{L^p}. \]

Linfty_norm 

The maximum absolute value of the function:

\[ E_K = \sup_K \max_c |e_c| \, w_c \]

and

\[ E = \max_K E_K = \sup_\Omega \max_c |e_c| \, w_c \]

or, for \(w \equiv 1\):

\[ E = \sup_\Omega \|e\|_\infty = \| e \|_{L^\infty}. \]

H1_seminorm 

L2_norm of the gradient:

\[ E_K = \sqrt{ \int_K \sum_c (\nabla e_c)^2 \, w_c } \]

and

\[ E = \sqrt{\sum_K E_K^2} = \sqrt{ \int_\Omega \sum_c (\nabla e_c)^2 \, w_c } \]

or, for \(w \equiv 1\):

\[ E = \| \nabla e \|_{L^2}. \]

Hdiv_seminorm 

L2_norm of the divergence of a vector field. The function \(f\) is expected to have \(c \geq \text{dim}\) components and the first dim will be used to compute the divergence:

\[ E_K = \sqrt{ \int_K \left( \sum_c \frac{\partial e_c}{\partial x_c} \, \sqrt{w_c} \right)^2 } \]

and

\[ E = \sqrt{\sum_K E_K^2} = \sqrt{ \int_\Omega \left( \sum_c \frac{\partial e_c}{\partial x_c} \, \sqrt{w_c} \right)^2 } \]

or, for \(w \equiv 1\):

\[ E = \| \nabla \cdot e \|_{L^2}. \]

H1_norm 

The square of this norm is the square of the L2_norm plus the square of the H1_seminorm:

\[ E_K = \sqrt{ \int_K \sum_c (e_c^2 + (\nabla e_c)^2) \, w_c } \]

and

\[ E = \sqrt{\sum_K E_K^2} = \sqrt{ \int_\Omega \sum_c (e_c^2 + (\nabla e_c)^2) \, w_c } \]

or, for \(w \equiv 1\):

\[ E = \left( \| e \|_{L^2}^2 + \| \nabla e \|_{L^2}^2 \right)^{1/2}. \]

W1p_seminorm 

Lp_norm of the gradient:

\[ E_K = \left( \int_K \sum_c |\nabla e_c|^p \, w_c \right)^{1/p} \]

and

\[ E = \left( \sum_K E_K^p \right)^{1/p} = \left( \int_\Omega \sum_c |\nabla e_c|^p \, w_c \right)^{1/p} \]

or, for \(w \equiv 1\):

\[ E = \| \nabla e \|_{L^p}. \]

W1p_norm 

The same as the H1_norm but using Lp:

\[ E_K = \left( \int_K \sum_c (|e_c|^p + |\nabla e_c|^p) \, w_c \right)^{1/p} \]

and

\[ E = \left( \sum_K E_K^p \right)^{1/p} = \left( \int_\Omega \sum_c (|e_c|^p + |\nabla e_c|^p) \, w_c \right)^{1/p} \]

or, for \(w \equiv 1\):

\[ E = \left( \| e \|_{L^p}^p + \| \nabla e \|_{L^p}^p \right)^{1/p}. \]

W1infty_seminorm 

Linfty_norm of the gradient:

\[ E_K = \sup_K \max_c |\nabla e_c| \, w_c \]

and

\[ E = \max_K E_K = \sup_\Omega \max_c |\nabla e_c| \, w_c \]

or, for \(w \equiv 1\):

\[ E = \| \nabla e \|_{L^\infty}. \]

W1infty_norm 

The sum of Linfty_norm and W1infty_seminorm:

\[ E_K = \sup_K \max_c |e_c| \, w_c + \sup_K \max_c |\nabla e_c| \, w_c. \]

The global norm is not implemented in compute_global_error(), because it is impossible to compute the sum of the global norms from the values \(E_K\). As a work-around, you can compute the global Linfty_norm and W1infty_seminorm separately and then add them to get (with \(w \equiv 1\)):

\[ E = \| e \|_{L^\infty} + \| \nabla e \|_{L^\infty}. \]

Definition at line 52 of file vector_tools_common.h.

Function Documentation

◆ interpolate_boundary_values() [1/6]

template<int dim, int spacedim, typename number >
void VectorTools::interpolate_boundary_values ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const std::map< types::boundary_id, const Function< spacedim, number > * > &  function_map,
std::map< types::global_dof_index, number > &  boundary_values,
const ComponentMask component_mask = {} 
)

Compute constraints on the solution that corresponds to the imposition of Dirichlet boundary conditions on parts of the boundary. This function creates 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, its boundary value will be overwritten if its index already exists in boundary_values. 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 key of this map corresponds 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. For an example of how to use this argument with a non-empty map, see the step-16 tutorial program.

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.

Note
This function applies the same ComponentMask to all parts of the boundary indicated by keys of the function_map argument. If you want to apply different component masks to parts of the boundary represented by different boundary indicators, this function needs to be called multiple times. For performance reasons, it might be reasonable to use the present function by grouping together all boundary indicators with the same ComponentMask. An alternative is to use one of the other functions with this name, which take only one boundary indicator with corresponding boundary function, to be called separately for every boundary indicator.
Mathematically, boundary conditions can only be applied to a part of the boundary that has a nonzero \((d-1)\)-dimensional measure; in other words, it must be the union of faces of a mesh, rather than a set of edges in 3d, or even just a few vertices. That is because applying boundary conditions on individual vertices (rather than on the entire face of which this vertex might be a part) would correspond to using Dirac delta functions as boundary values, and this generally leads to singular solutions that can not adequately be resolved with the finite element method. These considerations notwithstanding, people often do apply boundary conditions at individual vertices – in particular in solid mechanics, where one would then impose constraints on one or all components of the displacement at a vertex. This function does not support this operation: It works solely by looping over faces, checking whether the boundary indicator of the face is one of the ones of interest, and then considers all of the degrees of freedom on the face; it does not consider vertices (or, in 3d, edges) separately from the faces they are part of. But you can impose constraints on individual vertices by looping over all cells, over all vertices of each cell, and identifying whether this is the vertex you care about; then you use DoFAccessor::vertex_dof_index() to obtain the indices of the DoFs located on it. You can then entries for these degrees of freedom by hand to the std::map or AffineConstraints object you typically use to represent boundary value constraints.
When solving a partial differential equation with boundary conditions \(u|_{\partial\Omega}=g\) (on the entire boundary \(\partial\Omega\), or perhaps only on parts \(\Gamma\subset\partial\Omega\) of the boundary), then this boundary condition is in general not satisfiable exactly using finite elements in the form \(u_h|_{\partial\Omega}=g\). That is because the function \(g\) is generally not a polynomial, whereas \(u_h|_{\partial\Omega}\) is a polynomial on each face of the mesh that is located at the boundary. In other words, it is in general not possible to impose such boundary condition; what one can do, however, is to impose

\[ u_h|_{\partial\Omega}=I_h^{\partial\Omega} g, \]

where \(I_h^{\partial\Omega} g\) is a function that equals \(g\) at each node of the finite element space located on the boundary, and is piecewise polynomial in between. In other words, \(I_h^{\partial\Omega}\) is an interpolation operator and \(I_h^{\partial\Omega} g\) are the interpolated boundary values – thus the name. The use of \(I_h^{\partial\Omega} g\) instead of \(g\) as boundary values imposes an additional error (in the same spirit as using quadrature introduces an additional error compared to being able to compute the integrals of the weak form exactly). In most cases, this additional error is of the same order as the other error terms in the finite element method, though there are some subtle differences when measuring the error in the \(L^2\) norm. For some details, see [14] .
An alternative to using the interpolant,

\[ u_h|_{\partial\Omega}=I_h^{\partial\Omega} g \]

is to use the projection of the boundary values \(g\) onto the finite element space on the boundary:

\[ u_h|_{\partial\Omega}=\Pi_h^{\partial\Omega} g. \]

The projection is available using the project_boundary_values() function. Using the projection may have some theoretical advantages (see again [14]) but has the practical disadvantage that computing the projection is far more expensive than computing the interpolation because the latter can be done one face at a time whereas the projection requires the solution of a problem on the entire boundary. On the other hand, interpolation is only possible for "nodal" finite element spaces (such as FE_Q, but not FE_Q_Hierarchical), whereas the projection is always possible.

See the general documentation of this namespace for more information.

◆ interpolate_boundary_values() [2/6]

template<int dim, int spacedim, typename number >
void VectorTools::interpolate_boundary_values ( const hp::MappingCollection< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const std::map< types::boundary_id, const Function< spacedim, number > * > &  function_map,
std::map< types::global_dof_index, number > &  boundary_values,
const ComponentMask component_mask = {} 
)

Like the previous function, but take a mapping collection to go with DoFHandler objects with hp-capabilities.

◆ interpolate_boundary_values() [3/6]

template<int dim, int spacedim, typename number >
void VectorTools::interpolate_boundary_values ( const DoFHandler< dim, spacedim > &  dof,
const std::map< types::boundary_id, const Function< spacedim, number > * > &  function_map,
std::map< types::global_dof_index, number > &  boundary_values,
const ComponentMask component_mask = {} 
)

Like the previous functions but without Mapping argument, using mapping=MappingQ<dim,spacedim>(1) internally.

◆ interpolate_boundary_values() [4/6]

template<int dim, int spacedim, typename number >
void VectorTools::interpolate_boundary_values ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const types::boundary_id  boundary_indicator,
const Function< spacedim, number > &  boundary_function,
std::map< types::global_dof_index, number > &  boundary_values,
const ComponentMask component_mask = {} 
)

Take only one boundary indicator with corresponding boundary function.

See also
Glossary entry on boundary indicators

◆ interpolate_boundary_values() [5/6]

template<int dim, int spacedim, typename number >
void VectorTools::interpolate_boundary_values ( const hp::MappingCollection< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const types::boundary_id  boundary_indicator,
const Function< spacedim, number > &  boundary_function,
std::map< types::global_dof_index, number > &  boundary_values,
const ComponentMask component_mask = {} 
)

Like the previous function, but take a mapping collection to go with DoFHandler objects with hp-capabilities.

◆ interpolate_boundary_values() [6/6]

template<int dim, int spacedim, typename number >
void VectorTools::interpolate_boundary_values ( const DoFHandler< dim, spacedim > &  dof,
const types::boundary_id  boundary_indicator,
const Function< spacedim, number > &  boundary_function,
std::map< types::global_dof_index, number > &  boundary_values,
const ComponentMask component_mask = {} 
)

Like the previous functions but without Mapping argument, using mapping=MappingQ<dim,spacedim>(1) internally.

See also
Glossary entry on boundary indicators

◆ project_boundary_values() [1/4]

template<int dim, int spacedim, typename number >
void VectorTools::project_boundary_values ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const std::map< types::boundary_id, const Function< spacedim, number > * > &  boundary_functions,
const Quadrature< dim - 1 > &  q,
std::map< types::global_dof_index, number > &  boundary_values,
std::vector< unsigned int component_mapping = {} 
)

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.

In case this function is used with \(H_{div}\) conforming finite element space, the solution of a different problem is computed, namely: Find \(\vec{u}_h \in V_h \subset H(\text{div}; \Omega)\) so that

\begin{align*} \int_{\Gamma} (\vec{\varphi}_i \cdot \vec{n}) (\vec{u}_h \cdot \vec{n}) = \sum_{k \in {\cal K}} \int_{\Gamma_k} (\vec{\varphi}_i \cdot \vec{n}) (\vec{f}_k \cdot \vec{n}), \qquad \forall \vec{\varphi_i} \in V_h, \end{align*}

where \(\vec{n}\) is an outward normal vector.

This function throws an exception if used with \(H_\text{curl}\) conforming elements, so the project_boundary_values_curl_conforming_l2() should be used instead.

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.
Note
Using the projection rather than the interpolation of boundary values makes relatively little difference in practice. That said, it is far more computationally expensive to compute projections because the require the solution of a problem that couples all unknowns on the boundary, whereas interpolation works on one face at a time. On the other hand, interpolation is only possible for "nodal" finite element spaces (such as FE_Q, but not FE_Q_Hierarchical), whereas the projection is always possible. (For some more theoretical considerations, see the documentation of the first interpolate_boundary_values() function above.)

◆ project_boundary_values() [2/4]

template<int dim, int spacedim, typename number >
void VectorTools::project_boundary_values ( const DoFHandler< dim, spacedim > &  dof,
const std::map< types::boundary_id, const Function< spacedim, number > * > &  boundary_function,
const Quadrature< dim - 1 > &  q,
std::map< types::global_dof_index, number > &  boundary_values,
std::vector< unsigned int component_mapping = {} 
)

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

◆ project_boundary_values() [3/4]

template<int dim, int spacedim, typename number >
void VectorTools::project_boundary_values ( const hp::MappingCollection< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const std::map< types::boundary_id, const Function< spacedim, number > * > &  boundary_functions,
const hp::QCollection< dim - 1 > &  q,
std::map< types::global_dof_index, number > &  boundary_values,
std::vector< unsigned int component_mapping = {} 
)

Same as above, but with hp-capabilities.

◆ project_boundary_values() [4/4]

template<int dim, int spacedim, typename number >
void VectorTools::project_boundary_values ( const DoFHandler< dim, spacedim > &  dof,
const std::map< types::boundary_id, const Function< spacedim, number > * > &  boundary_function,
const hp::QCollection< dim - 1 > &  q,
std::map< types::global_dof_index, number > &  boundary_values,
std::vector< unsigned int component_mapping = {} 
)

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

◆ integrate_difference() [1/4]

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

Compute the cellwise 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.

To compute a global error norm of a finite element solution, use VectorTools::compute_global_error() with the output vector computed with this function.

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 others 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 ignored if a norm other than NormType::Lp_norm, NormType::W1p_norm, or NormType::W1p_seminorm 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, see VectorTools::compute_global_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>.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ integrate_difference() [2/4]

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

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

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ integrate_difference() [3/4]

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

Same as above for hp.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ integrate_difference() [4/4]

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

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

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ compute_global_error()

template<int dim, int spacedim, class InVector >
double VectorTools::compute_global_error ( const Triangulation< dim, spacedim > &  tria,
const InVector &  cellwise_error,
const NormType norm,
const double  exponent = 2. 
)

Take a Vector cellwise_error of errors on each cell with tria.n_active_cells() entries and return the global error as given by norm.

The cellwise_error vector is typically an output produced by VectorTools::integrate_difference() and you normally want to supply the same value for norm as you used in VectorTools::integrate_difference().

If the given Triangulation is a parallel::TriangulationBase, entries in cellwise_error that do not correspond to locally owned cells are assumed to be 0.0 and a parallel reduction using MPI is done to compute the global error.

Parameters
triaThe Triangulation with active cells corresponding with the entries in cellwise_error.
cellwise_errorVector of errors on each active cell.
normThe type of norm to compute.
exponentThe exponent \(p\) to use for \(L^p\)-norms and \(W^{1,p}\)-norms. The value is ignored if a norm other than NormType::Lp_norm, NormType::W1p_norm, or NormType::W1p_seminorm is chosen.
Note
Instantiated for type Vector<double> and Vector<float>.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ interpolate() [1/3]

template<int dim, int spacedim, typename VectorType >
void VectorTools::interpolate ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const Function< spacedim, typename VectorType::value_type > &  function,
VectorType &  vec,
const ComponentMask component_mask = {} 
)

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.

See the general documentation of this namespace for further information.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ interpolate() [2/3]

template<int dim, int spacedim, typename VectorType >
void VectorTools::interpolate ( const hp::MappingCollection< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const Function< spacedim, typename VectorType::value_type > &  function,
VectorType &  vec,
const ComponentMask component_mask = {} 
)

Same as above but in an hp-context.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ interpolate() [3/3]

template<int dim, int spacedim, typename VectorType >
void VectorTools::interpolate ( const DoFHandler< dim, spacedim > &  dof,
const Function< spacedim, typename VectorType::value_type > &  function,
VectorType &  vec,
const ComponentMask component_mask = {} 
)

Call the interpolate() function above with mapping=MappingQ<dim,spacedim>(1).

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ DEAL_II_CXX20_REQUIRES()

template<int dim, class InVector , class OutVector , int spacedim>
VectorTools::DEAL_II_CXX20_REQUIRES ( concepts::is_dealii_vector_type< InVector > &&concepts::is_writable_dealii_vector_type< OutVector >  ) const

Interpolate different finite element spaces. The interpolation of vector data_1 (which is assumed to be ghosted, see GlossGhostedVector) 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.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ interpolate_based_on_material_id()

template<int dim, int spacedim, typename VectorType >
void VectorTools::interpolate_based_on_material_id ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof_handler,
const std::map< types::material_id, const Function< spacedim, typename VectorType::value_type > * > &  function_map,
VectorType &  dst,
const ComponentMask component_mask = {} 
)

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

Parameters
[in]mappingThe mapping to use to determine the location of support points at which the functions are to be evaluated.
[in]dof_handlerDoFHandler initialized with Triangulation and FiniteElement objects and that defines the finite element space.
[in]function_mapA std::map reflecting the correspondence between material ids on those cells on which something should be interpolated, and the functions to be interpolated onto the finite element space.
[out]dstThe global finie element vector holding the output of the interpolated values.
[in]component_maskA mask of components that shall be interpolated.
Note
If the algorithm encounters a cell whose material id is not listed in the given function_map, then dst will not be updated in the respective degrees of freedom of the output vector. For example, if dst was initialized to zero, then those zeros which correspond to the missed material ids will still remain in dst 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 the order of cells is somewhat 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.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ interpolate_to_different_mesh() [1/3]

template<int dim, int spacedim, typename VectorType >
void VectorTools::interpolate_to_different_mesh ( const DoFHandler< dim, spacedim > &  dof_handler_1,
const VectorType &  u1,
const DoFHandler< dim, spacedim > &  dof_handler_2,
VectorType &  u2 
)

Compute the interpolation of a dof_handler_1 function u1 (i.e., a function defined on the mesh that underlies dof_handler_1, using the finite element associated with that DoFHandler) to a dof_handler_2 function u2 (i.e., a function defined on the mesh and finite element associated with dof_handler_2). This function assumes that dof_handler_1 and dof_handler_2 are defined on (possibly different) triangulations that have 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 AffineConstraints argument, see below, or make the field conforming yourself by calling the AffineConstraints::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). In practice, this is rarely the case because two triangulations, partitioned in their own ways, will not typically have corresponding cells owned by the same process, and implementing the interpolation procedure would require transfering data between processes in ways that are difficult to implement efficiently. However, some special cases can more easily be implemented, namely the case where one of the meshes is strictly coarser or finer than the other. For these cases, see the interpolate_to_coarser_mesh() and interpolate_to_finer_mesh().
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ interpolate_to_different_mesh() [2/3]

template<int dim, int spacedim, typename VectorType >
void VectorTools::interpolate_to_different_mesh ( const DoFHandler< dim, spacedim > &  dof_handler_1,
const VectorType &  u1,
const DoFHandler< dim, spacedim > &  dof_handler_2,
const AffineConstraints< typename VectorType::value_type > &  constraints,
VectorType &  u2 
)

This is a variation of the previous function that takes an additional constraint object as argument.

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.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ interpolate_to_different_mesh() [3/3]

template<int dim, int spacedim, typename VectorType >
void VectorTools::interpolate_to_different_mesh ( const InterGridMap< DoFHandler< dim, spacedim > > &  intergridmap,
const VectorType &  u1,
const AffineConstraints< typename VectorType::value_type > &  constraints,
VectorType &  u2 
)

The same function as above, but takes an InterGridMap object directly as a parameter. This function is 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.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ interpolate_to_coarser_mesh()

template<int dim, int spacedim, typename VectorType >
void VectorTools::interpolate_to_coarser_mesh ( const DoFHandler< dim, spacedim > &  dof_handler_fine,
const VectorType &  u_fine,
const DoFHandler< dim, spacedim > &  dof_handler_coarse,
const AffineConstraints< typename VectorType::value_type > &  constraints_coarse,
VectorType &  u_coarse 
)

This function addresses one of the limitations of the interpolate_to_different_mesh() function, namely that the latter only works on parallel triangulations in very specific cases where both triangulations involved happen to be partitioned in such a way that if a process owns a cell of one mesh, it also needs to own the corresponding parent of child cells of the other mesh. In practice, this is rarely the case.

This function does not have this restriction, and consequently also works in parallel, as long as the mesh we interpolate onto can be obtained from the mesh we interpolate off by coarsening. In other words, the target mesh is coarser than the source mesh.

The function takes an additional constraints argument that is used after interpolation to ensure that the output vector is conforming (that is, that all entries of the output vector conform to their constraints). constraints_coarse therefore needs to correspond to dof_handler_coarse .

The opposite operation, interpolation from a coarser to a finer mesh, is implemented in the interpolate_to_finer_mesh() function.

◆ interpolate_to_finer_mesh()

template<int dim, int spacedim, typename VectorType >
void VectorTools::interpolate_to_finer_mesh ( const DoFHandler< dim, spacedim > &  dof_handler_coarse,
const VectorType &  u_coarse,
const DoFHandler< dim, spacedim > &  dof_handler_fine,
const AffineConstraints< typename VectorType::value_type > &  constraints_fine,
VectorType &  u_fine 
)

This function addresses one of the limitations of the interpolate_to_different_mesh() function, namely that the latter only works on parallel triangulations in very specific cases where both triangulations involved happen to be partitioned in such a way that if a process owns a cell of one mesh, it also needs to own the corresponding parent of child cells of the other mesh. In practice, this is rarely the case.

This function does not have this restriction, and consequently also works in parallel, as long as the mesh we interpolate onto can be obtained from the mesh we interpolate off by refinement. In other words, the target mesh is finer than the source mesh.

The function takes an additional constraints argument that is used after interpolation to ensure that the output vector is conforming (that is, that all entries of the output vector conform to their constraints). constraints_finee therefore needs to correspond to dof_handler_fine .

The opposite operation, interpolation from a finer to a coarser mesh, is implemented in the interpolate_to_coarser_mesh() function.

◆ get_position_vector() [1/2]

template<int dim, int spacedim, typename VectorType >
void VectorTools::get_position_vector ( const DoFHandler< dim, spacedim > &  dh,
VectorType &  vector,
const ComponentMask mask = {} 
)

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-constructed mask is used, which is then interpreted as saying that 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.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ get_position_vector() [2/2]

template<int dim, int spacedim, typename VectorType >
void VectorTools::get_position_vector ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dh,
VectorType &  vector,
const ComponentMask mask = {} 
)

Like the above function but also taking mapping as argument. This will introduce an additional approximation between the true geometry specified by the manifold if the degree of the mapping is lower than the degree of the finite element in the DoFHandler dh, but more importantly it allows to fill location vectors for mappings that do not preserve vertex locations (like Eulerian mappings).

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ subtract_mean_value()

template<typename VectorType >
void VectorTools::subtract_mean_value ( VectorType &  v,
const std::vector< bool > &  p_select = {} 
)

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.
Warning
This function can only be used for distributed vector classes provided the boolean mask is empty, i.e. selecting the whole vector.
Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ add_constant()

template<typename VectorType , int dim, int spacedim = dim>
void VectorTools::add_constant ( VectorType &  solution,
const DoFHandler< dim, spacedim > &  dof_handler,
const unsigned int  component,
const typename VectorType::value_type  constant_adjustment 
)

Add the constant constant_adjustment to the specified component of the finite element function given by the coefficient vector solution defined by the given DoFHandler.

This operation is a common operation to compute a solution with mean pressure zero for a Stokes flow problem. Here, one can use VectorTools::compute_mean_value() to compute the value to subtract.

For a nodal finite element like FE_Q, this function will simply add the value constant_adjustment to each coefficient of the corresponding component. If you have the component in a separate block b, you could directly use solution.block(b) += constant_adjustment instead of calling this function (and this would be more efficient). For other finite element spaces like FE_DGP, the logic is more complicated and handled correctly by this function.

Note
Not all kinds of finite elements are supported and the selected component must not be part of a non-primitive element for this implementation to work.
In contrast to subtract_mean_value(), this function can adjust a single component of a distributed vector.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ compute_mean_value() [1/3]

template<int dim, typename Number , int spacedim>
Number VectorTools::compute_mean_value ( const hp::MappingCollection< dim, spacedim > &  mapping_collection,
const DoFHandler< dim, spacedim > &  dof,
const hp::QCollection< dim > &  q_collection,
const ReadVector< Number > &  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.

◆ compute_mean_value() [2/3]

template<int dim, typename Number , int spacedim>
Number VectorTools::compute_mean_value ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim > &  quadrature,
const ReadVector< Number > &  v,
const unsigned int  component 
)

Calls the other compute_mean_value() function, see above, for the non-hp case. That means, it requires a single FiniteElement, a single Quadrature, and a single Mapping object.

◆ compute_mean_value() [3/3]

template<int dim, typename Number , int spacedim>
Number VectorTools::compute_mean_value ( const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim > &  quadrature,
const ReadVector< Number > &  v,
const unsigned int  component 
)

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

◆ point_gradient() [1/6]

template<int dim, typename VectorType , int spacedim>
void VectorTools::point_gradient ( const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim, double > &  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.

This function is not particularly cheap. This is because it first needs to find which cell a given point is in, then find the point on the reference cell that matches the given evaluation point, and then evaluate the shape functions there. You probably do not want to use this function to evaluate the solution at many points. For this kind of application, the FEFieldFunction class offers at least some optimizations. On the other hand, if you want to evaluate many solutions at the same point, you may want to look at the VectorTools::create_point_source_vector() function.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the gradient of the finite element field either here or there, depending on which cell the point is found in. Since the gradient is, for most elements, discontinuous from one cell or the other, you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

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.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the gradient of the finite element field either here or there, depending on which cell the point is found in. Since the gradient is, for most elements, discontinuous from one cell or the other, you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ point_gradient() [2/6]

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, double > &  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.

This function is not particularly cheap. This is because it first needs to find which cell a given point is in, then find the point on the reference cell that matches the given evaluation point, and then evaluate the shape functions there. You probably do not want to use this function to evaluate the solution at many points. For this kind of application, the FEFieldFunction class offers at least some optimizations. On the other hand, if you want to evaluate many solutions at the same point, you may want to look at the VectorTools::create_point_source_vector() function.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the gradient of the finite element field either here or there, depending on which cell the point is found in. Since the gradient is, for most elements, discontinuous from one cell or the other, you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

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.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the gradient of the finite element field either here or there, depending on which cell the point is found in. Since the gradient is, for most elements, discontinuous from one cell or the other, you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ point_gradient() [3/6]

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, double > &  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.

This function is not particularly cheap. This is because it first needs to find which cell a given point is in, then find the point on the reference cell that matches the given evaluation point, and then evaluate the shape functions there. You probably do not want to use this function to evaluate the solution at many points. For this kind of application, the FEFieldFunction class offers at least some optimizations. On the other hand, if you want to evaluate many solutions at the same point, you may want to look at the VectorTools::create_point_source_vector() function.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the gradient of the finite element field either here or there, depending on which cell the point is found in. Since the gradient is, for most elements, discontinuous from one cell or the other, you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ point_gradient() [4/6]

template<int dim, typename VectorType , int spacedim>
void VectorTools::point_gradient ( const hp::MappingCollection< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim, double > &  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.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the gradient of the finite element field either here or there, depending on which cell the point is found in. Since the gradient is, for most elements, discontinuous from one cell or the other, you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ point_gradient() [5/6]

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, double > &  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.

This function is not particularly cheap. This is because it first needs to find which cell a given point is in, then find the point on the reference cell that matches the given evaluation point, and then evaluate the shape functions there. You probably do not want to use this function to evaluate the solution at many points. For this kind of application, the FEFieldFunction class offers at least some optimizations. On the other hand, if you want to evaluate many solutions at the same point, you may want to look at the VectorTools::create_point_source_vector() function.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the gradient of the finite element field either here or there, depending on which cell the point is found in. Since the gradient is, for most elements, discontinuous from one cell or the other, you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ point_gradient() [6/6]

template<int dim, typename VectorType , int spacedim>
Tensor< 1, spacedim, typename VectorType::value_type > VectorTools::point_gradient ( const hp::MappingCollection< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim, double > &  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.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the gradient of the finite element field either here or there, depending on which cell the point is found in. Since the gradient is, for most elements, discontinuous from one cell or the other, you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ create_point_source_vector() [1/6]

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof_handler,
const Point< spacedim, double > &  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) \varphi_i(x) dx\) where \(\varphi_i\) are the shape functions described by dof_handler and p is the point at which the delta function is located. Prior content of the given rhs_vector vector is deleted. This function is for the case of a scalar finite element.

This function is typically used in one of these two contexts:

  • Let's say you want to solve the same kind of problems many times over, with different values for right hand sides or coefficients, and then evaluate the solution at the same point every time. You could do this by calling VectorTools::point_value() after each solve, or you could realize that to evaluate the solution \(u_h\) at a point \(p\), you could rearrange operations like this:

    \begin{align*} u_h(p) &= \sum_j U_j \varphi_j(p) = \sum_j U_j F_j \\ &= U \cdot F \end{align*}

    with the vector as defined above. In other words, point evaluation can be achieved with just a single vector-vector product, and the vector \(F\) can be computed once and for all and reused for each solve, without having to go through the mesh every time to find out which cell (and where in the cell) the point \(p\) is located.
  • This function is also useful if you wanted to compute the Green's function for the problem you are solving. This is because the Green's function \(G(x,p)\) is defined by

    \begin{align*} L G(x,p) &= \delta(x-p) \end{align*}

    where \(L\) is the differential operator of your problem. The discrete version then requires computing the right hand side vector \(F_i = \int_\Omega \varphi_i(x) \delta(x-p)\), which is exactly the vector computed by the current function.

While maybe not relevant for documenting what this function does, it may be interesting to note that delta functions do not exist in reality, and consequently, using this function does not model any real situation. This is, because no real object is able to focus an infinite force density at an infinitesimally small part of the domain (rather, all real devices will spread out the force over a finite area); nor is it possible to measure values at individual points (but all measurements will somehow be averaged over small areas). Only if this area is so small that it cannot be resolved by any mesh does it make sense to model the situation in a way that uses a delta function with the same overall force or sensitivity. On the other hand, a situation that is probably more fruitfully simulated with a delta function is the electric potential of a point source; in this case, the solution is known to have a logarithmic singularity (in 2d) or a \(\frac{1}{r}\) singularity (in 3d), neither of which is bounded.

Mathematically, the use of delta functions typically leads to exact solutions to which the numerically obtained, approximate solution does not converge. This is because, taking the Laplace equation as an example, the error between exact and numerical solution can be bounded by the expression

\begin{align*} \| u-u_h \|_{L_2} \le C h \| \nabla u \|_{L_2} \end{align*}

but when using a delta function on the right hand side, the term \(\| \nabla u \|_{L_2} = |u|_{H^1}\) is not finite. This can be seen by using the a-priori bound for solutions of the Laplace equation \(-\Delta u = f\) that states that \(|u|_{H^1} \le \|f\|_{H^{-1}}\). When using a delta function as right hand side, \(f(x)=\delta(x-p)\), one would need to take the \(H^{-1}\) norm of a delta function, which however is not finite because \(\delta(\cdot-p) \not\in H^{-1}\).

The consequence of all of this is that the exact solution of the Laplace equation with a delta function on the right hand side – i.e., the Green's function – has a singularity at \(p\) that is so strong that it cannot be resolved by a finite element solution, and consequently finite element approximations do not converge towards the exact solution in any of the usual norms.

All of this is also the case for all of the other usual second-order partial differential equations in dimensions two or higher. (Because in dimension two and higher, \(H^1\) functions are not necessarily continuous, and consequently the delta function is not in the dual space \(H^{-1}\).)

◆ create_point_source_vector() [2/6]

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

Like the previous function, but for hp-objects.

◆ create_point_source_vector() [3/6]

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

Call the create_point_source_vector() function, see above, with an implied default \(Q_1\) mapping object.

Note that if your 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.

◆ create_point_source_vector() [4/6]

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof_handler,
const Point< spacedim, double > &  p,
const Point< dim, double > &  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 \varphi_i(x) dx\). Note here that \(\varphi_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 discussion of the first create_point_source_vector() variant for more on the use of delta functions.

◆ create_point_source_vector() [5/6]

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

Like the previous function, but for hp-objects.

◆ create_point_source_vector() [6/6]

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

Call the create_point_source_vector() function for vector-valued finite elements, see above, with an implied default \(Q_1\) mapping object.

Note that if your 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.

◆ point_difference() [1/2]

template<int dim, typename VectorType , int spacedim>
void VectorTools::point_difference ( const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Function< spacedim, typename VectorType::value_type > &  exact_solution,
Vector< typename VectorType::value_type > &  difference,
const Point< spacedim, double > &  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.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ point_difference() [2/2]

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, typename VectorType::value_type > &  exact_solution,
Vector< typename VectorType::value_type > &  difference,
const Point< spacedim, double > &  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.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ point_value() [1/6]

template<int dim, typename VectorType , int spacedim>
void VectorTools::point_value ( const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim, double > &  point,
Vector< typename VectorType::value_type > &  value 
)

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

This function uses a \(Q_1\)-mapping for the cell the point is evaluated in. If you need to evaluate using a different mapping (for example when using curved boundaries), use the point_difference() function that takes a mapping.

This function is not particularly cheap. This is because it first needs to find which cell a given point is in, then find the point on the reference cell that matches the given evaluation point, and then evaluate the shape functions there. You probably do not want to use this function to evaluate the solution at many points. For this kind of application, the FEFieldFunction class offers at least some optimizations. On the other hand, if you want to evaluate many solutions at the same point, you may want to look at the VectorTools::create_point_source_vector() function.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the value of the finite element field either here or there, depending on which cell the point is found in. This does not matter (to within the same tolerance) if the finite element field is continuous. On the other hand, if the finite element in use is not continuous, then you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

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.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the value of the finite element field either here or there, depending on which cell the point is found in. This does not matter (to within the same tolerance) if the finite element field is continuous. On the other hand, if the finite element in use is not continuous, then you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ point_value() [2/6]

template<int dim, typename VectorType , int spacedim>
VectorType::value_type VectorTools::point_value ( const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim, double > &  point 
)

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

This function uses a Q1-mapping for the cell the point is evaluated in. If you need to evaluate using a different mapping (for example when using curved boundaries), use the point_difference() function that takes a mapping.

This function is not particularly cheap. This is because it first needs to find which cell a given point is in, then find the point on the reference cell that matches the given evaluation point, and then evaluate the shape functions there. You probably do not want to use this function to evaluate the solution at many points. For this kind of application, the FEFieldFunction class offers at least some optimizations. On the other hand, if you want to evaluate many solutions at the same point, you may want to look at the VectorTools::create_point_source_vector() function.

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.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the value of the finite element field either here or there, depending on which cell the point is found in. This does not matter (to within the same tolerance) if the finite element field is continuous. On the other hand, if the finite element in use is not continuous, then you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

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.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the value of the finite element field either here or there, depending on which cell the point is found in. This does not matter (to within the same tolerance) if the finite element field is continuous. On the other hand, if the finite element in use is not continuous, then you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ point_value() [3/6]

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, double > &  point,
Vector< typename VectorType::value_type > &  value 
)

Evaluate a possibly vector-valued finite element function defined by the given DoFHandler and nodal vector fe_function at the given point 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 point value.

This function is not particularly cheap. This is because it first needs to find which cell a given point is in, then find the point on the reference cell that matches the given evaluation point, and then evaluate the shape functions there. You probably do not want to use this function to evaluate the solution at many points. For this kind of application, the FEFieldFunction class offers at least some optimizations. On the other hand, if you want to evaluate many solutions at the same point, you may want to look at the VectorTools::create_point_source_vector() function.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the value of the finite element field either here or there, depending on which cell the point is found in. This does not matter (to within the same tolerance) if the finite element field is continuous. On the other hand, if the finite element in use is not continuous, then you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ point_value() [4/6]

template<int dim, typename VectorType , int spacedim>
void VectorTools::point_value ( const hp::MappingCollection< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim, double > &  point,
Vector< 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.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the value of the finite element field either here or there, depending on which cell the point is found in. This does not matter (to within the same tolerance) if the finite element field is continuous. On the other hand, if the finite element in use is not continuous, then you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ point_value() [5/6]

template<int dim, typename VectorType , int spacedim>
VectorType::value_type VectorTools::point_value ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim, double > &  point 
)

Evaluate a scalar finite element function defined by the given DoFHandler and nodal vector fe_function at the given point 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.

This function is not particularly cheap. This is because it first needs to find which cell a given point is in, then find the point on the reference cell that matches the given evaluation point, and then evaluate the shape functions there. You probably do not want to use this function to evaluate the solution at many points. For this kind of application, the FEFieldFunction class offers at least some optimizations. On the other hand, if you want to evaluate many solutions at the same point, you may want to look at the VectorTools::create_point_source_vector() function.

Note
If the cell in which the point is found is not locally owned, an exception of type VectorTools::ExcPointNotAvailableHere is thrown.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the value of the finite element field either here or there, depending on which cell the point is found in. This does not matter (to within the same tolerance) if the finite element field is continuous. On the other hand, if the finite element in use is not continuous, then you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ point_value() [6/6]

template<int dim, typename VectorType , int spacedim>
VectorType::value_type VectorTools::point_value ( const hp::MappingCollection< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const VectorType &  fe_function,
const Point< spacedim, double > &  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.
This function needs to find the cell within which a point lies, and this can only be done up to a certain numerical tolerance of course. Consequently, for points that are on, or close to, the boundary of a cell, you may get the value of the finite element field either here or there, depending on which cell the point is found in. This does not matter (to within the same tolerance) if the finite element field is continuous. On the other hand, if the finite element in use is not continuous, then you will get unpredictable values for points on or close to the boundary of the cell, as one would expect when trying to evaluate point values of discontinuous functions.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ project() [1/7]

template<int dim, typename VectorType , int spacedim>
void VectorTools::project ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const AffineConstraints< typename VectorType::value_type > &  constraints,
const Quadrature< dim > &  quadrature,
const Function< spacedim, typename VectorType::value_type > &  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. In other words, given a function \(f(\mathbf x)\), the current function computes a finite element function \(f_h(\mathbf x)=\sum_j F_j \varphi_j(\mathbf x)\) characterized by the (output) vector of nodal values \(F\) that satisfies the equation

\begin{align*} (\varphi_i, f_h)_\Omega = (\varphi_i,f)_\Omega \end{align*}

for all test functions \(\varphi_i\). This requires solving a linear system involving the mass matrix since the equation above is equivalent to the linear system

\begin{align*} \sum_j (\varphi_i, \varphi_j)_\Omega F_j = (\varphi_i,f)_\Omega \end{align*}

which can also be written as \(MF = \Phi\) with \(M_{ij} = (\varphi_i, \varphi_j)_\Omega\) and \(\Phi_i = (\varphi_i,f)_\Omega\).

By default, no boundary values for \(f_h\) are needed nor imposed, but there are optional parameters to this function that allow imposing either zero boundary values or, in a first step, to project the boundary values of \(f\) onto the finite element space on the boundary of the mesh in a similar way to above, and then using these values as the imposed boundary values for \(f_h\). The ordering of arguments to this function is such that you need not give a second quadrature formula (of type Quadrature<dim-1> and used for the computation of the matrix and right hand side for the projection of boundary values) if you don't want to project to the boundary first, but that you must if you want to do so.

A MatrixFree implementation is used if the following conditions are met:

In this case, this function performs numerical quadrature using the given quadrature formula for integration of the right hand side \(\Phi_i\) and for the mass operator. In the case of hypercube cells, a QGauss(fe_degree+2) object is used for the mass operator. You should therefore make sure that the given quadrature formula is sufficiently accurate for creating the right-hand side.

Otherwise, only serial Triangulations are supported and the mass matrix is assembled using MatrixTools::create_mass_matrix. The given quadrature rule is then used for both the matrix and the right-hand side. You should therefore make sure that the given quadrature formula is also sufficient for creating the mass matrix. In particular, the degree of the quadrature formula must be sufficiently high to ensure that the mass matrix is invertible. For example, if you are using a FE_Q(k) element, then the integrand of the matrix entries \(M_{ij}\) is of polynomial degree \(2k\) in each variable, and you need a Gauss quadrature formula with \(k+1\) points in each coordinate direction to ensure that \(M\) is invertible.

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.

Parameters
[in]mappingThe mapping object to use.
[in]dofThe DoFHandler the describes the finite element space to project into and that corresponds to vec.
[in]constraintsConstraints to be used when assembling the mass matrix, typically needed when you have hanging nodes.
[in]quadratureThe quadrature formula to be used for assembling the mass matrix.
[in]functionThe function to project into the finite element space.
[out]vecThe output vector where the projected function will be stored in. This vector is required to be already initialized and must not have ghost elements.
[in]enforce_zero_boundaryIf true, vec will have zero boundary conditions.
[in]q_boundaryQuadrature rule to be used if project_to_boundary_first is true.
[in]project_to_boundary_firstIf true, perform a projection on the boundary before projecting the interior of the function.
Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ project() [2/7]

template<int dim, typename VectorType , int spacedim>
void VectorTools::project ( const DoFHandler< dim, spacedim > &  dof,
const AffineConstraints< typename VectorType::value_type > &  constraints,
const Quadrature< dim > &  quadrature,
const Function< spacedim, typename VectorType::value_type > &  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 
)

Call the project() function above, with mapping=MappingQ<dim>(1).

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ project() [3/7]

template<int dim, typename VectorType , int spacedim>
void VectorTools::project ( const hp::MappingCollection< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const AffineConstraints< typename VectorType::value_type > &  constraints,
const hp::QCollection< dim > &  quadrature,
const Function< spacedim, typename VectorType::value_type > &  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 with hp-capabilities.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ project() [4/7]

template<int dim, typename VectorType , int spacedim>
void VectorTools::project ( const DoFHandler< dim, spacedim > &  dof,
const AffineConstraints< typename VectorType::value_type > &  constraints,
const hp::QCollection< dim > &  quadrature,
const Function< spacedim, typename VectorType::value_type > &  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 
)

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

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ project() [5/7]

template<int dim, typename VectorType , int spacedim>
void VectorTools::project ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const AffineConstraints< typename VectorType::value_type > &  constraints,
const Quadrature< dim > &  quadrature,
const std::function< typename VectorType::value_type(const typename DoFHandler< dim, spacedim >::active_cell_iterator &, const unsigned int)> &  func,
VectorType &  vec_result 
)

The same as above for projection of scalar-valued quadrature data. The user provided function should return a value at the quadrature point based on the cell iterator and quadrature number and of course should be consistent with the provided quadrature object, which will be used to assemble the right-hand-side.

This function can be used with lambdas:

(mapping,
dof_handler,
constraints,
quadrature_formula,
[&] (const typename DoFHandler<dim>::active_cell_iterator & cell,
const unsigned int q) -> double
{
return qp_data.get_data(cell)[q]->density;
},
field);
typename ActiveSelector::active_cell_iterator active_cell_iterator
void project(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const AffineConstraints< typename VectorType::value_type > &constraints, const Quadrature< dim > &quadrature, const Function< spacedim, typename VectorType::value_type > &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)

where qp_data is a CellDataStorage object, which stores quadrature point data.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ project() [6/7]

template<int dim, typename VectorType >
void VectorTools::project ( std::shared_ptr< const MatrixFree< dim, typename VectorType::value_type, VectorizedArray< typename VectorType::value_type > > >  data,
const AffineConstraints< typename VectorType::value_type > &  constraints,
const unsigned int  n_q_points_1d,
const std::function< VectorizedArray< typename VectorType::value_type >(const unsigned int, const unsigned int)> &  func,
VectorType &  vec_result,
const unsigned int  fe_component = 0 
)

The same as above for projection of scalar-valued MatrixFree quadrature data. The user provided function func should return a VectorizedArray value at the quadrature point based on the cell number and quadrature number and should be consistent with the n_q_points_1d.

This function can be used with lambdas:

(matrix_free_data,
constraints,
3,
[&] (const unsigned int cell,
const unsigned int q) -> VectorizedArray<double>
{
return qp_data(cell,q);
},
field);

where qp_data is a an object of type Table<2, VectorizedArray<double> >, which stores quadrature point data.

fe_component allow to additionally specify which component of data to use in case it was constructed with an std::vector<const DoFHandler<dim>*>. It will be used internally in constructor of FEEvaluation object.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ project() [7/7]

template<int dim, typename VectorType >
void VectorTools::project ( std::shared_ptr< const MatrixFree< dim, typename VectorType::value_type, VectorizedArray< typename VectorType::value_type > > >  data,
const AffineConstraints< typename VectorType::value_type > &  constraints,
const std::function< VectorizedArray< typename VectorType::value_type >(const unsigned int, const unsigned int)> &  func,
VectorType &  vec_result,
const unsigned int  fe_component = 0 
)

Same as above but for n_q_points_1d = matrix_free.get_dof_handler().get_fe().degree+1.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ create_right_hand_side() [1/4]

template<int dim, int spacedim, typename VectorType >
void VectorTools::create_right_hand_side ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim > &  q,
const Function< spacedim, typename VectorType::value_type > &  rhs,
VectorType &  rhs_vector,
const AffineConstraints< typename VectorType::value_type > &  constraints = AffineConstraints< typename VectorType::value_type >() 
)

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.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ create_right_hand_side() [2/4]

template<int dim, int spacedim, typename VectorType >
void VectorTools::create_right_hand_side ( const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim > &  q,
const Function< spacedim, typename VectorType::value_type > &  rhs,
VectorType &  rhs_vector,
const AffineConstraints< typename VectorType::value_type > &  constraints = AffineConstraints< typename VectorType::value_type >() 
)

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

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ create_right_hand_side() [3/4]

template<int dim, int spacedim, typename VectorType >
void VectorTools::create_right_hand_side ( const hp::MappingCollection< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const hp::QCollection< dim > &  q,
const Function< spacedim, typename VectorType::value_type > &  rhs,
VectorType &  rhs_vector,
const AffineConstraints< typename VectorType::value_type > &  constraints = AffineConstraints< typename VectorType::value_type >() 
)

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

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ create_right_hand_side() [4/4]

template<int dim, int spacedim, typename VectorType >
void VectorTools::create_right_hand_side ( const DoFHandler< dim, spacedim > &  dof,
const hp::QCollection< dim > &  q,
const Function< spacedim, typename VectorType::value_type > &  rhs,
VectorType &  rhs_vector,
const AffineConstraints< typename VectorType::value_type > &  constraints = AffineConstraints< typename VectorType::value_type >() 
)

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

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ create_boundary_right_hand_side() [1/6]

template<int dim, int spacedim, typename VectorType >
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, typename VectorType::value_type > &  rhs,
VectorType &  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
Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ create_boundary_right_hand_side() [2/6]

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

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

See also
Glossary entry on boundary indicators
Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ create_boundary_right_hand_side() [3/6]

template<int dim, int spacedim, typename VectorType >
void VectorTools::create_boundary_right_hand_side ( const hp::MappingCollection< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const hp::QCollection< dim - 1 > &  q,
const Function< spacedim, typename VectorType::value_type > &  rhs,
VectorType &  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
Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ create_boundary_right_hand_side() [4/6]

template<int dim, int spacedim, typename VectorType >
void VectorTools::create_boundary_right_hand_side ( const DoFHandler< dim, spacedim > &  dof,
const hp::QCollection< dim - 1 > &  q,
const Function< spacedim, typename VectorType::value_type > &  rhs,
VectorType &  rhs_vector,
const std::set< types::boundary_id > &  boundary_ids = std::set< types::boundary_id >() 
)

Call 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
Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ create_boundary_right_hand_side() [5/6]

template<>
void VectorTools::create_boundary_right_hand_side ( const Mapping< 1, 1 > &  ,
const DoFHandler< 1, 1 > &  ,
const Quadrature< 0 > &  ,
const Function< 1 > &  ,
Vector< double > &  ,
const std::set< types::boundary_id > &   
)

Definition at line 27 of file vector_tools_rhs.cc.

◆ create_boundary_right_hand_side() [6/6]

template<>
void VectorTools::create_boundary_right_hand_side ( const Mapping< 1, 2 > &  ,
const DoFHandler< 1, 2 > &  ,
const Quadrature< 0 > &  ,
const Function< 2 > &  ,
Vector< double > &  ,
const std::set< types::boundary_id > &   
)

Definition at line 41 of file vector_tools_rhs.cc.

Variable Documentation

◆ dof_1

spacedim& VectorTools::dof_1

Definition at line 135 of file vector_tools_interpolate.h.

◆ dof_2

spacedim const DoFHandler<dim, spacedim>& VectorTools::dof_2

Definition at line 136 of file vector_tools_interpolate.h.

◆ transfer

spacedim const DoFHandler<dim, spacedim> const FullMatrix<double>& VectorTools::transfer

Definition at line 137 of file vector_tools_interpolate.h.

◆ data_1

spacedim const DoFHandler<dim, spacedim> const FullMatrix<double> const InVector& VectorTools::data_1

Definition at line 138 of file vector_tools_interpolate.h.

◆ data_2

spacedim const DoFHandler<dim, spacedim> const FullMatrix<double> const InVector OutVector& VectorTools::data_2

Definition at line 139 of file vector_tools_interpolate.h.