From: Peter Munch Date: Sun, 13 Dec 2020 14:45:09 +0000 (+0100) Subject: Add missing versions of VectorTools::interpolate_boundary_values for p::MappingCollection X-Git-Tag: v9.3.0-rc1~756^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F11367%2Fhead;p=dealii.git Add missing versions of VectorTools::interpolate_boundary_values for p::MappingCollection --- diff --git a/include/deal.II/numerics/vector_tools_boundary.h b/include/deal.II/numerics/vector_tools_boundary.h index 69cd87a0f7..6593be957b 100644 --- a/include/deal.II/numerics/vector_tools_boundary.h +++ b/include/deal.II/numerics/vector_tools_boundary.h @@ -139,6 +139,20 @@ namespace VectorTools std::map &boundary_values, const ComponentMask &component_mask = ComponentMask()); + /** + * Like the previous function, but take a mapping collection to go with + * DoFHandler objects with hp-capabilities. + */ + template + void + interpolate_boundary_values( + const hp::MappingCollection &mapping, + const DoFHandler & dof, + const types::boundary_id boundary_component, + const Function & boundary_function, + std::map & boundary_values, + const ComponentMask &component_mask = ComponentMask()); + /** * Call the other interpolate_boundary_values() function, see above, with * mapping=MappingQGeneric@(1). The same comments @@ -245,6 +259,20 @@ namespace VectorTools AffineConstraints &constraints, const ComponentMask & component_mask = ComponentMask()); + /** + * Like the previous function, but take a mapping collection to go with + * DoFHandler objects with hp-capabilities. + */ + template + void + interpolate_boundary_values( + const hp::MappingCollection &mapping, + const DoFHandler & dof, + const std::map *> + & function_map, + AffineConstraints &constraints, + 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 @@ -266,6 +294,20 @@ namespace VectorTools AffineConstraints & constraints, const ComponentMask & component_mask = ComponentMask()); + /** + * Like the previous function, but take a mapping collection to go with + * DoFHandler objects with hp-capabilities. + */ + template + void + interpolate_boundary_values( + const hp::MappingCollection &mapping, + const DoFHandler & dof, + const types::boundary_id boundary_component, + const Function & boundary_function, + AffineConstraints & constraints, + const ComponentMask &component_mask = ComponentMask()); + /** * Call the other interpolate_boundary_values() function, see above, with * mapping=MappingQGeneric@(1). The same comments diff --git a/include/deal.II/numerics/vector_tools_boundary.templates.h b/include/deal.II/numerics/vector_tools_boundary.templates.h index 7801f0d42e..cacf289443 100644 --- a/include/deal.II/numerics/vector_tools_boundary.templates.h +++ b/include/deal.II/numerics/vector_tools_boundary.templates.h @@ -440,6 +440,25 @@ namespace VectorTools + template + void + interpolate_boundary_values( + const hp::MappingCollection &mapping, + const DoFHandler & dof, + const types::boundary_id boundary_component, + const Function & boundary_function, + std::map & boundary_values, + const ComponentMask & component_mask) + { + std::map *> + function_map; + function_map[boundary_component] = &boundary_function; + interpolate_boundary_values( + mapping, dof, function_map, boundary_values, component_mask); + } + + + template void interpolate_boundary_values( @@ -530,6 +549,54 @@ namespace VectorTools + template + void + interpolate_boundary_values( + const hp::MappingCollection &mapping, + const DoFHandler & dof, + const std::map *> + & function_map, + AffineConstraints &constraints, + const ComponentMask & component_mask_) + { + std::map boundary_values; + interpolate_boundary_values( + mapping, dof, function_map, boundary_values, component_mask_); + typename std::map::const_iterator + boundary_value = boundary_values.begin(); + for (; boundary_value != boundary_values.end(); ++boundary_value) + { + if (constraints.can_store_line(boundary_value->first) && + !constraints.is_constrained(boundary_value->first)) + { + constraints.add_line(boundary_value->first); + constraints.set_inhomogeneity(boundary_value->first, + boundary_value->second); + } + } + } + + + + template + void + interpolate_boundary_values( + const hp::MappingCollection &mapping, + const DoFHandler & dof, + const types::boundary_id boundary_component, + const Function & boundary_function, + AffineConstraints & constraints, + const ComponentMask & component_mask) + { + std::map *> + function_map; + function_map[boundary_component] = &boundary_function; + interpolate_boundary_values( + mapping, dof, function_map, constraints, component_mask); + } + + + template void interpolate_boundary_values( diff --git a/source/numerics/vector_tools_boundary.inst.in b/source/numerics/vector_tools_boundary.inst.in index 468ceeb805..f2e24d69f3 100644 --- a/source/numerics/vector_tools_boundary.inst.in +++ b/source/numerics/vector_tools_boundary.inst.in @@ -40,6 +40,18 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS; template void interpolate_boundary_values( + const hp::MappingCollection + &, + const DoFHandler &, + const std::map *> &, + std::map &, + const ComponentMask &); + + template void + interpolate_boundary_values( + const hp::MappingCollection + &, const DoFHandler &, const types::boundary_id, const Function &, @@ -49,25 +61,13 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS; template void interpolate_boundary_values( const DoFHandler &, - const std::map *> &, + const types::boundary_id, + const Function &, std::map &, const ComponentMask &); - \} -#endif - } - -for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS; - number : REAL_AND_COMPLEX_SCALARS) - { -#if deal_II_dimension <= deal_II_space_dimension - namespace VectorTools - \{ template void interpolate_boundary_values( - const hp::MappingCollection - &, const DoFHandler &, const std::map *> &, @@ -104,6 +104,26 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS; AffineConstraints &, const ComponentMask &); + template void + interpolate_boundary_values( + const hp::MappingCollection + &, + const DoFHandler &, + const std::map *> &, + AffineConstraints &, + const ComponentMask &); + + template void + interpolate_boundary_values( + const hp::MappingCollection + &, + const DoFHandler &, + const types::boundary_id, + const Function &, + AffineConstraints &, + const ComponentMask &); + template void interpolate_boundary_values( const DoFHandler &,