Reference documentation for deal.II version 9.0.0
|
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 ::ExceptionBase & | ExcPointNotAvailableHere () |
Interpolation and projection | |
template<int dim, int spacedim, typename VectorType , template< int, int > class DoFHandlerType> | |
void | interpolate (const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof, const Function< spacedim, typename VectorType::value_type > &function, VectorType &vec, const ComponentMask &component_mask=ComponentMask()) |
template<int dim, int spacedim, typename VectorType , template< int, int > class DoFHandlerType> | |
void | interpolate (const DoFHandlerType< dim, spacedim > &dof, const Function< spacedim, typename VectorType::value_type > &function, VectorType &vec, const ComponentMask &component_mask=ComponentMask()) |
template<int dim, class InVector , class OutVector , int spacedim> | |
void | interpolate (const DoFHandler< dim, spacedim > &dof_1, const DoFHandler< dim, spacedim > &dof_2, const FullMatrix< double > &transfer, const InVector &data_1, OutVector &data_2) |
template<int dim, int spacedim, typename VectorType , template< int, int > class DoFHandlerType> | |
void | interpolate_based_on_material_id (const Mapping< dim, spacedim > &mapping, const DoFHandlerType< 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=ComponentMask()) |
template<int dim, int spacedim, typename VectorType , template< int, int > class DoFHandlerType> | |
void | interpolate_to_different_mesh (const DoFHandlerType< dim, spacedim > &dof1, const VectorType &u1, const DoFHandlerType< dim, spacedim > &dof2, VectorType &u2) |
template<int dim, int spacedim, typename VectorType , template< int, int > class DoFHandlerType> | |
void | interpolate_to_different_mesh (const DoFHandlerType< dim, spacedim > &dof1, const VectorType &u1, const DoFHandlerType< dim, spacedim > &dof2, const ConstraintMatrix &constraints, VectorType &u2) |
template<int dim, int spacedim, typename VectorType , template< int, int > class DoFHandlerType> | |
void | interpolate_to_different_mesh (const InterGridMap< DoFHandlerType< dim, spacedim > > &intergridmap, const VectorType &u1, const ConstraintMatrix &constraints, VectorType &u2) |
template<int dim, typename VectorType , int spacedim> | |
void | project (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const ConstraintMatrix &constraints, const Quadrature< dim > &quadrature, const Function< spacedim, 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 ConstraintMatrix &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 hp::DoFHandler< dim, spacedim > &dof, const ConstraintMatrix &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 hp::DoFHandler< dim, spacedim > &dof, const ConstraintMatrix &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 ConstraintMatrix &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 > > data, const ConstraintMatrix &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 > > data, const ConstraintMatrix &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) |
template<int dim, int spacedim, template< int, int > class DoFHandlerType, typename number > | |
void | interpolate_boundary_values (const Mapping< dim, spacedim > &mapping, const DoFHandlerType< 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=ComponentMask()) |
template<int dim, int spacedim, typename number > | |
void | interpolate_boundary_values (const hp::MappingCollection< dim, spacedim > &mapping, const hp::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=ComponentMask()) |
template<int dim, int spacedim, template< int, int > class DoFHandlerType, typename number > | |
void | interpolate_boundary_values (const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof, const types::boundary_id boundary_component, const Function< spacedim, number > &boundary_function, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask=ComponentMask()) |
template<int dim, int spacedim, template< int, int > class DoFHandlerType, typename number > | |
void | interpolate_boundary_values (const DoFHandlerType< dim, spacedim > &dof, const types::boundary_id boundary_component, const Function< spacedim, number > &boundary_function, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask=ComponentMask()) |
template<int dim, int spacedim, template< int, int > class DoFHandlerType, typename number > | |
void | interpolate_boundary_values (const DoFHandlerType< 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=ComponentMask()) |
template<int dim, int spacedim, template< int, int > class DoFHandlerType, typename number > | |
void | interpolate_boundary_values (const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > *> &function_map, ConstraintMatrix &constraints, const ComponentMask &component_mask=ComponentMask()) |
template<int dim, int spacedim, template< int, int > class DoFHandlerType, typename number > | |
void | interpolate_boundary_values (const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof, const types::boundary_id boundary_component, const Function< spacedim, number > &boundary_function, ConstraintMatrix &constraints, const ComponentMask &component_mask=ComponentMask()) |
template<int dim, int spacedim, template< int, int > class DoFHandlerType, typename number > | |
void | interpolate_boundary_values (const DoFHandlerType< dim, spacedim > &dof, const types::boundary_id boundary_component, const Function< spacedim, number > &boundary_function, ConstraintMatrix &constraints, const ComponentMask &component_mask=ComponentMask()) |
template<int dim, int spacedim, template< int, int > class DoFHandlerType, typename number > | |
void | interpolate_boundary_values (const DoFHandlerType< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > *> &function_map, ConstraintMatrix &constraints, const ComponentMask &component_mask=ComponentMask()) |
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=std::vector< unsigned int >()) |
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=std::vector< unsigned int >()) |
template<int dim, int spacedim, typename number > | |
void | project_boundary_values (const hp::MappingCollection< dim, spacedim > &mapping, const hp::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=std::vector< unsigned int >()) |
template<int dim, int spacedim, typename number > | |
void | project_boundary_values (const hp::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=std::vector< unsigned int >()) |
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, ConstraintMatrix &constraints, std::vector< unsigned int > component_mapping=std::vector< unsigned int >()) |
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, ConstraintMatrix &constraints, std::vector< unsigned int > component_mapping=std::vector< unsigned int >()) |
template<int dim> | |
void | project_boundary_values_curl_conforming (const DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, double > &boundary_function, const types::boundary_id boundary_component, ConstraintMatrix &constraints, const Mapping< dim > &mapping=StaticMappingQ1< dim >::mapping) |
template<int dim> | |
void | project_boundary_values_curl_conforming (const hp::DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, double > &boundary_function, const types::boundary_id boundary_component, ConstraintMatrix &constraints, const hp::MappingCollection< dim, dim > &mapping_collection=hp::StaticMappingQ1< dim >::mapping_collection) |
template<int dim> | |
void | project_boundary_values_curl_conforming_l2 (const DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, double > &boundary_function, const types::boundary_id boundary_component, ConstraintMatrix &constraints, const Mapping< dim > &mapping=StaticMappingQ1< dim >::mapping) |
template<int dim> | |
void | project_boundary_values_curl_conforming_l2 (const hp::DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, double > &boundary_function, const types::boundary_id boundary_component, ConstraintMatrix &constraints, const hp::MappingCollection< dim, dim > &mapping_collection=hp::StaticMappingQ1< dim >::mapping_collection) |
template<int dim> | |
void | project_boundary_values_div_conforming (const DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, double > &boundary_function, const types::boundary_id boundary_component, ConstraintMatrix &constraints, const Mapping< dim > &mapping=StaticMappingQ1< dim >::mapping) |
template<int dim> | |
void | project_boundary_values_div_conforming (const hp::DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, double > &boundary_function, const types::boundary_id boundary_component, ConstraintMatrix &constraints, const hp::MappingCollection< dim, dim > &mapping_collection=hp::StaticMappingQ1< dim >::mapping_collection) |
template<int dim, int spacedim, template< int, int > class DoFHandlerType> | |
void | compute_nonzero_normal_flux_constraints (const DoFHandlerType< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, typename FunctionMap< spacedim >::type &function_map, ConstraintMatrix &constraints, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim >::mapping) |
template<int dim, int spacedim, template< int, int > class DoFHandlerType> | |
void | compute_no_normal_flux_constraints (const DoFHandlerType< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, ConstraintMatrix &constraints, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim >::mapping) |
template<int dim, int spacedim, template< int, int > class DoFHandlerType> | |
void | compute_nonzero_tangential_flux_constraints (const DoFHandlerType< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, typename FunctionMap< spacedim >::type &function_map, ConstraintMatrix &constraints, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim >::mapping) |
template<int dim, int spacedim, template< int, int > class DoFHandlerType> | |
void | compute_normal_flux_constraints (const DoFHandlerType< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, ConstraintMatrix &constraints, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim >::mapping) |
Assembling of right hand sides | |
template<int dim, int spacedim, 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 ConstraintMatrix &constraints=ConstraintMatrix()) |
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 ConstraintMatrix &constraints=ConstraintMatrix()) |
template<int dim, int spacedim, typename VectorType > | |
void | create_right_hand_side (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim > &q, const Function< spacedim, typename VectorType::value_type > &rhs, VectorType &rhs_vector, const ConstraintMatrix &constraints=ConstraintMatrix()) |
template<int dim, int spacedim, typename VectorType > | |
void | create_right_hand_side (const hp::DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim > &q, const Function< spacedim, typename VectorType::value_type > &rhs, VectorType &rhs_vector, const ConstraintMatrix &constraints=ConstraintMatrix()) |
template<int dim, int spacedim> | |
void | create_point_source_vector (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim > &p, Vector< double > &rhs_vector) |
template<int dim, int spacedim> | |
void | create_point_source_vector (const DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim > &p, Vector< double > &rhs_vector) |
template<int dim, int spacedim> | |
void | create_point_source_vector (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim > &p, Vector< double > &rhs_vector) |
template<int dim, int spacedim> | |
void | create_point_source_vector (const hp::DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim > &p, Vector< double > &rhs_vector) |
template<int dim, int spacedim> | |
void | create_point_source_vector (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim > &p, const Point< dim > &direction, Vector< double > &rhs_vector) |
template<int dim, int spacedim> | |
void | create_point_source_vector (const DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim > &p, const Point< dim > &direction, Vector< double > &rhs_vector) |
template<int dim, int spacedim> | |
void | create_point_source_vector (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim > &p, const Point< dim > &direction, Vector< double > &rhs_vector) |
template<int dim, int spacedim> | |
void | create_point_source_vector (const hp::DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim > &p, const Point< dim > &direction, Vector< double > &rhs_vector) |
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 hp::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 hp::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 >()) |
Evaluation of functions and errors | |
template<int dim, class InVector , class OutVector , int spacedim> | |
void | integrate_difference (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim, double > &exact_solution, OutVector &difference, const Quadrature< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=nullptr, const double exponent=2.) |
template<int dim, class InVector , class OutVector , int spacedim> | |
void | integrate_difference (const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim, double > &exact_solution, OutVector &difference, const Quadrature< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=nullptr, const double exponent=2.) |
template<int dim, class InVector , class OutVector , int spacedim> | |
void | integrate_difference (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim, double > &exact_solution, OutVector &difference, const hp::QCollection< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=nullptr, const double exponent=2.) |
template<int dim, class InVector , class OutVector , int spacedim> | |
void | integrate_difference (const hp::DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim, double > &exact_solution, OutVector &difference, const hp::QCollection< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=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_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 > &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 > &point) |
template<int dim, typename VectorType , int spacedim> | |
void | point_value (const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, Vector< typename VectorType::value_type > &value) |
template<int dim, typename VectorType , int spacedim> | |
void | point_value (const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, Vector< 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 > &point) |
template<int dim, typename VectorType , int spacedim> | |
VectorType::value_type | point_value (const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point) |
template<int dim, typename VectorType , int spacedim> | |
void | point_value (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, Vector< typename VectorType::value_type > &value) |
template<int dim, typename VectorType , int spacedim> | |
void | point_value (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, Vector< 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 > &point) |
template<int dim, typename VectorType , int spacedim> | |
VectorType::value_type | point_value (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point) |
template<int dim, typename VectorType , int spacedim> | |
void | point_gradient (const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > &value) |
template<int dim, typename VectorType , int spacedim> | |
void | point_gradient (const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > &value) |
template<int dim, typename VectorType , int spacedim> | |
Tensor< 1, spacedim, typename VectorType::value_type > | point_gradient (const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point) |
template<int dim, typename VectorType , int spacedim> | |
Tensor< 1, spacedim, typename VectorType::value_type > | point_gradient (const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point) |
template<int dim, typename VectorType , int spacedim> | |
void | point_gradient (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > &value) |
template<int dim, typename VectorType , int spacedim> | |
void | point_gradient (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > &value) |
template<int dim, typename VectorType , int spacedim> | |
Tensor< 1, spacedim, typename VectorType::value_type > | point_gradient (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point) |
template<int dim, typename VectorType , int spacedim> | |
Tensor< 1, spacedim, typename VectorType::value_type > | point_gradient (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point) |
template<typename VectorType > | |
void | subtract_mean_value (VectorType &v, const std::vector< bool > &p_select=std::vector< bool >()) |
template<int dim, typename VectorType , int spacedim> | |
VectorType::value_type | compute_mean_value (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &quadrature, const VectorType &v, const unsigned int component) |
template<int dim, typename VectorType , int spacedim> | |
VectorType::value_type | compute_mean_value (const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &quadrature, const VectorType &v, const unsigned int component) |
template<int dim, int spacedim, template< int, int > class DoFHandlerType, typename VectorType > | |
void | get_position_vector (const DoFHandlerType< dim, spacedim > &dh, VectorType &vector, const ComponentMask &mask=ComponentMask()) |
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.
This collection of methods offers the following operations:
Interpolation: assign each degree of freedom in the vector to be the value of the function given as argument. This is identical to saying that the resulting finite element function (which is isomorphic to the output vector) has exact function values in all support points of trial functions. The support point of a trial function is the point where its value equals one, e.g. for linear trial functions the support points are four corners of an element. This function therefore relies on the assumption that a finite element is used for which the degrees of freedom are function values (Lagrange elements) rather than gradients, normal derivatives, second derivatives, etc (Hermite elements, quintic Argyris element, etc.).
It seems inevitable that some values of the vector to be created are set twice or even more than that. The reason is that we have to loop over all cells and get the function values for each of the trial functions located thereon. This applies also to the functions located on faces and corners which we thus visit more than once. While setting the value in the vector is not an expensive operation, the evaluation of the given function may be, taking into account that a virtual function has to be called.
Projection: compute the L2-projection of the given function onto the finite element space, i.e. if f is the function to be projected, compute fh in Vh such that (fh,vh)=(f,vh) for all discrete test functions vh. This is done through the solution of the linear system of equations M v = f where M is the mass matrix \(m_{ij} = \int_\Omega \phi_i(x) \phi_j(x) dx\) and \(f_i = \int_\Omega f(x) \phi_i(x) dx\). The solution vector \(v\) then is the nodal representation of the projection fh. The project() functions are used in the step-21 and step-23 tutorial programs.
In order to get proper results, it be may necessary to treat boundary conditions right. Below are listed some cases where this may be needed. If needed, this is done by L2-projection of the trace of the given function onto the finite element space restricted to the boundary of the domain, then taking this information and using it to eliminate the boundary nodes from the mass matrix of the whole domain, using the MatrixTools::apply_boundary_values() function. The projection of the trace of the function to the boundary is done with the VectorTools::project_boundary_values() (see below) function, which is called with a map of boundary functions std::map<types::boundary_id, const Function<spacedim,number>*> in which all boundary indicators from zero to numbers::internal_face_boundary_id-1 (numbers::internal_face_boundary_id is used for other purposes, see the Triangulation class documentation) point to the function to be projected. The projection to the boundary takes place using a second quadrature formula on the boundary given to the project() function. The first quadrature formula is used to compute the right hand side and for numerical quadrature of the mass matrix.
The projection of the boundary values first, then eliminating them from the global system of equations is not needed usually. It may be necessary if you want to enforce special restrictions on the boundary values of the projected function, for example in time dependent problems: you may want to project the initial values but need consistency with the boundary values for later times. Since the latter are projected onto the boundary in each time step, it is necessary that we also project the boundary values of the initial values, before projecting them to the whole domain.
Obviously, the results of the two schemes for projection are different. Usually, when projecting to the boundary first, the L2-norm of the difference between original function and projection over the whole domain will be larger (factors of five have been observed) while the L2-norm of the error integrated over the boundary should of course be less. The reverse should also hold if no projection to the boundary is performed.
The selection whether the projection to the boundary first is needed is done with the project_to_boundary_first
flag passed to the function. If false
is given, the additional quadrature formula for faces is ignored.
You should be aware of the fact that if no projection to the boundary is requested, a function with zero boundary values may not have zero boundary values after projection. There is a flag for this especially important case, which tells the function to enforce zero boundary values on the respective boundary parts. Since enforced zero boundary values could also have been reached through projection, but are more economically obtain using other methods, the project_to_boundary_first
flag is ignored if the enforce_zero_boundary
flag is set.
The solution of the linear system is presently done using a simple CG method without preconditioning and without multigrid. This is clearly not too efficient, but sufficient in many cases and simple to implement. This detail may change in the future.
Creation of right hand side vectors: The create_right_hand_side() function computes the vector \(f_i = \int_\Omega f(x) \phi_i(x) dx\). This is the same as what the MatrixCreator::create_*
functions which take a right hand side do, but without assembling a matrix.
Creation of right hand side vectors for point sources: The create_point_source_vector() function computes the vector \(f_i = \int_\Omega \delta(x-x_0) \phi_i(x) dx\).
Creation of boundary right hand side vectors: The create_boundary_right_hand_side() function computes the vector \(f_i = \int_{\partial\Omega} g(x) \phi_i(x) dx\). This is the right hand side contribution of boundary forces when having inhomogeneous Neumann boundary values in Laplace's equation or other second order operators. This function also takes an optional argument denoting over which parts of the boundary the integration shall extend. If the default argument is used, it is applied to all boundaries.
Interpolation of boundary values: The MatrixTools::apply_boundary_values() function takes a list of boundary nodes and their values. You can get such a list by interpolation of a boundary function using the interpolate_boundary_values() function. To use it, you have to specify a list of pairs of boundary indicators (of type types::boundary_id
; see the section in the documentation of the Triangulation class for more details) and the according functions denoting the Dirichlet boundary values of the nodes on boundary faces with this boundary indicator.
Usually, all other boundary conditions, such as inhomogeneous Neumann values or mixed boundary conditions are handled in the weak formulation. No attempt is made to include these into the process of matrix and vector assembly therefore.
Within this function, boundary values are interpolated, i.e. a node is given the point value of the boundary function. In some cases, it may be necessary to use the L2-projection of the boundary function or any other method. For this purpose we refer to the project_boundary_values() function below.
You should be aware that the boundary function may be evaluated at nodes on the interior of faces. These, however, need not be on the true boundary, but rather are on the approximation of the boundary represented by the mapping of the unit cell to the real cell. Since this mapping will in most cases not be the exact one at the face, the boundary function is evaluated at points which are not on the boundary and you should make sure that the returned values are reasonable in some sense anyway.
In 1d the situation is a bit different since there faces (i.e. vertices) have no boundary indicator. It is assumed that if the boundary indicator zero is given in the list of boundary functions, the left boundary point is to be interpolated while the right boundary point is associated with the boundary index 1 in the map. The respective boundary functions are then evaluated at the place of the respective boundary point.
Projection of boundary values: The project_boundary_values() function acts similar to the interpolate_boundary_values() function, apart from the fact that it does not get the nodal values of boundary nodes by interpolation but rather through the L2-projection of the trace of the function to the boundary.
The projection takes place on all boundary parts with boundary indicators listed in the map (std::map<types::boundary_id, const Function<spacedim,number>*>) of boundary functions. These boundary parts may or may not be continuous. For these boundary parts, the mass matrix is assembled using the MatrixTools::create_boundary_mass_matrix() function, as well as the appropriate right hand side. Then the resulting system of equations is solved using a simple CG method (without preconditioning), which is in most cases sufficient for the present purpose.
Computing errors: The function integrate_difference() performs the calculation of the error between a given (continuous) reference function and the finite element solution in different norms. The integration is performed using a given quadrature formula and assumes that the given finite element objects equals that used for the computation of the solution.
The result is stored in a vector (named difference
), where each entry equals the given norm of the difference on a cell. The order of entries is the same as a cell_iterator
takes when started with begin_active
and promoted with the ++
operator.
This data, one number per active cell, can be used to generate graphical output by directly passing it to the DataOut class through the DataOut::add_data_vector function. Alternatively, the global error can be computed using VectorTools::compute_global_error(). Finally, the output per cell from VectorTools::integrate_difference() can be interpolated to the nodal points of a finite element field using the DoFTools::distribute_cell_to_dof_vector function.
Presently, there is the possibility to compute the following values from the difference, on each cell: mean
, L1_norm
, L2_norm
, Linfty_norm
, H1_seminorm
and H1_norm
, see VectorTools::NormType. For the mean difference value, the reference function minus the numerical solution is computed, not the other way round.
The infinity norm of the difference on a given cell returns the maximum absolute value of the difference at the quadrature points given by the quadrature formula parameter. This will in some cases not be too good an approximation, since for example the Gauss quadrature formulae do not evaluate the difference at the end or corner points of the cells. You may want to choose a quadrature formula with more quadrature points or one with another distribution of the quadrature points in this case. You should also take into account the superconvergence properties of finite elements in some points: for example in 1D, the standard finite element method is a collocation method and should return the exact value at nodal points. Therefore, the trapezoidal rule should always return a vanishing L-infinity error. Conversely, in 2D the maximum L-infinity error should be located at the vertices or at the center of the cell, which would make it plausible to use the Simpson quadrature rule. On the other hand, there may be superconvergence at Gauss integration points. These examples are not intended as a rule of thumb, rather they are thought to illustrate that the use of the wrong quadrature formula may show a significantly wrong result and care should be taken to chose the right formula.
The H1 seminorm is the L2 norm of the gradient of the difference. The square of the full H1 norm is the sum of the square of seminorm and the square of the L2 norm.
To get the global L1 error, you have to sum up the entries in difference
, e.g. using Vector::l1_norm() function. For the global L2 difference, you have to sum up the squares of the entries and take the root of the sum, e.g. using Vector::l2_norm(). These two operations represent the l1 and l2 norms of the vectors, but you need not take the absolute value of each entry, since the cellwise norms are already positive.
To get the global mean difference, simply sum up the elements as above. To get the \(L_\infty\) norm, take the maximum of the vector elements, e.g. using the Vector::linfty_norm() function.
For the global H1 norm and seminorm, the same rule applies as for the L2 norm: compute the l2 norm of the cell error vector.
Note that, in the codimension one case, if you ask for a norm that requires the computation of a gradient, then the provided function is automatically projected along the curve, and the difference is only computed on the tangential part of the gradient, since no information is available on the normal component of the gradient anyway.
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.
Vector<float>, Vector<double>, BlockVector<float>, BlockVector<double>
; others can be generated in application code (see the section on Template instantiations in the manual).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 \[ 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 \[ 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 347 of file vector_tools.h.
void VectorTools::interpolate | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandlerType< dim, spacedim > & | dof, | ||
const Function< spacedim, typename VectorType::value_type > & | function, | ||
VectorType & | vec, | ||
const ComponentMask & | component_mask = ComponentMask() |
||
) |
Compute the interpolation of function
at the support points to the finite element space described by the Triangulation and FiniteElement object with which the given DoFHandler argument is initialized. It is assumed that the number of components of function
matches that of the finite element used by dof
.
Note that you may have to call hanging_nodes.distribute(vec)
with the hanging nodes from space dof
afterwards, to make the result continuous again.
The template argument DoFHandlerType
may either be of type DoFHandler or hp::DoFHandler.
See the general documentation of this namespace for further information.
mapping
argument should be replaced by a hp::MappingCollection in case of a hp::DoFHandler. void VectorTools::interpolate | ( | const DoFHandlerType< dim, spacedim > & | dof, |
const Function< spacedim, typename VectorType::value_type > & | function, | ||
VectorType & | vec, | ||
const ComponentMask & | component_mask = ComponentMask() |
||
) |
Call the interpolate()
function above with mapping=MappingQGeneric1<dim>@()
.
void VectorTools::interpolate | ( | const DoFHandler< dim, spacedim > & | dof_1, |
const DoFHandler< dim, spacedim > & | dof_2, | ||
const FullMatrix< double > & | transfer, | ||
const InVector & | data_1, | ||
OutVector & | data_2 | ||
) |
Interpolate different finite element spaces. The interpolation of vector data_1
is executed from the FE space represented by dof_1
to the vector data_2
on FE space dof_2
. The interpolation on each cell is represented by the matrix transfer
. Curved boundaries are neglected so far.
Note that you may have to call hanging_nodes.distribute(data_2)
with the hanging nodes from space dof_2
afterwards, to make the result continuous again.
void VectorTools::interpolate_based_on_material_id | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandlerType< 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 = ComponentMask() |
||
) |
This function is a kind of generalization or modification of the very first interpolate() function in the series. It interpolations a set of functions onto the finite element space given by the DoFHandler argument where the determination which function to use is made based on the material id (see GlossMaterialId) of each cell.
mapping | - The mapping to use to determine the location of support points at which the functions are to be evaluated. |
dof_handler | - DoFHandler initialized with Triangulation and FiniteElement objects, |
function_map | - std::map reflecting the correspondence between material ids and functions, |
dst | - global FE vector at the support points, |
component_mask | - mask of components that shall be interpolated |
function_map
, then dst
will not be updated in the respective degrees of freedom of the output vector For example, if dst
was successfully initialized to capture the degrees of freedom of the dof_handler
of the problem with all zeros in it, then those zeros which correspond to the missed material ids will still remain in dst
even after calling this function.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.void VectorTools::interpolate_to_different_mesh | ( | const DoFHandlerType< dim, spacedim > & | dof1, |
const VectorType & | u1, | ||
const DoFHandlerType< dim, spacedim > & | dof2, | ||
VectorType & | u2 | ||
) |
Compute the interpolation of a dof1-function
u1
to a dof2-function
u2
, where dof1
and dof2
represent different triangulations with a common coarse grid.
dof1 and dof2 need to have the same finite element discretization.
Note that for continuous elements on grids with hanging nodes (i.e. locally refined grids) this function does not give the expected output. Indeed, the resulting output vector does not necessarily respect continuity requirements at hanging nodes, due to local cellwise interpolation.
For this case (continuous elements on grids with hanging nodes), please use the interpolate_to_different_mesh function with an additional ConstraintMatrix argument, see below, or make the field conforming yourself by calling the ConstraintsMatrix::distribute
function of your hanging node constraints object.
void VectorTools::interpolate_to_different_mesh | ( | const DoFHandlerType< dim, spacedim > & | dof1, |
const VectorType & | u1, | ||
const DoFHandlerType< dim, spacedim > & | dof2, | ||
const ConstraintMatrix & | constraints, | ||
VectorType & | u2 | ||
) |
Compute the interpolation of a dof1-function
u1
to a dof2-function
u2
, where dof1
and dof2
represent different triangulations with a common coarse grid.
dof1 and dof2 need to have the same finite element discretization.
constraints
is a hanging node constraints object corresponding to dof2
. This object is particularly important when interpolating onto continuous elements on grids with hanging nodes (locally refined grids): Without it - due to cellwise interpolation - the resulting output vector does not necessarily respect continuity requirements at hanging nodes.
void VectorTools::interpolate_to_different_mesh | ( | const InterGridMap< DoFHandlerType< dim, spacedim > > & | intergridmap, |
const VectorType & | u1, | ||
const ConstraintMatrix & | constraints, | ||
VectorType & | u2 | ||
) |
The same function as above, but takes an InterGridMap object directly as a parameter. Useful for interpolating several vectors at the same time.
intergridmap
has to be initialized via InterGridMap::make_mapping pointing from a source DoFHandler to a destination DoFHandler.
void VectorTools::project | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandler< dim, spacedim > & | dof, | ||
const ConstraintMatrix & | 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.
By default, projection to the boundary and enforcement of zero boundary values are disabled. The ordering of arguments to this function is such that you need not give a second quadrature formula if you don't want to project to the boundary first, but that you must if you want to do so.
A MatrixFree implementation is used if the following conditions are met:
enforce_zero_boundary
is false,project_to_boundary_first
is false,In this case, this function performs numerical quadrature using the given quadrature formula for integration of the provided function while a QGauss(fe_degree+2) object is used for the mass operator. You should therefore make sure that the given quadrature formula is sufficient for creating the right-hand side.
Otherwise, only serial Triangulations are supported and the mass matrix is assembled exactly using MatrixTools::create_mass_matrix and the same quadrature rule as for the right-hand side. You should therefore make sure that the given quadrature formula is also sufficient for creating the mass matrix.
See the general documentation of this namespace for further information.
In 1d, the default value of the boundary quadrature formula is an invalid object since integration on the boundary doesn't happen in 1d.
[in] | mapping | The mapping object to use. |
[in] | dof | The DoFHandler the describes the finite element space to project into and that corresponds to vec . |
[in] | constraints | Constraints to be used when assembling the mass matrix, typically needed when you have hanging nodes. |
[in] | quadrature | The quadrature formula to be used for assembling the mass matrix. |
[in] | function | The function to project into the finite element space. |
[out] | vec | The 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_boundary | If true, vec will have zero boundary conditions. |
[in] | q_boundary | Quadrature rule to be used if project_to_boundary_first is true. |
[in] | project_to_boundary_first | If true, perform a projection on the boundary before projecting the interior of the function. |
void VectorTools::project | ( | const DoFHandler< dim, spacedim > & | dof, |
const ConstraintMatrix & | 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=MappingQGeneric<dim>(1)
.
void VectorTools::project | ( | const hp::MappingCollection< dim, spacedim > & | mapping, |
const hp::DoFHandler< dim, spacedim > & | dof, | ||
const ConstraintMatrix & | constraints, | ||
const hp::QCollection< dim > & | quadrature, | ||
const Function< spacedim, 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 for arguments of type hp::DoFHandler, hp::QCollection, and hp::MappingCollection.
void VectorTools::project | ( | const hp::DoFHandler< dim, spacedim > & | dof, |
const ConstraintMatrix & | 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.
void VectorTools::project | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandler< dim, spacedim > & | dof, | ||
const ConstraintMatrix & | 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:
where qp_data
is a CellDataStorage object, which stores quadrature point data.
void VectorTools::project | ( | std::shared_ptr< const MatrixFree< dim, typename VectorType::value_type > > | data, |
const ConstraintMatrix & | 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:
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.
void VectorTools::project | ( | std::shared_ptr< const MatrixFree< dim, typename VectorType::value_type > > | data, |
const ConstraintMatrix & | 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
.
void VectorTools::interpolate_boundary_values | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandlerType< 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 = ComponentMask() |
||
) |
Compute Dirichlet boundary conditions. This function makes up a map of degrees of freedom subject to Dirichlet boundary conditions and the corresponding values to be assigned to them, by interpolation around the boundary. For each degree of freedom at the boundary, if its index already exists in boundary_values
then its boundary value will be overwritten, otherwise a new entry with proper index and boundary value for this degree of freedom will be inserted into boundary_values
.
The parameter function_map
provides a list of boundary indicators to be handled by this function and corresponding boundary value functions. The keys of this map correspond to the number boundary_id
of the face. numbers::internal_face_boundary_id is an illegal value for this key since it is reserved for interior faces. For an example of how to use this argument with a non-empty map, see the step-16 tutorial program. (Note that the argument type is equal to the FunctionMap type.)
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)
.
function_map
must match that of the finite element used by dof
. In other words, for the example above, you need to provide a Function object that has 3 components (the two velocities and the pressure), even though you are only interested in the first two of them. interpolate_boundary_values() will then call this function to obtain a vector of 3 values at each interpolation point but only take the first two and discard the third. In other words, you are free to return whatever you like in the third component of the vector returned by Function::vector_value, but the Function object must state that it has 3 components.If the finite element used has shape functions that are non-zero in more than one component (in deal.II speak: they are non-primitive), then these components can presently not be used for interpolating boundary values. Thus, the elements in the component mask corresponding to the components of these non-primitive shape functions must be false
.
See the general documentation of this namespace for more information.
void VectorTools::interpolate_boundary_values | ( | const hp::MappingCollection< dim, spacedim > & | mapping, |
const hp::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 = ComponentMask() |
||
) |
Like the previous function, but take a mapping collection to go with the hp::DoFHandler object.
void VectorTools::interpolate_boundary_values | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandlerType< dim, spacedim > & | dof, | ||
const types::boundary_id | boundary_component, | ||
const Function< spacedim, number > & | boundary_function, | ||
std::map< types::global_dof_index, number > & | boundary_values, | ||
const ComponentMask & | component_mask = ComponentMask() |
||
) |
Same function as above, but taking only one pair of boundary indicator and corresponding boundary function. The same comments apply as for the previous function, in particular about the use of the component mask and the requires size of the function object.
void VectorTools::interpolate_boundary_values | ( | const DoFHandlerType< dim, spacedim > & | dof, |
const types::boundary_id | boundary_component, | ||
const Function< spacedim, number > & | boundary_function, | ||
std::map< types::global_dof_index, number > & | boundary_values, | ||
const ComponentMask & | component_mask = ComponentMask() |
||
) |
Call the other interpolate_boundary_values() function, see above, with mapping=MappingQGeneric<dim,spacedim>(1)
. The same comments apply as for the previous function, in particular about the use of the component mask and the requires size of the function object.
void VectorTools::interpolate_boundary_values | ( | const DoFHandlerType< 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 = ComponentMask() |
||
) |
Call the other interpolate_boundary_values() function, see above, with mapping=MappingQGeneric<dim,spacedim>(1)
. The same comments apply as for the previous function, in particular about the use of the component mask and the requires size of the function object.
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 = std::vector< unsigned int >() |
||
) |
Project a function or a set of functions to the boundary of the domain. In other words, compute the solution of the following problem: Find \(u_h \in V_h\) (where \(V_h\) is the finite element space represented by the DoFHandler argument of this function) so that
\begin{align*} \int_{\Gamma} \varphi_i u_h = \sum_{k \in {\cal K}} \int_{\Gamma_k} \varphi_i f_k, \qquad \forall \varphi_i \in V_h \end{align*}
where \(\Gamma = \bigcup_{k \in {\cal K}} \Gamma_k\), \(\Gamma_k \subset \partial\Omega\), \(\cal K\) is the set of indices and \(f_k\) the corresponding boundary functions represented in the function map argument boundary_values
to this function, and the integrals are evaluated by quadrature. This problem has a non-unique solution in the interior, but it is well defined for the degrees of freedom on the part of the boundary, \(\Gamma\), for which we do the integration. The values of \(u_h|_\Gamma\), i.e., the nodal values of the degrees of freedom of this function along the boundary, are then what is computed by this function.
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 exception if used with \(H_{curl}\) conforming elements, so the project_boundary_values_curl_conforming() should be used instead.
[in] | mapping | The mapping that will be used in the transformations necessary to integrate along the boundary. |
[in] | dof | The DoFHandler that describes the finite element space and the numbering of degrees of freedom. |
[in] | boundary_functions | A 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] | q | The face quadrature used in the integration necessary to compute the mass matrix and right hand side of the projection. |
[out] | boundary_values | The 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_mapping | It 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. |
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 = std::vector< unsigned int >() |
||
) |
Call the project_boundary_values() function, see above, with mapping=MappingQGeneric<dim,spacedim>(1)
.
void VectorTools::project_boundary_values | ( | const hp::MappingCollection< dim, spacedim > & | mapping, |
const hp::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 = std::vector< unsigned int >() |
||
) |
Same as above, but for objects of type hp::DoFHandler
void VectorTools::project_boundary_values | ( | const hp::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 = std::vector< unsigned int >() |
||
) |
Call the project_boundary_values() function, see above, with mapping=MappingQGeneric<dim,spacedim>(1)
.
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 ConstraintMatrix & | constraints = ConstraintMatrix() |
||
) |
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.
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 ConstraintMatrix & | constraints = ConstraintMatrix() |
||
) |
Call the create_right_hand_side() function, see above, with mapping=MappingQGeneric<dim>(1)
.
void VectorTools::create_right_hand_side | ( | const hp::MappingCollection< dim, spacedim > & | mapping, |
const hp::DoFHandler< dim, spacedim > & | dof, | ||
const hp::QCollection< dim > & | q, | ||
const Function< spacedim, typename VectorType::value_type > & | rhs, | ||
VectorType & | rhs_vector, | ||
const ConstraintMatrix & | constraints = ConstraintMatrix() |
||
) |
Like the previous set of functions, but for hp objects.
void VectorTools::create_right_hand_side | ( | const hp::DoFHandler< dim, spacedim > & | dof, |
const hp::QCollection< dim > & | q, | ||
const Function< spacedim, typename VectorType::value_type > & | rhs, | ||
VectorType & | rhs_vector, | ||
const ConstraintMatrix & | constraints = ConstraintMatrix() |
||
) |
Like the previous set of functions, but for hp objects.
void VectorTools::create_point_source_vector | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandler< dim, spacedim > & | dof_handler, | ||
const Point< spacedim > & | p, | ||
Vector< double > & | rhs_vector | ||
) |
Create a right hand side vector for a point source at point p
. In other words, it creates a vector \(F\) so that \(F_i = \int_\Omega \delta(x-p) \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.
It is worth noting 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. 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. 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}\).)
void VectorTools::create_point_source_vector | ( | const DoFHandler< dim, spacedim > & | dof_handler, |
const Point< spacedim > & | p, | ||
Vector< double > & | rhs_vector | ||
) |
Call the create_point_source_vector() function, see above, with an implied default \(Q_1\) mapping object.
void VectorTools::create_point_source_vector | ( | const hp::MappingCollection< dim, spacedim > & | mapping, |
const hp::DoFHandler< dim, spacedim > & | dof_handler, | ||
const Point< spacedim > & | p, | ||
Vector< double > & | rhs_vector | ||
) |
Like the previous set of functions, but for hp objects.
void VectorTools::create_point_source_vector | ( | const hp::DoFHandler< dim, spacedim > & | dof_handler, |
const Point< spacedim > & | p, | ||
Vector< double > & | rhs_vector | ||
) |
Like the previous set of functions, but for hp objects. The function uses an implied default \(Q_1\) mapping object. Note that if your hp::DoFHandler uses any active fe index other than zero, then you need to call the function above that provides a mapping object for each active fe index.
void VectorTools::create_point_source_vector | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandler< dim, spacedim > & | dof_handler, | ||
const Point< spacedim > & | p, | ||
const Point< dim > & | direction, | ||
Vector< double > & | rhs_vector | ||
) |
Create a right hand side vector for a point source at point p
. This variation of the function is meant for vector-valued problems with exactly dim components (it will also work for problems with more than dim components, and in this case simply consider only the first dim components of the shape functions). It computes a right hand side that corresponds to a forcing function that is equal to a delta function times a given direction. In other words, it creates a vector \(F\) so that \(F_i = \int_\Omega [\mathbf d \delta(x-p)] \cdot \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.
void VectorTools::create_point_source_vector | ( | const DoFHandler< dim, spacedim > & | dof_handler, |
const Point< spacedim > & | p, | ||
const Point< dim > & | 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.
void VectorTools::create_point_source_vector | ( | const hp::MappingCollection< dim, spacedim > & | mapping, |
const hp::DoFHandler< dim, spacedim > & | dof_handler, | ||
const Point< spacedim > & | p, | ||
const Point< dim > & | direction, | ||
Vector< double > & | rhs_vector | ||
) |
Like the previous set of functions, but for hp objects.
void VectorTools::create_point_source_vector | ( | const hp::DoFHandler< dim, spacedim > & | dof_handler, |
const Point< spacedim > & | p, | ||
const Point< dim > & | direction, | ||
Vector< double > & | rhs_vector | ||
) |
Like the previous set of functions, but for hp objects. The function uses an implied default \(Q_1\) mapping object. Note that if your hp::DoFHandler uses any active fe index other than zero, then you need to call the function above that provides a mapping object for each active fe index.
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.
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=MappingQGeneric<dim>(1)
.
void VectorTools::create_boundary_right_hand_side | ( | const hp::MappingCollection< dim, spacedim > & | mapping, |
const hp::DoFHandler< dim, spacedim > & | dof, | ||
const hp::QCollection< dim-1 > & | q, | ||
const Function< spacedim, 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.
void VectorTools::create_boundary_right_hand_side | ( | const hp::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.
void VectorTools::integrate_difference | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandler< dim, spacedim > & | dof, | ||
const InVector & | fe_function, | ||
const Function< spacedim, double > & | exact_solution, | ||
OutVector & | difference, | ||
const Quadrature< dim > & | q, | ||
const NormType & | norm, | ||
const Function< spacedim, double > * | weight = 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.
[in] | mapping | The mapping that is used when integrating the difference \(u-u_h\). |
[in] | dof | The DoFHandler object that describes the finite element space in which the solution vector lives. |
[in] | fe_function | A 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_solution | The exact solution that is used to compute the error. |
[out] | difference | The vector of values \(d_K\) computed as above. |
[in] | q | The quadrature formula used to approximate the integral shown above. Note that some quadrature formulas are more useful than other in integrating \(u-u_h\). For example, it is known that the \(Q_1\) approximation \(u_h\) to the exact solution \(u\) of a Laplace equation is particularly accurate (in fact, superconvergent, i.e. accurate to higher order) at the 4 Gauss points of a cell in 2d (or 8 points in 3d) that correspond to a QGauss(2) object. Consequently, because a QGauss(2) formula only evaluates the two solutions at these particular points, choosing this quadrature formula may indicate an error far smaller than it actually is. |
[in] | norm | The 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] | weight | The 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] | exponent | This 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.
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>.
void VectorTools::integrate_difference | ( | const DoFHandler< dim, spacedim > & | dof, |
const InVector & | fe_function, | ||
const Function< spacedim, double > & | exact_solution, | ||
OutVector & | difference, | ||
const Quadrature< dim > & | q, | ||
const NormType & | norm, | ||
const Function< spacedim, double > * | weight = nullptr , |
||
const double | exponent = 2. |
||
) |
Call the integrate_difference() function, see above, with mapping=MappingQGeneric<dim>(1)
.
void VectorTools::integrate_difference | ( | const hp::MappingCollection< dim, spacedim > & | mapping, |
const hp::DoFHandler< dim, spacedim > & | dof, | ||
const InVector & | fe_function, | ||
const Function< spacedim, double > & | exact_solution, | ||
OutVector & | difference, | ||
const hp::QCollection< dim > & | q, | ||
const NormType & | norm, | ||
const Function< spacedim, double > * | weight = nullptr , |
||
const double | exponent = 2. |
||
) |
Same as above for hp.
void VectorTools::integrate_difference | ( | const hp::DoFHandler< dim, spacedim > & | dof, |
const InVector & | fe_function, | ||
const Function< spacedim, double > & | exact_solution, | ||
OutVector & | difference, | ||
const hp::QCollection< dim > & | q, | ||
const NormType & | norm, | ||
const Function< spacedim, double > * | weight = nullptr , |
||
const double | exponent = 2. |
||
) |
Call the integrate_difference() function, see above, with mapping=MappingQGeneric<dim>(1)
.
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::Triangulation, 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.
tria | The Triangulation with active cells corresponding with the entries in cellwise_error . |
cellwise_error | Vector of errors on each active cell. |
norm | The type of norm to compute. |
exponent | The 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. |
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 > & | 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.
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 > & | 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.
void VectorTools::point_value | ( | const DoFHandler< dim, spacedim > & | dof, |
const VectorType & | fe_function, | ||
const Point< spacedim > & | point, | ||
Vector< 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) value of this function through the last argument.
This is a wrapper function using a Q1-mapping for cell boundaries to call the other point_difference() function.
void VectorTools::point_value | ( | const hp::DoFHandler< dim, spacedim > & | dof, |
const VectorType & | fe_function, | ||
const Point< spacedim > & | point, | ||
Vector< typename VectorType::value_type > & | value | ||
) |
Same as above for hp.
VectorType::value_type VectorTools::point_value | ( | const DoFHandler< dim, spacedim > & | dof, |
const VectorType & | fe_function, | ||
const Point< spacedim > & | point | ||
) |
Evaluate a scalar finite element function defined by the given DoFHandler and nodal vector at the given point, and return the value of this function.
Compared with the other function of the same name, this is a wrapper function using a Q1-mapping for cells.
This function is used in the "Possibilities for extensions" part of the results section of step-3.
VectorType::value_type VectorTools::point_value | ( | const hp::DoFHandler< dim, spacedim > & | dof, |
const VectorType & | fe_function, | ||
const Point< spacedim > & | point | ||
) |
Same as above for hp.
void VectorTools::point_value | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandler< dim, spacedim > & | dof, | ||
const VectorType & | fe_function, | ||
const Point< spacedim > & | point, | ||
Vector< 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) value of this function through the last argument.
Compared with the other function of the same name, this function uses an arbitrary mapping to evaluate the difference.
void VectorTools::point_value | ( | const hp::MappingCollection< dim, spacedim > & | mapping, |
const hp::DoFHandler< dim, spacedim > & | dof, | ||
const VectorType & | fe_function, | ||
const Point< spacedim > & | point, | ||
Vector< typename VectorType::value_type > & | value | ||
) |
Same as above for hp.
VectorType::value_type VectorTools::point_value | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandler< dim, spacedim > & | dof, | ||
const VectorType & | fe_function, | ||
const Point< spacedim > & | point | ||
) |
Evaluate a scalar finite element function defined by the given DoFHandler and nodal vector at the given point, and return the value of this function.
Compared with the other function of the same name, this function uses an arbitrary mapping to evaluate the difference.
VectorType::value_type VectorTools::point_value | ( | const hp::MappingCollection< dim, spacedim > & | mapping, |
const hp::DoFHandler< dim, spacedim > & | dof, | ||
const VectorType & | fe_function, | ||
const Point< spacedim > & | point | ||
) |
Same as above for hp.
void VectorTools::point_gradient | ( | const DoFHandler< dim, spacedim > & | dof, |
const VectorType & | fe_function, | ||
const Point< spacedim > & | point, | ||
std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > & | value | ||
) |
Evaluate a possibly vector-valued finite element function defined by the given DoFHandler and nodal vector at the given point, and return the (vector) gradient of this function through the last argument.
This is a wrapper function using a Q1-mapping for cell boundaries to call the other point_gradient() function.
void VectorTools::point_gradient | ( | const hp::DoFHandler< dim, spacedim > & | dof, |
const VectorType & | fe_function, | ||
const Point< spacedim > & | point, | ||
std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > & | value | ||
) |
Same as above for hp.
Tensor<1, spacedim, typename VectorType::value_type> VectorTools::point_gradient | ( | const DoFHandler< dim, spacedim > & | dof, |
const VectorType & | fe_function, | ||
const Point< spacedim > & | point | ||
) |
Evaluate a scalar finite element function defined by the given DoFHandler and nodal vector at the given point, and return the gradient of this function.
Compared with the other function of the same name, this is a wrapper function using a Q1-mapping for cells.
Tensor<1, spacedim, typename VectorType::value_type> VectorTools::point_gradient | ( | const hp::DoFHandler< dim, spacedim > & | dof, |
const VectorType & | fe_function, | ||
const Point< spacedim > & | point | ||
) |
Same as above for hp.
void VectorTools::point_gradient | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandler< dim, spacedim > & | dof, | ||
const VectorType & | fe_function, | ||
const Point< spacedim > & | point, | ||
std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > & | value | ||
) |
Evaluate a possibly vector-valued finite element function defined by the given DoFHandler and nodal vector at the given point, and return the gradients of this function through the last argument.
Compared with the other function of the same name, this function uses an arbitrary mapping for evaluation.
void VectorTools::point_gradient | ( | const hp::MappingCollection< dim, spacedim > & | mapping, |
const hp::DoFHandler< dim, spacedim > & | dof, | ||
const VectorType & | fe_function, | ||
const Point< spacedim > & | point, | ||
std::vector< Tensor< 1, spacedim, typename VectorType::value_type > > & | value | ||
) |
Same as above for hp.
Tensor<1, spacedim, typename VectorType::value_type> VectorTools::point_gradient | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandler< dim, spacedim > & | dof, | ||
const VectorType & | fe_function, | ||
const Point< spacedim > & | point | ||
) |
Evaluate a scalar finite element function defined by the given DoFHandler and nodal vector at the given point, and return the gradient of this function.
Compared with the other function of the same name, this function uses an arbitrary mapping for evaluation.
Tensor<1, spacedim, typename VectorType::value_type> VectorTools::point_gradient | ( | const hp::MappingCollection< dim, spacedim > & | mapping, |
const hp::DoFHandler< dim, spacedim > & | dof, | ||
const VectorType & | fe_function, | ||
const Point< spacedim > & | point | ||
) |
Same as above for hp.
void VectorTools::subtract_mean_value | ( | VectorType & | v, |
const std::vector< bool > & | p_select = std::vector< bool >() |
||
) |
Mean value operations Subtract the (algebraic) mean value from a vector.
This function is most frequently used as a mean-value filter for Stokes: The pressure in Stokes' equations with only Dirichlet boundaries for the velocities is only determined up to a constant. This function allows to subtract the mean value of the pressure. It is usually called in a preconditioner and generates updates with mean value zero. The mean value is computed as the mean value of the degrees of freedom values as given by the input vector; they are not weighted by the area of cells, i.e. the mean is computed as \(\sum_i v_i\), rather than as \(\int_\Omega v(x) = \int_\Omega \sum_i v_i \phi_i(x)\). The latter can be obtained from the VectorTools::compute_mean_function, however.
Apart from the vector v
to operate on, this function takes a boolean mask p_select
that has a true entry for every element of the vector for which the mean value shall be computed and later subtracted. The argument is used to denote which components of the solution vector correspond to the pressure, and avoid touching all other components of the vector, such as the velocity components. (Note, however, that the mask is not a GlossComponentMask operating on the vector components of the finite element the solution vector v
may be associated with; rather, it is a mask on the entire vector, without reference to what the vector elements mean.)
The boolean mask p_select
has an empty vector as default value, which will be interpreted as selecting all vector elements, hence, subtracting the algebraic mean value on the whole vector. This allows to call this function without a boolean mask if the whole vector should be processed.
VectorType::value_type VectorTools::compute_mean_value | ( | const Mapping< dim, spacedim > & | mapping, |
const DoFHandler< dim, spacedim > & | dof, | ||
const Quadrature< dim > & | quadrature, | ||
const VectorType & | v, | ||
const unsigned int | component | ||
) |
Compute the mean value of one component of the solution.
This function integrates the chosen component over the whole domain and returns the result, i.e. it computes \(\frac{1}{|\Omega|}\int_\Omega [u_h(x)]_c \; dx\) where \(c\) is the vector component and \(u_h\) is the function representation of the nodal vector given as fourth argument. The integral is evaluated numerically using the quadrature formula given as third argument.
This function is used in the "Possibilities for extensions" part of the results section of step-3.
VectorType::value_type VectorTools::compute_mean_value | ( | const DoFHandler< dim, spacedim > & | dof, |
const Quadrature< dim > & | quadrature, | ||
const VectorType & | v, | ||
const unsigned int | component | ||
) |
Call the other compute_mean_value() function, see above, with mapping=MappingQGeneric<dim>(1)
.
void VectorTools::get_position_vector | ( | const DoFHandlerType< dim, spacedim > & | dh, |
VectorType & | vector, | ||
const ComponentMask & | mask = ComponentMask() |
||
) |
Geometrical interpolation Given a DoFHandler containing at least a spacedim vector field, this function interpolates the Triangulation at the support points of a FE_Q() finite element of the same degree as the degree of the required components.
Curved manifold are respected, and the resulting VectorType will be geometrically consistent. The resulting map is guaranteed to be interpolatory at the support points of a FE_Q() finite element of the same degree as the degree of the required components.
If the underlying finite element is an FE_Q(1)^spacedim, then the resulting VectorType
is a finite element field representation of the vertices of the Triangulation.
The optional ComponentMask argument can be used to specify what components of the FiniteElement to use to describe the geometry. If no mask is specified at construction time, then a default one is used, i.e., the first spacedim components of the FiniteElement are assumed to represent the geometry of the problem.
This function is only implemented for FiniteElements where the specified components are primitive.