]> https://gitweb.dealii.org/ - dealii.git/commitdiff
doc: give example of Utilities::MPI::Partitioner 6095/head
authorDenis Davydov <davydden@gmail.com>
Fri, 23 Mar 2018 20:26:04 +0000 (21:26 +0100)
committerDenis Davydov <davydden@gmail.com>
Sun, 25 Mar 2018 10:35:11 +0000 (12:35 +0200)
doc/doxygen/images/partitioner.png [new file with mode: 0644]
include/deal.II/base/partitioner.h
tests/mpi/parallel_partitioner_03b.cc [new file with mode: 0644]
tests/mpi/parallel_partitioner_03b.mpirun=4.output [new file with mode: 0644]

diff --git a/doc/doxygen/images/partitioner.png b/doc/doxygen/images/partitioner.png
new file mode 100644 (file)
index 0000000..a4c3534
Binary files /dev/null and b/doc/doxygen/images/partitioner.png differ
index 7dd443a9de617b60f13413f7ef8950434f40688d..1429531f26b76321a4393712a377f03c72690691 100644 (file)
@@ -53,6 +53,14 @@ namespace Utilities
      * information is gathered once when constructing the partitioner, which
      * obviates subsequent global communication steps when exchanging data.
      *
+     * The figure below gives an example of index space $[0,74)$ being split into
+     * four processes.
+     * @image html partitioner.png
+     * Here process 0 will import 5 DoFs from process 1 (first pair of import targets),
+     * which corresponds to the first 3 elements of its import indices. Whereas
+     * process 2 will import 3 DoFs from process 0, which corresponds to the first
+     * two elements of its import indices.
+     *
      * The partitioner includes a mechanism for converting global to local and
      * local to global indices. Internally, this class stores vector elements
      * using the convention as follows: The local range is associated with
