31 explicit operator bool() const noexcept {
return flag; }
41#if defined(AMREX_USE_GPU)
53template <
typename T,
bool SINGLE_WORD>
struct BlockStatus {};
56struct BlockStatus<T, true>
62 void operator=(Data<U>
const&) =
delete;
63 void operator=(Data<U> &&) =
delete;
68 void write (
char a_status, T a_value) {
69#if defined(AMREX_USE_CUDA)
70 volatile uint64_t tmp;
71 reinterpret_cast<STVA<T> volatile&
>(tmp).status = a_status;
72 reinterpret_cast<STVA<T> volatile&
>(tmp).value = a_value;
73 reinterpret_cast<uint64_t&
>(d.s) = tmp;
76 tmp.s = {a_status, a_value};
77 static_assert(
sizeof(
unsigned long long) ==
sizeof(uint64_t),
78 "HIP/SYCL: unsigned long long must be 64 bits");
79 Gpu::Atomic::Exch(
reinterpret_cast<unsigned long long*
>(&d),
80 reinterpret_cast<unsigned long long&
>(tmp));
85 T get_aggregate()
const {
return d.s.value; }
88 STVA<T> read ()
volatile {
89#if defined(AMREX_USE_CUDA)
90 volatile uint64_t tmp =
reinterpret_cast<uint64_t volatile&
>(d);
91 return {
reinterpret_cast<STVA<T> volatile&
>(tmp).status,
92 reinterpret_cast<STVA<T> volatile&
>(tmp).value };
94 static_assert(
sizeof(
unsigned long long) ==
sizeof(uint64_t),
95 "HIP/SYCL: unsigned long long must be 64 bits");
96 unsigned long long tmp = Gpu::Atomic::Add
97 (
reinterpret_cast<unsigned long long*
>(
const_cast<Data<T>*
>(&d)), 0ull);
98 return (*
reinterpret_cast<Data<T>*
>(&tmp)).s;
103 void set_status (
char a_status) { d.s.status = a_status; }
106 STVA<T> wait ()
volatile {
109#if defined(AMREX_USE_SYCL)
110 sycl::atomic_fence(sycl::memory_order::acq_rel, sycl::memory_scope::work_group);
112 __threadfence_block();
115 }
while (r.status ==
'x');
121struct BlockStatus<T, false>
128 void write (
char a_status, T a_value) {
129 if (a_status ==
'a') {
134#if defined(AMREX_USE_SYCL)
135 sycl::atomic_fence(sycl::memory_order::acq_rel, sycl::memory_scope::device);
143 T get_aggregate()
const {
return aggregate; }
146 STVA<T> read ()
volatile {
147#if defined(AMREX_USE_SYCL)
148 constexpr auto mo = sycl::memory_order::relaxed;
149 constexpr auto ms = sycl::memory_scope::device;
150 constexpr auto as = sycl::access::address_space::global_space;
154 }
else if (status ==
'a') {
155#if defined(AMREX_USE_SYCL)
156 sycl::atomic_ref<T,mo,ms,as> ar{
const_cast<T&
>(aggregate)};
157 return {
'a', ar.load()};
159 return {
'a', aggregate};
162#if defined(AMREX_USE_SYCL)
163 sycl::atomic_ref<T,mo,ms,as> ar{
const_cast<T&
>(
inclusive)};
164 return {
'p', ar.load()};
172 void set_status (
char a_status) { status = a_status; }
175 STVA<T> wait ()
volatile {
179#if defined(AMREX_USE_SYCL)
180 sycl::atomic_fence(sycl::memory_order::acq_rel, sycl::memory_scope::device);
184 }
while (r.status ==
'x');
192#if defined(AMREX_USE_SYCL)
194#ifndef AMREX_SYCL_NO_MULTIPASS_SCAN
195template <
typename T, std::
integral N,
typename FIN,
typename FOUT,
typename TYPE>
196T PrefixSum_mp (N n, FIN
const& fin, FOUT
const& fout, TYPE, RetSum a_ret_sum)
198 if (n <= 0) {
return 0; }
199 constexpr int nwarps_per_block = 8;
201 constexpr int nchunks = 12;
202 constexpr int nelms_per_block = nthreads * nchunks;
204 std::numeric_limits<int>::max())*nelms_per_block);
205 int nblocks = (
static_cast<Long>(n) + nelms_per_block - 1) / nelms_per_block;
209 std::size_t nbytes_blockresult =
Arena::align(
sizeof(T)*n);
210 std::size_t nbytes_blocksum =
Arena::align(
sizeof(T)*nblocks);
215 T* blockresult_p = (T*)dp;
216 T* blocksum_p = (T*)(dp + nbytes_blockresult);
217 T* totalsum_p = (T*)(dp + nbytes_blockresult + nbytes_blocksum);
219 amrex::launch<nthreads>(nblocks, sm, stream,
222 sycl::sub_group
const& sg = gh.item->get_sub_group();
223 int lane = sg.get_local_id()[0];
224 int warp = sg.get_group_id()[0];
225 int nwarps = sg.get_group_range()[0];
227 int threadIdxx = gh.item->get_local_id(0);
228 int blockIdxx = gh.item->get_group_linear_id();
229 int blockDimx = gh.item->get_local_range(0);
231 T* shared = (T*)(gh.local);
235 N ibegin =
static_cast<N
>(nelms_per_block) * blockIdxx;
236 N iend =
amrex::min(
static_cast<N
>(ibegin+nelms_per_block), n);
241 T sum_prev_chunk = 0;
242 for (
int ichunk = 0; ichunk < nchunks; ++ichunk) {
243 N
offset = ibegin + ichunk*blockDimx;
244 if (
offset >= iend) {
break; }
251 T s = sycl::shift_group_right(sg,
x, i);
252 if (lane >= i) {
x += s; }
262 gh.item->barrier(sycl::access::fence_space::local_space);
267 T
y = (lane < nwarps) ? shared[lane] : 0;
269 T s = sycl::shift_group_right(sg,
y, i);
270 if (lane >= i) {
y += s; }
273 if (lane < nwarps) { shared2[lane] =
y; }
276 gh.item->barrier(sycl::access::fence_space::local_space);
283 T sum_prev_warp = (warp == 0) ? 0 : shared2[warp-1];
284 T tmp_out = sum_prev_warp + sum_prev_chunk +
285 (std::is_same_v<std::decay_t<TYPE>,Type::Inclusive> ?
x :
x-x0);
286 sum_prev_chunk += shared2[nwarps-1];
289 blockresult_p[
offset] = tmp_out;
294 if (threadIdxx == 0) {
295 blocksum_p[blockIdxx] = sum_prev_chunk;
299 amrex::launch<nthreads>(1, sm, stream,
302 sycl::sub_group
const& sg = gh.item->get_sub_group();
303 int lane = sg.get_local_id()[0];
304 int warp = sg.get_group_id()[0];
305 int nwarps = sg.get_group_range()[0];
307 int threadIdxx = gh.item->get_local_id(0);
308 int blockDimx = gh.item->get_local_range(0);
310 T* shared = (T*)(gh.local);
313 T sum_prev_chunk = 0;
318 T s = sycl::shift_group_right(sg,
x, i);
319 if (lane >= i) {
x += s; }
329 gh.item->barrier(sycl::access::fence_space::local_space);
334 T
y = (lane < nwarps) ? shared[lane] : 0;
336 T s = sycl::shift_group_right(sg,
y, i);
337 if (lane >= i) {
y += s; }
340 if (lane < nwarps) { shared2[lane] =
y; }
343 gh.item->barrier(sycl::access::fence_space::local_space);
350 T sum_prev_warp = (warp == 0) ? 0 : shared2[warp-1];
351 T tmp_out = sum_prev_warp + sum_prev_chunk +
x;
352 sum_prev_chunk += shared2[nwarps-1];
355 blocksum_p[
offset] = tmp_out;
360 if (threadIdxx == 0) {
361 *totalsum_p = sum_prev_chunk;
365 amrex::launch<nthreads>(nblocks, 0, stream,
368 int threadIdxx = gh.item->get_local_id(0);
369 int blockIdxx = gh.item->get_group_linear_id();
370 int blockDimx = gh.item->get_local_range(0);
373 N ibegin =
static_cast<N
>(nelms_per_block) * blockIdxx;
374 N iend =
amrex::min(
static_cast<N
>(ibegin+nelms_per_block), n);
375 T prev_sum = (blockIdxx == 0) ? 0 : blocksum_p[blockIdxx-1];
397template <
typename T, std::
integral N,
typename FIN,
typename FOUT,
typename TYPE>
398requires (std::same_as<std::decay_t<TYPE>,Type::Inclusive> ||
399 std::same_as<std::decay_t<TYPE>,Type::Exclusive>)
400T
PrefixSum (N n, FIN && fin, FOUT && fout, TYPE type, RetSum a_ret_sum =
retSum)
402 if (n <= 0) {
return 0; }
403 constexpr int nwarps_per_block = 8;
405 constexpr int nchunks = 12;
406 constexpr int nelms_per_block = nthreads * nchunks;
408 std::numeric_limits<int>::max())*nelms_per_block);
409 int nblocks = (
static_cast<Long>(n) + nelms_per_block - 1) / nelms_per_block;
411#ifndef AMREX_SYCL_NO_MULTIPASS_SCAN
413 return PrefixSum_mp<T>(n, std::forward<FIN>(fin), std::forward<FOUT>(fout), type, a_ret_sum);
420 using BlockStatusT = std::conditional_t<
sizeof(detail::STVA<T>) <= 8,
421 detail::BlockStatus<T,true>, detail::BlockStatus<T,false> >;
423 std::size_t nbytes_blockstatus =
Arena::align(
sizeof(BlockStatusT)*nblocks);
424 std::size_t nbytes_blockid =
Arena::align(
sizeof(
unsigned int));
430 unsigned int*
AMREX_RESTRICT virtual_block_id_p = (
unsigned int*)(dp + nbytes_blockstatus);
431 T*
AMREX_RESTRICT totalsum_p = (T*)(dp + nbytes_blockstatus + nbytes_blockid);
434 BlockStatusT& block_status = block_status_p[i];
435 block_status.set_status(
'x');
437 *virtual_block_id_p = 0;
442 amrex::launch<nthreads>(nblocks, sm, stream,
445 sycl::sub_group
const& sg = gh.item->get_sub_group();
446 int lane = sg.get_local_id()[0];
447 int warp = sg.get_group_id()[0];
448 int nwarps = sg.get_group_range()[0];
450 int threadIdxx = gh.item->get_local_id(0);
451 int blockDimx = gh.item->get_local_range(0);
452 int gridDimx = gh.item->get_group_range(0);
454 T* shared = (T*)(gh.local);
460 int virtual_block_id = 0;
462 int& virtual_block_id_shared = *((
int*)(shared2+nwarps));
463 if (threadIdxx == 0) {
465 virtual_block_id_shared = bid;
467 gh.item->barrier(sycl::access::fence_space::local_space);
468 virtual_block_id = virtual_block_id_shared;
472 N ibegin =
static_cast<N
>(nelms_per_block) * virtual_block_id;
473 N iend =
amrex::min(
static_cast<N
>(ibegin+nelms_per_block), n);
474 BlockStatusT& block_status = block_status_p[virtual_block_id];
485 T sum_prev_chunk = 0;
487 for (
int ichunk = 0; ichunk < nchunks; ++ichunk) {
488 N
offset = ibegin + ichunk*blockDimx;
489 if (
offset >= iend) {
break; }
493 if (std::is_same_v<std::decay_t<TYPE>,Type::Exclusive> &&
offset == n-1) {
499 T s = sycl::shift_group_right(sg,
x, i);
500 if (lane >= i) {
x += s; }
510 gh.item->barrier(sycl::access::fence_space::local_space);
515 T
y = (lane < nwarps) ? shared[lane] : 0;
517 T s = sycl::shift_group_right(sg,
y, i);
518 if (lane >= i) {
y += s; }
521 if (lane < nwarps) { shared2[lane] =
y; }
524 gh.item->barrier(sycl::access::fence_space::local_space);
531 T sum_prev_warp = (warp == 0) ? 0 : shared2[warp-1];
532 tmp_out[ichunk] = sum_prev_warp + sum_prev_chunk +
533 (std::is_same_v<std::decay_t<TYPE>,Type::Inclusive> ?
x :
x-x0);
534 sum_prev_chunk += shared2[nwarps-1];
538 if (threadIdxx == 0 && gridDimx > 1) {
539 block_status.write((virtual_block_id == 0) ?
'p' :
'a',
543 if (virtual_block_id == 0) {
544 for (
int ichunk = 0; ichunk < nchunks; ++ichunk) {
545 N
offset = ibegin + ichunk*blockDimx + threadIdxx;
546 if (
offset >= iend) {
break; }
547 fout(
offset, tmp_out[ichunk]);
549 *totalsum_p += tmp_out[ichunk];
552 }
else if (virtual_block_id > 0) {
555 T exclusive_prefix = 0;
556 BlockStatusT
volatile* pbs = block_status_p;
559 int iblock = iblock0-lane;
560 detail::STVA<T> stva{
'p', 0};
562 stva = pbs[iblock].wait();
568 unsigned status_bf = (stva.status ==
'p') ? (0x1u << lane) : 0;
570 status_bf |= sycl::permute_group_by_xor(sg, status_bf, i);
573 bool stop_lookback = status_bf & 0x1u;
574 if (stop_lookback ==
false) {
575 if (status_bf != 0) {
577 if (lane > 0) {
x = 0; }
578 unsigned int bit_mask = 0x1u;
581 if (i == lane) {
x =
y; }
582 if (status_bf & bit_mask) {
583 stop_lookback =
true;
590 x += sycl::shift_group_left(sg,
x,i);
594 if (lane == 0) { exclusive_prefix +=
x; }
595 if (stop_lookback) {
break; }
599 block_status.write(
'p', block_status.get_aggregate() + exclusive_prefix);
600 shared[0] = exclusive_prefix;
604 gh.item->barrier(sycl::access::fence_space::local_space);
606 T exclusive_prefix = shared[0];
608 for (
int ichunk = 0; ichunk < nchunks; ++ichunk) {
609 N
offset = ibegin + ichunk*blockDimx + threadIdxx;
610 if (
offset >= iend) {
break; }
611 T t = tmp_out[ichunk] + exclusive_prefix;
637#elif defined(AMREX_USE_HIP)
639template <
typename T, std::
integral N,
typename FIN,
typename FOUT,
typename TYPE>
640requires (std::same_as<std::decay_t<TYPE>,Type::Inclusive> ||
641 std::same_as<std::decay_t<TYPE>,Type::Exclusive>)
642T
PrefixSum (N n, FIN
const& fin, FOUT
const& fout, TYPE, RetSum a_ret_sum =
retSum)
644 if (n <= 0) {
return 0; }
645 constexpr int nwarps_per_block = 4;
647 constexpr int nelms_per_thread =
sizeof(T) >= 8 ? 8 : 16;
648 constexpr int nelms_per_block = nthreads * nelms_per_thread;
650 std::numeric_limits<int>::max())*nelms_per_block);
651 int nblocks = (n + nelms_per_block - 1) / nelms_per_block;
655 using ScanTileState = rocprim::detail::lookback_scan_state<T>;
656 using OrderedBlockId = rocprim::detail::ordered_block_id<unsigned int>;
658#if (defined(HIP_VERSION_MAJOR) && (HIP_VERSION_MAJOR < 6)) || \
659 (defined(HIP_VERSION_MAJOR) && (HIP_VERSION_MAJOR == 6) && \
660 defined(HIP_VERSION_MINOR) && (HIP_VERSION_MINOR == 0))
662 std::size_t nbytes_tile_state = rocprim::detail::align_size
663 (ScanTileState::get_storage_size(nblocks));
664 std::size_t nbytes_block_id = OrderedBlockId::get_storage_size();
666 auto dp = (
char*)(
The_Arena()->
alloc(nbytes_tile_state+nbytes_block_id));
668 ScanTileState tile_state = ScanTileState::create(dp, nblocks);
672 std::size_t nbytes_tile_state;
673 AMREX_HIP_SAFE_CALL(ScanTileState::get_storage_size(nblocks, stream, nbytes_tile_state));
674 nbytes_tile_state = rocprim::detail::align_size(nbytes_tile_state);
676 std::size_t nbytes_block_id = OrderedBlockId::get_storage_size();
678 auto dp = (
char*)(
The_Arena()->
alloc(nbytes_tile_state+nbytes_block_id));
680 ScanTileState tile_state;
681 AMREX_HIP_SAFE_CALL(ScanTileState::create(tile_state, dp, nblocks, stream));
685 auto ordered_block_id = OrderedBlockId::create
686 (
reinterpret_cast<OrderedBlockId::id_type*
>(dp + nbytes_tile_state));
689 amrex::launch<nthreads>((nblocks+nthreads-1)/nthreads, 0, stream, [=]
AMREX_GPU_DEVICE ()
691 auto& scan_tile_state =
const_cast<ScanTileState&
>(tile_state);
692 auto& scan_bid =
const_cast<OrderedBlockId&
>(ordered_block_id);
693 const unsigned int gid = blockIdx.x*nthreads + threadIdx.x;
694 if (gid == 0) { scan_bid.reset(); }
695 scan_tile_state.initialize_prefix(gid, nblocks);
700 amrex::launch_global<nthreads> <<<nblocks, nthreads, sm, stream>>> (
703 using BlockLoad = rocprim::block_load<T, nthreads, nelms_per_thread,
704 rocprim::block_load_method::block_load_transpose>;
705 using BlockScan = rocprim::block_scan<T, nthreads,
706 rocprim::block_scan_algorithm::using_warp_scan>;
707 using BlockExchange = rocprim::block_exchange<T, nthreads, nelms_per_thread>;
708 using LookbackScanPrefixOp = rocprim::detail::lookback_scan_prefix_op
709 <T, rocprim::plus<T>, ScanTileState>;
711 __shared__
struct TempStorage {
712 typename OrderedBlockId::storage_type ordered_bid;
714 typename BlockLoad::storage_type load;
715 typename BlockExchange::storage_type exchange;
716 typename BlockScan::storage_type scan;
721 auto& scan_tile_state =
const_cast<ScanTileState&
>(tile_state);
722 auto& scan_bid =
const_cast<OrderedBlockId&
>(ordered_block_id);
724 auto const virtual_block_id = scan_bid.get(threadIdx.x, temp_storage.ordered_bid);
727 N ibegin =
static_cast<N
>(nelms_per_block) * virtual_block_id;
728 N iend =
amrex::min(
static_cast<N
>(ibegin+nelms_per_block), n);
730 auto input_begin = rocprim::make_transform_iterator(
731 rocprim::make_counting_iterator(N(0)),
732 [&] (N i) -> T {
return fin(i+ibegin); });
734 T data[nelms_per_thread];
735 if (
static_cast<int>(iend-ibegin) == nelms_per_block) {
736 BlockLoad().load(input_begin, data, temp_storage.load);
739 BlockLoad().load(input_begin, data, iend-ibegin, 0, temp_storage.load);
744 constexpr bool is_exclusive = std::is_same_v<std::decay_t<TYPE>,Type::Exclusive>;
746 if (virtual_block_id == 0) {
749 BlockScan().exclusive_scan(data, data, T{0}, block_agg, temp_storage.scan);
751 BlockScan().inclusive_scan(data, data, block_agg, temp_storage.scan);
753 if (threadIdx.x == 0) {
755 scan_tile_state.set_complete(0, block_agg);
756 }
else if (nblocks == 1 && totalsum_p) {
757 *totalsum_p = block_agg;
761 T last = data[nelms_per_thread-1];
763 LookbackScanPrefixOp prefix_op(virtual_block_id, rocprim::plus<T>(), scan_tile_state);
765 BlockScan().exclusive_scan(data, data, temp_storage.scan, prefix_op,
768 BlockScan().inclusive_scan(data, data, temp_storage.scan, prefix_op,
772 if (iend == n && threadIdx.x == nthreads-1) {
773 T tsum = data[nelms_per_thread-1];
782 BlockExchange().blocked_to_striped(data, data, temp_storage.exchange);
784 for (
int i = 0; i < nelms_per_thread; ++i) {
785 N
offset = ibegin + i*nthreads + threadIdx.x;
799 T ret = (a_ret_sum) ? *totalsum_p : T(0);
805#elif defined(AMREX_USE_CUDA)
807template <
typename T, std::
integral N,
typename FIN,
typename FOUT,
typename TYPE>
808requires (std::same_as<std::decay_t<TYPE>,Type::Inclusive> ||
809 std::same_as<std::decay_t<TYPE>,Type::Exclusive>)
812 if (n <= 0) {
return 0; }
813 constexpr int nwarps_per_block = 8;
815 constexpr int nelms_per_thread =
sizeof(T) >= 8 ? 4 : 8;
816 constexpr int nelms_per_block = nthreads * nelms_per_thread;
818 std::numeric_limits<int>::max())*nelms_per_block);
819 int nblocks = (n + nelms_per_block - 1) / nelms_per_block;
823 using ScanTileState = cub::ScanTileState<T>;
824 std::size_t tile_state_size = 0;
825 ScanTileState::AllocationSize(nblocks, tile_state_size);
827 std::size_t nbytes_tile_state =
Arena::align(tile_state_size);
828 auto tile_state_p = (
char*)(
The_Arena()->
alloc(nbytes_tile_state));
830 ScanTileState tile_state;
831 tile_state.Init(nblocks, tile_state_p, tile_state_size);
835 amrex::launch<nthreads>((nblocks+nthreads-1)/nthreads, 0, stream, [=]
AMREX_GPU_DEVICE ()
837 const_cast<ScanTileState&
>(tile_state).InitializeStatus(nblocks);
843 amrex::launch_global<nthreads> <<<nblocks, nthreads, sm, stream>>> (
846 using BlockLoad = cub::BlockLoad<T, nthreads, nelms_per_thread, cub::BLOCK_LOAD_WARP_TRANSPOSE>;
847 using BlockScan = cub::BlockScan<T, nthreads, cub::BLOCK_SCAN_WARP_SCANS>;
848 using BlockExchange = cub::BlockExchange<T, nthreads, nelms_per_thread>;
850#ifdef AMREX_CUDA_CCCL_VER_GE_2_8
851 using Sum = cuda::std::plus<T>;
853 using Sum = cub::Sum;
855 using TilePrefixCallbackOp = cub::TilePrefixCallbackOp<T, Sum, ScanTileState>;
857 __shared__
union TempStorage
859 typename BlockLoad::TempStorage load;
860 typename BlockExchange::TempStorage exchange;
862 typename BlockScan::TempStorage scan;
863 typename TilePrefixCallbackOp::TempStorage prefix;
868 auto& scan_tile_state =
const_cast<ScanTileState&
>(tile_state);
870 int virtual_block_id = blockIdx.x;
873 N ibegin =
static_cast<N
>(nelms_per_block) * virtual_block_id;
874 N iend =
amrex::min(
static_cast<N
>(ibegin+nelms_per_block), n);
876 auto input_lambda = [&] (N i) -> T {
return fin(i+ibegin); };
877#ifdef AMREX_CUDA_CCCL_VER_GE_2_8
878 thrust::transform_iterator<
decltype(input_lambda),thrust::counting_iterator<N> >
879 input_begin(thrust::counting_iterator<N>(0), input_lambda);
881 cub::TransformInputIterator<T,
decltype(input_lambda),cub::CountingInputIterator<N> >
882 input_begin(cub::CountingInputIterator<N>(0), input_lambda);
885 T data[nelms_per_thread];
886 if (
static_cast<int>(iend-ibegin) == nelms_per_block) {
887 BlockLoad(temp_storage.load).Load(input_begin, data);
889 BlockLoad(temp_storage.load).Load(input_begin, data, iend-ibegin, 0);
894 constexpr bool is_exclusive = std::is_same_v<std::decay_t<TYPE>,
Type::Exclusive>;
896 if (virtual_block_id == 0) {
899 BlockScan(temp_storage.scan_storeage.scan).ExclusiveSum(data, data, block_agg);
901 BlockScan(temp_storage.scan_storeage.scan).InclusiveSum(data, data, block_agg);
903 if (threadIdx.x == 0) {
905 scan_tile_state.SetInclusive(0, block_agg);
906 }
else if (nblocks == 1 && totalsum_p) {
907 *totalsum_p = block_agg;
911 T last = data[nelms_per_thread-1];
913 TilePrefixCallbackOp prefix_op(scan_tile_state, temp_storage.scan_storeage.prefix,
914 Sum{}, virtual_block_id);
916 BlockScan(temp_storage.scan_storeage.scan).ExclusiveSum(data, data, prefix_op);
918 BlockScan(temp_storage.scan_storeage.scan).InclusiveSum(data, data, prefix_op);
921 if (iend == n && threadIdx.x == nthreads-1) {
922 T tsum = data[nelms_per_thread-1];
931 BlockExchange(temp_storage.exchange).BlockedToStriped(data);
933 for (
int i = 0; i < nelms_per_thread; ++i) {
934 N
offset = ibegin + i*nthreads + threadIdx.x;
948 T ret = (a_ret_sum) ? *totalsum_p : T(0);
967template <std::
integral N,
typename T >
970 if (n <= 0) {
return 0; }
971#if defined(AMREX_USE_CUDA)
972 void* d_temp =
nullptr;
973 std::size_t temp_bytes = 0;
987#elif defined(AMREX_USE_HIP)
988 void* d_temp =
nullptr;
989 std::size_t temp_bytes = 0;
1003#elif defined(AMREX_USE_SYCL) && defined(AMREX_USE_ONEDPL)
1004 auto policy = oneapi::dpl::execution::make_device_policy(Gpu::Device::streamQueue());
1005 std::inclusive_scan(policy, in, in+n, out, std::plus<T>(), T(0));
1014 if (
static_cast<Long>(n) <=
static_cast<Long>(std::numeric_limits<int>::max())) {
1015 return PrefixSum<T>(
static_cast<int>(n),
1020 return PrefixSum<T>(n,
1039template <std::
integral N,
typename T >
1042 if (n <= 0) {
return 0; }
1043#if defined(AMREX_USE_CUDA)
1048 void* d_temp =
nullptr;
1049 std::size_t temp_bytes = 0;
1062 return in_last+out_last;
1063#elif defined(AMREX_USE_HIP)
1068 void* d_temp =
nullptr;
1069 std::size_t temp_bytes = 0;
1082 return in_last+out_last;
1083#elif defined(AMREX_USE_SYCL) && defined(AMREX_USE_ONEDPL)
1088 auto policy = oneapi::dpl::execution::make_device_policy(Gpu::Device::streamQueue());
1089 std::exclusive_scan(policy, in, in+n, out, T(0), std::plus<T>());
1096 return in_last+out_last;
1098 if (
static_cast<Long>(n) <=
static_cast<Long>(std::numeric_limits<int>::max())) {
1099 return PrefixSum<T>(
static_cast<int>(n),
1104 return PrefixSum<T>(n,
1114template <
typename T, std::
integral N,
typename FIN,
typename FOUT,
typename TYPE>
1115requires (std::same_as<std::decay_t<TYPE>,Type::Inclusive> ||
1116 std::same_as<std::decay_t<TYPE>,Type::Exclusive>)
1117T
PrefixSum (N n, FIN
const& fin, FOUT
const& fout, TYPE, RetSum =
retSum)
1119 if (n <= 0) {
return 0; }
1121 for (N i = 0; i < n; ++i) {
1134template <std::
integral N,
typename T >
1137#if (__cplusplus >= 201703L) && (!defined(_GLIBCXX_RELEASE) || _GLIBCXX_RELEASE >= 10)
1139 std::inclusive_scan(in, in+n, out);
1141 std::partial_sum(in, in+n, out);
1143 return (n > 0) ? out[n-1] : T(0);
1147template <std::
integral N,
typename T >
1150 if (n <= 0) {
return 0; }
1152 auto in_last = in[n-1];
1153#if (__cplusplus >= 201703L) && (!defined(_GLIBCXX_RELEASE) || _GLIBCXX_RELEASE >= 10)
1155 std::exclusive_scan(in, in+n, out, 0);
1158 std::partial_sum(in, in+n-1, out+1);
1160 return in_last + out[n-1];