]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Move is_globally_ascending into an own function
authorDaniel Arndt <d.arndt@math.uni-goettingen.de>
Sun, 11 Sep 2016 10:56:29 +0000 (12:56 +0200)
committerDaniel Arndt <d.arndt@math.uni-goettingen.de>
Sun, 11 Sep 2016 21:34:49 +0000 (23:34 +0200)
include/deal.II/base/index_set.h
source/base/index_set.cc

index 8d2bc790bcaad30180a8527c3a16c65d7c5c1a8e..2ec7a302709786fe9383798272e3b0c16361a619 100644 (file)
@@ -218,13 +218,14 @@ public:
   bool is_contiguous () const;
 
   /**
-   * Return whether the index set stored by this object defines a linear
-   * range, i.e., each index is contained in exactly one IndexSet,
+   * Return whether the IndexSets are ascending with respect to MPI process
+   * number, i.e., each index is contained in exactly one IndexSet,
    * the first indices are contained in the IndexSet of the first MPI process,
    * the second indices are contained in the IndexSet of the second MPI process
    * and so on.
+   * In case there is only one MPI process, this is always true.
    */
-  bool is_linear() const;
+  bool is_globally_ascending() const;
 
   /**
    * Return the number of elements stored in this index set.
index b3ff1ab4b3c7129ffda761868d6c1f3b1978bc97..24d5c8ca9340886360391a08e8a5a58a2d4c6b0e 100644 (file)
@@ -512,56 +512,11 @@ IndexSet::make_trilinos_map (const MPI_Comm &communicator,
     }
 #endif
 
-  // Find out if the IndexSet is linear
-  // If there is only one process, the IndexSet is always linear
-  const unsigned int n_ranks = Utilities::MPI::n_mpi_processes(communicator);
-  bool is_linear = (n_ranks==1);
-  const bool all_contiguous = (Utilities::MPI::min (is_contiguous() ? 1 : 0, communicator) == 1);
-  // overlapping IndexSets or non_contiguous IndexSets can't be linear
-  if ((all_contiguous) && (!overlapping) && (!is_linear))
-    {
-      is_linear = true;
-      // we know that there is only one interval
-      types::global_dof_index local_dofs[2];
-      local_dofs[0] = (n_elements()>0) ? *(begin_intervals()->begin())
-                      : numbers::invalid_dof_index ;
-      local_dofs[1] = (n_elements()>0) ? begin_intervals()->last()
-                      : numbers::invalid_dof_index;
-
-      const unsigned int my_rank = Utilities::MPI::this_mpi_process(communicator);
-      // first gather all information on process 0
-      const unsigned int gather_size = (my_rank==0)?2*n_ranks:1;
-      std::vector<types::global_dof_index> global_dofs(gather_size);
-
-      MPI_Gather(local_dofs, 2, DEAL_II_DOF_INDEX_MPI_TYPE,
-                 &(global_dofs[0]), 2, DEAL_II_DOF_INDEX_MPI_TYPE, 0,
-                 communicator);
-      if (my_rank == 0)
-        {
-          // find out if the received std::vector is linear
-          types::global_dof_index old_dof = global_dofs[0], new_dof = 0;
-          types::global_dof_index index = 0;
-          while (global_dofs[index] == numbers::invalid_dof_index)
-            ++index;
-          old_dof = global_dofs[index++];
-          for (; index<global_dofs.size(); ++index)
-            {
-              new_dof = global_dofs[index];
-              if (new_dof == numbers::invalid_dof_index)
-                new_dof = old_dof;
-              else if (new_dof<=old_dof)
-                {
-                  is_linear = false;
-                  break;
-                }
-            }
-        }
-    }
-
-  // now broadcast the result
-  MPI_Bcast(&is_linear, 1, MPI::BOOL, 0, communicator);
+  // Find out if the IndexSet is ascending.
+  // Overlapping IndexSets are never ascending.
+  const bool ascending = overlapping ? false : is_globally_ascending();
 
-  if (is_linear)
+  if (ascending)
     return Epetra_Map (TrilinosWrappers::types::int_type(size()),
                        TrilinosWrappers::types::int_type(n_elements()),
                        0,
@@ -591,9 +546,66 @@ IndexSet::make_trilinos_map (const MPI_Comm &communicator,
                         );
     }
 }
+#endif
 
 
-#endif
+
+bool
+IndexSet::is_globally_ascending () const
+{
+  // If there is only one process, the IndexSet is always ascending.
+  const unsigned int n_ranks = Utilities::MPI::n_mpi_processes(communicator);
+  if (n_ranks==1)
+    return true;
+  // / Non-contiguous IndexSets can't be linear.
+  const bool all_contiguous = (Utilities::MPI::min (is_contiguous() ? 1 : 0, communicator) == 1);
+  if (!all_contiguous)
+    return false;
+
+  bool is_globally_ascending = true;
+  // we know that there is only one interval
+  types::global_dof_index local_dofs[2];
+  local_dofs[0] = (n_elements()>0) ? *(begin_intervals()->begin())
+                  : numbers::invalid_dof_index ;
+  local_dofs[1] = (n_elements()>0) ? begin_intervals()->last()
+                  : numbers::invalid_dof_index;
+
+  const unsigned int my_rank = Utilities::MPI::this_mpi_process(communicator);
+  // first gather all information on process 0
+  const unsigned int gather_size = (my_rank==0)?2*n_ranks:1;
+  std::vector<types::global_dof_index> global_dofs(gather_size);
+
+  MPI_Gather(local_dofs, 2, DEAL_II_DOF_INDEX_MPI_TYPE,
+             &(global_dofs[0]), 2, DEAL_II_DOF_INDEX_MPI_TYPE, 0,
+             communicator);
+  if (my_rank == 0)
+    {
+      // find out if the received std::vector is linear
+      types::global_dof_index old_dof = global_dofs[0], new_dof = 0;
+      types::global_dof_index index = 0;
+      while (global_dofs[index] == numbers::invalid_dof_index)
+        ++index;
+      old_dof = global_dofs[index++];
+      for (; index<global_dofs.size(); ++index)
+        {
+          new_dof = global_dofs[index];
+          if (new_dof == numbers::invalid_dof_index)
+            new_dof = old_dof;
+          else if (new_dof<old_dof)
+            {
+              is_globally_ascending = false;
+              break;
+            }
+        }
+    }
+
+  // now broadcast the result
+  int is_ascending = is_globally_ascending ? 1 : 0;
+  MPI_Bcast(&is_ascending, 1, MPI_INT, 0, communicator);
+
+  return (is_ascending==1);
+}
+
 
 
 std::size_t

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.