diff --git a/tests/mpi/parallel_partitioner_03b.cc b/tests/mpi/parallel_partitioner_03b.cc
new file mode 100644 (file)
index 0000000..80f2187
--- /dev/null
@@ -0,0 +1,144 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2012 - 2017 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// a variation of parallel_partitioner_03, which is used to draw a figure in
+// MPI::Partitioner documentation.
+
+#include "../tests.h"
+#include <deal.II/base/utilities.h>
+#include <deal.II/base/index_set.h>
+#include <deal.II/base/partitioner.h>
+#include <iostream>
+#include <vector>
+
+
+void test ()
+{
+  unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+  unsigned int numproc = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD);
+
+  if (myid==0) deallog << "numproc=" << numproc << std::endl;
+
+  const unsigned int set = 20;
+  AssertIndexRange (numproc, set-2);
+  const unsigned int local_size = set - myid;
+  unsigned int global_size = 0;
+  unsigned int my_start = 0;
+  for (unsigned int i=0; i<numproc; ++i)
+    {
+      global_size += set - i;
+      if (i<myid)
+        my_start += set - i;
+    }
+  // each processor owns some indices and all
+  // are ghosting elements from three
+  // processors (the second). some entries
+  // are right around the border between two
+  // processors
+  IndexSet local_owned(global_size);
+  local_owned.add_range(my_start, my_start + local_size);
+  IndexSet local_relevant(global_size);
+  local_relevant = local_owned;
+
+  if (myid==0)
+    {
+      std::vector<unsigned int>ghost_indices = {1, 2, 13, set-2, set-1, set, set+1, 2*set,
+                                                2*set+1, 2*set+3
+                                               };
+      local_relevant.add_indices (ghost_indices.begin(), ghost_indices.end());
+    }
+  else if (myid==1)
+    {
+      std::vector<unsigned int>ghost_indices = {1, 2, 13, set-2, set-1, set, set+1, 2*set,
+                                                2*set+1, 2*set+3, 3*set
+                                               };
+      local_relevant.add_indices (ghost_indices.begin(), ghost_indices.end());
+    }
+  else if (myid==2)
+    {
+      std::vector<unsigned int> ghost_indices = {          set-2, set-1, set, set+1, 2*set,
+                                                           2*set+1, 2*set+3
+                                                };
+      local_relevant.add_indices (ghost_indices.begin(), ghost_indices.end());
+
+    }
+  else
+    {
+      std::vector<unsigned int>ghost_indices = {1, 2, 13, 2*set+3
+                                               };
+      local_relevant.add_indices (ghost_indices.begin(), ghost_indices.end());
+    }
+
+  Utilities::MPI::Partitioner v(local_owned, local_relevant, MPI_COMM_WORLD);
+
+  // write the info on ghost processors and import indices to file
+  {
+    std::ofstream file((std::string("dat.") + Utilities::int_to_string(myid)).c_str());
+    file << "**** proc " << myid << std::endl;
+    file  << "owned:" << std::endl;
+    local_owned.print(file);
+    file <<"relevant:" << std::endl;
+    local_relevant.print(file);
+    file << "ghost targets: ";
+    for (unsigned int i=0; i<v.ghost_targets().size(); ++i)
+      file << "[" << v.ghost_targets()[i].first << "/"
+           << v.ghost_targets()[i].second << "] ";
+    file << std::endl;
+    file << "import targets: ";
+    for (unsigned int i=0; i<v.import_targets().size(); ++i)
+      file << "[" << v.import_targets()[i].first << "/"
+           << v.import_targets()[i].second << "] ";
+    file << std::endl;
+    file << "import indices:" << std::endl;
+    for (unsigned int i=0; i<v.import_indices().size(); ++i)
+      file << "[" << v.import_indices()[i].first << "/"
+           << v.import_indices()[i].second << ")" << std::endl;
+    file << "****" << std::endl;
+  }
+
+  MPI_Barrier(MPI_COMM_WORLD);
+
+  if (myid==0)
+    {
+      for (unsigned int i=0; i<numproc; ++i)
+        {
+          cat_file((std::string("dat.") + Utilities::int_to_string(i)).c_str());
+        }
+
+    }
+
+}
+
+
+
+int main (int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, testing_max_num_threads());
+
+  unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+  deallog.push(Utilities::int_to_string(myid));
+
+  if (myid == 0)
+    {
+      initlog();
+      deallog << std::setprecision(4);
+
+      test();
+    }
+  else
+    test();
+
+}
diff --git a/tests/mpi/parallel_partitioner_03b.mpirun=4.output b/tests/mpi/parallel_partitioner_03b.mpirun=4.output
new file mode 100644 (file)
index 0000000..cbcb676
--- /dev/null
@@ -0,0 +1,56 @@
+
+DEAL:0::numproc=4
+**** proc 0
+owned:
+{[0,19]}
+relevant:
+{[0,21], [40,41], 43}
+ghost targets: [1/2] [2/3] 
+import targets: [1/5] [2/2] [3/3] 
+import indices:
+[1/3)
+[13/14)
+[18/20)
+[18/20)
+[1/3)
+[13/14)
+****
+
+**** proc 1
+owned:
+{[20,38]}
+relevant:
+{[1,2], 13, [18,38], [40,41], 43, 60}
+ghost targets: [0/5] [2/3] [3/1] 
+import targets: [0/2] [2/2] 
+import indices:
+[0/2)
+[0/2)
+****
+
+**** proc 2
+owned:
+{[39,56]}
+relevant:
+{[18,21], [39,56]}
+ghost targets: [0/2] [1/2] 
+import targets: [0/3] [1/3] [3/1] 
+import indices:
+[1/3)
+[4/5)
+[1/3)
+[4/5)
+[4/5)
+****
+
+**** proc 3
+owned:
+{[57,73]}
+relevant:
+{[1,2], 13, 43, [57,73]}
+ghost targets: [0/3] [2/1] 
+import targets: [1/1] 
+import indices:
+[3/4)
+****
+

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.