]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Change get_active_fe_indices() to return the result. 14207/head
authorMarc Fehling <mafehling.git@gmail.com>
Sun, 21 Aug 2022 02:24:23 +0000 (20:24 -0600)
committerMarc Fehling <mafehling.git@gmail.com>
Sun, 21 Aug 2022 03:40:07 +0000 (21:40 -0600)
doc/news/changes/incompatibilities/20220820Fehling [new file with mode: 0644]
include/deal.II/dofs/dof_accessor.templates.h
include/deal.II/dofs/dof_handler.h
include/deal.II/dofs/dof_tools.h
source/dofs/dof_handler.cc
source/dofs/dof_tools.cc
tests/hp/get_active_fe_indices.cc

diff --git a/doc/news/changes/incompatibilities/20220820Fehling b/doc/news/changes/incompatibilities/20220820Fehling
new file mode 100644 (file)
index 0000000..654e48a
--- /dev/null
@@ -0,0 +1,5 @@
+Changed: The function DoFHandler::get_active_fe_indices() now returns
+the result vector. The previous version that writes into the argument
+as well as DoFTools::get_active_fe_indices() have been deprecated.
+<br>
+(Marc Fehling, 08/20/2022)
index 982da7e0a07d107064e03bec13869c160e9b8e0e..cbad821178958f2e112ec11a959d450d170165e2 100644 (file)
@@ -2333,15 +2333,10 @@ namespace internal
                  accessor.dof_handler->hp_cell_future_fe_indices.size(),
                ExcMessage("DoFHandler not initialized"));
 
-        // TODO
-        using active_fe_index_type = unsigned short int;
-        static const active_fe_index_type invalid_active_fe_index =
-          static_cast<active_fe_index_type>(-1);
-
         accessor.dof_handler
           ->hp_cell_future_fe_indices[accessor.level()]
                                      [accessor.present_index] =
-          invalid_active_fe_index;
+          DoFHandler<dim, spacedim>::invalid_active_fe_index;
       }
     };
   } // namespace DoFCellAccessorImplementation
index 2c9417b7be6b5ecb4b557cc3cc33dad4fd74aa18..957ce517ea2650f52bb33ce14617e42325d517f3 100644 (file)
@@ -573,16 +573,26 @@ public:
 
   /**
    * Go through the triangulation and set the active FE indices of all
-   * active cells to the values given in @p active_fe_indices.
+   * locally owned cells to the values given in @p active_fe_indices.
    */
   void
   set_active_fe_indices(const std::vector<unsigned int> &active_fe_indices);
 
+  /**
+   * Go through the triangulation and return a vector of active FE indices of
+   * all locally relevant cells.
+   */
+  std::vector<unsigned int>
+  get_active_fe_indices() const;
+
   /**
    * Go through the triangulation and store the active FE indices of all
-   * active cells to the vector @p active_fe_indices. This vector is
-   * resized, if necessary.
+   * locally relevant cells to the vector @p active_fe_indices. This vector
+   * is resized, if necessary.
+   *
+   * @deprecated Use the function that returns the result vector.
    */
+  DEAL_II_DEPRECATED_EARLY
   void
   get_active_fe_indices(std::vector<unsigned int> &active_fe_indices) const;
 
index 24320c89cca47ab35b204bdfb8e7766f51257070..7d9c9838034df9a4fbc473fa0c13cd14d8512ee7 100644 (file)
@@ -2133,9 +2133,11 @@ namespace DoFTools
    * For DoFHandler objects without hp-capabilities given as first argument, the
    * returned vector will consist of only zeros, indicating that all cells use
    * the same finite element. In hp-mode, the values may be different, though.
+   *
+   * @deprecated Use DoFHandler::get_active_fe_indices() instead.
    */
   template <int dim, int spacedim>
-  void
+  DEAL_II_DEPRECATED_EARLY void
   get_active_fe_indices(const DoFHandler<dim, spacedim> &dof_handler,
                         std::vector<unsigned int> &      active_fe_indices);
 
