DiFfRG
Loading...
Searching...
No Matches
integrator_3D_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_3d(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 uint len_x = gridDim.x * blockDim.x;
26 uint len_y = gridDim.y * blockDim.y;
27 uint idx_x = (blockIdx.x * blockDim.x) + threadIdx.x;
28 uint idx_y = (blockIdx.y * blockDim.y) + threadIdx.y;
29 uint idx_z = (blockIdx.z * blockDim.z) + threadIdx.z;
30 uint idx = idx_z * len_x * len_y + idx_y * len_x + idx_x;
31
32 const ctype cos = 2 * (ang_quadrature_p[idx_y] - (ctype)0.5);
33 const ctype phi = 2 * (ctype)M_PI * ang_quadrature_p[idx_z];
34 const ctype weight =
35 2 * (ctype)M_PI * ang_quadrature_w[idx_z] * 2 * ang_quadrature_w[idx_y] * x_quadrature_w[idx_x] * x_extent;
36 const ctype q = k * sqrt(x_quadrature_p[idx_x] * x_extent);
37 const ctype int_element = (powr<1>(q) / (ctype)2 * powr<2>(k)) // x = p^2 / k^2 integral
38 / powr<3>(2 * (ctype)M_PI); // fourier factor
39
40 NT res = KERNEL::kernel(q, cos, phi, k, t...);
41
42 dest[idx] = int_element * res * weight;
43 }
44
45 template <typename NT, typename KERNEL> class Integrator3DGPU
46 {
47 public:
51 using ctype = typename get_type::ctype<NT>;
52
53 Integrator3DGPU(QuadratureProvider &quadrature_provider, const std::array<uint, 3> grid_sizes, const ctype x_extent,
54 const JSONValue &json)
55 : Integrator3DGPU(quadrature_provider, grid_sizes, x_extent, json.get_uint("/integration/cudathreadsperblock"))
56 {
57 }
58
59 Integrator3DGPU(QuadratureProvider &quadrature_provider, std::array<uint, 3> grid_sizes, const ctype x_extent,
60 const uint max_block_size = 256)
62 pool(rmm::mr::get_current_device_resource(), (device_data_size / 256 + 1) * 256)
63 {
64 if (grid_sizes[2] != grid_sizes[1])
65 throw std::runtime_error("Grid sizes must be currently equal for all angular dimensions!");
66
67 ptr_x_quadrature_p = quadrature_provider.get_device_points<ctype>(grid_sizes[0]);
68 ptr_x_quadrature_w = quadrature_provider.get_device_weights<ctype>(grid_sizes[0]);
69 ptr_ang_quadrature_p = quadrature_provider.get_device_points<ctype>(grid_sizes[1]);
71
72 block_sizes = {max_block_size, max_block_size, max_block_size};
73 // choose block sizes such that the size is both as close to max_block_size as possible and the individual sizes
74 // are as close to each other as possible
75 uint optimize_dim = 2;
76 while (block_sizes[0] * block_sizes[1] * block_sizes[2] > max_block_size || block_sizes[0] > grid_sizes[0] ||
77 block_sizes[1] > grid_sizes[1] || block_sizes[2] > grid_sizes[2]) {
78 if (block_sizes[optimize_dim] > 1) block_sizes[optimize_dim]--;
79 while (grid_sizes[optimize_dim] % block_sizes[optimize_dim] != 0)
80 block_sizes[optimize_dim]--;
81 optimize_dim = (optimize_dim + 2) % 3;
82 }
83
84 uint blocks1 = grid_sizes[0] / block_sizes[0];
85 uint threads1 = block_sizes[0];
86 uint blocks2 = grid_sizes[1] / block_sizes[1];
87 uint threads2 = block_sizes[1];
88 uint blocks3 = grid_sizes[2] / block_sizes[2];
89 uint threads3 = block_sizes[2];
90
91 num_blocks = dim3(blocks1, blocks2, blocks3);
92 threads_per_block = dim3(threads1, threads2, threads3);
93 }
94
105
106 template <typename... T> NT get(const ctype k, const T &...t) const
107 {
108 const auto cuda_stream = cuda_stream_pool.get_stream();
109 rmm::device_uvector<NT> device_data(device_data_size, cuda_stream, &pool);
112 x_extent, k, t...);
113 check_cuda();
114 return KERNEL::constant(k, t...) + thrust::reduce(thrust::cuda::par.on(cuda_stream.value()), device_data.begin(),
115 device_data.end(), NT(0.), thrust::plus<NT>());
116 }
117
128 template <typename... T> std::future<NT> request(const ctype k, const T &...t) const
129 {
130 const auto cuda_stream = cuda_stream_pool.get_stream();
131 std::shared_ptr<rmm::device_uvector<NT>> device_data =
132 std::make_shared<rmm::device_uvector<NT>>(device_data_size, cuda_stream, &pool);
135 x_extent, k, t...);
136 check_cuda();
137 const NT constant = KERNEL::constant(k, t...);
138
139 return std::async(std::launch::deferred, [=, this]() {
140 return constant + thrust::reduce(thrust::cuda::par.on(cuda_stream.value()), (*device_data).begin(),
141 (*device_data).end(), NT(0.), thrust::plus<NT>());
142 });
143 }
144
145 private:
146 const std::array<uint, 3> grid_sizes;
147 std::array<uint, 3> block_sizes;
148
150
155
157
160
161 using PoolMR = rmm::mr::pool_memory_resource<rmm::mr::device_memory_resource>;
162 mutable PoolMR pool;
163 const rmm::cuda_stream_pool cuda_stream_pool;
164 };
165} // namespace DiFfRG
166
167#else
168
169#ifdef USE_CUDA
170
171#else
172
174
175namespace DiFfRG
176{
177 template <typename NT, typename KERNEL> class Integrator3DGPU : public Integrator3DTBB<NT, KERNEL>
178 {
179 public:
180 using ctype = typename get_type::ctype<NT>;
181
182 Integrator3DGPU(QuadratureProvider &quadrature_provider, std::array<uint, 3> grid_sizes, const ctype x_extent,
183 const uint max_block_size = 256)
184 : Integrator3DTBB<NT, KERNEL>(quadrature_provider, grid_sizes, x_extent)
185 {
186 }
187
188 Integrator3DGPU(QuadratureProvider &quadrature_provider, const std::array<uint, 3> grid_sizes, const ctype x_extent,
189 const JSONValue &)
190 : Integrator3DTBB<NT, KERNEL>(quadrature_provider, grid_sizes, x_extent)
191 {
192 }
193 };
194} // namespace DiFfRG
195
196#endif
197
198#endif
Definition integrator_3D_gpu.hh:46
dim3 num_blocks
Definition integrator_3D_gpu.hh:158
NT get(const ctype k, const T &...t) const
Definition integrator_3D_gpu.hh:106
const ctype * ptr_x_quadrature_w
Definition integrator_3D_gpu.hh:152
dim3 threads_per_block
Definition integrator_3D_gpu.hh:159
const ctype * ptr_x_quadrature_p
Definition integrator_3D_gpu.hh:151
Integrator3DGPU(QuadratureProvider &quadrature_provider, std::array< uint, 3 > grid_sizes, const ctype x_extent, const uint max_block_size=256)
Definition integrator_3D_gpu.hh:59
Integrator3DGPU(QuadratureProvider &quadrature_provider, const std::array< uint, 3 > grid_sizes, const ctype x_extent, const JSONValue &json)
Definition integrator_3D_gpu.hh:53
const rmm::cuda_stream_pool cuda_stream_pool
Definition integrator_3D_gpu.hh:163
const uint device_data_size
Definition integrator_3D_gpu.hh:149
const ctype * ptr_ang_quadrature_p
Definition integrator_3D_gpu.hh:153
std::future< NT > request(const ctype k, const T &...t) const
Request a future for the integral of the kernel.
Definition integrator_3D_gpu.hh:128
const ctype * ptr_ang_quadrature_w
Definition integrator_3D_gpu.hh:154
std::array< uint, 3 > block_sizes
Definition integrator_3D_gpu.hh:147
rmm::mr::pool_memory_resource< rmm::mr::device_memory_resource > PoolMR
Definition integrator_3D_gpu.hh:161
Integrator3DGPU(const Integrator3DGPU &other)
Definition integrator_3D_gpu.hh:95
const ctype x_extent
Definition integrator_3D_gpu.hh:156
PoolMR pool
Definition integrator_3D_gpu.hh:162
const std::array< uint, 3 > grid_sizes
Definition integrator_3D_gpu.hh:146
typename get_type::ctype< NT > ctype
Numerical type to be used for integration tasks e.g. the argument or possible jacobians.
Definition integrator_3D_gpu.hh:51
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
__global__ void gridreduce_3d(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_3D_gpu.hh:21
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