From: Eldar Date: Tue, 4 Sep 2018 21:12:15 +0000 (-0500) Subject: Simplified p2p comm computation X-Git-Tag: v9.1.0-rc1~749^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=fe0d62c0f64c188065d590738bfac833d2d97908;p=dealii.git Simplified p2p comm computation --- diff --git a/doc/news/changes/minor/20180904EldarKhattatov b/doc/news/changes/minor/20180904EldarKhattatov new file mode 100644 index 0000000000..0434ff1d61 --- /dev/null +++ b/doc/news/changes/minor/20180904EldarKhattatov @@ -0,0 +1,6 @@ +Added: New function +Utilities::MPI::compute_n_point_to_point_communications() +is implemented to be used for computing number of processes +in an MPI universe to expect communication from. +
+(Eldar Khattatov, 2018/09/04) diff --git a/include/deal.II/base/mpi.h b/include/deal.II/base/mpi.h index f5966c2bc0..d2482be792 100644 --- a/include/deal.II/base/mpi.h +++ b/include/deal.II/base/mpi.h @@ -117,6 +117,30 @@ namespace Utilities const MPI_Comm & mpi_comm, const std::vector &destinations); + /** + * Simplified (for efficiency) version of the + * compute_point_to_point_communication_pattern() + * which only computes the number of processes in an MPI universe to expect + * communication from. + * + * @param mpi_comm A + * @ref GlossMPICommunicator "communicator" + * that describes the processors that are going to communicate with each + * other. + * + * @param destinations The list of processors the current process wants to + * send information to. This list need not be sorted in any way. If it + * contains duplicate entries that means that multiple messages are + * intended for a given destination. + * + * @return A number of processors that want to send something to the current + * processor. + */ + unsigned int + compute_n_point_to_point_communications( + const MPI_Comm & mpi_comm, + const std::vector &destinations); + /** * Given a * @ref GlossMPICommunicator "communicator", diff --git a/include/deal.II/fe/fe_tools_extrapolate.templates.h b/include/deal.II/fe/fe_tools_extrapolate.templates.h index 01492c650f..6f073f5ffd 100644 --- a/include/deal.II/fe/fe_tools_extrapolate.templates.h +++ b/include/deal.II/fe/fe_tools_extrapolate.templates.h @@ -1170,14 +1170,14 @@ namespace FETools Assert(destinations.size() == cells_to_send.size(), ExcInternalError()); - std::vector senders = - Utilities::MPI::compute_point_to_point_communication_pattern( - communicator, destinations); + const unsigned int n_senders = + Utilities::MPI::compute_n_point_to_point_communications(communicator, + destinations); // receive data std::vector receive; CellData cell_data; - for (unsigned int index = 0; index < senders.size(); ++index) + for (unsigned int index = 0; index < n_senders; ++index) { MPI_Status status; int len; diff --git a/source/base/mpi.cc b/source/base/mpi.cc index b7fa8393c2..2ecf9fcd5c 100644 --- a/source/base/mpi.cc +++ b/source/base/mpi.cc @@ -293,6 +293,68 @@ namespace Utilities } + + unsigned int + compute_n_point_to_point_communications( + const MPI_Comm & mpi_comm, + const std::vector &destinations) + { + const unsigned int myid = Utilities::MPI::this_mpi_process(mpi_comm); + const unsigned int n_procs = Utilities::MPI::n_mpi_processes(mpi_comm); + + for (unsigned int i = 0; i < destinations.size(); ++i) + { + Assert(destinations[i] < n_procs, + ExcIndexRange(destinations[i], 0, n_procs)); + Assert(destinations[i] != myid, + ExcMessage( + "There is no point in communicating with ourselves.")); + } + + // Calculate the number of messages to send to each process + std::vector dest_vector(n_procs); + for (const auto &el : destinations) + ++dest_vector[el]; + +# if DEAL_II_MPI_VERSION_GTE(2, 2) + // Find out how many processes will send to this one + // MPI_Reduce_scatter(_block) does exactly this + unsigned int n_recv_from; + const int ierr = MPI_Reduce_scatter_block( + &dest_vector[0], &n_recv_from, 1, MPI_UNSIGNED, MPI_SUM, mpi_comm); + + AssertThrowMPI(ierr); + + return n_recv_from; +# else + // Find out how many processes will send to this one + // by reducing with sum and then scattering the + // results over all processes + std::vector buffer(dest_vector.size()); + unsigned int n_recv_from; + + MPI_Reduce(&dest_vector[0], + &buffer[0], + dest_vector.size(), + MPI_UNSIGNED, + MPI_SUM, + 0, + mpi_comm); + MPI_Scatter(&buffer[0], + 1, + MPI_UNSIGNED, + &n_recv_from, + 1, + MPI_UNSIGNED, + 0, + mpi_comm); + + return n_recv_from; +# endif + } + + + namespace { // custom MIP_Op for calculate_collective_mpi_min_max_avg diff --git a/source/distributed/tria.cc b/source/distributed/tria.cc index 8e16f9c3da..d3cb89ba98 100644 --- a/source/distributed/tria.cc +++ b/source/distributed/tria.cc @@ -4262,14 +4262,14 @@ namespace parallel // collect the neighbors // that are going to send stuff to us - const std::vector senders = - Utilities::MPI::compute_point_to_point_communication_pattern( + const unsigned int n_senders = + Utilities::MPI::compute_n_point_to_point_communications( this->get_communicator(), destinations); // receive ghostcelldata std::vector receive; CommunicateLocallyMovedVertices::CellInfo cellinfo; - for (unsigned int i = 0; i < senders.size(); ++i) + for (unsigned int i = 0; i < n_senders; ++i) { MPI_Status status; int len; @@ -4329,7 +4329,7 @@ namespace parallel // check all msgs got sent and received Assert(Utilities::MPI::sum(needs_to_get_cells.size(), this->get_communicator()) == - Utilities::MPI::sum(senders.size(), this->get_communicator()), + Utilities::MPI::sum(n_senders, this->get_communicator()), ExcInternalError()); } diff --git a/source/lac/sparsity_tools.cc b/source/lac/sparsity_tools.cc index d060d882cd..0d36b69f9d 100644 --- a/source/lac/sparsity_tools.cc +++ b/source/lac/sparsity_tools.cc @@ -987,9 +987,8 @@ namespace SparsityTools send_to.push_back(it->first); num_receive = - Utilities::MPI::compute_point_to_point_communication_pattern(mpi_comm, - send_to) - .size(); + Utilities::MPI::compute_n_point_to_point_communications(mpi_comm, + send_to); } std::vector requests(send_data.size()); @@ -1134,9 +1133,8 @@ namespace SparsityTools send_to.push_back(it->first); num_receive = - Utilities::MPI::compute_point_to_point_communication_pattern(mpi_comm, - send_to) - .size(); + Utilities::MPI::compute_n_point_to_point_communications(mpi_comm, + send_to); } std::vector requests(send_data.size()); diff --git a/tests/mpi/point_to_point_pattern_01.cc b/tests/mpi/point_to_point_pattern_01.cc index caad9858ad..d654185656 100644 --- a/tests/mpi/point_to_point_pattern_01.cc +++ b/tests/mpi/point_to_point_pattern_01.cc @@ -17,7 +17,13 @@ // Check if point-to-point communication // works correctly with MPI_Reduce_scatter_block -// implementation +// implementation. +// Also check that the simplified point-to-point communication +// computation agrees with the original function. +// Namely, the size of +// compute_point_to_point_communication_pattern() +// returned vector should agree with the number returned +// by compute_point_to_point_communications() #include @@ -76,6 +82,14 @@ test_mpi() std::vector origins = Utilities::MPI::compute_point_to_point_communication_pattern(MPI_COMM_WORLD, destinations); + + const unsigned int n_origins = + Utilities::MPI::compute_n_point_to_point_communications(MPI_COMM_WORLD, + destinations); + + if (origins.size() != n_origins) + deallog << "Size mismatch!" << std::endl; + std::sort(origins.begin(), origins.end()); if (myid == 0)