Sacado Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Sacado_DynamicArrayTraits.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Sacado Package
4 //
5 // Copyright 2006 NTESS and the Sacado contributors.
6 // SPDX-License-Identifier: LGPL-2.1-or-later
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef SACADO_DYNAMICARRAYTRAITS_HPP
11 #define SACADO_DYNAMICARRAYTRAITS_HPP
12 
13 #include <new>
14 #include <cstring>
15 #include <stdint.h>
16 
17 #include "Sacado_Traits.hpp"
18 #if defined(HAVE_SACADO_KOKKOS)
19 #include "Kokkos_Core.hpp"
20 #endif
21 
22 namespace Sacado {
23 
24  template <typename ExecSpace>
25  void createGlobalMemoryPool(const ExecSpace& space
26  , const size_t min_total_alloc_size
27  , const uint32_t min_block_alloc_size
28  , const uint32_t max_block_alloc_size
29  , const uint32_t min_superblock_size
30  ) {}
31 
32  template <typename ExecSpace>
33  void destroyGlobalMemoryPool(const ExecSpace& space) {}
34 
35 #if 0 && defined(HAVE_SACADO_KOKKOS) && defined(KOKKOS_ENABLE_OPENMP)
36  namespace Impl {
37  extern const Kokkos::MemoryPool<Kokkos::OpenMP>* global_sacado_openmp_memory_pool;
38  }
39 
40  inline void
41  createGlobalMemoryPool(const ExecSpace& space
42  , const size_t min_total_alloc_size
43  , const uint32_t min_block_alloc_size
44  , const uint32_t max_block_alloc_size
45  , const uint32_t min_superblock_size
46  )
47  {
48  typedef Kokkos::MemoryPool<Kokkos::OpenMP> pool_t;
49  Impl::global_sacado_openmp_memory_pool =
50  new pool_t(typename Kokkos::OpenMP::memory_space(),
51  min_total_alloc_size,
52  min_block_alloc_size,
53  max_block_alloc_size,
54  min_superblock_size);
55  }
56 
57  inline void destroyGlobalMemoryPool(const Kokkos::OpenMP& space)
58  {
59  delete Impl::global_sacado_openmp_memory_pool;
60  }
61 #endif
62 
63 #if defined(HAVE_SACADO_KOKKOS) && defined(SACADO_KOKKOS_USE_MEMORY_POOL) && !defined(SACADO_DISABLE_CUDA_IN_KOKKOS) && defined(KOKKOS_ENABLE_CUDA) && defined(__CUDACC__)
64 
65  namespace Impl {
66 
67  extern const Kokkos::MemoryPool<Kokkos::Cuda>* global_sacado_cuda_memory_pool_host;
68  extern const Kokkos::MemoryPool<Kokkos::Cuda>* global_sacado_cuda_memory_pool_device;
69 #ifdef KOKKOS_ENABLE_CUDA_RELOCATABLE_DEVICE_CODE
70  extern __device__ const Kokkos::MemoryPool<Kokkos::Cuda>* global_sacado_cuda_memory_pool_on_device;
71 #else
72  __device__ const Kokkos::MemoryPool<Kokkos::Cuda>* global_sacado_cuda_memory_pool_on_device = 0;
73 #endif
74 
75  struct SetMemoryPoolPtr {
76  Kokkos::MemoryPool<Kokkos::Cuda>* pool_device;
77  __device__ inline void operator()(int) const {
78  global_sacado_cuda_memory_pool_on_device = pool_device;
79  };
80  };
81 
82  }
83 
84  // For some reason we get memory errors if these functions are defined in
85  // Sacado_DynamicArrayTraits.cpp
86  inline void
87  createGlobalMemoryPool(const Kokkos::Cuda& space
88  , const size_t min_total_alloc_size
89  , const uint32_t min_block_alloc_size
90  , const uint32_t max_block_alloc_size
91  , const uint32_t min_superblock_size
92  )
93  {
94  typedef Kokkos::MemoryPool<Kokkos::Cuda> pool_t;
95  pool_t* pool =
96  new pool_t(typename Kokkos::Cuda::memory_space(),
97  min_total_alloc_size,
98  min_block_alloc_size,
99  max_block_alloc_size,
100  min_superblock_size);
101  Impl::SetMemoryPoolPtr f;
102  KOKKOS_IMPL_CUDA_SAFE_CALL( cudaMalloc( &f.pool_device, sizeof(pool_t) ) );
103  KOKKOS_IMPL_CUDA_SAFE_CALL( cudaMemcpy( f.pool_device, pool,
104  sizeof(pool_t),
105  cudaMemcpyHostToDevice ) );
106  Kokkos::parallel_for(Kokkos::RangePolicy<Kokkos::Cuda>(0,1),f);
107  Impl::global_sacado_cuda_memory_pool_host = pool;
108  Impl::global_sacado_cuda_memory_pool_device = f.pool_device;
109  }
110 
111  inline void destroyGlobalMemoryPool(const Kokkos::Cuda& space)
112  {
113  KOKKOS_IMPL_CUDA_SAFE_CALL( cudaFree( (void*) Impl::global_sacado_cuda_memory_pool_device ) );
114  delete Impl::global_sacado_cuda_memory_pool_host;
115  }
116 
117 #endif
118 
119 #if !defined(SACADO_DISABLE_CUDA_IN_KOKKOS) && defined(KOKKOS_ENABLE_CUDA) && defined(__CUDACC__)
120 
121  namespace Impl {
122 
123  // Compute warp lane/thread index
124  __device__ inline int warpLane(const int warp_size = 32) {
125  return ( threadIdx.x + (threadIdx.y + threadIdx.z*blockDim.y)*blockDim.x ) % warp_size;
126  }
127 
128  // Reduce y across the warp and broadcast to all lanes
129  template <typename T>
130  __device__ inline T warpReduce(T y, const int warp_size = 32) {
131  for (int i=1; i<warp_size; i*=2) {
132  y += Kokkos::shfl_down(y, i, warp_size);
133  }
134  y = Kokkos::shfl(y, 0, warp_size);
135  return y;
136  }
137 
138  // Non-inclusive plus-scan up the warp, replacing the first entry with 0
139  template <typename T>
140  __device__ inline int warpScan(T y, const int warp_size = 32) {
141  const int lane = warpLane();
142  y = Kokkos::shfl_up(y, 1, warp_size);
143  if (lane == 0)
144  y = T(0);
145  for (int i=1; i<warp_size; i*=2) {
146  T t = Kokkos::shfl_up(y, i, warp_size);
147  if (lane > i)
148  y += t;
149  }
150  return y;
151  }
152 
153  template <typename T>
154  __device__ inline T warpBcast(T y, int id, const int warp_size = 32) {
155  return Kokkos::shfl(y, id, warp_size);
156  }
157 
158  }
159 
160 #endif
161 
162  namespace Impl {
163 
164  template <typename T>
166  static T* ds_alloc(const int sz) {
167 #if defined( CUDA_VERSION ) && ( 6000 <= CUDA_VERSION ) && defined(KOKKOS_ENABLE_CUDA_UVM) && !defined( __CUDA_ARCH__ )
168  T* m = 0;
169  if (sz > 0)
170  KOKKOS_IMPL_CUDA_SAFE_CALL( cudaMallocManaged( (void**) &m, sz*sizeof(T), cudaMemAttachGlobal ) );
171 #elif defined(HAVE_SACADO_KOKKOS) && defined(SACADO_KOKKOS_USE_MEMORY_POOL) && !defined(SACADO_DISABLE_CUDA_IN_KOKKOS) && defined(__CUDA_ARCH__)
172  // This code assumes all threads enter ds_alloc, even those with sz == 0
173  T* m = 0;
174  const int total_sz = warpReduce(sz);
175  const int lane = warpLane();
176  if (total_sz > 0 && lane == 0) {
177  m = static_cast<T*>(global_sacado_cuda_memory_pool_on_device->allocate(total_sz*sizeof(T)));
178  if (m == 0)
179  Kokkos::abort("Allocation failed. Kokkos memory pool is out of memory");
180  }
181  m = warpBcast(m,0);
182  m += warpScan(sz);
183 #elif 0 && defined(HAVE_SACADO_KOKKOS) && defined(SACADO_KOKKOS_USE_MEMORY_POOL) && defined(KOKKOS_ENABLE_OPENMP)
184  T* m = 0;
185  if (sz > 0) {
186  if (global_sacado_openmp_memory_pool != 0) {
187  m = static_cast<T*>(global_sacado_openmp_memory_pool->allocate(sz*sizeof(T)));
188  if (m == 0)
189  Kokkos::abort("Allocation failed. Kokkos memory pool is out of memory");
190  }
191  else
192  m = static_cast<T* >(operator new(sz*sizeof(T)));
193  }
194 #else
195  T* m = 0;
196  if (sz > 0) {
197  m = static_cast<T* >(operator new(sz*sizeof(T)));
198 #if defined(HAVE_SACADO_KOKKOS)
199  if (m == 0)
200  Kokkos::abort("Allocation failed.");
201 #endif
202  }
203 #endif
204  return m;
205  }
206 
207  template <typename T>
209  static void ds_free(T* m, int sz) {
210 #if defined( CUDA_VERSION ) && ( 6000 <= CUDA_VERSION ) && defined(KOKKOS_ENABLE_CUDA_UVM) && !defined( __CUDA_ARCH__ )
211  if (sz > 0)
212  KOKKOS_IMPL_CUDA_SAFE_CALL( cudaFree(m) );
213 #elif defined(HAVE_SACADO_KOKKOS) && defined(SACADO_KOKKOS_USE_MEMORY_POOL) && !defined(SACADO_DISABLE_CUDA_IN_KOKKOS) && defined(__CUDA_ARCH__)
214  const int total_sz = warpReduce(sz);
215  const int lane = warpLane();
216  if (total_sz > 0 && lane == 0) {
217  global_sacado_cuda_memory_pool_on_device->deallocate((void*) m, total_sz*sizeof(T));
218  }
219 #elif 0 && defined(HAVE_SACADO_KOKKOS) && defined(SACADO_KOKKOS_USE_MEMORY_POOL) && defined(KOKKOS_ENABLE_OPENMP)
220  if (sz > 0) {
221  if (global_sacado_openmp_memory_pool != 0)
222  global_sacado_openmp_memory_pool->deallocate((void*) m, sz*sizeof(T));
223  else
224  operator delete((void*) m);
225  }
226 #else
227  if (sz > 0)
228  operator delete((void*) m);
229 #endif
230  }
231 
232  }
233 
238  struct ds_array {
239 
242  static T* get(int sz) {
243  T* m = Impl::ds_alloc<T>(sz);
244  T* p = m;
245  for (int i=0; i<sz; ++i)
246  new (p++) T();
247  return m;
248  }
249 
252  static T* get_and_fill(int sz) {
253  T* m = Impl::ds_alloc<T>(sz);
254  T* p = m;
255  for (int i=0; i<sz; ++i)
256  new (p++) T(0.0);
257  return m;
258  }
259 
265  static T* get_and_fill(const T* src, int sz) {
266  T* m = Impl::ds_alloc<T>(sz);
267  T* p = m;
268  for (int i=0; i<sz; ++i)
269  new (p++) T(*(src++));
270  return m;
271  }
272 
278  static T* strided_get_and_fill(const T* src, int stride, int sz) {
279  T* m = Impl::ds_alloc<T>(sz);
280  T* p = m;
281  for (int i=0; i<sz; ++i) {
282  new (p++) T(*(src));
283  src += stride;
284  }
285  return m;
286  }
287 
290  static void copy(const T* src, T* dest, int sz) {
291  for (int i=0; i<sz; ++i)
292  *(dest++) = *(src++);
293  }
294 
297  static void strided_copy(const T* src, int src_stride,
298  T* dest, int dest_stride, int sz) {
299  for (int i=0; i<sz; ++i) {
300  *(dest) = *(src);
301  dest += dest_stride;
302  src += src_stride;
303  }
304  }
305 
308  static void zero(T* dest, int sz) {
309  for (int i=0; i<sz; ++i)
310  *(dest++) = T(0.);
311  }
312 
315  static void strided_zero(T* dest, int stride, int sz) {
316  for (int i=0; i<sz; ++i) {
317  *(dest) = T(0.);
318  dest += stride;
319  }
320  }
321 
324  static void destroy_and_release(T* m, int sz) {
325  T* e = m+sz;
326  for (T* b = m; b!=e; b++)
327  b->~T();
328  Impl::ds_free(m, sz);
329  }
330  };
331 
332 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL_DFAD) && !defined(SACADO_DISABLE_CUDA_IN_KOKKOS) && defined(__CUDA_ARCH__)
333 
334  namespace Impl {
335 
336  template <typename T>
338  static T* ds_strided_alloc(const int sz) {
339  T* m = 0;
340  // Only do strided memory allocations when we are doing hierarchical
341  // parallelism with a vector dimension of 32. The limitation on the
342  // memory pool allowing only a single thread in a warp to allocate
343  // makes it too difficult to do otherwise.
344  if (blockDim.x == 32) {
345  //const int lane = warpLane();
346  const int lane = threadIdx.x;
347  if (sz > 0 && lane == 0) {
348 #if defined(HAVE_SACADO_KOKKOS) && defined(SACADO_KOKKOS_USE_MEMORY_POOL)
349  m = static_cast<T*>(global_sacado_cuda_memory_pool_on_device->allocate(sz*sizeof(T)));
350  if (m == 0)
351  Kokkos::abort("Allocation failed. Kokkos memory pool is out of memory");
352 #else
353  m = static_cast<T* >(operator new(sz*sizeof(T)));
354 #if defined(HAVE_SACADO_KOKKOS)
355  if (m == 0)
356  Kokkos::abort("Allocation failed.");
357 #endif
358 #endif
359  }
360  m = warpBcast(m,0,blockDim.x);
361  }
362  else {
363  if (sz > 0) {
364  m = static_cast<T* >(operator new(sz*sizeof(T)));
365 #if defined(HAVE_SACADO_KOKKOS)
366  if (m == 0)
367  Kokkos::abort("Allocation failed.");
368 #endif
369  }
370  }
371 
372  return m;
373  }
374 
375  template <typename T>
377  static void ds_strided_free(T* m, int sz) {
378  if (blockDim.x == 32) {
379  // const int lane = warpLane();
380  const int lane = threadIdx.x;
381  if (sz > 0 && lane == 0) {
382 #if defined(HAVE_SACADO_KOKKOS) && defined(SACADO_KOKKOS_USE_MEMORY_POOL)
383  global_sacado_cuda_memory_pool_on_device->deallocate((void*) m, sz*sizeof(T));
384 #else
385  operator delete((void*) m);
386 #endif
387  }
388  }
389  else {
390  if (sz > 0)
391  operator delete((void*) m);
392  }
393 
394  }
395 
396  }
397 
402  template <typename T>
403  struct ds_array<T,true> {
404 
407  static T* get(int sz) {
408  T* m = Impl::ds_strided_alloc<T>(sz);
409  return m;
410  }
411 
414  static T* get_and_fill(int sz) {
415  T* m = Impl::ds_strided_alloc<T>(sz);
416  for (int i=threadIdx.x; i<sz; i+=blockDim.x)
417  m[i] = 0.0;
418  return m;
419  }
420 
426  static T* get_and_fill(const T* src, int sz) {
427  T* m = Impl::ds_strided_alloc<T>(sz);
428  for (int i=threadIdx.x; i<sz; i+=blockDim.x)
429  m[i] = src[i];
430  return m;
431  }
432 
438  static T* strided_get_and_fill(const T* src, int stride, int sz) {
439  T* m = Impl::ds_strided_alloc<T>(sz);
440  for (int i=threadIdx.x; i<sz; i+=blockDim.x)
441  m[i] = src[i*stride];
442  return m;
443  }
444 
447  static void copy(const T* src, T* dest, int sz) {
448  for (int i=threadIdx.x; i<sz; i+=blockDim.x)
449  dest[i] = src[i];
450  }
451 
454  static void strided_copy(const T* src, int src_stride,
455  T* dest, int dest_stride, int sz) {
456  for (int i=threadIdx.x; i<sz; i+=blockDim.x) {
457  dest[i*dest_stride] = src[i*src_stride];
458  }
459  }
460 
463  static void zero(T* dest, int sz) {
464  for (int i=threadIdx.x; i<sz; i+=blockDim.x)
465  dest[i] = T(0.);
466  }
467 
470  static void strided_zero(T* dest, int stride, int sz) {
471  for (int i=threadIdx.x; i<sz; i+=blockDim.x) {
472  dest[i*stride] = T(0.);
473  }
474  }
475 
478  static void destroy_and_release(T* m, int sz) {
479  Impl::ds_strided_free(m, sz);
480  }
481  };
482 
483 #elif defined(SACADO_VIEW_CUDA_HIERARCHICAL_DFAD_STRIDED) && !defined(SACADO_DISABLE_CUDA_IN_KOKKOS) && defined(__CUDA_ARCH__)
484 
485  namespace Impl {
486 
487  template <typename T>
489  static T* ds_strided_alloc(const int sz) {
490  T* m = 0;
491  // Only do strided memory allocations when we are doing hierarchical
492  // parallelism with a vector dimension of 32. The limitation on the
493  // memory pool allowing only a single thread in a warp to allocate
494  // makes it too difficult to do otherwise.
495  if (blockDim.x == 32) {
496  // const int total_sz = warpReduce(sz);
497  // const int lane = warpLane();
498  const int total_sz = warpReduce(sz, blockDim.x);
499  const int lane = threadIdx.x;
500  if (total_sz > 0 && lane == 0) {
501 #if defined(HAVE_SACADO_KOKKOS) && defined(SACADO_KOKKOS_USE_MEMORY_POOL)
502  m = static_cast<T*>(global_sacado_cuda_memory_pool_on_device->allocate(total_sz*sizeof(T)));
503  if (m == 0)
504  Kokkos::abort("Allocation failed. Kokkos memory pool is out of memory");
505 #else
506  m = static_cast<T* >(operator new(total_sz*sizeof(T)));
507 #if defined(HAVE_SACADO_KOKKOS)
508  if (m == 0)
509  Kokkos::abort("Allocation failed.");
510 #endif
511 #endif
512  }
513  m = warpBcast(m,0,blockDim.x);
514  m += lane;
515  }
516  else {
517  if (sz > 0) {
518  m = static_cast<T* >(operator new(sz*sizeof(T)));
519 #if defined(HAVE_SACADO_KOKKOS)
520  if (m == 0)
521  Kokkos::abort("Allocation failed.");
522 #endif
523  }
524  }
525 
526  return m;
527  }
528 
529  template <typename T>
531  static void ds_strided_free(T* m, int sz) {
532  if (blockDim.x == 32) {
533  // const int total_sz = warpReduce(sz);
534  // const int lane = warpLane();
535  const int total_sz = warpReduce(sz, blockDim.x);
536  const int lane = threadIdx.x;
537  if (total_sz > 0 && lane == 0) {
538 #if defined(HAVE_SACADO_KOKKOS) && defined(SACADO_KOKKOS_USE_MEMORY_POOL)
539  global_sacado_cuda_memory_pool_on_device->deallocate((void*) m, total_sz*sizeof(T));
540 #else
541  operator delete((void*) m);
542 #endif
543  }
544  }
545  else {
546  if (sz > 0)
547  operator delete((void*) m);
548  }
549  }
550  }
551 
556  template <typename T>
557  struct ds_array<T,true> {
558 
561  static T* get(int sz) {
562  T* m = Impl::ds_strided_alloc<T>(sz);
563  return m;
564  }
565 
568  static T* get_and_fill(int sz) {
569  T* m = Impl::ds_strided_alloc<T>(sz);
570  for (int i=0; i<sz; ++i)
571  m[i*blockDim.x] = 0.0;
572  return m;
573  }
574 
580  static T* get_and_fill(const T* src, int sz) {
581  T* m = Impl::ds_strided_alloc<T>(sz);
582  for (int i=0; i<sz; ++i)
583  m[i*blockDim.x] = src[i*blockDim.x];
584  return m;
585  }
586 
592  static T* strided_get_and_fill(const T* src, int stride, int sz) {
593  T* m = Impl::ds_strided_alloc<T>(sz);
594  for (int i=0; i<sz; ++i)
595  m[i*blockDim.x] = src[i*stride];
596  return m;
597  }
598 
601  static void copy(const T* src, T* dest, int sz) {
602  for (int i=0; i<sz; ++i)
603  dest[i*blockDim.x] = src[i*blockDim.x];
604  }
605 
608  static void strided_copy(const T* src, int src_stride,
609  T* dest, int dest_stride, int sz) {
610  for (int i=0; i<sz; ++i) {
611  *(dest) = *(src);
612  dest += dest_stride;
613  src += src_stride;
614  }
615  }
616 
619  static void zero(T* dest, int sz) {
620  for (int i=0; i<sz; ++i)
621  dest[i*blockDim.x] = T(0.);
622  }
623 
626  static void strided_zero(T* dest, int stride, int sz) {
627  for (int i=0; i<sz; ++i) {
628  *(dest) = T(0.);
629  dest += stride;
630  }
631  }
632 
635  static void destroy_and_release(T* m, int sz) {
636  Impl::ds_strided_free(m, sz);
637  }
638  };
639 
640 #else
641 
646  template <typename T>
647  struct ds_array<T,true> {
648 
651  static T* get(int sz) {
652  T* m = Impl::ds_alloc<T>(sz);
653  return m;
654  }
655 
658  static T* get_and_fill(int sz) {
659  T* m = Impl::ds_alloc<T>(sz);
660 #if defined(__CUDACC__ ) || defined(__HIPCC__ )
661  for (int i=0; i<sz; ++i)
662  m[i] = 0.0;
663 #else
664  if (sz > 0)
665  std::memset(m,0,sz*sizeof(T));
666 #endif
667  return m;
668  }
669 
675  static T* get_and_fill(const T* src, int sz) {
676  T* m = Impl::ds_alloc<T>(sz);
677  for (int i=0; i<sz; ++i)
678  m[i] = src[i];
679  return m;
680  }
681 
687  static T* strided_get_and_fill(const T* src, int stride, int sz) {
688  T* m = Impl::ds_alloc<T>(sz);
689  for (int i=0; i<sz; ++i)
690  m[i] = src[i*stride];
691  return m;
692  }
693 
696  static void copy(const T* src, T* dest, int sz) {
697  if (sz > 0 && dest != NULL && src != NULL)
698 #if defined( __CUDACC__) || defined(__HIPCC__ )
699  for (int i=0; i<sz; ++i)
700  dest[i] = src[i];
701 #else
702  std::memcpy(dest,src,sz*sizeof(T));
703 #endif
704  }
705 
708  static void strided_copy(const T* src, int src_stride,
709  T* dest, int dest_stride, int sz) {
710  for (int i=0; i<sz; ++i) {
711  *(dest) = *(src);
712  dest += dest_stride;
713  src += src_stride;
714  }
715  }
716 
719  static void zero(T* dest, int sz) {
720  if (sz > 0 && dest != NULL)
721 #if defined(__CUDACC__ ) || defined(__HIPCC__ )
722  for (int i=0; i<sz; ++i)
723  dest[i] = T(0.);
724 #else
725  std::memset(dest,0,sz*sizeof(T));
726 #endif
727  }
728 
731  static void strided_zero(T* dest, int stride, int sz) {
732  for (int i=0; i<sz; ++i) {
733  *(dest) = T(0.);
734  dest += stride;
735  }
736  }
737 
740  static void destroy_and_release(T* m, int sz) {
741  Impl::ds_free(m, sz);
742  }
743  };
744 
745 #endif
746 
747 } // namespace Sacado
748 
749 #endif // SACADO_DYNAMICARRAY_HPP
const char * p
static SACADO_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.
void f()
static SACADO_INLINE_FUNCTION void strided_zero(T *dest, int stride, int sz)
Zero out array dest of length sz.
static SACADO_INLINE_FUNCTION void copy(const T *src, T *dest, int sz)
Copy array from src to dest of length sz.
void createGlobalMemoryPool(const ExecSpace &space, const size_t min_total_alloc_size, const uint32_t min_block_alloc_size, const uint32_t max_block_alloc_size, const uint32_t min_superblock_size)
static SACADO_INLINE_FUNCTION void zero(T *dest, int sz)
Zero out array dest of length sz.
static SACADO_INLINE_FUNCTION void destroy_and_release(T *m, int sz)
Destroy array elements and release memory.
static SACADO_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 SACADO_INLINE_FUNCTION void zero(T *dest, int sz)
Zero out array dest of length sz.
static SACADO_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.
expr true
static SACADO_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.
#define T
Definition: Sacado_rad.hpp:553
static SACADO_INLINE_FUNCTION void ds_free(T *m, int sz)
static SACADO_INLINE_FUNCTION void destroy_and_release(T *m, int sz)
Destroy array elements and release memory.
int value
static SACADO_INLINE_FUNCTION T * ds_alloc(const int sz)
static SACADO_INLINE_FUNCTION void copy(const T *src, T *dest, int sz)
Copy array from src to dest of length sz.
static SACADO_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.
void destroyGlobalMemoryPool(const ExecSpace &space)
#define SACADO_INLINE_FUNCTION
static SACADO_INLINE_FUNCTION T * get_and_fill(int sz)
Get memory for new array of length sz and fill with zeros.
static SACADO_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 SACADO_INLINE_FUNCTION void strided_zero(T *dest, int stride, int sz)
Zero out array dest of length sz.
static SACADO_INLINE_FUNCTION T * get_and_fill(int sz)
Get memory for new array of length sz and fill with zeros.
const double y