DiFfRG
Loading...
Searching...
No Matches
integrator_4D_gpu.hh
Go to the documentation of this file.
1#pragma once
2
3#ifdef __CUDACC__
4
5// standard library
6#include <future>
7
8// external libraries
9#include <rmm/cuda_stream_pool.hpp>
10#include <rmm/device_uvector.hpp>
11#include <rmm/mr/device/pool_memory_resource.hpp>
12#include <thrust/reduce.h>
13
14// DiFfRG
17
18namespace DiFfRG
19{
20 template <typename ctype, typename NT, typename KERNEL, typename... T>
21 __global__ void gridreduce_4d(NT *dest, const ctype *x_quadrature_p, const ctype *x_quadrature_w,
22 const ctype *ang_quadrature_p, const ctype *ang_quadrature_w, const ctype x_extent,
23 const ctype k, T... t)
24 {
25 const uint len_x = gridDim.x * blockDim.x;
26 const uint len_y = gridDim.y * blockDim.y;
27 const uint idx_x = (blockIdx.x * blockDim.x) + threadIdx.x;
28 const uint idx_y = (blockIdx.y * blockDim.y) + threadIdx.y;
29 const uint idx_z = (blockIdx.z * blockDim.z) + threadIdx.z;
30 const uint idx = idx_z * len_x * len_y + idx_y * len_x + idx_x;
31
32 const ctype q = k * sqrt(x_quadrature_p[idx_x] * x_extent);
33 const ctype cos1 = 2 * (ang_quadrature_p[idx_y] - (ctype)0.5);
34 const ctype phi = 2 * (ctype)M_PI * ang_quadrature_p[idx_z];
35 const ctype int_element = (powr<2>(q) * (ctype)0.5 * powr<2>(k)) // x = p^2 / k^2 integral
36 * sqrt(1. - powr<2>(cos1)) // cos1 integral jacobian
37 / powr<4>(2 * (ctype)M_PI); // fourier factor
38 const ctype weight = 2 * (ctype)M_PI * ang_quadrature_w[idx_z] // phi weight
39 * 2 * ang_quadrature_w[idx_y] // cos1 weight
40 * x_quadrature_w[idx_x] * x_extent; // x weight
41 NT mres = 0.;
42 for (uint idx_int = 0; idx_int < len_y; ++idx_int) {
43 const ctype cos2 = 2 * (ang_quadrature_p[idx_int] - (ctype)0.5);
44 mres += KERNEL::kernel(q, cos1, cos2, phi, k, t...) // kernel
45 * int_element * weight // other weights and integration elements
46 * 2 * ang_quadrature_w[idx_int]; // cos2 weight
47 }
48 dest[idx] = mres;
49 }
50
62 template <typename NT, typename KERNEL> class Integrator4DGPU
63 {
64 using PoolMR = rmm::mr::pool_memory_resource<rmm::mr::device_memory_resource>;
65
66 public:
70 using ctype = typename get_type::ctype<NT>;
71
73 const JSONValue &json)
74 : Integrator4DGPU(quadrature_provider, grid_sizes, x_extent, json.get_uint("/integration/cudathreadsperblock"))
75 {
76 }
77
79 const uint max_block_size = 256)
82 {
83 cudaGetDeviceCount(&n_devices);
84 if (n_devices == 0) throw std::runtime_error("No CUDA devices found!");
85
86 for (int device = 0; device < n_devices; ++device) {
87 const rmm::cuda_device_id device_id(device);
88 pool.emplace_back(
89 std::make_shared<PoolMR>(rmm::mr::get_per_device_resource(device_id), (device_data_size / 256 + 1) * 256));
90 }
91
92 if (grid_sizes[2] != grid_sizes[1] || grid_sizes[3] != grid_sizes[1])
93 throw std::runtime_error("Grid sizes must be currently equal for all angular dimensions!");
94
95 block_sizes = {max_block_size, max_block_size, max_block_size};
96 // choose block sizes such that the size is both as close to max_block_size as possible and the individual sizes
97 // are as close to each other as possible
98 uint optimize_dim = 2;
99 while (block_sizes[0] * block_sizes[1] * block_sizes[2] > max_block_size || block_sizes[0] > grid_sizes[0] ||
100 block_sizes[1] > grid_sizes[1] || block_sizes[2] > grid_sizes[2]) {
101 if (block_sizes[optimize_dim] > 1) block_sizes[optimize_dim]--;
102 while (grid_sizes[optimize_dim] % block_sizes[optimize_dim] != 0)
103 block_sizes[optimize_dim]--;
104 optimize_dim = (optimize_dim + 2) % 3;
105 }
106
107 uint blocks1 = grid_sizes[0] / block_sizes[0];
108 uint threads1 = block_sizes[0];
109 uint blocks2 = grid_sizes[1] / block_sizes[1];
110 uint threads2 = block_sizes[1];
111 uint blocks3 = grid_sizes[2] / block_sizes[2];
112 uint threads3 = block_sizes[2];
113
114 num_blocks = dim3(blocks1, blocks2, blocks3);
115 threads_per_block = dim3(threads1, threads2, threads3);
116
117 cudaSetDevice(0);
118 }
119
131
142 template <typename... T> NT get(const ctype k, const T &...t) const
143 {
144 const int m_device = 0; // evaluations++ % n_devices;
145 cudaSetDevice(0);
146
151
152 const auto cuda_stream = cuda_stream_pool.get_stream();
153 rmm::device_uvector<NT> device_data(device_data_size, cuda_stream, pool[m_device].get());
156 x_extent, k, t...);
157 check_cuda();
158 return KERNEL::constant(k, t...) + thrust::reduce(thrust::cuda::par.on(cuda_stream.value()), device_data.begin(),
159 device_data.end(), NT(0.), thrust::plus<NT>());
160
161 cudaSetDevice(0);
162 }
163
174 template <typename... T> std::future<NT> request(const ctype k, const T &...t) const
175 {
176 const int m_device = 0; // evaluations++ % n_devices;
177 cudaSetDevice(0);
178
183
184 const auto cuda_stream = cuda_stream_pool.get_stream();
185 std::shared_ptr<rmm::device_uvector<NT>> device_data =
186 std::make_shared<rmm::device_uvector<NT>>(device_data_size, cuda_stream, pool[m_device].get());
189 x_extent, k, t...);
190 check_cuda();
191 const NT constant = KERNEL::constant(k, t...);
192
193 return std::async(std::launch::deferred, [=, this]() {
194 cudaSetDevice(0);
195
196 return constant + thrust::reduce(thrust::cuda::par.on(cuda_stream.value()), (*device_data).begin(),
197 (*device_data).end(), NT(0.), thrust::plus<NT>());
198 });
199
200 cudaSetDevice(0);
201 }
202
203 private:
205 const std::array<uint, 4> grid_sizes;
206 std::array<uint, 3> block_sizes;
207
209
214
216
219
221 mutable std::vector<std::shared_ptr<PoolMR>> pool;
222 const rmm::cuda_stream_pool cuda_stream_pool;
223
224 mutable std::atomic_ullong evaluations;
225 };
226} // namespace DiFfRG
227
228#else
229
230#ifdef USE_CUDA
231
232namespace DiFfRG
233{
234 template <typename NT, typename KERNEL> class Integrator4DGPU;
235}
236
237#else
238
240
241namespace DiFfRG
242{
243 template <typename NT, typename KERNEL> class Integrator4DGPU : public Integrator4DTBB<NT, KERNEL>
244 {
245 public:
246 using ctype = typename get_type::ctype<NT>;
247
248 Integrator4DGPU(QuadratureProvider &quadrature_provider, std::array<uint, 4> grid_sizes, const ctype x_extent,
249 const uint max_block_size = 256)
250 : Integrator4DTBB<NT, KERNEL>(quadrature_provider, grid_sizes, x_extent)
251 {
252 }
253
254 Integrator4DGPU(QuadratureProvider &quadrature_provider, const std::array<uint, 4> grid_sizes, const ctype x_extent,
255 const JSONValue &)
256 : Integrator4DTBB<NT, KERNEL>(quadrature_provider, grid_sizes, x_extent)
257 {
258 }
259 };
260} // namespace DiFfRG
261
262#endif
263
264#endif
GPU integrator for the integration of a 4D function with three angles with CUDA. Calculates.
Definition integrator_4D_gpu.hh:63
rmm::mr::pool_memory_resource< rmm::mr::device_memory_resource > PoolMR
Definition integrator_4D_gpu.hh:64
const rmm::cuda_stream_pool cuda_stream_pool
Definition integrator_4D_gpu.hh:222
const ctype x_extent
Definition integrator_4D_gpu.hh:215
Integrator4DGPU(QuadratureProvider &quadrature_provider, std::array< uint, 4 > grid_sizes, const ctype x_extent, const uint max_block_size=256)
Definition integrator_4D_gpu.hh:78
Integrator4DGPU(QuadratureProvider &quadrature_provider, const std::array< uint, 4 > grid_sizes, const ctype x_extent, const JSONValue &json)
Definition integrator_4D_gpu.hh:72
dim3 threads_per_block
Definition integrator_4D_gpu.hh:218
const ctype * ptr_x_quadrature_w
Definition integrator_4D_gpu.hh:211
const ctype * ptr_ang_quadrature_w
Definition integrator_4D_gpu.hh:213
dim3 num_blocks
Definition integrator_4D_gpu.hh:217
std::atomic_ullong evaluations
Definition integrator_4D_gpu.hh:224
const std::array< uint, 4 > grid_sizes
Definition integrator_4D_gpu.hh:205
int n_devices
Definition integrator_4D_gpu.hh:220
typename get_type::ctype< NT > ctype
Numerical type to be used for integration tasks e.g. the argument or possible jacobians.
Definition integrator_4D_gpu.hh:70
const ctype * ptr_ang_quadrature_p
Definition integrator_4D_gpu.hh:212
Integrator4DGPU(const Integrator4DGPU &other)
Definition integrator_4D_gpu.hh:120
QuadratureProvider & quadrature_provider
Definition integrator_4D_gpu.hh:204
const ctype * ptr_x_quadrature_p
Definition integrator_4D_gpu.hh:210
std::future< NT > request(const ctype k, const T &...t) const
Request a future for the integral of the kernel.
Definition integrator_4D_gpu.hh:174
std::array< uint, 3 > block_sizes
Definition integrator_4D_gpu.hh:206
NT get(const ctype k, const T &...t) const
Get the integral of the kernel.
Definition integrator_4D_gpu.hh:142
const uint device_data_size
Definition integrator_4D_gpu.hh:208
std::vector< std::shared_ptr< PoolMR > > pool
Definition integrator_4D_gpu.hh:221
A wrapper around the boost json value class.
Definition json.hh:19
A class that provides quadrature points and weights, in host and device memory. The quadrature points...
Definition quadrature_provider.hh:139
const NT * get_device_weights(const size_t order, const int device=0, const QuadratureType type=QuadratureType::legendre)
Get the device-side quadrature weights for a quadrature of size quadrature_size.
Definition quadrature_provider.hh:211
const NT * get_device_points(const size_t order, const int device=0, const QuadratureType type=QuadratureType::legendre)
Get the device-side quadrature points for a quadrature of size quadrature_size.
Definition quadrature_provider.hh:198
typename internal::_ctype< CT >::value ctype
Definition types.hh:106
Definition complex_math.hh:14
constexpr __forceinline__ __host__ __device__ NumberType powr(const NumberType x)
A compile-time evaluatable power function for whole number exponents.
Definition math.hh:45
void check_cuda(std::string prefix="")
Check if a CUDA error occurred and print an error message if it did.
unsigned int uint
Definition utils.hh:22
__global__ void gridreduce_4d(NT *dest, const ctype *x_quadrature_p, const ctype *x_quadrature_w, const ctype *ang_quadrature_p, const ctype *ang_quadrature_w, const ctype x_extent, const ctype k, T... t)
Definition integrator_4D_gpu.hh:21