Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

spgemm handle: check that A,B,C graphs never change #1742

Merged
merged 7 commits into from
Apr 5, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
43 changes: 43 additions & 0 deletions common/src/KokkosKernels_SimpleUtils.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -380,6 +380,49 @@ KOKKOS_FORCEINLINE_FUNCTION Value xorshiftHash(Value v) {
: static_cast<Value>(x * 2685821657736338717ULL - 1);
}

struct ViewHashFunctor {
ViewHashFunctor(const uint8_t *data_) : data(data_) {}

KOKKOS_INLINE_FUNCTION void operator()(size_t i, uint32_t &lhash) const {
// Compute a hash/digest of both the index i, and data[i]. Then add that to
// overall hash.
uint32_t x = uint32_t(i);
x ^= x << 13;
x ^= x >> 17;
x ^= x << 5;
x ^= uint32_t(data[i]);
brian-kelley marked this conversation as resolved.
Show resolved Hide resolved
x ^= x << 13;
x ^= x >> 17;
x ^= x << 5;
lhash += x;
}

const uint8_t *data;
};

/// \brief Compute a hash of a view.
/// \param v: the view to hash. Must be contiguous, and its element type must
/// not contain any padding bytes.
template <typename View>
uint32_t hashView(const View &v) {
assert(v.span_is_contiguous());
// Note: This type trait is supposed to be part of C++17,
// but it's not defined on Intel 19 (with GCC 7.2.0 standard library).
// So just check if it's available before using.
#ifdef __cpp_lib_has_unique_object_representations
static_assert(std::has_unique_object_representations<
typename View::non_const_value_type>::value,
"KokkosKernels::Impl::hashView: the view's element type must "
"not have any padding bytes.");
#endif
size_t nbytes = v.span() * sizeof(typename View::value_type);
brian-kelley marked this conversation as resolved.
Show resolved Hide resolved
uint32_t h;
Kokkos::parallel_reduce(
Kokkos::RangePolicy<typename View::execution_space, size_t>(0, nbytes),
ViewHashFunctor(reinterpret_cast<const uint8_t *>(v.data())), h);
return h;
}

template <typename V>
struct SequentialFillFunctor {
using size_type = typename V::size_type;
Expand Down
54 changes: 54 additions & 0 deletions sparse/src/KokkosSparse_spgemm_handle.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -776,6 +776,60 @@ class SPGEMMHandle {
}

bool get_compression_step() { return is_compression_single_step; }

private:
// An SpGEMM handle can be reused for multiple products C = A*B, but only if
// the sparsity patterns of A and B do not change. Enforce this (in debug
// builds only) by recording hashes of the graphs, and then checking they
// match in later calls.
bool computedInputHashes = false;
uint32_t a_graph_hash = 0U;
uint32_t b_graph_hash = 0U;

public:
template <typename a_rowptrs_t, typename a_entries_t, typename b_rowptrs_t,
typename b_entries_t>
bool checkMatrixIdentitiesSymbolic(const a_rowptrs_t &a_rowptrsIn,
const a_entries_t &a_entriesIn,
const b_rowptrs_t &b_rowptrsIn,
const b_entries_t &b_entriesIn) {
#ifndef NDEBUG
// If this is the first symbolic call, assign the handle's CRS pointers to
// check against later
if (!computedInputHashes) {
a_graph_hash = KokkosKernels::Impl::hashView(a_rowptrsIn) ^
KokkosKernels::Impl::hashView(a_entriesIn);
b_graph_hash = KokkosKernels::Impl::hashView(b_rowptrsIn) ^
KokkosKernels::Impl::hashView(b_entriesIn);
computedInputHashes = true;
} else {
if (a_graph_hash != (KokkosKernels::Impl::hashView(a_rowptrsIn) ^
KokkosKernels::Impl::hashView(a_entriesIn)))
return false;
if (b_graph_hash != (KokkosKernels::Impl::hashView(b_rowptrsIn) ^
KokkosKernels::Impl::hashView(b_entriesIn)))
return false;
}
#endif
return true;
}

template <typename a_rowptrs_t, typename a_entries_t, typename b_rowptrs_t,
typename b_entries_t>
bool checkMatrixIdentitiesNumeric(const a_rowptrs_t &a_rowptrsIn,
const a_entries_t &a_entriesIn,
const b_rowptrs_t &b_rowptrsIn,
const b_entries_t &b_entriesIn) {
#ifndef NDEBUG
if (a_graph_hash != (KokkosKernels::Impl::hashView(a_rowptrsIn) ^
KokkosKernels::Impl::hashView(a_entriesIn)))
return false;
if (b_graph_hash != (KokkosKernels::Impl::hashView(b_rowptrsIn) ^
KokkosKernels::Impl::hashView(b_entriesIn)))
return false;
#endif
return true;
}
};

