]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Deprecate a DoFTools function that returns by argument, rather than directly.
authorWolfgang Bangerth <bangerth@colostate.edu>
Wed, 31 Mar 2021 03:25:49 +0000 (21:25 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Mon, 5 Apr 2021 16:20:22 +0000 (10:20 -0600)
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.

include/deal.II/dofs/dof_tools.h
source/dofs/dof_tools.cc
source/dofs/dof_tools.inst.in

index 9b526fdb45bbf7587a52f759d6b4ff09d8de093b..7b1a460e243524932ad056f57d369ea540a04634 100644 (file)
@@ -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 <int dim, int spacedim>
-  void
+  IndexSet
+  extract_boundary_dofs(const DoFHandler<dim, spacedim> &   dof_handler,
+                        const ComponentMask &               component_mask,
+                        const std::set<types::boundary_id> &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 <int dim, int spacedim>
+  DEAL_II_DEPRECATED_EARLY void
   extract_boundary_dofs(const DoFHandler<dim, spacedim> &   dof_handler,
                         const ComponentMask &               component_mask,
                         IndexSet &                          selected_dofs,
index 55f29fad2123e390e5c6ea24a12e4415f7c704c8..6f8a63e8e2404e90bad83a44766ff52e449553ae 100644 (file)
@@ -641,6 +641,19 @@ namespace DoFTools
                         const ComponentMask &               component_mask,
                         IndexSet &                          selected_dofs,
                         const std::set<types::boundary_id> &boundary_ids)
+  {
+    // Simply forward to the other function
+    selected_dofs =
+      extract_boundary_dofs(dof_handler, component_mask, boundary_ids);
+  }
+
+
+
+  template <int dim, int spacedim>
+  IndexSet
+  extract_boundary_dofs(const DoFHandler<dim, spacedim> &   dof_handler,
+                        const ComponentMask &               component_mask,
+                        const std::set<types::boundary_id> &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;
   }
 
 
index eca9cc1feeeb2e33c8515992214f3e4c346c69ec..36162aa3debda467725ce8ae0052d816d6cd5831 100644 (file)
@@ -188,6 +188,12 @@ for (deal_II_dimension : DIMENSIONS)
       IndexSet &,
       const std::set<types::boundary_id> &);
 
+    template IndexSet
+    DoFTools::extract_boundary_dofs<deal_II_dimension, deal_II_dimension>(
+      const DoFHandler<deal_II_dimension> &,
+      const ComponentMask &,
+      const std::set<types::boundary_id> &);
+
     template void DoFTools::extract_dofs_with_support_on_boundary<
       deal_II_dimension,
       deal_II_dimension>(const DoFHandler<deal_II_dimension> &,

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.