From: Doug Shi-Dong Date: Thu, 30 Apr 2020 00:09:18 +0000 (-0400) Subject: Add U::MPI::create_evenly_distributed_partitioning X-Git-Tag: v9.2.0-rc1~133^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=7d3205f73a34fbed2ee927c5fd007603b9320d8d;p=dealii.git Add U::MPI::create_evenly_distributed_partitioning --- diff --git a/doc/news/changes/minor/20200429DougShiDong b/doc/news/changes/minor/20200429DougShiDong new file mode 100644 index 0000000000..ad7514835b --- /dev/null +++ b/doc/news/changes/minor/20200429DougShiDong @@ -0,0 +1,6 @@ +New: Add Utilities::create_evenly_distributed_partitioning(my_partition_id, n_partitions, total_size) +to create a one-to-one evenly distributed ascending partitioning from total size. +Add Utilities::MPI::create_evenly_distributed_partitioning(comm, total_size) to use processor ID and +number of MPI processes to determine partitioning. +
+(Doug Shi-Dong, 2020/04/29) diff --git a/include/deal.II/base/mpi.h b/include/deal.II/base/mpi.h index 68bdb8e12b..5880b1889b 100644 --- a/include/deal.II/base/mpi.h +++ b/include/deal.II/base/mpi.h @@ -104,6 +104,23 @@ class IndexSet; namespace Utilities { + /** + * Given the total number of elements @p total_size, create an evenly + * distributed 1:1 partitioning of the elements for across @p n_partitions. + * The local sizes will be equal to the @p total_size divided by the number + * of partitions plus the remainder being divided amongst the first + * processes. Each process will store a contiguous subset of indices, and the + * index set on process p+1 starts at the index one larger than the last one + * stored on process p. + * For example, a @p total_size of 11 with 3 processes will result + * in the IndexSets { [0,4), [4,8), [8,11)] }, and this function will + * return the @p my_partition_id 's IndexSet. + */ + IndexSet + create_evenly_distributed_partitioning(const unsigned int my_partition_id, + const unsigned int n_partitions, + const IndexSet::size_type total_size); + /** * A namespace for utility functions that abstract certain operations using * the Message Passing Interface (MPI) or provide fallback operations in @@ -431,6 +448,18 @@ namespace Utilities create_ascending_partitioning(const MPI_Comm & comm, const IndexSet::size_type local_size); + /** + * Given the total number of elements @p total_size, create an evenly + * distributed 1:1 partitioning of the elements across the + * MPI communicator @p comm. + * Uses @p comm to determine number of partitions and processor ID to call the + * @p create_evenly_distributed_partitioning() function above. + */ + IndexSet + create_evenly_distributed_partitioning( + const MPI_Comm & comm, + const IndexSet::size_type total_size); + #ifdef DEAL_II_WITH_MPI /** * Calculate mean and standard deviation across the MPI communicator @p comm diff --git a/source/base/mpi.cc b/source/base/mpi.cc index 9295a2c0c8..3e91985641 100644 --- a/source/base/mpi.cc +++ b/source/base/mpi.cc @@ -68,6 +68,24 @@ DEAL_II_NAMESPACE_OPEN namespace Utilities { + IndexSet + create_evenly_distributed_partitioning(const unsigned int my_partition_id, + const unsigned int n_partitions, + const IndexSet::size_type total_size) + { + const unsigned int remain = total_size % n_partitions; + + const IndexSet::size_type min_size = total_size / n_partitions; + + const IndexSet::size_type begin = + min_size * my_partition_id + std::min(my_partition_id, remain); + const IndexSet::size_type end = + min_size * (my_partition_id + 1) + std::min(my_partition_id + 1, remain); + IndexSet result(total_size); + result.add_range(begin, end); + return result; + } + namespace MPI { MinMaxAvg @@ -243,6 +261,18 @@ namespace Utilities return res; } + IndexSet + create_evenly_distributed_partitioning(const MPI_Comm & comm, + const IndexSet::size_type total_size) + { + const unsigned int this_proc = this_mpi_process(comm); + const unsigned int n_proc = n_mpi_processes(comm); + + return Utilities::create_evenly_distributed_partitioning(this_proc, + n_proc, + total_size); + } + /** @@ -664,6 +694,13 @@ namespace Utilities return std::vector(1, complete_index_set(local_size)); } + IndexSet + create_evenly_distributed_partitioning(const MPI_Comm & /*comm*/, + const IndexSet::size_type total_size) + { + return complete_index_set(total_size); + } + MPI_Comm diff --git a/tests/base/mpi_evenly_distributed_partitioning_01.cc b/tests/base/mpi_evenly_distributed_partitioning_01.cc new file mode 100644 index 0000000000..b64ace6ceb --- /dev/null +++ b/tests/base/mpi_evenly_distributed_partitioning_01.cc @@ -0,0 +1,97 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2017 - 2020 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.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + +// Test Utilities::MPI::create_evenly_distributed_partitioning() + +#include + +#include + +#include "../tests.h" + +void +test() +{ + MPI_Comm comm(MPI_COMM_WORLD); + const auto my_proc = Utilities::MPI::this_mpi_process(comm); + const auto n_proc = Utilities::MPI::n_mpi_processes(comm); + + const auto total_size = 17; + const IndexSet index_set = + Utilities::MPI::create_evenly_distributed_partitioning(comm, total_size); + + std::vector vec(n_proc); + for (unsigned int i_proc = 0; i_proc < n_proc; ++i_proc) + { + vec[i_proc] = + Utilities::create_evenly_distributed_partitioning(i_proc, + n_proc, + total_size); + } + // correct number: + AssertThrow(index_set == vec[my_proc], ExcInternalError()); + + // ascending one-to-one + AssertThrow(vec[my_proc].is_ascending_and_one_to_one(comm), + ExcInternalError()); + + // same size: + const auto size = vec[0].size(); + for (unsigned int i = 1; i < vec.size(); ++i) + AssertThrow(size == vec[i].size(), ExcInternalError()); + + // Local sizes differ by one or less + auto local_size_1 = vec[0].n_elements(); + for (unsigned int i = 1; i < vec.size(); ++i) + { + auto local_size_2 = vec[i].n_elements(); + AssertThrow(std::abs((int)local_size_1 - (int)local_size_2) <= 1, + ExcInternalError()); + } + + // 1:1 + for (unsigned int i = 0; i < vec.size(); ++i) + for (unsigned int j = i + 1; j < vec.size(); ++j) + { + const IndexSet intersection = vec[i] & vec[j]; + AssertThrow(intersection.is_empty(), ExcInternalError()); + } + + // union of all + IndexSet all(size); + for (unsigned int i = 0; i < vec.size(); ++i) + all.add_indices(vec[i]); + + AssertThrow(all == complete_index_set(size), ExcInternalError()); + + // finally print: + deallog << "Total size=" << size << std::endl; + for (unsigned int p = 0; p < vec.size(); ++p) + { + deallog << "p=" << p << ":" << std::endl; + vec[p].print(deallog.get_file_stream()); + } +} + +int +main(int argc, char *argv[]) +{ + Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1); + + MPILogInitAll log; + + test(); + return 0; +} diff --git a/tests/base/mpi_evenly_distributed_partitioning_01.mpirun=1.output b/tests/base/mpi_evenly_distributed_partitioning_01.mpirun=1.output new file mode 100644 index 0000000000..e29cd194c0 --- /dev/null +++ b/tests/base/mpi_evenly_distributed_partitioning_01.mpirun=1.output @@ -0,0 +1,4 @@ + +DEAL:0::Total size=17 +DEAL:0::p=0: +{[0,16]} diff --git a/tests/base/mpi_evenly_distributed_partitioning_01.mpirun=2.output b/tests/base/mpi_evenly_distributed_partitioning_01.mpirun=2.output new file mode 100644 index 0000000000..d5dc835a11 --- /dev/null +++ b/tests/base/mpi_evenly_distributed_partitioning_01.mpirun=2.output @@ -0,0 +1,13 @@ + +DEAL:0::Total size=17 +DEAL:0::p=0: +{[0,8]} +DEAL:0::p=1: +{[9,16]} + +DEAL:1::Total size=17 +DEAL:1::p=0: +{[0,8]} +DEAL:1::p=1: +{[9,16]} + diff --git a/tests/base/mpi_evenly_distributed_partitioning_01.mpirun=4.output b/tests/base/mpi_evenly_distributed_partitioning_01.mpirun=4.output new file mode 100644 index 0000000000..808d3928a1 --- /dev/null +++ b/tests/base/mpi_evenly_distributed_partitioning_01.mpirun=4.output @@ -0,0 +1,43 @@ + +DEAL:0::Total size=17 +DEAL:0::p=0: +{[0,4]} +DEAL:0::p=1: +{[5,8]} +DEAL:0::p=2: +{[9,12]} +DEAL:0::p=3: +{[13,16]} + +DEAL:1::Total size=17 +DEAL:1::p=0: +{[0,4]} +DEAL:1::p=1: +{[5,8]} +DEAL:1::p=2: +{[9,12]} +DEAL:1::p=3: +{[13,16]} + + +DEAL:2::Total size=17 +DEAL:2::p=0: +{[0,4]} +DEAL:2::p=1: +{[5,8]} +DEAL:2::p=2: +{[9,12]} +DEAL:2::p=3: +{[13,16]} + + +DEAL:3::Total size=17 +DEAL:3::p=0: +{[0,4]} +DEAL:3::p=1: +{[5,8]} +DEAL:3::p=2: +{[9,12]} +DEAL:3::p=3: +{[13,16]} +