--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 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.
+//
+// ---------------------------------------------------------------------
+
+
+// check the parallel partitioner with an additional index set to describe a
+// larger set for which the partitioner provides us with some renumbering
+#include "../tests.h"
+#include <deal.II/base/utilities.h>
+#include <deal.II/base/index_set.h>
+#include <deal.II/base/partitioner.h>
+#include <fstream>
+#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);
+
+ const unsigned int set = 200;
+ AssertIndexRange (numproc, set-2);
+ const unsigned int local_size = set - myid;
+ types::global_dof_index global_size = 0;
+ types::global_dof_index 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;
+ types::global_dof_index ghost_indices [10] = { 1, 2, 13, set-2, set-1, set, set+1, 2*set,
+ 2*set+1, 2*set+3
+ };
+ local_relevant.add_indices (&ghost_indices[0], &ghost_indices[0]+10);
+ types::global_dof_index before_start = myid > 0 ? my_start - set/4 : 0;
+ types::global_dof_index after_end = myid < numproc-1 ? my_start + local_size + set/3 : my_start+local_size;
+ if (before_start < my_start)
+ local_relevant.add_range(before_start, my_start);
+ if (after_end > my_start + local_size)
+ local_relevant.add_range(my_start + local_size, after_end);
+
+ Utilities::MPI::Partitioner v(local_owned, local_relevant, MPI_COMM_WORLD);
+
+ std::vector<types::global_dof_index> restricted_indices;
+ for (types::global_dof_index i=before_start; i<my_start; ++i)
+ if (i%2==0)
+ restricted_indices.push_back(i);
+ for (types::global_dof_index i=my_start+local_size; i<after_end; ++i)
+ if (((i/4)%3)==0)
+ restricted_indices.push_back(i);
+ restricted_indices.push_back(13);
+ restricted_indices.push_back(2*set);
+ restricted_indices.push_back(2*set+1);
+ IndexSet restricted_set(global_size);
+ restricted_set.add_indices (restricted_indices.begin(), restricted_indices.end());
+ Utilities::MPI::Partitioner w(local_owned, MPI_COMM_WORLD);
+ w.set_ghost_indices(restricted_set, v.ghost_indices());
+
+ IndexSet restricted_set2(global_size);
+ restricted_set2.add_index(2);
+ if (before_start < my_start)
+ restricted_set2.add_range(before_start, my_start);
+ Utilities::MPI::Partitioner x(local_owned, MPI_COMM_WORLD);
+ x.set_ghost_indices(restricted_set2, v.ghost_indices());
+
+ // print the additional ghost index number to the log stream
+ deallog << "Ghost subset in " << v.n_ghost_indices() << " indices: ";
+ for (unsigned int i=0; i<v.ghost_indices_within_larger_ghost_set().size(); ++i)
+ deallog << "[" << v.ghost_indices_within_larger_ghost_set()[i].first << ", "
+ << v.ghost_indices_within_larger_ghost_set()[i].second << ") ";
+ deallog << std::endl;
+
+ deallog << "Ghost subset in " << w.n_ghost_indices() << " indices: ";
+ for (unsigned int i=0; i<w.ghost_indices_within_larger_ghost_set().size(); ++i)
+ deallog << "[" << w.ghost_indices_within_larger_ghost_set()[i].first << ", "
+ << w.ghost_indices_within_larger_ghost_set()[i].second << ") ";
+ deallog << std::endl;
+
+ deallog << "Ghost subset in " << x.n_ghost_indices() << " indices: ";
+ for (unsigned int i=0; i<x.ghost_indices_within_larger_ghost_set().size(); ++i)
+ deallog << "[" << x.ghost_indices_within_larger_ghost_set()[i].first << ", "
+ << x.ghost_indices_within_larger_ghost_set()[i].second << ") ";
+ deallog << std::endl;
+}
+
+
+
+int main (int argc, char **argv)
+{
+ Utilities::MPI::MPI_InitFinalize mpi(argc, argv);
+ MPILogInitAll log;
+ test();
+}
--- /dev/null
+
+DEAL:0::Ghost subset in 69 indices: [200, 269)
+DEAL:0::Ghost subset in 24 indices: [4, 8) [16, 20) [28, 32) [40, 44) [52, 56) [64, 66) [66, 68)
+DEAL:0::Ghost subset in 0 indices:
+
+DEAL:1::Ghost subset in 119 indices: [199, 318)
+DEAL:1::Ghost subset in 49 indices: [2, 4) [5, 6) [7, 8) [9, 10) [11, 12) [13, 14) [15, 16) [17, 18) [19, 20) [21, 22) [23, 24) [25, 26) [27, 28) [29, 30) [31, 32) [33, 34) [35, 36) [37, 38) [39, 40) [41, 42) [43, 44) [45, 46) [47, 48) [49, 50) [51, 52) [53, 56) [62, 66) [74, 78) [86, 90) [98, 102) [110, 114)
+DEAL:1::Ghost subset in 51 indices: [1, 2) [3, 53)
+
+
+DEAL:2::Ghost subset in 123 indices: [198, 321)
+DEAL:2::Ghost subset in 49 indices: [2, 3) [8, 9) [10, 11) [12, 13) [14, 15) [16, 17) [18, 19) [20, 21) [22, 23) [24, 25) [26, 27) [28, 29) [30, 31) [32, 33) [34, 35) [36, 37) [38, 39) [40, 41) [42, 43) [44, 45) [46, 47) [48, 49) [50, 51) [52, 53) [54, 55) [56, 57) [60, 64) [72, 76) [84, 88) [96, 100) [108, 112) [120, 123)
+DEAL:2::Ghost subset in 51 indices: [1, 2) [7, 57)
+
+
+DEAL:3::Ghost subset in 60 indices: [197, 257)
+DEAL:3::Ghost subset in 28 indices: [2, 3) [7, 9) [11, 12) [13, 14) [15, 16) [17, 18) [19, 20) [21, 22) [23, 24) [25, 26) [27, 28) [29, 30) [31, 32) [33, 34) [35, 36) [37, 38) [39, 40) [41, 42) [43, 44) [45, 46) [47, 48) [49, 50) [51, 52) [53, 54) [55, 56) [57, 58) [59, 60)
+DEAL:3::Ghost subset in 51 indices: [1, 2) [10, 60)
+