inline SPGEMMAlgorithm StringToSPGEMMAlgorithm(std::string &name) {
Expand Down
19 changes: 18 additions & 1 deletion sparse/src/KokkosSparse_spgemm_numeric.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -236,7 +236,24 @@ void spgemm_numeric(KernelHandle *handle,
return;
}

auto algo = tmp_handle.get_spgemm_handle()->get_algorithm_type();
auto spgemmHandle = tmp_handle.get_spgemm_handle();

if (!spgemmHandle) {
throw std::invalid_argument(
"KokkosSparse::spgemm_numeric: the given KernelHandle does not have "
"an SpGEMM handle associated with it.");
}

if (!spgemmHandle->checkMatrixIdentitiesNumeric(const_a_r, const_a_l,
const_b_r, const_b_l)) {
throw std::invalid_argument(
"KokkosSparse::spgemm_numeric: once used, an spgemm handle cannot be "
"reused for a product with a different sparsity pattern.\n"
"The rowptrs and entries of A and B must be identical to those "
"passed to the first spgemm_symbolic and spgemm_numeric calls.");
}

auto algo = spgemmHandle->get_algorithm_type();

if (algo == SPGEMM_DEBUG || algo == SPGEMM_SERIAL) {
// Never call a TPL if serial/debug is requested (this is needed for
Expand Down
19 changes: 18 additions & 1 deletion sparse/src/KokkosSparse_spgemm_symbolic.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -151,7 +151,24 @@ void spgemm_symbolic(KernelHandle *handle,
"rows. May use KokkosSparse::sort_crs_matrix to sort it.");
#endif

auto algo = tmp_handle.get_spgemm_handle()->get_algorithm_type();
auto spgemmHandle = tmp_handle.get_spgemm_handle();

if (!spgemmHandle) {
throw std::invalid_argument(
"KokkosSparse::spgemm_symbolic: the given KernelHandle does not have "
"an SpGEMM handle associated with it.");
}

if (!spgemmHandle->checkMatrixIdentitiesSymbolic(const_a_r, const_a_l,
const_b_r, const_b_l)) {
throw std::invalid_argument(
"KokkosSparse::spgemm_symbolic: once used, an spgemm handle cannot be "
"reused for a product with a different sparsity pattern.\n"
"The rowptrs and entries of A and B must be identical to those "
"passed to the first spgemm_symbolic call.");
}

auto algo = spgemmHandle->get_algorithm_type();

if (algo == SPGEMM_DEBUG || algo == SPGEMM_SERIAL) {
// Never call a TPL if serial/debug is requested (this is needed for
Expand Down
66 changes: 62 additions & 4 deletions sparse/unit_test/Test_Sparse_spgemm.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -78,12 +78,13 @@ void run_spgemm_noreuse(crsMat_t A, crsMat_t B, crsMat_t &C) {
}

template <typename crsMat_t, typename device>
int run_spgemm(crsMat_t A, crsMat_t B,
int run_spgemm(crsMat_t &A, crsMat_t &B,
KokkosSparse::SPGEMMAlgorithm spgemm_algorithm, crsMat_t &C,
bool testReuse) {
typedef typename crsMat_t::size_type size_type;
typedef typename crsMat_t::ordinal_type lno_t;
typedef typename crsMat_t::value_type scalar_t;
typedef typename crsMat_t::values_type::non_const_type scalar_view_t;

typedef KokkosKernels::Experimental::KokkosKernelsHandle<
size_type, lno_t, scalar_t, typename device::execution_space,
Expand Down Expand Up @@ -113,7 +114,14 @@ int run_spgemm(crsMat_t A, crsMat_t B,
EXPECT_TRUE(sh->is_numeric_called());

if (testReuse) {
// Give A and B completely new random values, and re-run just numeric
// Give A and B completely new random values (changing both the pointer
brian-kelley marked this conversation as resolved.
Show resolved Hide resolved
// and contents), and re-run just numeric.
A.values = scalar_view_t(
Kokkos::view_alloc(Kokkos::WithoutInitializing, "new A values"),
A.nnz());
B.values = scalar_view_t(
Kokkos::view_alloc(Kokkos::WithoutInitializing, "new B values"),
B.nnz());
randomize_matrix_values(A.values);
randomize_matrix_values(B.values);
KokkosSparse::spgemm_numeric(kh, A, false, B, false, C);
Expand All @@ -127,7 +135,7 @@ int run_spgemm(crsMat_t A, crsMat_t B,
}

template <typename crsMat_t, typename device>
int run_spgemm_old_interface(crsMat_t A, crsMat_t B,
int run_spgemm_old_interface(crsMat_t &A, crsMat_t &B,
KokkosSparse::SPGEMMAlgorithm spgemm_algorithm,
crsMat_t &result, bool testReuse) {
typedef typename crsMat_t::StaticCrsGraphType graph_t;
Expand Down Expand Up @@ -188,7 +196,14 @@ int run_spgemm_old_interface(crsMat_t A, crsMat_t B,
EXPECT_TRUE(sh->is_numeric_called());

if (testReuse) {
// Give A and B completely new random values, and re-run just numeric
// Give A and B completely new random values (changing both the pointer
// and contents), and re-run just numeric.
A.values = scalar_view_t(
Kokkos::view_alloc(Kokkos::WithoutInitializing, "new A values"),
A.nnz());
B.values = scalar_view_t(
Kokkos::view_alloc(Kokkos::WithoutInitializing, "new B values"),
B.nnz());
randomize_matrix_values(A.values);
randomize_matrix_values(B.values);
KokkosSparse::Experimental::spgemm_numeric(
Expand Down Expand Up @@ -468,6 +483,48 @@ void test_issue402() {
<< "SpGEMM still has issue 402 bug; C=AA' is incorrect!\n";
}

template <typename scalar_t, typename lno_t, typename size_type,
typename device>
void test_issue1738() {
// Make sure that std::invalid_argument is thrown if you:
// - call numeric where an input matrix's entries have changed.
// - try to reuse an spgemm handle by calling symbolic with new input
// matrices
// This check is only enabled in debug builds.
#ifndef NDEBUG
using crsMat_t = CrsMatrix<scalar_t, lno_t, device, void, size_type>;
using KernelHandle = KokkosKernels::Experimental::KokkosKernelsHandle<
size_type, lno_t, scalar_t, typename device::execution_space,
typename device::memory_space, typename device::memory_space>;
crsMat_t A1 = KokkosSparse::Impl::kk_generate_diag_matrix<crsMat_t>(100);
crsMat_t B1 = KokkosSparse::Impl::kk_generate_diag_matrix<crsMat_t>(100);
crsMat_t A2 = KokkosSparse::Impl::kk_generate_diag_matrix<crsMat_t>(50);
crsMat_t B2 = KokkosSparse::Impl::kk_generate_diag_matrix<crsMat_t>(50);
{
KernelHandle kh;
kh.create_spgemm_handle();
crsMat_t C1;
KokkosSparse::spgemm_symbolic(kh, A1, false, B1, false, C1);
KokkosSparse::spgemm_numeric(kh, A1, false, B1, false, C1);
crsMat_t C2;
EXPECT_THROW(KokkosSparse::spgemm_symbolic(kh, A2, false, B2, false, C2),
std::invalid_argument);
}
{
KernelHandle kh;
kh.create_spgemm_handle();
crsMat_t C1;
KokkosSparse::spgemm_symbolic(kh, A1, false, B1, false, C1);
// Note: A1 is a 100x100 diagonal matrix, so the first entry in the first
// row is 0. Change it to a 1 and make sure spgemm_numeric notices that it
// changed.
Kokkos::deep_copy(Kokkos::subview(A1.graph.entries, 0), 1);
EXPECT_THROW(KokkosSparse::spgemm_numeric(kh, A1, false, B1, false, C1),
std::invalid_argument);
}
#endif
}

#define KOKKOSKERNELS_EXECUTE_TEST(SCALAR, ORDINAL, OFFSET, DEVICE) \
TEST_F(TestCategory, \
sparse##_##spgemm##_##SCALAR##_##ORDINAL##_##OFFSET##_##DEVICE) { \
Expand Down Expand Up @@ -513,6 +570,7 @@ void test_issue402() {
test_spgemm_symbolic<SCALAR, ORDINAL, OFFSET, DEVICE>(true, false); \
test_spgemm_symbolic<SCALAR, ORDINAL, OFFSET, DEVICE>(false, false); \
test_issue402<SCALAR, ORDINAL, OFFSET, DEVICE>(); \
test_issue1738<SCALAR, ORDINAL, OFFSET, DEVICE>(); \
}

// test_spgemm<SCALAR,ORDINAL,OFFSET,DEVICE>(50000, 50000 * 30, 100, 10);
Expand Down