]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Deprecate a function in DoFTools.
authorWolfgang Bangerth <bangerth@colostate.edu>
Wed, 31 Mar 2021 03:18:42 +0000 (21:18 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Mon, 5 Apr 2021 16:20:22 +0000 (10:20 -0600)
While there also update the documentation of its replacement a bit.

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

index 5e4926f8221816aa9185aed18c7cda473f6f8a3e..9b526fdb45bbf7587a52f759d6b4ff09d8de093b 100644 (file)
@@ -1405,9 +1405,9 @@ namespace DoFTools
    * specified components of the solution. The function returns its results in
    * the last non-default-valued parameter which contains @p true if a degree
    * of freedom is at the boundary and belongs to one of the selected
-   * components, and @p false otherwise. The function is used in step-15.
+   * components, and @p false otherwise.
    *
-   * By specifying the @p boundary_id variable, you can select which boundary
+   * By specifying the @p boundary_ids variable, you can select which boundary
    * indicators the faces have to have on which the degrees of freedom are
    * located that shall be extracted. If it is an empty list, then all
    * boundary indicators are accepted.
@@ -1424,13 +1424,15 @@ namespace DoFTools
    * component mask is used that corresponds to the first non-zero components.
    * Elements in the mask corresponding to later components are ignored.
    *
-   * @note This function will not work for DoFHandler objects that are built
+   * @deprecated This function will not work for DoFHandler objects that are built
    * on a parallel::distributed::Triangulation object. The reasons is that the
    * output argument @p selected_dofs has to have a length equal to <i>all</i>
    * global degrees of freedom. Consequently, this does not scale to very
-   * large problems. If you need the functionality of this function for
+   * large problems, and this is also why the function is deprecated. If you
+   * need the functionality of this function for
    * parallel triangulations, then you need to use the other
-   * DoFTools::extract_boundary_dofs function.
+   * DoFTools::extract_boundary_dofs() function that returns its information
+   * via an IndexSet object.
    *
    * @param[in] dof_handler The object that describes which degrees of freedom
    * live on which cell.
@@ -1453,24 +1455,31 @@ namespace DoFTools
    * @ref GlossBoundaryIndicator "Glossary entry on boundary indicators"
    */
   template <int dim, int spacedim>
-  void
+  DEAL_II_DEPRECATED_EARLY void
   extract_boundary_dofs(const DoFHandler<dim, spacedim> &   dof_handler,
                         const ComponentMask &               component_mask,
                         std::vector<bool> &                 selected_dofs,
-                        const std::set<types::boundary_id> &boundary_ids =
-                          std::set<types::boundary_id>());
+                        const std::set<types::boundary_id> &boundary_ids = {});
 
   /**
-   * This function does the same as the previous one but it returns its result
-   * as an IndexSet rather than a std::vector@<bool@>. Thus, it can also be
-   * called for DoFHandler objects that are defined on
-   * parallel::distributed::Triangulation objects.
+   * Extract all degrees of freedom which are at the boundary and belong to
+   * specified components of the solution. The function returns its results in
+   * the form of an IndexSet that contains those entries that correspond to
+   * these selected degrees of freedom, i.e., which are at the boundary and
+   * belong to one of the selected components.
+   *
+   * By specifying the @p boundary_ids variable, you can select which boundary
+   * indicators the faces have to have on which the degrees of freedom are
+   * located that shall be extracted. If it is an empty list (the default), then
+   * all boundary indicators are accepted.
    *
-   * @note If the DoFHandler object is indeed defined on a
-   * parallel::distributed::Triangulation, then the @p selected_dofs index set
+   * @note If the DoFHandler object is defined on a
+   * parallel Triangulation object, then the computed index set
    * will contain only those degrees of freedom on the boundary that belong to
    * the locally relevant set (see
-   * @ref GlossLocallyRelevantDof "locally relevant DoFs").
+   * @ref GlossLocallyRelevantDof "locally relevant DoFs"), i.e., the function
+   * only considers faces of locally owned and ghost cells, but not of
+   * artificial cells.
    *
    * @param[in] dof_handler The object that describes which degrees of freedom
    * live on which cell.
@@ -1495,8 +1504,7 @@ namespace DoFTools
   extract_boundary_dofs(const DoFHandler<dim, spacedim> &   dof_handler,
                         const ComponentMask &               component_mask,
                         IndexSet &                          selected_dofs,
-                        const std::set<types::boundary_id> &boundary_ids =
-                          std::set<types::boundary_id>());
+                        const std::set<types::boundary_id> &boundary_ids = {});
 
   /**
    * This function is similar to the extract_boundary_dofs() function but it
index 30fc0849fe9acbe84558f9d5254c35893542c261..55f29fad2123e390e5c6ea24a12e4415f7c704c8 100644 (file)
@@ -634,6 +634,7 @@ namespace DoFTools
   }
 
 
+
   template <int dim, int spacedim>
   void
   extract_boundary_dofs(const DoFHandler<dim, spacedim> &   dof_handler,
@@ -675,7 +676,6 @@ namespace DoFTools
     // boundary line is also part of a boundary face which we will be
     // visiting sooner or later
     for (const auto &cell : dof_handler.active_cell_iterators())
-
       // only work on cells that are either locally owned or at least ghost
       // cells
       if (cell->is_artificial() == false)

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.