From: Wolfgang Bangerth Date: Wed, 31 Mar 2021 03:25:49 +0000 (-0600) Subject: Deprecate a DoFTools function that returns by argument, rather than directly. X-Git-Tag: v9.3.0-rc1~243^2~5 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=fd1b86e3381fe05dd840e6b0aee0a1cf214d086a;p=dealii.git Deprecate a DoFTools function that returns by argument, rather than directly. We want to move away from these kinds of functions, so deprecate it and replace it by a function that does the right thing instead. --- diff --git a/include/deal.II/dofs/dof_tools.h b/include/deal.II/dofs/dof_tools.h index 9b526fdb45..7b1a460e24 100644 --- a/include/deal.II/dofs/dof_tools.h +++ b/include/deal.II/dofs/dof_tools.h @@ -1486,21 +1486,33 @@ namespace DoFTools * @param[in] component_mask A mask denoting the vector components of the * finite element that should be considered (see also * @ref GlossComponentMask). - * @param[out] selected_dofs The IndexSet object that is returned and that - * will contain the indices of degrees of freedom that are located on the - * boundary (and correspond to the selected vector components and boundary - * indicators, depending on the values of the @p component_mask and @p - * boundary_ids arguments). * @param[in] boundary_ids If empty, this function extracts the indices of the * degrees of freedom for all parts of the boundary. If it is a non- empty * list, then the function only considers boundary faces with the boundary * indicators listed in this argument. + * @return The IndexSet object that + * will contain the indices of degrees of freedom that are located on the + * boundary (and correspond to the selected vector components and boundary + * indicators, depending on the values of the @p component_mask and @p + * boundary_ids arguments). * * @see * @ref GlossBoundaryIndicator "Glossary entry on boundary indicators" */ template - void + IndexSet + extract_boundary_dofs(const DoFHandler & dof_handler, + const ComponentMask & component_mask, + const std::set &boundary_ids = {}); + + /** + * The same as the previous function, except that it returns its information + * via the third argument. + * + * @deprecated Use the previous function instead. + */ + template + DEAL_II_DEPRECATED_EARLY void extract_boundary_dofs(const DoFHandler & dof_handler, const ComponentMask & component_mask, IndexSet & selected_dofs, diff --git a/source/dofs/dof_tools.cc b/source/dofs/dof_tools.cc index 55f29fad21..6f8a63e8e2 100644 --- a/source/dofs/dof_tools.cc +++ b/source/dofs/dof_tools.cc @@ -641,6 +641,19 @@ namespace DoFTools const ComponentMask & component_mask, IndexSet & selected_dofs, const std::set &boundary_ids) + { + // Simply forward to the other function + selected_dofs = + extract_boundary_dofs(dof_handler, component_mask, boundary_ids); + } + + + + template + IndexSet + extract_boundary_dofs(const DoFHandler & dof_handler, + const ComponentMask & component_mask, + const std::set &boundary_ids) { Assert(component_mask.represents_n_components( dof_handler.get_fe_collection().n_components()), @@ -649,9 +662,7 @@ namespace DoFTools boundary_ids.end(), ExcInvalidBoundaryIndicator()); - // first reset output argument - selected_dofs.clear(); - selected_dofs.set_size(dof_handler.n_dofs()); + IndexSet selected_dofs(dof_handler.n_dofs()); // let's see whether we have to check for certain boundary indicators // or whether we can accept all @@ -760,6 +771,8 @@ namespace DoFTools } } } + + return selected_dofs; } diff --git a/source/dofs/dof_tools.inst.in b/source/dofs/dof_tools.inst.in index eca9cc1fee..36162aa3de 100644 --- a/source/dofs/dof_tools.inst.in +++ b/source/dofs/dof_tools.inst.in @@ -188,6 +188,12 @@ for (deal_II_dimension : DIMENSIONS) IndexSet &, const std::set &); + template IndexSet + DoFTools::extract_boundary_dofs( + const DoFHandler &, + const ComponentMask &, + const std::set &); + template void DoFTools::extract_dofs_with_support_on_boundary< deal_II_dimension, deal_II_dimension>(const DoFHandler &,