DiFfRG
Loading...
Searching...
No Matches
integrator_2Dpq0_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_2dpq0(NT *dest, const ctype *x_quadrature_p, const ctype *x_quadrature_w,
22 const ctype *ang_quadrature_p, const ctype *ang_quadrature_w,
23 const ctype *q0_quadrature_p, const ctype *q0_quadrature_w, const ctype x_extent,
24 const ctype q0_extent, const ctype k, T... t)
25 {
26 uint len_x = gridDim.x * blockDim.x;
27 uint len_y = gridDim.y * blockDim.y;
28 uint idx_x = (blockIdx.x * blockDim.x) + threadIdx.x;
29 uint idx_y = (blockIdx.y * blockDim.y) + threadIdx.y;
30 uint idx_z = (blockIdx.z * blockDim.z) + threadIdx.z;
31 uint idx = idx_z * len_x * len_y + idx_y * len_x + idx_x;
32
33 const ctype cos = 2 * (ang_quadrature_p[idx_y] - (ctype)0.5);
34 const ctype q0 = q0_quadrature_p[idx_z] * q0_extent;
35 const ctype weight =
36 q0_quadrature_w[idx_z] * q0_extent * 2 * ang_quadrature_w[idx_y] * x_quadrature_w[idx_x] * x_extent;
37 const ctype q = k * sqrt(x_quadrature_p[idx_x] * x_extent);
38 const ctype S_2d = 2 * (ctype)M_PI;
39 const ctype int_element = (1 / (ctype)2 * powr<2>(k)) // x = p^2 / k^2 integral
40 * (S_2d) * (1 / (ctype)2) // S_2 and divide the cos integral out
41 / powr<3>(2 * (ctype)M_PI); // fourier factor
42
43 NT res = KERNEL::kernel(q, cos, q0, k, t...) + KERNEL::kernel(q, cos, -q0, k, t...);
44
45 dest[idx] = int_element * res * weight;
46 }
47
58 template <typename NT, typename KERNEL> class Integrator2Dpq0GPU
59 {
60 public:
64 using ctype = typename get_type::ctype<NT>;
65
75 Integrator2Dpq0GPU(QuadratureProvider &quadrature_provider, const std::array<uint, 3> grid_sizes,
76 const ctype x_extent, const ctype q0_extent, const JSONValue &json)
77 : Integrator2Dpq0GPU(quadrature_provider, grid_sizes, x_extent, q0_extent,
78 json.get_uint("/integration/cudathreadsperblock"))
79 {
80 }
81
91 Integrator2Dpq0GPU(QuadratureProvider &quadrature_provider, std::array<uint, 3> grid_sizes, const ctype x_extent,
92 const ctype q0_extent, const uint max_block_size = 256)
94 q0_extent(q0_extent), pool(rmm::mr::get_current_device_resource(), (device_data_size / 256 + 1) * 256)
95 {
96 ptr_x_quadrature_p = quadrature_provider.get_device_points(grid_sizes[0]);
97 ptr_x_quadrature_w = quadrature_provider.get_device_weights(grid_sizes[0]);
98 ptr_ang_quadrature_p = quadrature_provider.get_device_points(grid_sizes[1]);
99 ptr_ang_quadrature_w = quadrature_provider.get_device_weights(grid_sizes[1]);
100 ptr_q0_quadrature_p = quadrature_provider.get_device_points(grid_sizes[2]);
101 ptr_q0_quadrature_w = quadrature_provider.get_device_weights(grid_sizes[2]);
102
103 block_sizes = {max_block_size, max_block_size, max_block_size};
104 // choose block sizes such that the size is both as close to max_block_size as possible and the individual sizes
105 // are as close to each other as possible
106 uint optimize_dim = 2;
107 while (block_sizes[0] * block_sizes[1] * block_sizes[2] > max_block_size || block_sizes[0] > grid_sizes[0] ||
108 block_sizes[1] > grid_sizes[1] || block_sizes[2] > grid_sizes[2]) {
109 if (block_sizes[optimize_dim] > 1) block_sizes[optimize_dim]--;
110 while (grid_sizes[optimize_dim] % block_sizes[optimize_dim] != 0)
111 block_sizes[optimize_dim]--;
112 optimize_dim = (optimize_dim + 2) % 3;
113 }
114
115 uint blocks1 = grid_sizes[0] / block_sizes[0];
116 uint threads1 = block_sizes[0];
117 uint blocks2 = grid_sizes[1] / block_sizes[1];
118 uint threads2 = block_sizes[1];
119 uint blocks3 = grid_sizes[2] / block_sizes[2];
120 uint threads3 = block_sizes[2];
121
122 num_blocks = dim3(blocks1, blocks2, blocks3);
123 threads_per_block = dim3(threads1, threads2, threads3);
124 }
125
138
149 template <typename... T> NT get(const ctype k, const T &...t) const
150 {
151 const auto cuda_stream = cuda_stream_pool.get_stream();
152 rmm::device_uvector<NT> device_data(device_data_size, cuda_stream, &pool);
156 check_cuda();
157 return KERNEL::constant(k, t...) + thrust::reduce(thrust::cuda::par.on(cuda_stream.value()), device_data.begin(),
158 device_data.end(), NT(0.), thrust::plus<NT>());
159 }
160
171 template <typename... T> std::future<NT> request(const ctype k, const T &...t) const
172 {
173 const auto cuda_stream = cuda_stream_pool.get_stream();
174 std::shared_ptr<rmm::device_uvector<NT>> device_data =
175 std::make_shared<rmm::device_uvector<NT>>(device_data_size, cuda_stream, &pool);
179 check_cuda();
180 const NT constant = KERNEL::constant(k, t...);
181
182 return std::async(std::launch::deferred, [=, this]() {
183 return constant + thrust::reduce(thrust::cuda::par.on(cuda_stream.value()), (*device_data).begin(),
184 (*device_data).end(), NT(0.), thrust::plus<NT>());
185 });
186 }
187
188 private:
189 const std::array<uint, 3> grid_sizes;
190 std::array<uint, 3> block_sizes;
191
193
200
203
206
207 using PoolMR = rmm::mr::pool_memory_resource<rmm::mr::device_memory_resource>;
208 mutable PoolMR pool;
209 const rmm::cuda_stream_pool cuda_stream_pool;
210 };
211} // namespace DiFfRG
212
213#else
214
215#ifdef USE_CUDA
216
217#else
218
220
221namespace DiFfRG
222{
223 template <typename NT, typename KERNEL> class Integrator2Dpq0GPU : public Integrator2Dpq0TBB<NT, KERNEL>
224 {
225 public:
226 using ctype = typename get_type::ctype<NT>;
227
228 Integrator2Dpq0GPU(QuadratureProvider &quadrature_provider, std::array<uint, 3> grid_sizes, const ctype x_extent,
229 const ctype q0_extent, const uint max_block_size = 256)
230 : Integrator2Dpq0TBB<NT, KERNEL>(quadrature_provider, grid_sizes, x_extent, q0_extent)
231 {
232 }
233
234 Integrator2Dpq0GPU(QuadratureProvider &quadrature_provider, const std::array<uint, 3> grid_sizes,
235 const ctype x_extent, const ctype q0_extent, const JSONValue &)
236 : Integrator2Dpq0TBB<NT, KERNEL>(quadrature_provider, grid_sizes, x_extent, q0_extent)
237 {
238 }
239 };
240} // namespace DiFfRG
241
242#endif
243
244#endif
Integrator for 2+1D integrals over p, q0 and an angle on the GPU. Calculates.
Definition integrator_2Dpq0_gpu.hh:59
const ctype * ptr_ang_quadrature_w
Definition integrator_2Dpq0_gpu.hh:197
NT get(const ctype k, const T &...t) const
Get the integral of the kernel.
Definition integrator_2Dpq0_gpu.hh:149
const ctype x_extent
Definition integrator_2Dpq0_gpu.hh:201
PoolMR pool
Definition integrator_2Dpq0_gpu.hh:208
const ctype q0_extent
Definition integrator_2Dpq0_gpu.hh:202
Integrator2Dpq0GPU(const Integrator2Dpq0GPU &other)
Definition integrator_2Dpq0_gpu.hh:126
Integrator2Dpq0GPU(QuadratureProvider &quadrature_provider, std::array< uint, 3 > grid_sizes, const ctype x_extent, const ctype q0_extent, const uint max_block_size=256)
Construct a new Integrator2Dpq0TBB object.
Definition integrator_2Dpq0_gpu.hh:91
std::future< NT > request(const ctype k, const T &...t) const
Request a future for the integral of the kernel.
Definition integrator_2Dpq0_gpu.hh:171
dim3 threads_per_block
Definition integrator_2Dpq0_gpu.hh:205
const uint device_data_size
Definition integrator_2Dpq0_gpu.hh:192
rmm::mr::pool_memory_resource< rmm::mr::device_memory_resource > PoolMR
Definition integrator_2Dpq0_gpu.hh:207
std::array< uint, 3 > block_sizes
Definition integrator_2Dpq0_gpu.hh:190
Integrator2Dpq0GPU(QuadratureProvider &quadrature_provider, const std::array< uint, 3 > grid_sizes, const ctype x_extent, const ctype q0_extent, const JSONValue &json)
Construct a new Integrator2Dpq0TBB object.
Definition integrator_2Dpq0_gpu.hh:75
const rmm::cuda_stream_pool cuda_stream_pool
Definition integrator_2Dpq0_gpu.hh:209
const ctype * ptr_x_quadrature_w
Definition integrator_2Dpq0_gpu.hh:195
dim3 num_blocks
Definition integrator_2Dpq0_gpu.hh:204
const std::array< uint, 3 > grid_sizes
Definition integrator_2Dpq0_gpu.hh:189
const ctype * ptr_x_quadrature_p
Definition integrator_2Dpq0_gpu.hh:194
const ctype * ptr_q0_quadrature_w
Definition integrator_2Dpq0_gpu.hh:199
const ctype * ptr_ang_quadrature_p
Definition integrator_2Dpq0_gpu.hh:196
typename get_type::ctype< NT > ctype
Numerical type to be used for integration tasks e.g. the argument or possible jacobians.
Definition integrator_2Dpq0_gpu.hh:64
const ctype * ptr_q0_quadrature_p
Definition integrator_2Dpq0_gpu.hh:198
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_2dpq0(NT *dest, const ctype *x_quadrature_p, const ctype *x_quadrature_w, const ctype *ang_quadrature_p, const ctype *ang_quadrature_w, const ctype *q0_quadrature_p, const ctype *q0_quadrature_w, const ctype x_extent, const ctype q0_extent, const ctype k, T... t)
Definition integrator_2Dpq0_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