Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_CSR.H
Go to the documentation of this file.
1#ifndef AMREX_CSR_H_
2#define AMREX_CSR_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_Gpu.H>
6#include <AMReX_INT.H>
7#include <AMReX_OpenMP.H>
8
9#if defined(AMREX_USE_CUDA)
10#include <cub/cub.cuh> // for Clang
11#endif
12
13#include <algorithm>
14#include <climits>
15#include <type_traits>
16
17namespace amrex {
18
32template <typename T>
33struct CsrView {
34 using U = std::conditional_t<std::is_const_v<T>, Long const, Long>;
35 T* AMREX_RESTRICT mat = nullptr;
38 Long nnz = 0;
40};
41
48template <typename T, template <typename> class V>
49struct CSR {
50 V<T> mat;
51 V<Long> col_index;
52 V<Long> row_offset;
53 Long nnz = 0;
54
56 [[nodiscard]] Long nrows () const {
57 return row_offset.empty() ? Long(0) : Long(row_offset.size())-1;
58 }
59
69 void resize (Long num_rows, Long num_non_zeros) {
70 mat.resize(num_non_zeros);
71 col_index.resize(num_non_zeros);
72 row_offset.resize(num_rows+1);
73 nnz = num_non_zeros;
74 }
75
78 return CsrView<T>{mat.data(), col_index.data(), row_offset.data(),
79 nnz, Long(row_offset.size())-1};
80 }
81
83 [[nodiscard]] CsrView<T const> view () const {
84 return CsrView<T>{mat.data(), col_index.data(), row_offset.data(),
85 nnz, Long(row_offset.size())-1};
86 }
87
89 [[nodiscard]] CsrView<T const> const_view () const {
90 return CsrView<T const>{mat.data(), col_index.data(), row_offset.data(),
91 nnz, Long(row_offset.size())-1};
92 }
93
97 void sort ();
98
103};
104
105template <typename C, typename T, template<typename> class AD, template<typename> class AS,
106 std::enable_if_t<std::is_same_v<C,Gpu::HostToDevice> ||
107 std::is_same_v<C,Gpu::DeviceToHost> ||
108 std::is_same_v<C,Gpu::DeviceToDevice>, int> = 0>
116void duplicateCSR (C c, CSR<T,AD>& dst, CSR<T,AS> const& src)
117{
118 dst.mat.resize(src.mat.size());
119 dst.col_index.resize(src.col_index.size());
120 dst.row_offset.resize(src.row_offset.size());
122 src.mat.begin(),
123 src.mat.end(),
124 dst.mat.begin());
126 src.col_index.begin(),
127 src.col_index.end(),
128 dst.col_index.begin());
130 src.row_offset.begin(),
131 src.row_offset.end(),
132 dst.row_offset.begin());
133 dst.nnz = src.nnz;
134}
135
136template <typename T, template <typename> class V>
138{
139 if (nnz <= 0) { return; }
140
141#ifdef AMREX_USE_GPU
142
143#if defined(AMREX_USE_CUDA) || defined(AMREX_USE_HIP)
144
145 // The function is synchronous. If that is no longer the case, we might
146 // need to update SpMatrix::define.
147
148 constexpr int nthreads = 256;
149 constexpr int nwarps_per_block = nthreads / Gpu::Device::warp_size;
150
151 AMREX_ALWAYS_ASSERT((nrows()+nwarps_per_block-1) < Long(std::numeric_limits<int>::max()));
152
153 auto nr = int(nrows());
154 int nblocks = (nr + nwarps_per_block-1) / nwarps_per_block;
155 auto const& stream = Gpu::gpuStream();
156
157 auto* pmat = mat.data();
158 auto* pcol = col_index.data();
159 auto* prow = row_offset.data();
160
161 Gpu::Buffer<int> needs_fallback({0});
162 auto* d_needs_fallback = needs_fallback.data();
163
164 amrex::launch_global<nthreads><<<nblocks, nthreads, 0, stream>>>
165 ([=] AMREX_GPU_DEVICE () noexcept
166 {
167 int wid = int(threadIdx.x)/Gpu::Device::warp_size;
168 int r = int(blockIdx.x)*nwarps_per_block + wid;
169 if (r >= nr) return;
170
171 Long const b = prow[r];
172 Long const e = prow[r+1];
173 auto const len = int(e - b);
174
175 if (len <= 1) return;
176
177 int lane = int(threadIdx.x) - wid * Gpu::Device::warp_size;
178
179 bool sorted = true;
180 for (Long i = lane + 1; i < len; i += Gpu::Device::warp_size) {
181 sorted = sorted && (pcol[b+i-1] <= pcol[b+i]);
182 }
183#if defined(AMREX_USE_CUDA)
184 if (__all_sync(0xffffffff, sorted)) { return; }
185#else
186 if (__all(sorted)) { return; }
187#endif
188
189 constexpr int ITEMS_PER_THREAD = AMREX_HIP_OR_CUDA(2,4);
190 constexpr int ITEMS_PER_WARP = Gpu::Device::warp_size * ITEMS_PER_THREAD;
191
192 if (len <= ITEMS_PER_WARP)
193 {
194#if defined(AMREX_USE_CUDA)
195 using WarpSort = cub::WarpMergeSort<Long, ITEMS_PER_THREAD, Gpu::Device::warp_size, T>;
196 __shared__ typename WarpSort::TempStorage temp_storage[nwarps_per_block];
197#elif defined(AMREX_USE_HIP)
198 using WarpSort = rocprim::warp_sort<Long, Gpu::Device::warp_size, T>;
199 __shared__ typename WarpSort::storage_type temp_storage[nwarps_per_block];
200#endif
201
202 Long keys[ITEMS_PER_THREAD];
203 T values[ITEMS_PER_THREAD];
204
205 #pragma unroll
206 for (int i = 0; i < ITEMS_PER_THREAD; ++i) {
207 int idx = lane * ITEMS_PER_THREAD + i;
208 if (idx < len) {
209 keys[i] = pcol[b + idx];
210 values[i] = pmat[b + idx];
211 } else {
212 keys[i] = std::numeric_limits<Long>::max();
213 values[i] = T{};
214 }
215 }
216
218 WarpSort{}.sort(keys, values, temp_storage[wid]),
219 WarpSort(temp_storage[wid]).Sort(
220 keys, values, [](Long x, Long y) {return x < y;}));
221
222 #pragma unroll
223 for (int i = 0; i < ITEMS_PER_THREAD; ++i) {
224 int idx = lane * ITEMS_PER_THREAD + i;
225 if (idx < len) {
226 pcol[b + idx] = keys[i];
227 pmat[b + idx] = values[i];
228 }
229 }
230 } else {
231 if (lane == 0) {
232 Gpu::Atomic::AddNoRet(d_needs_fallback, 1);
233 }
234 }
235 });
236
237 auto* h_needs_fallback = needs_fallback.copyToHost();
238
239 if (*h_needs_fallback)
240 {
241 V<Long> col_index_out(col_index.size());
242 V<T> mat_out(mat.size());
243 auto* d_col_out = col_index_out.data();
244 auto* d_val_out = mat_out.data();
245
246 std::size_t temp_bytes = 0;
247
249 rocprim::segmented_radix_sort_pairs,
250 cub::DeviceSegmentedRadixSort::SortPairs)
251 (nullptr, temp_bytes, pcol, d_col_out, pmat, d_val_out,
252 nnz, nr, prow, prow+1, 0, int(sizeof(Long)*CHAR_BIT),
253 stream));
254
255 auto* d_temp = (void*) The_Arena()->alloc(temp_bytes);
256
258 rocprim::segmented_radix_sort_pairs,
259 cub::DeviceSegmentedRadixSort::SortPairs)
260 (d_temp, temp_bytes, pcol, d_col_out, pmat, d_val_out,
261 nnz, nr, prow, prow+1, 0, int(sizeof(Long)*CHAR_BIT),
262 stream));
263
264 std::swap(col_index, col_index_out);
265 std::swap(mat, mat_out);
266
268 The_Arena()->free(d_temp);
269 }
270
271 // let's test both by print matrix out to see if it's sorted.
272
274
275#elif defined(AMREX_USE_SYCL)
276
277 // xxxxx TODO SYCL: Let's not worry about performance for now.
279 duplicateCSR(Gpu::deviceToHost, h_csr, *this);
281 h_csr.sort_on_host();
282 duplicateCSR(Gpu::hostToDevice, *this, h_csr);
284
285#endif
286
287#else
288
289 sort_on_host();
290
291#endif
292}
293
294template <typename T, template <typename> class V>
296{
297 if (nnz <= 0) { return; }
298
299 constexpr int SMALL = 128;
300
301 Long nr = nrows();
302
303#ifdef AMREX_USE_OMP
304#pragma omp parallel
305#endif
306 {
307 V<Long> lcols;
308 V<T > lvals;
309 V<int > perm;
310
311 Long scols[SMALL];
312 T svals[SMALL];
313
314#ifdef AMREX_USE_OMP
315#pragma omp for
316#endif
317 for (Long r = 0; r < nr; ++r) {
318 Long const b = row_offset[r ];
319 Long const e = row_offset[r+1];
320 auto const len = int(e - b);
321
322 if (len <= 1) { continue; }
323
324 bool sorted = true;
325 for (int i = 1; i < len; ++i) {
326 if (col_index[b+i-1] > col_index[b+i]) {
327 sorted = false;
328 break;
329 }
330 }
331 if (sorted) { continue; }
332
333 if (len <= SMALL) {
334 // Insertion sort using arrays on stack
335 for (int i = 0; i < len; ++i) {
336 scols[i] = col_index[b+i];
337 svals[i] = mat [b+i];
338 }
339 for (int i = 1; i < len; ++i) {
340 auto c = scols[i];
341 auto v = svals[i];
342 auto j = i;
343 while (j > 0 && scols[j-1] > c) {
344 scols[j] = scols[j-1];
345 svals[j] = svals[j-1];
346 --j;
347 }
348 scols[j] = c;
349 svals[j] = v;
350 }
351 for (int i = 0; i < len; ++i) {
352 col_index[b+i] = scols[i];
353 mat [b+i] = svals[i];
354 }
355 } else {
356 lcols.resize(len);
357 lvals.resize(len);
358 perm.resize(len);
359
360 for (int i = 0; i < len; ++i) {
361 lcols[i] = col_index[b+i];
362 lvals[i] = mat [b+i];
363 perm [i] = i;
364 }
365
366 std::sort(perm.begin(), perm.end(),
367 [&] (int i0, int i1) {
368 return lcols[i0] < lcols[i1];
369 });
370
371 for (int out = 0; out < len; ++out) {
372 auto const in = perm[out];
373 col_index[b+out] = lcols[in];
374 mat [b+out] = lvals[in];
375 }
376 }
377 }
378 }
379}
380
381}
382
383#endif
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
#define AMREX_RESTRICT
Definition AMReX_Extension.H:32
#define AMREX_HIP_OR_CUDA(a, b)
Definition AMReX_GpuControl.H:21
#define AMREX_GPU_SAFE_CALL(call)
Definition AMReX_GpuError.H:63
#define AMREX_GPU_ERROR_CHECK()
Definition AMReX_GpuError.H:151
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
virtual void free(void *pt)=0
A pure virtual function for deleting the arena pointed to by pt.
virtual void * alloc(std::size_t sz)=0
Definition AMReX_GpuBuffer.H:18
T const * data() const noexcept
Definition AMReX_GpuBuffer.H:45
static constexpr int warp_size
Definition AMReX_GpuDevice.H:236
amrex_long Long
Definition AMReX_INT.H:30
Arena * The_Arena()
Definition AMReX_Arena.cpp:805
__host__ __device__ AMREX_FORCE_INLINE void AddNoRet(T *sum, T value) noexcept
Definition AMReX_GpuAtomic.H:284
void copyAsync(HostToDevice, InIter begin, InIter end, OutIter result) noexcept
A host-to-device copy routine. Note this is just a wrapper around memcpy, so it assumes contiguous st...
Definition AMReX_GpuContainers.H:228
static constexpr DeviceToHost deviceToHost
Definition AMReX_GpuContainers.H:106
static constexpr HostToDevice hostToDevice
Definition AMReX_GpuContainers.H:105
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
gpuStream_t gpuStream() noexcept
Definition AMReX_GpuDevice.H:291
Definition AMReX_Amr.cpp:49
void duplicateCSR(C c, CSR< T, AD > &dst, CSR< T, AS > const &src)
Copy CSR buffers between memory spaces asynchronously.
Definition AMReX_CSR.H:116
const int[]
Definition AMReX_BLProfiler.cpp:1664
Owning CSR container backed by AMReX resizable vectors.
Definition AMReX_CSR.H:49
V< Long > row_offset
Definition AMReX_CSR.H:52
Long nrows() const
Number of logical rows represented by the CSR offset array.
Definition AMReX_CSR.H:56
Long nnz
Definition AMReX_CSR.H:53
void sort()
Sort each row by column index. Uses GPU acceleration when possible.
Definition AMReX_CSR.H:137
CsrView< T > view()
Mutable view of the underlying buffers.
Definition AMReX_CSR.H:77
void sort_on_host()
Host-only fallback that sorts column indices row by row.
Definition AMReX_CSR.H:295
CsrView< T const > view() const
Const view of the underlying buffers.
Definition AMReX_CSR.H:83
CsrView< T const > const_view() const
Convenience alias for view() const.
Definition AMReX_CSR.H:89
void resize(Long num_rows, Long num_non_zeros)
Resize the storage to accommodate num_rows and num_non_zeros entries.
Definition AMReX_CSR.H:69
V< Long > col_index
Definition AMReX_CSR.H:51
V< T > mat
Definition AMReX_CSR.H:50
Lightweight non-owning CSR view that can point to host or device buffers.
Definition AMReX_CSR.H:33
std::conditional_t< std::is_const_v< T >, Long const, Long > U
Definition AMReX_CSR.H:34
T *__restrict__ mat
Definition AMReX_CSR.H:35
Long nrows
Definition AMReX_CSR.H:39
Long nnz
Definition AMReX_CSR.H:38
U *__restrict__ row_offset
Definition AMReX_CSR.H:37
U *__restrict__ col_index
Definition AMReX_CSR.H:36