]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Use std::pair<cell_iterator, CellStatus> for DistributedTriangulationBase::local_cell... 11913/head
authorPasquale Africa <pasqualeclaudio.africa@polimi.it>
Mon, 15 Mar 2021 14:20:59 +0000 (14:20 +0000)
committerPasquale Africa <pasqualeclaudio.africa@polimi.it>
Tue, 16 Mar 2021 08:29:38 +0000 (08:29 +0000)
include/deal.II/distributed/fully_distributed_tria.h
include/deal.II/distributed/shared_tria.h
include/deal.II/distributed/tria.h
include/deal.II/distributed/tria_base.h
source/distributed/fully_distributed_tria.cc
source/distributed/tria.cc
source/distributed/tria_base.cc

index 9f2ad34d6cb1456709ec4473e0a7ced60194ada4..4996109c63d353590a7b5ea3831c04dad7990e36 100644 (file)
@@ -266,7 +266,7 @@ namespace parallel
 
       /**
        * Go through all locally owned cells and store the relations between
-       * cells and their status in the private member local_cell_relations.
+       * cells and their CellStatus in the private member local_cell_relations.
        *
        * The stored vector will be ordered by the occurrence of cells.
        */
index 0d4efafe5e2442e11bd0121c8c9c3e26aad01462..5c03d90be9ea6bd567a1e14942d4c10bbf39454b 100644 (file)
@@ -30,7 +30,6 @@
 #include <functional>
 #include <list>
 #include <set>
-#include <tuple>
 #include <utility>
 #include <vector>
 
index dda7535d94fb1f6ebf15e101175bf0bd427c838f..8ba5939840bd32980017fc6a4b80ac75c28c2097 100644 (file)
@@ -33,7 +33,6 @@
 #include <functional>
 #include <list>
 #include <set>
-#include <tuple>
 #include <type_traits>
 #include <utility>
 #include <vector>
@@ -695,8 +694,8 @@ namespace parallel
       typename dealii::internal::p4est::types<dim>::ghost *parallel_ghost;
 
       /**
-       * Go through all p4est trees and store the relations between the status
-       * of locally owned quadrants and cells in the private member
+       * Go through all p4est trees and store the relations between a deal.II
+       * cell and its current CellStatus in the private member
        * local_cell_relations.
        *
        * The stored vector will be ordered by the occurrence of quadrants in
index b9849799acd67d733392ae9592002abb72b2188b..235939c50d210bd26e70a296b2f5e5259a0d5b3d 100644 (file)
@@ -30,7 +30,6 @@
 #include <functional>
 #include <list>
 #include <set>
-#include <tuple>
 #include <utility>
 #include <vector>
 
@@ -689,9 +688,10 @@ namespace parallel
                        const unsigned int n_attached_deserialize_variable);
 
     /**
-     * Go through all cells and store the relations between locally
-     * owned quadrants and cells in the private member
-     * local_cell_relations.
+     * Go through all cells and store the relations between a deal.II cell and
+     * its current CellStatus in the private member local_cell_relations.
+     * For an extensive description of CellStatus, see the documentation
+     * for the member function register_data_attach().
      *
      * The stored vector will be ordered by the occurrence of quadrants.
      */
@@ -704,12 +704,12 @@ namespace parallel
      * description of the latter, see the documentation for the member
      * function register_data_attach().
      */
