From 35cdd4e18e27e305d5ac56ca0433e71b3f7e482f Mon Sep 17 00:00:00 2001 From: Ian Rose Date: Wed, 22 Apr 2015 16:20:57 -0700 Subject: [PATCH] Add collective MPI sums to tensorial objects Fix indentation indentation Move functions to mpi.h, address other comments Try to avoid variable length array declaration Add entry to changelog Fix slash --- doc/news/changes.h | 8 +- include/deal.II/base/mpi.h | 69 +++++++++++ tests/mpi/collective_tensor.cc | 120 +++++++++++++++++++ tests/mpi/collective_tensor.mpirun=1.output | 7 ++ tests/mpi/collective_tensor.mpirun=10.output | 7 ++ tests/mpi/collective_tensor.mpirun=4.output | 7 ++ 6 files changed, 217 insertions(+), 1 deletion(-) create mode 100644 tests/mpi/collective_tensor.cc create mode 100644 tests/mpi/collective_tensor.mpirun=1.output create mode 100644 tests/mpi/collective_tensor.mpirun=10.output create mode 100644 tests/mpi/collective_tensor.mpirun=4.output diff --git a/doc/news/changes.h b/doc/news/changes.h index 7c10dafcb5..7d33088137 100644 --- a/doc/news/changes.h +++ b/doc/news/changes.h @@ -468,13 +468,19 @@ inconvenience this causes.

