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

Draft: Desul ordered atomic policies + litmus tests #1616

Draft
wants to merge 14 commits into
base: develop
Choose a base branch
from
Draft
204 changes: 129 additions & 75 deletions include/RAJA/policy/desul/atomic.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -12,10 +12,9 @@

#if defined(RAJA_ENABLE_DESUL_ATOMICS)

#include "RAJA/util/macros.hpp"

#include "RAJA/policy/atomic_builtin.hpp"

#include "RAJA/policy/desul/policy.hpp"
#include "RAJA/util/macros.hpp"
#include "desul/atomics.hpp"

// Default desul options for RAJA
Expand All @@ -26,153 +25,208 @@ using raja_default_desul_scope = desul::MemoryScopeDevice;
namespace RAJA
{

namespace detail
{
template <typename T>
struct DesulAtomicPolicy {
using memory_order = raja_default_desul_order;
using memory_scope = raja_default_desul_scope;
};

template <typename OrderingPolicy, typename ScopePolicy>
struct DesulAtomicPolicy<detail_atomic_t<OrderingPolicy, ScopePolicy>> {
using memory_order = OrderingPolicy;
using memory_scope = ScopePolicy;
};

} // namespace detail

RAJA_SUPPRESS_HD_WARN
template <typename AtomicPolicy, typename T>
RAJA_HOST_DEVICE
RAJA_INLINE T
atomicAdd(AtomicPolicy, T volatile *acc, T value) {
return desul::atomic_fetch_add(const_cast<T*>(acc),
RAJA_HOST_DEVICE RAJA_INLINE T atomicAdd(AtomicPolicy, T volatile *acc, T value)
{
using desul_order =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_order;
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Does this mean that AtomicPolicy has to be detail_atomic_t<...> instead of DesulAtomicPolicy<...>

using desul_scope =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_scope;
return desul::atomic_fetch_add(const_cast<T *>(acc),
value,
raja_default_desul_order{},
raja_default_desul_scope{});
desul_order{},
desul_scope{});
}

RAJA_SUPPRESS_HD_WARN
template <typename AtomicPolicy, typename T>
RAJA_HOST_DEVICE
RAJA_INLINE T
atomicSub(AtomicPolicy, T volatile *acc, T value) {
return desul::atomic_fetch_sub(const_cast<T*>(acc),
RAJA_HOST_DEVICE RAJA_INLINE T atomicSub(AtomicPolicy, T volatile *acc, T value)
{
using desul_order =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_order;
using desul_scope =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_scope;
return desul::atomic_fetch_sub(const_cast<T *>(acc),
value,
raja_default_desul_order{},
raja_default_desul_scope{});
desul_order{},
desul_scope{});
}

RAJA_SUPPRESS_HD_WARN
template <typename AtomicPolicy, typename T>
RAJA_HOST_DEVICE
RAJA_INLINE T atomicMin(AtomicPolicy, T volatile *acc, T value)
RAJA_HOST_DEVICE RAJA_INLINE T atomicMin(AtomicPolicy, T volatile *acc, T value)
{
return desul::atomic_fetch_min(const_cast<T*>(acc),
using desul_order =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_order;
using desul_scope =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_scope;
return desul::atomic_fetch_min(const_cast<T *>(acc),
value,
raja_default_desul_order{},
raja_default_desul_scope{});
desul_order{},
desul_scope{});
}

RAJA_SUPPRESS_HD_WARN
template <typename AtomicPolicy, typename T>
RAJA_HOST_DEVICE
RAJA_INLINE T atomicMax(AtomicPolicy, T volatile *acc, T value)
RAJA_HOST_DEVICE RAJA_INLINE T atomicMax(AtomicPolicy, T volatile *acc, T value)
{
return desul::atomic_fetch_max(const_cast<T*>(acc),
using desul_order =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_order;
using desul_scope =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_scope;
return desul::atomic_fetch_max(const_cast<T *>(acc),
value,
raja_default_desul_order{},
raja_default_desul_scope{});
desul_order{},
desul_scope{});
}

RAJA_SUPPRESS_HD_WARN
template <typename AtomicPolicy, typename T>
RAJA_HOST_DEVICE
RAJA_INLINE T atomicInc(AtomicPolicy, T volatile *acc)
RAJA_HOST_DEVICE RAJA_INLINE T atomicInc(AtomicPolicy, T volatile *acc)
{
return desul::atomic_fetch_inc(const_cast<T*>(acc),
raja_default_desul_order{},
raja_default_desul_scope{});
using desul_order =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_order;
using desul_scope =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_scope;
return desul::atomic_fetch_inc(const_cast<T *>(acc),
desul_order{},
desul_scope{});
}

RAJA_SUPPRESS_HD_WARN
template <typename AtomicPolicy, typename T>
RAJA_HOST_DEVICE
RAJA_INLINE T atomicInc(AtomicPolicy, T volatile *acc, T val)
RAJA_HOST_DEVICE RAJA_INLINE T atomicInc(AtomicPolicy, T volatile *acc, T val)
{
using desul_order =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_order;
using desul_scope =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_scope;
// See:
// http://docs.nvidia.com/cuda/cuda-c-programming-guide/index.html#atomicinc
return desul::atomic_fetch_inc_mod(const_cast<T*>(acc),
val,
raja_default_desul_order{},
raja_default_desul_scope{});
return desul::atomic_fetch_inc_mod(const_cast<T *>(acc),
val,
desul_order{},
desul_scope{});
}

RAJA_SUPPRESS_HD_WARN
template <typename AtomicPolicy, typename T>
RAJA_HOST_DEVICE
RAJA_INLINE T atomicDec(AtomicPolicy, T volatile *acc)
RAJA_HOST_DEVICE RAJA_INLINE T atomicDec(AtomicPolicy, T volatile *acc)
{
return desul::atomic_fetch_dec(const_cast<T*>(acc),
raja_default_desul_order{},
raja_default_desul_scope{});
using desul_order =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_order;
using desul_scope =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_scope;
return desul::atomic_fetch_dec(const_cast<T *>(acc),
desul_order{},
desul_scope{});
}

RAJA_SUPPRESS_HD_WARN
template <typename AtomicPolicy, typename T>
RAJA_HOST_DEVICE
RAJA_INLINE T atomicDec(AtomicPolicy, T volatile *acc, T val)
RAJA_HOST_DEVICE RAJA_INLINE T atomicDec(AtomicPolicy, T volatile *acc, T val)
{
using desul_order =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_order;
using desul_scope =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_scope;
// See:
// http://docs.nvidia.com/cuda/cuda-c-programming-guide/index.html#atomicdec
return desul::atomic_fetch_dec_mod(const_cast<T*>(acc),
val,
raja_default_desul_order{},
raja_default_desul_scope{});
return desul::atomic_fetch_dec_mod(const_cast<T *>(acc),
val,
desul_order{},
desul_scope{});
}

RAJA_SUPPRESS_HD_WARN
template <typename AtomicPolicy, typename T>
RAJA_HOST_DEVICE
RAJA_INLINE T atomicAnd(AtomicPolicy, T volatile *acc, T value)
RAJA_HOST_DEVICE RAJA_INLINE T atomicAnd(AtomicPolicy, T volatile *acc, T value)
{
return desul::atomic_fetch_and(const_cast<T*>(acc),
using desul_order =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_order;
using desul_scope =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_scope;
return desul::atomic_fetch_and(const_cast<T *>(acc),
value,
raja_default_desul_order{},
raja_default_desul_scope{});
desul_order{},
desul_scope{});
}

RAJA_SUPPRESS_HD_WARN
template <typename AtomicPolicy, typename T>
RAJA_HOST_DEVICE
RAJA_INLINE T atomicOr(AtomicPolicy, T volatile *acc, T value)
RAJA_HOST_DEVICE RAJA_INLINE T atomicOr(AtomicPolicy, T volatile *acc, T value)
{
return desul::atomic_fetch_or(const_cast<T*>(acc),
using desul_order =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_order;
using desul_scope =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_scope;
return desul::atomic_fetch_or(const_cast<T *>(acc),
value,
raja_default_desul_order{},
raja_default_desul_scope{});
desul_order{},
desul_scope{});
}

RAJA_SUPPRESS_HD_WARN
template <typename AtomicPolicy, typename T>
RAJA_HOST_DEVICE
RAJA_INLINE T atomicXor(AtomicPolicy, T volatile *acc, T value)
RAJA_HOST_DEVICE RAJA_INLINE T atomicXor(AtomicPolicy, T volatile *acc, T value)
{
return desul::atomic_fetch_xor(const_cast<T*>(acc),
using desul_order =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_order;
using desul_scope =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_scope;
return desul::atomic_fetch_xor(const_cast<T *>(acc),
value,
raja_default_desul_order{},
raja_default_desul_scope{});
desul_order{},
desul_scope{});
}

RAJA_SUPPRESS_HD_WARN
template <typename AtomicPolicy, typename T>
RAJA_HOST_DEVICE
RAJA_INLINE T atomicExchange(AtomicPolicy, T volatile *acc, T value)
RAJA_HOST_DEVICE RAJA_INLINE T atomicExchange(AtomicPolicy,
T volatile *acc,
T value)
{
return desul::atomic_exchange(const_cast<T*>(acc),
using desul_order =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_order;
using desul_scope =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_scope;
return desul::atomic_exchange(const_cast<T *>(acc),
value,
raja_default_desul_order{},
raja_default_desul_scope{});
desul_order{},
desul_scope{});
}

RAJA_SUPPRESS_HD_WARN
template <typename AtomicPolicy, typename T>
RAJA_HOST_DEVICE
RAJA_INLINE T atomicCAS(AtomicPolicy, T volatile *acc, T compare, T value)
RAJA_HOST_DEVICE RAJA_INLINE T
atomicCAS(AtomicPolicy, T volatile *acc, T compare, T value)
{
return desul::atomic_compare_exchange(const_cast<T*>(acc),
compare,
value,
raja_default_desul_order{},
raja_default_desul_scope{});
using desul_order =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_order;
using desul_scope =
typename detail::DesulAtomicPolicy<AtomicPolicy>::memory_scope;
return desul::atomic_compare_exchange(
const_cast<T *>(acc), compare, value, desul_order{}, desul_scope{});
}

} // namespace RAJA

#endif // RAJA_ENABLE_DESUL_ATOMICS
#endif // guard
#endif // guard
62 changes: 62 additions & 0 deletions include/RAJA/policy/desul/policy.hpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,62 @@
//~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~//
// Copyright (c) 2016-24, Lawrence Livermore National Security, LLC
// and RAJA project contributors. See the RAJA/LICENSE file for details.
//
// SPDX-License-Identifier: (BSD-3-Clause)
//~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~//

#ifndef RAJA_policy_desul_HPP
#define RAJA_policy_desul_HPP

#include "RAJA/config.hpp"

#if defined(RAJA_ENABLE_DESUL_ATOMICS)

#include "desul/atomics.hpp"

namespace RAJA
{

// Policy to perform an atomic operation with a given memory ordering.
template <typename OrderingPolicy, typename ScopePolicy>
struct detail_atomic_t {
};

using atomic_seq_cst =
detail_atomic_t<desul::MemoryOrderSeqCst, desul::MemoryScopeDevice>;
using atomic_acq_rel =
detail_atomic_t<desul::MemoryOrderAcqRel, desul::MemoryScopeDevice>;
using atomic_acquire =
detail_atomic_t<desul::MemoryOrderAcquire, desul::MemoryScopeDevice>;
using atomic_release =
detail_atomic_t<desul::MemoryOrderRelease, desul::MemoryScopeDevice>;
using atomic_relaxed =
detail_atomic_t<desul::MemoryOrderRelaxed, desul::MemoryScopeDevice>;

using atomic_seq_cst_block =
detail_atomic_t<desul::MemoryOrderSeqCst, desul::MemoryScopeCore>;
using atomic_acq_rel_block =
detail_atomic_t<desul::MemoryOrderAcqRel, desul::MemoryScopeCore>;
using atomic_acquire_block =
detail_atomic_t<desul::MemoryOrderAcquire, desul::MemoryScopeCore>;
using atomic_release_block =
detail_atomic_t<desul::MemoryOrderRelease, desul::MemoryScopeCore>;
using atomic_relaxed_block =
detail_atomic_t<desul::MemoryOrderRelaxed, desul::MemoryScopeCore>;

using atomic_seq_cst_sys =
detail_atomic_t<desul::MemoryOrderSeqCst, desul::MemoryScopeSystem>;
using atomic_acq_rel_sys =
detail_atomic_t<desul::MemoryOrderAcqRel, desul::MemoryScopeSystem>;
using atomic_acquire_sys =
detail_atomic_t<desul::MemoryOrderAcquire, desul::MemoryScopeSystem>;
using atomic_release_sys =
detail_atomic_t<desul::MemoryOrderRelease, desul::MemoryScopeSystem>;
using atomic_relaxed_sys =
detail_atomic_t<desul::MemoryOrderRelaxed, desul::MemoryScopeSystem>;

} // namespace RAJA

#endif // RAJA_ENABLE_DESUL_ATOMICS

#endif
35 changes: 35 additions & 0 deletions test/functional/forall/atomic-basic/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -28,3 +28,38 @@ foreach( ATOMIC_BACKEND ${FORALL_ATOMIC_BACKENDS} )
target_include_directories(test-forall-atomic-basic-unsigned-${ATOMIC_BACKEND}.exe
PRIVATE ${CMAKE_CURRENT_SOURCE_DIR}/tests)
endforeach()


set(ENABLE_LITMUS_TESTS OFF)

if(RAJA_ENABLE_DESUL_ATOMICS AND RAJA_ENABLE_CUDA)
set(LITMUS_BACKEND "Cuda")
set(ENABLE_LITMUS_TESTS ON)
endif()

if(RAJA_ENABLE_DESUL_ATOMICS AND RAJA_ENABLE_HIP)
set(LITMUS_BACKEND "Hip")
set(ENABLE_LITMUS_TESTS ON)
endif()

set(FORALL_LITMUS_TESTS
mp # Message Passing
sb # Store Buffer
lb # Load Buffer
store # Store
read # Read
write2x2 # 2+2 write
)

if (ENABLE_LITMUS_TESTS)
foreach ( LITMUS_TEST ${FORALL_LITMUS_TESTS} )
raja_add_test( NAME test-forall-atomic-litmus-${LITMUS_BACKEND}-${LITMUS_TEST}
SOURCES test-forall-atomic-litmus-${LITMUS_TEST}.cpp)
target_include_directories(test-forall-atomic-litmus-${LITMUS_BACKEND}-${LITMUS_TEST}.exe
PRIVATE ${CMAKE_CURRENT_SOURCE_DIR}/tests)
endforeach()
endif()

unset(FORALL_LITMUS_TESTS)
unset(ENABLE_LITMUS_TESTS)
unset(LITMUS_BACKEND)