10#include "Tpetra_ConfigDefs.hpp"
13#ifdef HAVE_TPETRACORE_MPI
14#include "Teuchos_DefaultMpiComm.hpp"
16#include "Teuchos_DefaultSerialComm.hpp"
21#ifdef HAVE_TPETRACORE_MPI
22std::string getMpiErrorString(
const int errCode) {
25 char errString[MPI_MAX_ERROR_STRING + 1];
26 int errStringLen = MPI_MAX_ERROR_STRING;
27 (void)MPI_Error_string(errCode, errString, &errStringLen);
32 if (errString[errStringLen - 1] !=
'\0') {
33 errString[errStringLen] =
'\0';
35 return std::string(errString);
41std::shared_ptr<CommRequest>
43 return std::shared_ptr<CommRequest>(
new CommRequest());
46#ifdef HAVE_TPETRACORE_MPI
50iallreduceRaw(
const void* sendbuf,
53 MPI_Datatype mpiDatatype,
54 const Teuchos::EReductionType op,
56 MPI_Op rawOp = ::Teuchos::Details::getMpiOpForEReductionType(op);
57 MPI_Request req = MPI_REQUEST_NULL;
58 int err = MPI_SUCCESS;
59 if (sendbuf == recvbuf) {
63 err = MPI_Iallreduce(MPI_IN_PLACE, recvbuf, count, mpiDatatype,
66 err = MPI_Iallreduce(sendbuf, recvbuf, count, mpiDatatype,
69 TEUCHOS_TEST_FOR_EXCEPTION(err != MPI_SUCCESS, std::runtime_error,
70 "MPI_Iallreduce failed with the following error: "
71 << getMpiErrorString(err));
76void allreduceRaw(
const void* sendbuf,
79 MPI_Datatype mpiDatatype,
80 const Teuchos::EReductionType op,
82 MPI_Op rawOp = ::Teuchos::Details::getMpiOpForEReductionType(op);
83 int err = MPI_SUCCESS;
84 if (sendbuf == recvbuf) {
85 err = MPI_Allreduce(MPI_IN_PLACE, recvbuf,
86 count, mpiDatatype, rawOp, comm);
89 (void)MPI_Allreduce(
const_cast<void*
>(sendbuf), recvbuf,
90 count, mpiDatatype, rawOp, comm);
92 TEUCHOS_TEST_FOR_EXCEPTION(err != MPI_SUCCESS, std::runtime_error,
93 "MPI_Allreduce failed with the following error: "
94 << getMpiErrorString(err));
101template <
class ValueType>
102std::shared_ptr<CommRequest>
103iallreduce(
const ValueType localValue,
104 ValueType& globalValue,
105 const ::Teuchos::EReductionType op,
106 const ::Teuchos::Comm<int>& comm) {
108 Kokkos::View<ValueType*, Kokkos::HostSpace> localView(
109 Kokkos::ViewAllocateWithoutInitializing(
"localValue"), 1);
110 localView(0) = localValue;
111 Kokkos::View<ValueType*, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>>
112 globalView(&globalValue, 1);
113 return ::Tpetra::Details::iallreduce<decltype(localView), decltype(globalView)>(localView, globalView, op, comm);
119template std::shared_ptr<Tpetra::Details::CommRequest> Tpetra::Details::iallreduce(
const int,
121 const ::Teuchos::EReductionType,
122 const ::Teuchos::Comm<int>&);
123template std::shared_ptr<Tpetra::Details::CommRequest> Tpetra::Details::iallreduce(
const Tpetra::global_size_t,
125 const ::Teuchos::EReductionType,
126 const ::Teuchos::Comm<int>&);
Declaration of Tpetra::Details::iallreduce.
Struct that holds views of the contents of a CrsMatrix.
Implementation details of Tpetra.
Namespace Tpetra contains the class and methods constituting the Tpetra library.