-    using cell_relation_t = typename std::tuple<CellStatus, cell_iterator>;
+    using cell_relation_t = typename std::pair<cell_iterator, CellStatus>;
 
     /**
-     * Vector of tuples, which each contain a deal.II cell
+     * Vector of pair, each containing a deal.II cell iterator
      * and their relation after refinement. To update its contents, use the
-     * compute_cell_relations member function.
+     * update_cell_relations() member function.
      *
      * The size of this vector is assumed to be equal to the number of locally
      * owned quadrants in the parallel_forest object.
index 5a964891fb437437a08b839fbb7da4a484ee028b..c49843ccd66b2b9de744d4bec238f71285140c84 100644 (file)
@@ -419,7 +419,7 @@ namespace parallel
             continue;
 
           this->local_cell_relations[cell_id] =
-            std::make_tuple(Triangulation<dim, spacedim>::CELL_PERSIST, cell);
+            std::make_pair(cell, Triangulation<dim, spacedim>::CELL_PERSIST);
           ++cell_id;
         }
     }
index ceb9d84b5fcaccfdb16229bff81c4cde92581a81..0f8d574f3a16b2e1f57e63fc3ac9e399fedce33d 100644 (file)
@@ -934,14 +934,14 @@ namespace
   }
 
   template <int dim, int spacedim>
-  using cell_relation_t = typename std::tuple<
-    typename dealii::Triangulation<dim, spacedim>::CellStatus,
-    typename dealii::Triangulation<dim, spacedim>::cell_iterator>;
+  using cell_relation_t = typename std::pair<
+    typename dealii::Triangulation<dim, spacedim>::cell_iterator,
+    typename dealii::Triangulation<dim, spacedim>::CellStatus>;
 
   /**
-   * Adds a tuple of a p4est quadrant, @p status and @p dealii_cell
+   * Adds a pair of a @p dealii_cell and its @p status
    * to the vector containing all relations @p cell_rel.
-   * The tuple will be inserted in the position corresponding to the one
+   * The pair will be inserted in the position corresponding to the one
    * of the p4est quadrant in the underlying p4est sc_array. The position
    * will be determined from @p idx, which is the position of the quadrant
    * in its corresponding @p tree. The p4est quadrant will be deduced from
@@ -962,7 +962,7 @@ namespace
     Assert(local_quadrant_index < cell_rel.size(), ExcInternalError());
 
     // store relation
-    cell_rel[local_quadrant_index] = std::make_tuple(status, dealii_cell);
+    cell_rel[local_quadrant_index] = std::make_pair(dealii_cell, status);
   }
 
 
@@ -1032,7 +1032,7 @@ namespace
     else if (!p4est_has_children && !dealii_cell->has_children())
       {
         // this active cell didn't change
-        // save tuple into corresponding position
+        // save pair into corresponding position
         add_single_cell_relation<dim, spacedim>(
           cell_rel,
           tree,
@@ -1455,7 +1455,7 @@ namespace parallel
         {
           (void)cell_rel;
           Assert(
-            (std::get<0>(cell_rel) == // cell_status
+            (cell_rel.second == // cell_status
              parallel::distributed::Triangulation<dim, spacedim>::CELL_PERSIST),
             ExcInternalError());
         }
@@ -3391,8 +3391,8 @@ namespace parallel
       // in the same order p4est will encounter them during repartitioning.
       for (const auto &cell_rel : this->local_cell_relations)
         {
-          const auto &cell_status = std::get<0>(cell_rel);
-          const auto &cell_it     = std::get<1>(cell_rel);
+          const auto &cell_it     = cell_rel.first;
+          const auto &cell_status = cell_rel.second;
 
           switch (cell_status)
             {
index 515081dac294d4b13a63f7ade0e76baabd7d4710..203a1365f8c813a9889dc4b351098854cf3e6e62 100644 (file)
@@ -745,7 +745,7 @@ namespace parallel
           {
             (void)cell_rel;
             Assert(
-              (std::get<0>(cell_rel) == // cell_status
+              (cell_rel.second == // cell_status
                parallel::DistributedTriangulationBase<dim,
                                                       spacedim>::CELL_PERSIST),
               ExcInternalError());
@@ -819,7 +819,7 @@ namespace parallel
       {
         // reset all cell_status entries after coarsening/refinement
         for (auto &cell_rel : local_cell_relations)
-          std::get<0>(cell_rel) =
+          cell_rel.second =
             parallel::DistributedTriangulationBase<dim, spacedim>::CELL_PERSIST;
       }
   }
@@ -893,8 +893,8 @@ namespace parallel
       auto data_cell_variable_it = packed_variable_size_data.begin();
       for (; cell_rel_it != cell_relations.cend(); ++cell_rel_it)
         {
-          const auto &cell_status = std::get<0>(*cell_rel_it);
-          const auto &dealii_cell = std::get<1>(*cell_rel_it);
+          const auto &dealii_cell = cell_rel_it->first;
+          const auto &cell_status = cell_rel_it->second;
 
           // Assertions about the tree structure.
           switch (cell_status)
@@ -1197,7 +1197,7 @@ namespace parallel
     for (; cell_rel_it != cell_relations.end();
          ++cell_rel_it, dest_fixed_it += sizes_fixed_cumulative.back())
       {
-        std::get<0>(*cell_rel_it) = // cell_status
+        cell_rel_it->second = // cell_status
           Utilities::unpack<typename parallel::DistributedTriangulationBase<
             dim,
             spacedim>::CellStatus>(dest_fixed_it,
@@ -1298,8 +1298,8 @@ namespace parallel
     auto dest_sizes_it = dest_sizes_variable.cbegin();
     for (; cell_rel_it != cell_relations.end(); ++cell_rel_it)
       {
-        const auto &cell_status = std::get<0>(*cell_rel_it);
-        const auto &dealii_cell = std::get<1>(*cell_rel_it);
+        const auto &dealii_cell = cell_rel_it->first;
+        const auto &cell_status = cell_rel_it->second;
 
         if (callback_variable_transfer)
           {

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.