10#ifndef MCGIDI_VECTOR_HPP
11#define MCGIDI_VECTOR_HPP
16#ifdef HAVE_OPENMP_TARGET
17 #ifdef USE_OPENMP_NO_GPU
28#define MCGIDI_SWAP(a,b,type) {type ttttttttt=a;a=b;b=ttttttttt;}
30#if defined(__CUDACC__) && !defined(__CUDA_ARCH__)
32#include <cuda_runtime.h>
33#include <cuda_runtime_api.h>
37#include <hip/hip_version.h>
38#include <hip/hip_runtime.h>
39#include <hip/hip_runtime_api.h>
40#include <hip/hip_common.h>
57 std::size_t _capacity;
69 if( s == 0 ){ _data =
nullptr;
return;}
70 switch ((
int)_mem_type){
72 _data =
new T [_capacity];
77#if defined(__CUDACC__) && !defined(__CUDA_ARCH__)
78 cudaMallocManaged(&ptr, _capacity*
sizeof(T), cudaMemAttachGlobal);
79#elif defined(__HIP__) && !defined(__HIP_DEVICE_COMPILE__)
80 hipMallocManaged(&ptr, _capacity*
sizeof(T), hipMemAttachGlobal);
82 _data =
new(ptr) T[_capacity];
86 _data =
new T [_capacity];
92 if( s == 0 ){ _data =
nullptr;
return;}
93 switch ( (
int) _mem_type){
95 _data =
new T [_capacity];
100#if defined(__CUDACC__) && !defined(__CUDA_ARCH__)
101 cudaMallocManaged(&ptr, _capacity*
sizeof(T), cudaMemAttachGlobal);
102#elif defined(__HIP__) && !defined(__HIP_DEVICE_COMPILE__)
103 hipMallocManaged(&ptr, _capacity*
sizeof(T), hipMemAttachGlobal);
105 _data =
new(ptr) T[_capacity];
109 _data =
new T [_capacity];
112 for (std::size_t ii = 0; ii < _capacity; ++ii)
117 : _data(0), _capacity(aa._capacity), _size(aa._size), _mem_type(aa._mem_type)
119 if( _capacity == 0 ){ _data =
nullptr;
return; }
121 switch ( (
int) _mem_type){
123 _data =
new T [_capacity];
128#if defined(__CUDACC__) && !defined(__CUDA_ARCH__)
129 cudaMallocManaged(&ptr, _capacity*
sizeof(T), cudaMemAttachGlobal);
130#elif defined(__HIP__) && !defined(__HIP_DEVICE_COMPILE__)
131 hipMallocManaged(&ptr, _capacity*
sizeof(T), hipMemAttachGlobal);
133 _data =
new(ptr) T[_capacity];
137 _data =
new T [_capacity];
141 for (std::size_t ii=0; ii<_size; ++ii)
142 _data[ii] = aa._data[ii];
148 if( _capacity == 0 ){ _data =
nullptr;
return;}
150 switch ( (
int) _mem_type){
152 _data =
new T [_capacity];
157#if defined(__CUDACC__) && !defined(__CUDA_ARCH__)
158 cudaMallocManaged(&ptr, _capacity*
sizeof(T), cudaMemAttachGlobal);
159#elif defined(__HIP__) && !defined(__HIP_DEVICE_COMPILE__)
160 hipMallocManaged(&ptr, _capacity*
sizeof(T), hipMemAttachGlobal);
162 _data =
new(ptr) T[_capacity];
166 _data =
new T [_capacity];
170 for (std::size_t ii=0; ii<_size; ++ii)
175 switch ( (
int) _mem_type){
180 for (std::size_t i=0; i < _size; ++i)
182#if defined(__CUDACC__) && !defined(__CUDA_ARCH__)
184#elif defined(__HIP__) && !defined(__HIP_DEVICE_COMPILE__)
206 MCGIDI_SWAP(_capacity, other._capacity, std::size_t);
236 assert( _size < _capacity );
237 _data[_size] = dataElem;
267 return _data[_size-1];
272 return _data[_size-1];
277 if (s == _capacity)
return;
278 assert( _capacity == 0 );
280 _mem_type = mem_flag;
281 if( s == 0 ){ _data =
nullptr;
return;}
282 switch ( (
int) _mem_type){
284 if (address ==
nullptr || *address ==
nullptr) _data =
new T [_capacity];
286 _data =
new(*address) T [_capacity];
287 *address +=
sizeof(T) * _capacity;
293#if defined(__CUDACC__) && !defined(__CUDA_ARCH__)
294 cudaMallocManaged(&ptr, _capacity*
sizeof(T), cudaMemAttachGlobal);
295#elif defined(__HIP__) && !defined(__HIP_DEVICE_COMPILE__)
296 hipMallocManaged(&ptr, _capacity*
sizeof(T), hipMemAttachGlobal);
298 _data =
new(ptr) T[_capacity];
302 if (address ==
nullptr || *address ==
nullptr) _data =
new T [_capacity];
304 _data =
new(*address) T [_capacity];
305 *address +=
sizeof(T) * _capacity;
313 if (_capacity != 0) {
314 assert( _capacity >= s);
318 assert( _capacity == 0 );
321 _mem_type = mem_flag;
322 if( s == 0 ){ _data =
nullptr;
return;}
323 switch ( (
int) _mem_type){
325 if (address ==
nullptr || *address ==
nullptr) {
326 _data =
new T [_capacity];
329 _data =
new(*address) T [_capacity];
330 std::size_t delta =
sizeof(T) * _capacity;
331 std::size_t sub = delta % 8;
332 if (sub != 0) delta += (8-sub);
339#if defined(__CUDACC__) && !defined(__CUDA_ARCH__)
340 cudaMallocManaged(&ptr, _capacity*
sizeof(T), cudaMemAttachGlobal);
341#elif defined(__HIP__) && !defined(__HIP_DEVICE_COMPILE__)
342 hipMallocManaged(&ptr, _capacity*
sizeof(T), hipMemAttachGlobal);
344 _data =
new(ptr) T[_capacity];
348 if (address ==
nullptr || *address ==
nullptr) _data =
new T [_capacity];
350 _data =
new(*address) T [_capacity];
351 std::size_t delta =
sizeof(T) * _capacity;
352 std::size_t sub = delta % 8;
353 if (sub != 0) delta += (8-sub);
362 assert( _capacity == 0 );
365 _mem_type = mem_flag;
366 if( s == 0 ){ _data =
nullptr;
return;}
367 switch ( (
int) _mem_type){
369 if (address ==
nullptr || *address ==
nullptr) _data =
new T [_capacity];
371 _data =
new(*address) T [_capacity];
372 std::size_t delta =
sizeof(T) * _capacity;
373 std::size_t sub = delta % 8;
374 if (sub != 0) delta += (8-sub);
381#if defined(__CUDACC__) && !defined(__CUDA_ARCH__)
382 cudaMallocManaged(&ptr, _capacity*
sizeof(T), cudaMemAttachGlobal);
383#elif defined(__HIP__) && !defined(__HIP_DEVICE_COMPILE__)
384 hipMallocManaged(&ptr, _capacity*
sizeof(T), hipMemAttachGlobal);
386 _data =
new(ptr) T[_capacity];
390 if (address ==
nullptr || *address ==
nullptr) _data =
new T [_capacity];
392 _data =
new(*address) T [_capacity];
393 std::size_t delta =
sizeof(T) * _capacity;
394 std::size_t sub = delta % 8;
395 if (sub != 0) delta += (8-sub);
397 *address +=
sizeof(T) * _capacity;
401 for (std::size_t ii = 0; ii < _capacity; ++ii)
407 return ( _size == 0 );
412 assert( NewEnd <= _size );
429 assert( _size + listSize < _capacity );
431 for( std::size_t i = _size; i < _size + listSize; i++ )
433 _data[i] = list[ i-_size ];
441 if (_size+inc > _capacity)
442 {
MCGIDI_PRINTF(
"inc too much (size %d, inc %d cap %d)\n", _size, inc, _capacity); abort(); }
443 assert(_size+inc <= _capacity);
447 {pos = _size; _size = _size + inc;}
454 std::size_t delta =
sizeof(T) * _size;
455 std::size_t sub = delta % 8;
456 if (sub != 0) delta += (8-sub);
#define MCGIDI_SWAP(a, b, type)
int MCGIDI_VectorSizeType
LUPI_HOST_DEVICE void push_back(const T &dataElem)
LUPI_HOST_DEVICE std::size_t internalSize() const
LUPI_HOST_DEVICE const_iterator begin() const
LUPI_HOST Vector< T > & operator=(const std::vector< T > &aa)
LUPI_HOST_DEVICE std::size_t size() const
LUPI_HOST_DEVICE void pop_back()
LUPI_HOST_DEVICE iterator begin()
LUPI_HOST_DEVICE void forceCreate(std::size_t a_size, T *a_data)
LUPI_HOST_DEVICE T & operator[](std::size_t index)
LUPI_HOST_DEVICE Vector(const Vector< T > &aa)
LUPI_HOST_DEVICE const_iterator end() const
LUPI_HOST_DEVICE Vector(std::size_t s, const T &d, bool mem_flag=CPU_MEM)
LUPI_HOST_DEVICE void eraseEnd(std::size_t NewEnd)
LUPI_HOST_DEVICE void swap(Vector< T > &other)
Needed for copy-swap idiom.
LUPI_HOST Vector(const std::vector< T > &aa)
LUPI_HOST_DEVICE void clear()
LUPI_HOST_DEVICE ~Vector()
LUPI_HOST_DEVICE Vector()
LUPI_HOST_DEVICE bool empty() const
LUPI_HOST_DEVICE std::size_t atomic_Index_Inc(std::size_t inc)
LUPI_HOST_DEVICE void appendList(std::size_t listSize, T *list)
LUPI_HOST_DEVICE Vector(std::size_t s, bool mem_flag=CPU_MEM)
LUPI_HOST_DEVICE int get_mem_type()
LUPI_HOST_DEVICE T & back() const
LUPI_HOST_DEVICE std::size_t capacity() const
LUPI_HOST_DEVICE T & back()
LUPI_HOST_DEVICE void resize(std::size_t s, const T &d, char **address=nullptr, bool mem_flag=CPU_MEM)
LUPI_HOST_DEVICE void resize(std::size_t s, char **address=nullptr, bool mem_flag=CPU_MEM)
LUPI_HOST_DEVICE iterator end()
LUPI_HOST_DEVICE Vector< T > & operator=(const Vector< T > &aa)
Implement assignment using copy-swap idiom.
LUPI_HOST_DEVICE const T & operator[](std::size_t index) const
LUPI_HOST_DEVICE void reserve(std::size_t s, char **address=nullptr, bool mem_flag=CPU_MEM)
Simple C++ string class, useful as replacement for std::string if this cannot be used,...