From: Denis Davydov Date: Thu, 9 Jul 2015 11:10:54 +0000 (+0200) Subject: implement Utilities::MPI::sum for Vector X-Git-Tag: v8.3.0-rc1~23^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F1093%2Fhead;p=dealii.git implement Utilities::MPI::sum for Vector --- diff --git a/include/deal.II/base/mpi.h b/include/deal.II/base/mpi.h index 54ed4017cf..5ffd0b4a71 100644 --- a/include/deal.II/base/mpi.h +++ b/include/deal.II/base/mpi.h @@ -57,6 +57,10 @@ DEAL_II_NAMESPACE_OPEN template class Tensor; template class SymmetricTensor; +//Forward type declaration to allow MPI sums over Vector type +template class Vector; + + namespace Utilities { /** @@ -166,6 +170,19 @@ namespace Utilities const MPI_Comm &mpi_communicator, std::vector &sums); + /** + * Like the previous function, but take the sums over the elements of a + * Vector. + * + * Input and output vectors may be the same. + */ + template + inline + void sum (const Vector &values, + const MPI_Comm &mpi_communicator, + Vector &sums); + + /** * Perform an MPI sum of the entries of a symmetric tensor. * @@ -547,6 +564,37 @@ namespace Utilities } } + template + inline + void all_reduce (const MPI_Op &mpi_op, + const Vector &values, + const MPI_Comm &mpi_communicator, + Vector &output) + { +#ifdef DEAL_II_WITH_MPI + if (job_supports_mpi()) + { + if (values.begin() != output.begin()) + output.reinit (values.size()); + + MPI_Allreduce ((values.begin() != output.begin() + ? + const_cast(static_cast(values.begin())) + : + MPI_IN_PLACE), + output.begin(), values.size(), internal::mpi_type_id((T *)0), mpi_op, + mpi_communicator); + } + else +#endif + { + (void)mpi_op; + (void)mpi_communicator; + output = values; + } + } + + } @@ -579,6 +627,16 @@ namespace Utilities internal::all_reduce(MPI_SUM, values, mpi_communicator, sums); } + template + inline + void sum (const Vector &values, + const MPI_Comm &mpi_communicator, + Vector &sums) + { + internal::all_reduce(MPI_SUM, values, mpi_communicator, sums); + } + + template inline Tensor diff --git a/tests/mpi/collective_02_dealii_vector.cc b/tests/mpi/collective_02_dealii_vector.cc new file mode 100644 index 0000000000..9ab0d66b72 --- /dev/null +++ b/tests/mpi/collective_02_dealii_vector.cc @@ -0,0 +1,88 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2009 - 2015 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 Utilities::MPI::sum() for dealii::Vector + +#include "../tests.h" +#include +#include +#include +#include + +void test() +{ + unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD); + const unsigned int numprocs = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD); + + { + Vector values(2); + values[0] = 1; + values[1] = 2; + Vector sums(2); + Utilities::MPI::sum (values, + MPI_COMM_WORLD, + sums); + Assert (sums[0] == numprocs, ExcInternalError()); + Assert (sums[1] == 2*numprocs, ExcInternalError()); + + if (myid==0) + deallog << sums[0] << ' ' << sums[1] << std::endl; + } + + { + Vector values(2); + values[0] = 1.5; + values[1] = 2.5; + Vector sums; + Utilities::MPI::sum (values, + MPI_COMM_WORLD, + sums); + Assert (sums[0] == 1.5 * numprocs, ExcInternalError()); + Assert (sums[1] == 2.5 * numprocs, ExcInternalError()); + + if (myid==0) + deallog << sums[0] << ' ' << sums[1] << std::endl; + } + +} + + +int main(int argc, char *argv[]) +{ +#ifdef DEAL_II_WITH_MPI + Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, numbers::invalid_unsigned_int); +#else + (void)argc; + (void)argv; + compile_time_error; + +#endif + + if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0) + { + std::ofstream logfile("output"); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + deallog.push("mpi"); + test(); + deallog.pop(); + } + else + test(); +} diff --git a/tests/mpi/collective_02_dealii_vector.mpirun=1.output b/tests/mpi/collective_02_dealii_vector.mpirun=1.output new file mode 100644 index 0000000000..af1578854f --- /dev/null +++ b/tests/mpi/collective_02_dealii_vector.mpirun=1.output @@ -0,0 +1,3 @@ + +DEAL:mpi::1 2 +DEAL:mpi::1.50000 2.50000 diff --git a/tests/mpi/collective_02_dealii_vector.mpirun=4.output b/tests/mpi/collective_02_dealii_vector.mpirun=4.output new file mode 100644 index 0000000000..8819d0d87b --- /dev/null +++ b/tests/mpi/collective_02_dealii_vector.mpirun=4.output @@ -0,0 +1,3 @@ + +DEAL:mpi::4 8 +DEAL:mpi::6.00000 10.0000 diff --git a/tests/mpi/collective_02_dealii_vector_in_place.cc b/tests/mpi/collective_02_dealii_vector_in_place.cc new file mode 100644 index 0000000000..c868f6d1aa --- /dev/null +++ b/tests/mpi/collective_02_dealii_vector_in_place.cc @@ -0,0 +1,86 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2009 - 2015 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 Utilities::MPI::sum() for dealii::Vector, but with input=output + +#include "../tests.h" +#include +#include +#include +#include + +void test() +{ + unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD); + const unsigned int numprocs = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD); + + { + Vector sums(2); + sums[0] = 1; + sums[1] = 2; + Utilities::MPI::sum (sums, + MPI_COMM_WORLD, + sums); + Assert (sums[0] == numprocs, ExcInternalError()); + Assert (sums[1] == 2*numprocs, ExcInternalError()); + + if (myid==0) + deallog << sums[0] << ' ' << sums[1] << std::endl; + } + + { + Vector sums(2); + sums[0] = 1.5; + sums[1] = 2.5; + Utilities::MPI::sum (sums, + MPI_COMM_WORLD, + sums); + Assert (sums[0] == 1.5 * numprocs, ExcInternalError()); + Assert (sums[1] == 2.5 * numprocs, ExcInternalError()); + + if (myid==0) + deallog << sums[0] << ' ' << sums[1] << std::endl; + } + +} + + +int main(int argc, char *argv[]) +{ +#ifdef DEAL_II_WITH_MPI + Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, numbers::invalid_unsigned_int); +#else + (void)argc; + (void)argv; + compile_time_error; + +#endif + + if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0) + { + std::ofstream logfile("output"); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + deallog.push("mpi"); + test(); + deallog.pop(); + } + else + test(); +} diff --git a/tests/mpi/collective_02_dealii_vector_in_place.mpirun=1.output b/tests/mpi/collective_02_dealii_vector_in_place.mpirun=1.output new file mode 100644 index 0000000000..af1578854f --- /dev/null +++ b/tests/mpi/collective_02_dealii_vector_in_place.mpirun=1.output @@ -0,0 +1,3 @@ + +DEAL:mpi::1 2 +DEAL:mpi::1.50000 2.50000 diff --git a/tests/mpi/collective_02_dealii_vector_in_place.mpirun=4.output b/tests/mpi/collective_02_dealii_vector_in_place.mpirun=4.output new file mode 100644 index 0000000000..8819d0d87b --- /dev/null +++ b/tests/mpi/collective_02_dealii_vector_in_place.mpirun=4.output @@ -0,0 +1,3 @@ + +DEAL:mpi::4 8 +DEAL:mpi::6.00000 10.0000