Specific improvements

    +
  1. New: There are now MPI sum functions for Tensors and SymmetricTensors + in the Utilities::MPI namespace. +
    + (Ian Rose, 2015/04/24) +
  2. +
  3. Fixed: project_boundary_values_curl_conforming_l2() produced incorrect results for non-uniform grids in 2D. Adjustment to the way 2D tangents to edges are computed fixes this.
    (Ross Kynch, 2015/04/23)
  4. - +
  5. Fixed: The TimerOutput class reported abnormally large cpu time when run with more than one process with MPI. Now this is fixed.
    diff --git a/include/deal.II/base/mpi.h b/include/deal.II/base/mpi.h index 34bbc4d462..09550080bd 100644 --- a/include/deal.II/base/mpi.h +++ b/include/deal.II/base/mpi.h @@ -50,6 +50,10 @@ namespace MPI DEAL_II_NAMESPACE_OPEN +//Forward type declarations to allow MPI sums over tensorial types +template class Tensor; +template class SymmetricTensor; + namespace Utilities { /** @@ -160,6 +164,28 @@ namespace Utilities const MPI_Comm &mpi_communicator, std::vector &sums); + /** + * Perform an MPI sum of the entries of a symmetric tensor. + * + * @relates SymmetricTensor + */ + template + inline + SymmetricTensor + sum (const SymmetricTensor &local, + const MPI_Comm &mpi_communicator); + + /** + * Perform an MPI sum of the entries of a tensor. + * + * @relates Tensor + */ + template + inline + Tensor + sum (const Tensor &local, + const MPI_Comm &mpi_communicator); + /** * Return the maximum over all processors of the value @p t. This function * is collective over all processors given in the communicator. If deal.II @@ -468,6 +494,49 @@ namespace Utilities } } + template + inline + Tensor + sum (const Tensor &local, + const MPI_Comm &mpi_communicator) + { + const unsigned int n_entries = Tensor::n_independent_components; + Number entries[ Tensor::n_independent_components ]; + + for (unsigned int i=0; i< n_entries; ++i) + entries[i] = local[ local.unrolled_to_component_indices(i) ]; + + Number global_entries[ Tensor::n_independent_components ]; + dealii::Utilities::MPI::sum( entries, mpi_communicator, global_entries ); + + Tensor global; + for (unsigned int i=0; i< n_entries; ++i) + global[ global.unrolled_to_component_indices(i) ] = global_entries[i]; + + return global; + } + + template + inline + SymmetricTensor + sum (const SymmetricTensor &local, + const MPI_Comm &mpi_communicator) + { + const unsigned int n_entries = SymmetricTensor::n_independent_components; + Number entries[ SymmetricTensor::n_independent_components ]; + + for (unsigned int i=0; i< n_entries; ++i) + entries[i] = local[ local.unrolled_to_component_indices(i) ]; + + Number global_entries[ SymmetricTensor::n_independent_components ]; + dealii::Utilities::MPI::sum( entries, mpi_communicator, global_entries ); + + SymmetricTensor global; + for (unsigned int i=0; i< n_entries; ++i) + global[ global.unrolled_to_component_indices(i) ] = global_entries[i]; + + return global; + } template inline diff --git a/tests/mpi/collective_tensor.cc b/tests/mpi/collective_tensor.cc new file mode 100644 index 0000000000..53f17460bd --- /dev/null +++ b/tests/mpi/collective_tensor.cc @@ -0,0 +1,120 @@ +// --------------------------------------------------------------------- +// +// 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 tensor objects + +#include "../tests.h" +#include +#include +#include +#include + +void test() +{ + Assert( Utilities::MPI::job_supports_mpi(), ExcInternalError()); + + unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD); + const unsigned int numprocs = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD); + + if (myid==0) + deallog << "Running on " << numprocs << " CPU(s)." << std::endl; + + Tensor<1,1> tensor_1; + + tensor_1[0] = 1.0; + + Tensor<1,1> result_1 = Utilities::MPI::sum( tensor_1, MPI_COMM_WORLD ); + + if (myid==0) + deallog << "1D tensor: " << result_1 << std::endl; + + Tensor<1,2> tensor_2; + + tensor_2[0] = 1.0; + tensor_2[1] = 2.0; + + Tensor<1,2> result_2 = Utilities::MPI::sum( tensor_2, MPI_COMM_WORLD ); + + if (myid==0) + deallog << "2D tensor: " << result_2 << std::endl; + + Tensor<1,3> tensor_3; + + tensor_3[0] = 1.0; + tensor_3[1] = 2.0; + tensor_3[2] = 3.0; + + Tensor<1,3> result_3 = Utilities::MPI::sum( tensor_3, MPI_COMM_WORLD ); + + if (myid==0) + deallog << "3D tensor: " << result_3 << std::endl; + + Tensor<2,3> rank_2_tensor; + + rank_2_tensor[0][0] = 1.0; + rank_2_tensor[0][1] = 2.0; + rank_2_tensor[0][2] = 3.0; + rank_2_tensor[1][0] = 4.0; + rank_2_tensor[1][1] = 5.0; + rank_2_tensor[1][2] = 6.0; + rank_2_tensor[2][0] = 7.0; + rank_2_tensor[2][1] = 8.0; + rank_2_tensor[2][2] = 9.0; + + Tensor<2,3> result_rank_2 = Utilities::MPI::sum( rank_2_tensor, MPI_COMM_WORLD ); + + if (myid==0) + deallog << "Rank 2 tensor: " << result_rank_2 << std::endl; + + SymmetricTensor<2,2> symmetric_tensor; + + symmetric_tensor[0][0] = 1.0; + symmetric_tensor[0][1] = 2.0; + symmetric_tensor[1][1] = 3.0; + + SymmetricTensor<2,2> result_symmetric = Utilities::MPI::sum( symmetric_tensor, MPI_COMM_WORLD ); + + if (myid==0) + deallog << "Symmetric tensor: " << result_symmetric << 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_tensor.mpirun=1.output b/tests/mpi/collective_tensor.mpirun=1.output new file mode 100644 index 0000000000..6fd92b294a --- /dev/null +++ b/tests/mpi/collective_tensor.mpirun=1.output @@ -0,0 +1,7 @@ + +DEAL:mpi::Running on 1 CPU(s). +DEAL:mpi::1D tensor: 1.00000 +DEAL:mpi::2D tensor: 1.00000 2.00000 +DEAL:mpi::3D tensor: 1.00000 2.00000 3.00000 +DEAL:mpi::Rank 2 tensor: 1.00000 2.00000 3.00000 4.00000 5.00000 6.00000 7.00000 8.00000 9.00000 +DEAL:mpi::Symmetric tensor: 1.00000 2.00000 2.00000 3.00000 diff --git a/tests/mpi/collective_tensor.mpirun=10.output b/tests/mpi/collective_tensor.mpirun=10.output new file mode 100644 index 0000000000..236824d14c --- /dev/null +++ b/tests/mpi/collective_tensor.mpirun=10.output @@ -0,0 +1,7 @@ + +DEAL:mpi::Running on 10 CPU(s). +DEAL:mpi::1D tensor: 10.0000 +DEAL:mpi::2D tensor: 10.0000 20.0000 +DEAL:mpi::3D tensor: 10.0000 20.0000 30.0000 +DEAL:mpi::Rank 2 tensor: 10.0000 20.0000 30.0000 40.0000 50.0000 60.0000 70.0000 80.0000 90.0000 +DEAL:mpi::Symmetric tensor: 10.0000 20.0000 20.0000 30.0000 diff --git a/tests/mpi/collective_tensor.mpirun=4.output b/tests/mpi/collective_tensor.mpirun=4.output new file mode 100644 index 0000000000..3e5d55c1c9 --- /dev/null +++ b/tests/mpi/collective_tensor.mpirun=4.output @@ -0,0 +1,7 @@ + +DEAL:mpi::Running on 4 CPU(s). +DEAL:mpi::1D tensor: 4.00000 +DEAL:mpi::2D tensor: 4.00000 8.00000 +DEAL:mpi::3D tensor: 4.00000 8.00000 12.0000 +DEAL:mpi::Rank 2 tensor: 4.00000 8.00000 12.0000 16.0000 20.0000 24.0000 28.0000 32.0000 36.0000 +DEAL:mpi::Symmetric tensor: 4.00000 8.00000 8.00000 12.0000 -- 2.39.5