index bd2b7bc949c5c66242ee8d1f62e849cc3a170a91..4bff5729694bbeefef603dbb047148b5fa20c4a2 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 1998 - 2021 by the deal.II authors
+// Copyright (C) 1998 - 2022 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -2645,11 +2645,11 @@ DoFHandler<dim, spacedim>::set_active_fe_indices(
 
 
 template <int dim, int spacedim>
-void
-DoFHandler<dim, spacedim>::get_active_fe_indices(
-  std::vector<unsigned int> &active_fe_indices) const
+std::vector<unsigned int>
+DoFHandler<dim, spacedim>::get_active_fe_indices() const
 {
-  active_fe_indices.resize(this->get_triangulation().n_active_cells());
+  std::vector<unsigned int> active_fe_indices(
+    this->get_triangulation().n_active_cells(), invalid_active_fe_index);
 
   // we could try to extract the values directly, since they are
   // stored as protected data of this object, but for simplicity we
@@ -2657,6 +2657,18 @@ DoFHandler<dim, spacedim>::get_active_fe_indices(
   for (const auto &cell : this->active_cell_iterators())
     if (!cell->is_artificial())
       active_fe_indices[cell->active_cell_index()] = cell->active_fe_index();
+
+  return active_fe_indices;
+}
+
+
+
+template <int dim, int spacedim>
+void
+DoFHandler<dim, spacedim>::get_active_fe_indices(
+  std::vector<unsigned int> &active_fe_indices) const
+{
+  active_fe_indices = get_active_fe_indices();
 }
 
 
@@ -3026,7 +3038,7 @@ DoFHandler<dim, spacedim>::prepare_for_serialization_of_active_fe_indices()
   // active FE indices since ownership of cells may change.
 
   // Gather all current active FE indices
-  get_active_fe_indices(active_fe_index_transfer->active_fe_indices);
+  active_fe_index_transfer->active_fe_indices = get_active_fe_indices();
 
   // Attach to transfer object
   active_fe_index_transfer->cell_data_transfer->prepare_for_serialization(
index 5785431fd8654b215e9fc15ffc8344c722382d4e..2283ff3f2573772349efa368e7d57f21a568d6fc 100644 (file)
@@ -1373,8 +1373,7 @@ namespace DoFTools
     AssertDimension(active_fe_indices.size(),
                     dof_handler.get_triangulation().n_active_cells());
 
-    for (const auto &cell : dof_handler.active_cell_iterators())
-      active_fe_indices[cell->active_cell_index()] = cell->active_fe_index();
+    active_fe_indices = dof_handler.get_active_fe_indices();
   }
 
   template <int dim, int spacedim>
index e67830f8683e06d77a226172a2a23628bd775bbb..e4319dba3764c9f402f13ff19ef6037068872e6e 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2005 - 2020 by the deal.II authors
+// Copyright (C) 2005 - 2022 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -16,7 +16,7 @@
 
 
 // distribute different finite elements randomly across the domain, then use
-// DoFTools::get_active_fe_indices()
+// DoFHandler::get_active_fe_indices()
 
 
 #include <deal.II/dofs/dof_accessor.h>
@@ -54,16 +54,13 @@ test()
 
   DoFHandler<dim> dof_handler(tria);
 
-  for (typename DoFHandler<dim>::active_cell_iterator cell =
-         dof_handler.begin_active();
-       cell != dof_handler.end();
-       ++cell)
+  for (const auto &cell : dof_handler.active_cell_iterators())
     cell->set_active_fe_index(Testing::rand() % fe_collection.size());
 
   dof_handler.distribute_dofs(fe_collection);
 
-  std::vector<unsigned int> active_fe_indices(tria.n_active_cells());
-  DoFTools::get_active_fe_indices(dof_handler, active_fe_indices);
+  std::vector<unsigned int> active_fe_indices =
+    dof_handler.get_active_fe_indices();
   for (unsigned int i = 0; i < tria.n_active_cells(); ++i)
     deallog << active_fe_indices[i] << std::endl;
 }

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.