30 #ifndef SACADO_DYNAMICARRAYTRAITS_HPP 31 #define SACADO_DYNAMICARRAYTRAITS_HPP 37 #if defined(HAVE_SACADO_KOKKOSCORE) 38 #include "Kokkos_Core.hpp" 46 template <typename T, bool isScalar = IsScalarType<T>::value>
51 #if defined(__CUDACC__) && defined( CUDA_VERSION ) && ( 6000 <= CUDA_VERSION ) && defined(KOKKOS_USE_CUDA_UVM) && !defined( __CUDA_ARCH__ ) 53 CUDA_SAFE_CALL( cudaMallocManaged( (
void**) &m, sz*
sizeof(
T), cudaMemAttachGlobal ) );
55 T* m =
static_cast<T*
>(
operator new(sz*
sizeof(
T)));
63 #if defined(__CUDACC__) && defined( CUDA_VERSION ) && ( 6000 <= CUDA_VERSION ) && defined(KOKKOS_USE_CUDA_UVM) && !defined( __CUDA_ARCH__ ) 64 CUDA_SAFE_CALL( cudaFree(m) );
66 operator delete((
void*) m);
73 static T*
get(
int sz) {
77 for (
int i=0; i<sz; ++i)
90 for (
int i=0; i<sz; ++i)
106 for (
int i=0; i<sz; ++i)
107 new (p++)
T(*(src++));
122 for (
int i=0; i<sz; ++i) {
133 static void copy(
const T* src,
T* dest,
int sz) {
134 for (
int i=0; i<sz; ++i)
135 *(dest++) = *(src++);
141 T* dest,
int dest_stride,
int sz) {
142 for (
int i=0; i<sz; ++i) {
152 for (
int i=0; i<sz; ++i)
159 for (
int i=0; i<sz; ++i) {
169 for (
T* b = m; b!=e; b++)
179 template <
typename T>
184 #if defined( CUDA_VERSION ) && ( 6000 <= CUDA_VERSION ) && defined(KOKKOS_USE_CUDA_UVM) && !defined( __CUDA_ARCH__ ) 186 CUDA_SAFE_CALL( cudaMallocManaged( (
void**) &m, sz*
sizeof(
T), cudaMemAttachGlobal ) );
188 T* m =
static_cast<T*
>(
operator new(sz*
sizeof(
T)));
196 #if defined( CUDA_VERSION ) && ( 6000 <= CUDA_VERSION ) && defined(KOKKOS_USE_CUDA_UVM) && !defined( __CUDA_ARCH__ ) 197 CUDA_SAFE_CALL( cudaFree(m) );
199 operator delete((
void*) m);
206 static T*
get(
int sz) {
220 for (
int i=0; i<sz; ++i)
223 std::memset(m,0,sz*
sizeof(
T));
238 for (
int i=0; i<sz; ++i)
253 for (
int i=0; i<sz; ++i)
254 m[i] = src[i*stride];
262 static void copy(
const T* src,
T* dest,
int sz) {
265 for (
int i=0; i<sz; ++i)
268 std::memcpy(dest,src,sz*
sizeof(
T));
275 T* dest,
int dest_stride,
int sz) {
276 for (
int i=0; i<sz; ++i) {
288 for (
int i=0; i<sz; ++i)
291 std::memset(dest,0,sz*
sizeof(
T));
298 for (
int i=0; i<sz; ++i) {
313 #endif // SACADO_DYNAMICARRAY_HPP static KOKKOS_INLINE_FUNCTION T * strided_get_and_fill(const T *src, int stride, int sz)
Get memory for new array of length sz and fill with entries from src.
static KOKKOS_INLINE_FUNCTION void destroy_and_release(T *m, int sz)
Destroy array elements and release memory.
static KOKKOS_INLINE_FUNCTION void copy(const T *src, T *dest, int sz)
Copy array from src to dest of length sz.
static KOKKOS_INLINE_FUNCTION T * my_alloc(const int sz)
static KOKKOS_INLINE_FUNCTION void my_free(T *m, int sz)
static KOKKOS_INLINE_FUNCTION void strided_zero(T *dest, int stride, int sz)
Zero out array dest of length sz.
static KOKKOS_INLINE_FUNCTION void destroy_and_release(T *m, int sz)
Destroy array elements and release memory.
static KOKKOS_INLINE_FUNCTION void strided_copy(const T *src, int src_stride, T *dest, int dest_stride, int sz)
Copy array from src to dest of length sz.
static KOKKOS_INLINE_FUNCTION T * my_alloc(const int sz)
#define KOKKOS_INLINE_FUNCTION
static KOKKOS_INLINE_FUNCTION void my_free(T *m, int sz)
static KOKKOS_INLINE_FUNCTION T * get_and_fill(int sz)
Get memory for new array of length sz and fill with zeros.
static KOKKOS_INLINE_FUNCTION void zero(T *dest, int sz)
Zero out array dest of length sz.
static KOKKOS_INLINE_FUNCTION void strided_copy(const T *src, int src_stride, T *dest, int dest_stride, int sz)
Copy array from src to dest of length sz.
static KOKKOS_INLINE_FUNCTION void copy(const T *src, T *dest, int sz)
Copy array from src to dest of length sz.
static KOKKOS_INLINE_FUNCTION T * strided_get_and_fill(const T *src, int stride, int sz)
Get memory for new array of length sz and fill with entries from src.
static KOKKOS_INLINE_FUNCTION void strided_zero(T *dest, int stride, int sz)
Zero out array dest of length sz.
static KOKKOS_INLINE_FUNCTION T * get_and_fill(const T *src, int sz)
Get memory for new array of length sz and fill with entries from src.
static KOKKOS_INLINE_FUNCTION T * get_and_fill(int sz)
Get memory for new array of length sz and fill with zeros.
static KOKKOS_INLINE_FUNCTION T * get_and_fill(const T *src, int sz)
Get memory for new array of length sz and fill with entries from src.
Dynamic array allocation class that works for any type.
static KOKKOS_INLINE_FUNCTION void zero(T *dest, int sz)
Zero out array dest of length sz.