Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
MCGIDI_vector.hpp
Go to the documentation of this file.
1/*
2# <<BEGIN-copyright>>
3# Copyright 2019, Lawrence Livermore National Security, LLC.
4# This file is part of the gidiplus package (https://github.com/LLNL/gidiplus).
5# gidiplus is licensed under the MIT license (see https://opensource.org/licenses/MIT).
6# SPDX-License-Identifier: MIT
7# <<END-copyright>>
8*/
9
10#ifndef MCGIDI_VECTOR_HPP
11#define MCGIDI_VECTOR_HPP
12
13#define CPU_MEM false
14#define UVM_MEM true
15
16#ifdef HAVE_OPENMP_TARGET
17 #ifdef USE_OPENMP_NO_GPU
18 #define VAR_MEM false
19 #else
20 #define VAR_MEM true
21 #endif
22#else
23 #define VAR_MEM false
24#endif
25
27
28#define MCGIDI_SWAP(a,b,type) {type ttttttttt=a;a=b;b=ttttttttt;}
29
30#if defined(__CUDACC__) && !defined(__CUDA_ARCH__)
31#include <cuda.h>
32#include <cuda_runtime.h>
33#include <cuda_runtime_api.h>
34#endif
35
36#if defined(__HIP__)
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>
41#endif
42
43#include <string.h>
44#include <stdio.h>
45#include "cassert"
46#include <algorithm>
47#include <LUPI_declareMacro.hpp>
48#include <vector>
49
50namespace MCGIDI {
51
52template <class T>
53class Vector
54{
55 private:
56 T* _data;
57 std::size_t _capacity;
58 std::size_t _size;
59 bool _mem_type;
60
61 public:
62 typedef T* iterator;
63 typedef T* const_iterator;
64
65 LUPI_HOST_DEVICE Vector() : _data(0), _capacity(0), _size(0), _mem_type(CPU_MEM) {};
66 LUPI_HOST_DEVICE Vector( std::size_t s, bool mem_flag = CPU_MEM ) : _data(0), _capacity(s), _size(s), _mem_type(mem_flag)
67 {
68
69 if( s == 0 ){ _data = nullptr; return;}
70 switch ((int)_mem_type){
71 case CPU_MEM:
72 _data = new T [_capacity];
73 break;
74 case UVM_MEM:
75 {
76 void *ptr = nullptr;
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);
81#endif
82 _data = new(ptr) T[_capacity];
83 break;
84 }
85 default:
86 _data = new T [_capacity];
87 break;
88 }
89 }
90 LUPI_HOST_DEVICE Vector( std::size_t s, const T& d, bool mem_flag = CPU_MEM ) : _data(0), _capacity(s), _size(s), _mem_type(mem_flag)
91 {
92 if( s == 0 ){ _data = nullptr; return;}
93 switch ( (int) _mem_type){
94 case CPU_MEM:
95 _data = new T [_capacity];
96 break;
97 case UVM_MEM:
98 {
99 void *ptr = nullptr;
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);
104#endif
105 _data = new(ptr) T[_capacity];
106 break;
107 }
108 default:
109 _data = new T [_capacity];
110 break;
111 }
112 for (std::size_t ii = 0; ii < _capacity; ++ii)
113 _data[ii] = d;
114 }
115
117 : _data(0), _capacity(aa._capacity), _size(aa._size), _mem_type(aa._mem_type)
118 {
119 if( _capacity == 0 ){ _data = nullptr; return; }
120
121 switch ( (int) _mem_type){
122 case CPU_MEM:
123 _data = new T [_capacity];
124 break;
125 case UVM_MEM:
126 {
127 void *ptr = nullptr;
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);
132#endif
133 _data = new(ptr) T[_capacity];
134 break;
135 }
136 default:
137 _data = new T [_capacity];
138 break;
139 }
140
141 for (std::size_t ii=0; ii<_size; ++ii)
142 _data[ii] = aa._data[ii];
143 }
144
145 LUPI_HOST Vector(const std::vector<T>& aa )
146 : _data(0), _capacity(aa.size()), _size(aa.size()), _mem_type(CPU_MEM)
147 {
148 if( _capacity == 0 ){ _data = nullptr; return;}
149
150 switch ( (int) _mem_type){
151 case CPU_MEM:
152 _data = new T [_capacity];
153 break;
154 case UVM_MEM:
155 {
156 void *ptr = nullptr;
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);
161#endif
162 _data = new(ptr) T[_capacity];
163 break;
164 }
165 default:
166 _data = new T [_capacity];
167 break;
168 }
169
170 for (std::size_t ii=0; ii<_size; ++ii)
171 _data[ii] = aa[ii];
172 }
173
175 switch ( (int) _mem_type){
176 case CPU_MEM:
177 delete[] _data;
178 break;
179 case UVM_MEM:
180 for (std::size_t i=0; i < _size; ++i)
181 _data[i].~T();
182#if defined(__CUDACC__) && !defined(__CUDA_ARCH__)
183 cudaFree(_data);
184#elif defined(__HIP__) && !defined(__HIP_DEVICE_COMPILE__)
185 hipFree(_data);
186#endif
187 break;
188 default:
189 delete[] _data;
190 break;
191 }
192 }
193
194 LUPI_HOST_DEVICE iterator begin() { return _data; }
195
196 LUPI_HOST_DEVICE const_iterator begin() const { return _data; }
197
198 LUPI_HOST_DEVICE iterator end() { return _data + _size; }
199
200 LUPI_HOST_DEVICE const_iterator end() const { return _data + _size; }
201
202 /// Needed for copy-swap idiom
204 {
205 MCGIDI_SWAP(_data, other._data, T*);
206 MCGIDI_SWAP(_capacity, other._capacity, std::size_t);
207 MCGIDI_SWAP(_size, other._size, std::size_t);
208 MCGIDI_SWAP(_mem_type, other._mem_type, bool);
209 }
210
211 /// Implement assignment using copy-swap idiom
213 {
214 if (&aa != this)
215 {
216 Vector<T> temp(aa);
217 this->swap(temp);
218 }
219 return *this;
220 }
221
222 LUPI_HOST Vector<T>& operator=(const std::vector<T>& aa)
223 {
224 Vector<T> temp(aa);
225 this->swap(temp);
226 return *this;
227 }
228
230 {
231 return _mem_type;
232 }
233
234 LUPI_HOST_DEVICE void push_back( const T& dataElem )
235 {
236 assert( _size < _capacity );
237 _data[_size] = dataElem;
238 _size++;
239 }
240
241 LUPI_HOST_DEVICE const T& operator[]( std::size_t index ) const
242 {
243 // assert( index < _capacity );
244 // assert( index >= 0); comment out pointless assertion size_t type is >= 0 by definition
245 return _data[index];
246 }
247
248 LUPI_HOST_DEVICE T& operator[]( std::size_t index )
249 {
250 // assert( index < _capacity );
251 // assert( index >= 0); comment out pointless assertion size_t type is >= 0 by definition
252 return _data[index];
253 }
254
255 LUPI_HOST_DEVICE std::size_t capacity() const
256 {
257 return _capacity;
258 }
259
260 LUPI_HOST_DEVICE std::size_t size() const
261 {
262 return _size;
263 }
264
266 {
267 return _data[_size-1];
268 }
269
271 {
272 return _data[_size-1];
273 }
274
275 LUPI_HOST_DEVICE void reserve( std::size_t s, char ** address = nullptr, bool mem_flag = CPU_MEM )
276 {
277 if (s == _capacity) return;
278 assert( _capacity == 0 );
279 _capacity = s;
280 _mem_type = mem_flag;
281 if( s == 0 ){ _data = nullptr; return;}
282 switch ( (int) _mem_type){
283 case CPU_MEM:
284 if (address == nullptr || *address == nullptr) _data = new T [_capacity];
285 else {
286 _data = new(*address) T [_capacity];
287 *address += sizeof(T) * _capacity;
288 }
289 break;
290 case UVM_MEM:
291 {
292 void *ptr = nullptr;
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);
297#endif
298 _data = new(ptr) T[_capacity];
299 break;
300 }
301 default:
302 if (address == nullptr || *address == nullptr) _data = new T [_capacity];
303 else {
304 _data = new(*address) T [_capacity];
305 *address += sizeof(T) * _capacity;
306 }
307 break;
308 }
309 }
310
311 LUPI_HOST_DEVICE void resize( std::size_t s, char ** address = nullptr, bool mem_flag = CPU_MEM )
312 {
313 if (_capacity != 0) {
314 assert( _capacity >= s);
315 _size = s;
316 return;
317 }
318 assert( _capacity == 0 );
319 _capacity = s;
320 _size = s;
321 _mem_type = mem_flag;
322 if( s == 0 ){ _data = nullptr; return;}
323 switch ( (int) _mem_type){
324 case CPU_MEM:
325 if (address == nullptr || *address == nullptr) {
326 _data = new T [_capacity];
327 }
328 else {
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);
333 *address += delta;
334 }
335 break;
336 case UVM_MEM:
337 {
338 void *ptr = nullptr;
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);
343#endif
344 _data = new(ptr) T[_capacity];
345 break;
346 }
347 default:
348 if (address == nullptr || *address == nullptr) _data = new T [_capacity];
349 else {
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);
354 *address += delta;
355 }
356 break;
357 }
358 }
359
360 LUPI_HOST_DEVICE void resize( std::size_t s, const T& d, char ** address = nullptr, bool mem_flag = CPU_MEM )
361 {
362 assert( _capacity == 0 );
363 _capacity = s;
364 _size = s;
365 _mem_type = mem_flag;
366 if( s == 0 ){ _data = nullptr; return;}
367 switch ( (int) _mem_type){
368 case CPU_MEM:
369 if (address == nullptr || *address == nullptr) _data = new T [_capacity];
370 else {
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);
375 *address += delta;
376 }
377 break;
378 case UVM_MEM:
379 {
380 void *ptr = nullptr;
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);
385#endif
386 _data = new(ptr) T[_capacity];
387 break;
388 }
389 default:
390 if (address == nullptr || *address == nullptr) _data = new T [_capacity];
391 else {
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);
396 *address += delta;
397 *address += sizeof(T) * _capacity;
398 }
399 break;
400 }
401 for (std::size_t ii = 0; ii < _capacity; ++ii)
402 _data[ii] = d;
403 }
404
406 {
407 return ( _size == 0 );
408 }
409
410 LUPI_HOST_DEVICE void eraseEnd( std::size_t NewEnd )
411 {
412 assert( NewEnd <= _size );
413 _size = NewEnd;
414 }
415
417 {
418 assert(_size > 0);
419 _size--;
420 }
421
423 {
424 _size = 0;
425 }
426
427 LUPI_HOST_DEVICE void appendList( std::size_t listSize, T* list )
428 {
429 assert( _size + listSize < _capacity );
430
431 for( std::size_t i = _size; i < _size + listSize; i++ )
432 {
433 _data[i] = list[ i-_size ];
434 }
435
436 }
437
438 //Atomically retrieve an availible index then increment that index some amount
439 LUPI_HOST_DEVICE std::size_t atomic_Index_Inc( std::size_t inc )
440 {
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);
444 std::size_t pos;
445
446// #include "mc_omp_atomic_capture.hh"
447 {pos = _size; _size = _size + inc;}
448
449 return pos;
450 }
451
452 // This will not work for a vector of base classes.
453 LUPI_HOST_DEVICE std::size_t internalSize() const {
454 std::size_t delta = sizeof(T) * _size;
455 std::size_t sub = delta % 8;
456 if (sub != 0) delta += (8-sub);
457 return delta;
458 }
459
460 LUPI_HOST_DEVICE void forceCreate(std::size_t a_size, T* a_data) {
461 _capacity = a_size;
462 _size = a_size;
463 _data = a_data;
464 }
465};
466
467}
468#endif
#define LUPI_HOST_DEVICE
#define LUPI_HOST
#define MCGIDI_PRINTF
Definition MCGIDI.hpp:40
#define MCGIDI_SWAP(a, b, type)
int MCGIDI_VectorSizeType
#define CPU_MEM
#define UVM_MEM
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,...
Definition MCGIDI.hpp:43