Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_Details_iallreduce.cpp
1 // @HEADER
2 // *****************************************************************************
3 // Tpetra: Templated Linear Algebra Services Package
4 //
5 // Copyright 2008 NTESS and the Tpetra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
11 
12 #ifdef HAVE_TPETRACORE_MPI
13 # include "Teuchos_DefaultMpiComm.hpp" // only needs to be in .cpp file
14 #endif // HAVE_TPETRACORE_MPI
15 #include "Teuchos_DefaultSerialComm.hpp" // only needs to be in .cpp file
16 
17 namespace Tpetra {
18 namespace Details {
19 
20 #ifdef HAVE_TPETRACORE_MPI
21 std::string getMpiErrorString (const int errCode) {
22  // Space for storing the error string returned by MPI.
23  // Leave room for null termination, since I don't know if MPI does this.
24  char errString [MPI_MAX_ERROR_STRING+1];
25  int errStringLen = MPI_MAX_ERROR_STRING; // output argument
26  (void) MPI_Error_string (errCode, errString, &errStringLen);
27  // errStringLen on output is the number of characters written.
28  // I'm not sure (the MPI 3.0 Standard doesn't say) if this
29  // includes the '\0', so I'll make sure. We reserved space for
30  // the extra '\0' if needed.
31  if (errString[errStringLen-1] != '\0') {
32  errString[errStringLen] = '\0';
33  }
34  return std::string (errString); // This copies the original string.
35 }
36 #endif // HAVE_TPETRACORE_MPI
37 
38 namespace Impl {
39 
40 std::shared_ptr<CommRequest>
41 emptyCommRequest ()
42 {
43  return std::shared_ptr<CommRequest> (new CommRequest());
44 }
45 
46 #ifdef HAVE_TPETRACORE_MPI
47 
48 #if MPI_VERSION >= 3
49 MPI_Request
50 iallreduceRaw (const void* sendbuf,
51  void* recvbuf,
52  const int count,
53  MPI_Datatype mpiDatatype,
54  const Teuchos::EReductionType op,
55  MPI_Comm comm)
56 {
57  MPI_Op rawOp = ::Teuchos::Details::getMpiOpForEReductionType (op);
58  MPI_Request req = MPI_REQUEST_NULL;
59  int err = MPI_SUCCESS;
60  if (sendbuf == recvbuf) {
61  // Fix for #850. This only works if rawComm is an
62  // intracommunicator. Intercommunicators don't have an in-place
63  // option for collectives.
64  err = MPI_Iallreduce (MPI_IN_PLACE, recvbuf, count, mpiDatatype,
65  rawOp, comm, &req);
66  }
67  else {
68  err = MPI_Iallreduce (sendbuf, recvbuf, count, mpiDatatype,
69  rawOp, comm, &req);
70  }
71  TEUCHOS_TEST_FOR_EXCEPTION
72  (err != MPI_SUCCESS, std::runtime_error,
73  "MPI_Iallreduce failed with the following error: "
74  << getMpiErrorString (err));
75  return req;
76 }
77 #endif //MPI >= 3
78 
79 void
80 allreduceRaw (const void* sendbuf,
81  void* recvbuf,
82  const int count,
83  MPI_Datatype mpiDatatype,
84  const Teuchos::EReductionType op,
85  MPI_Comm comm)
86 {
87  MPI_Op rawOp = ::Teuchos::Details::getMpiOpForEReductionType (op);
88  int err = MPI_SUCCESS;
89  if (sendbuf == recvbuf) {
90  err = MPI_Allreduce (MPI_IN_PLACE, recvbuf,
91  count, mpiDatatype, rawOp, comm);
92  }
93  else {
94  // OpenMPI 1.6.5 insists on void*, not const void*, for sendbuf.
95  (void) MPI_Allreduce (const_cast<void*> (sendbuf), recvbuf,
96  count, mpiDatatype, rawOp, comm);
97  }
98  TEUCHOS_TEST_FOR_EXCEPTION
99  (err != MPI_SUCCESS, std::runtime_error,
100  "MPI_Allreduce failed with the following error: "
101  << getMpiErrorString (err));
102 }
103 
104 #endif // HAVE_TPETRACORE_MPI
105 
106 } // namespace Impl
107 
108 std::shared_ptr<CommRequest>
109 iallreduce (const int localValue,
110  int& globalValue,
111  const ::Teuchos::EReductionType op,
112  const ::Teuchos::Comm<int>& comm)
113 {
114  //Input: needs to be an owning view containing localValue
115  Kokkos::View<int*, Kokkos::HostSpace> localView(
116  Kokkos::ViewAllocateWithoutInitializing("localValue"), 1);
117  localView(0) = localValue;
118  Kokkos::View<int*, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>>
119  globalView(&globalValue, 1);
120  return ::Tpetra::Details::iallreduce<decltype(localView), decltype(globalView)>(localView, globalView, op, comm);
121 }
122 
123 } // namespace Details
124 } // namespace Tpetra
125 
Declaration of Tpetra::Details::iallreduce.