Skip to content

Commit

Permalink
Add tests for Ready and Synchronous non-blocking send (isend)
Browse files Browse the repository at this point in the history
  • Loading branch information
dssgabriel committed Apr 11, 2024
1 parent c7e53ed commit 8f9ceed
Show file tree
Hide file tree
Showing 3 changed files with 190 additions and 0 deletions.
2 changes: 2 additions & 0 deletions unit_tests/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -44,6 +44,8 @@ target_link_libraries(test-mpi MPI::MPI_CXX)
# Kokkos Comm tests
add_executable(test-main test_main.cpp
test_isendrecv.cpp
test_irsendrecv.cpp
test_issendrecv.cpp
test_reduce.cpp
test_sendrecv.cpp
test_rsendrecv.cpp
Expand Down
94 changes: 94 additions & 0 deletions unit_tests/test_irsendrecv.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,94 @@
//@HEADER
// ************************************************************************
//
// Kokkos v. 4.0
// Copyright (2022) National Technology & Engineering
// Solutions of Sandia, LLC (NTESS).
//
// Under the terms of Contract DE-NA0003525 with NTESS,
// the U.S. Government retains certain rights in this software.
//
// Part of Kokkos, under the Apache License v2.0 with LLVM Exceptions.
// See https://kokkos.org/LICENSE for license information.
// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
//
//@HEADER

#include <gtest/gtest.h>

#include "KokkosComm.hpp"

template <typename T>
class IrsendRecv : public testing::Test {
public:
using Scalar = T;
};

using ScalarTypes =
::testing::Types<float, double, Kokkos::complex<float>,
Kokkos::complex<double>, int, unsigned, int64_t, size_t>;
TYPED_TEST_SUITE(IrsendRecv, ScalarTypes);

TYPED_TEST(IrsendRecv, 1D_contig) {
Kokkos::View<typename TestFixture::Scalar *> a("a", 1000);

int rank, size;
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
MPI_Comm_size(MPI_COMM_WORLD, &size);
if (size < 2) {
GTEST_SKIP() << "Requires >= 2 ranks (" << size << " provided)";
}

if (0 == rank) {
int dst = 1;
Kokkos::parallel_for(
a.extent(0), KOKKOS_LAMBDA(const int i) { a(i) = i; });
KokkosComm::Req req = KokkosComm::isend<KokkosComm::Mode::Ready>(
Kokkos::DefaultExecutionSpace(), a, dst, 0, MPI_COMM_WORLD);
req.wait();
} else if (1 == rank) {
int src = 0;
KokkosComm::recv(Kokkos::DefaultExecutionSpace(), a, src, 0,
MPI_COMM_WORLD);
int errs;
Kokkos::parallel_reduce(
a.extent(0),
KOKKOS_LAMBDA(const int &i, int &lsum) {
lsum += a(i) != typename TestFixture::Scalar(i);
},
errs);
ASSERT_EQ(errs, 0);
}
}

TYPED_TEST(IrsendRecv, 1D_noncontig) {
// this is C-style layout, i.e. b(0,0) is next to b(0,1)
Kokkos::View<typename TestFixture::Scalar **, Kokkos::LayoutRight> b("a", 10,
10);
auto a =
Kokkos::subview(b, Kokkos::ALL, 2); // take column 2 (non-contiguous)

int rank;
MPI_Comm_rank(MPI_COMM_WORLD, &rank);

if (0 == rank) {
int dst = 1;
Kokkos::parallel_for(
a.extent(0), KOKKOS_LAMBDA(const int i) { a(i) = i; });
KokkosComm::Req req = KokkosComm::isend<KokkosComm::Mode::Ready>(
Kokkos::DefaultExecutionSpace(), a, dst, 0, MPI_COMM_WORLD);
req.wait();
} else if (1 == rank) {
int src = 0;
KokkosComm::recv(Kokkos::DefaultExecutionSpace(), a, src, 0,
MPI_COMM_WORLD);
int errs;
Kokkos::parallel_reduce(
a.extent(0),
KOKKOS_LAMBDA(const int &i, int &lsum) {
lsum += a(i) != typename TestFixture::Scalar(i);
},
errs);
ASSERT_EQ(errs, 0);
}
}
94 changes: 94 additions & 0 deletions unit_tests/test_issendrecv.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,94 @@
//@HEADER
// ************************************************************************
//
// Kokkos v. 4.0
// Copyright (2022) National Technology & Engineering
// Solutions of Sandia, LLC (NTESS).
//
// Under the terms of Contract DE-NA0003525 with NTESS,
// the U.S. Government retains certain rights in this software.
//
// Part of Kokkos, under the Apache License v2.0 with LLVM Exceptions.
// See https://kokkos.org/LICENSE for license information.
// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
//
//@HEADER

#include <gtest/gtest.h>

#include "KokkosComm.hpp"

template <typename T>
class IssendRecv : public testing::Test {
public:
using Scalar = T;
};

using ScalarTypes =
::testing::Types<float, double, Kokkos::complex<float>,
Kokkos::complex<double>, int, unsigned, int64_t, size_t>;
TYPED_TEST_SUITE(IssendRecv, ScalarTypes);

TYPED_TEST(IssendRecv, 1D_contig) {
Kokkos::View<typename TestFixture::Scalar *> a("a", 1000);

int rank, size;
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
MPI_Comm_size(MPI_COMM_WORLD, &size);
if (size < 2) {
GTEST_SKIP() << "Requires >= 2 ranks (" << size << " provided)";
}

if (0 == rank) {
int dst = 1;
Kokkos::parallel_for(
a.extent(0), KOKKOS_LAMBDA(const int i) { a(i) = i; });
KokkosComm::Req req = KokkosComm::isend<KokkosComm::Mode::Synchronous>(
Kokkos::DefaultExecutionSpace(), a, dst, 0, MPI_COMM_WORLD);
req.wait();
} else if (1 == rank) {
int src = 0;
KokkosComm::recv(Kokkos::DefaultExecutionSpace(), a, src, 0,
MPI_COMM_WORLD);
int errs;
Kokkos::parallel_reduce(
a.extent(0),
KOKKOS_LAMBDA(const int &i, int &lsum) {
lsum += a(i) != typename TestFixture::Scalar(i);
},
errs);
ASSERT_EQ(errs, 0);
}
}

TYPED_TEST(IssendRecv, 1D_noncontig) {
// this is C-style layout, i.e. b(0,0) is next to b(0,1)
Kokkos::View<typename TestFixture::Scalar **, Kokkos::LayoutRight> b("a", 10,
10);
auto a =
Kokkos::subview(b, Kokkos::ALL, 2); // take column 2 (non-contiguous)

int rank;
MPI_Comm_rank(MPI_COMM_WORLD, &rank);

if (0 == rank) {
int dst = 1;
Kokkos::parallel_for(
a.extent(0), KOKKOS_LAMBDA(const int i) { a(i) = i; });
KokkosComm::Req req = KokkosComm::isend<KokkosComm::Mode::Synchronous>(
Kokkos::DefaultExecutionSpace(), a, dst, 0, MPI_COMM_WORLD);
req.wait();
} else if (1 == rank) {
int src = 0;
KokkosComm::recv(Kokkos::DefaultExecutionSpace(), a, src, 0,
MPI_COMM_WORLD);
int errs;
Kokkos::parallel_reduce(
a.extent(0),
KOKKOS_LAMBDA(const int &i, int &lsum) {
lsum += a(i) != typename TestFixture::Scalar(i);
},
errs);
ASSERT_EQ(errs, 0);
}
}

0 comments on commit 8f9ceed

Please sign in to comment.