DiFfRG
Loading...
Searching...
No Matches
integrator_3Dpq0_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_3dpq0(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 phi = 2 * (ctype)M_PI * ang_quadrature_p[idx_z / len_y];
35 const ctype q0 = q0_quadrature_p[idx_z % len_y] * q0_extent;
36 const ctype weight = q0_quadrature_w[idx_z % len_y] * q0_extent * 2 * (ctype)M_PI *
37 ang_quadrature_w[idx_z / len_y] * 2 * ang_quadrature_w[idx_y] * x_quadrature_w[idx_x] *
38 x_extent;
39 const ctype q = k * sqrt(x_quadrature_p[idx_x] * x_extent);
40 const ctype int_element = (powr<1>(q) / (ctype)2 * powr<2>(k)) // x = p^2 / k^2 integral
41 / powr<4>(2 * (ctype)M_PI); // fourier factor
42
43 NT res = KERNEL::kernel(q, cos, phi, q0, k, t...) + KERNEL::kernel(q, cos, phi, -q0, k, t...);
44
45 dest[idx] = int_element * res * weight;
46 }
47
48 template <typename NT, typename KERNEL> class Integrator3Dpq0GPU
49 {
50 public:
54 using ctype = typename get_type::ctype<NT>;
55
56 Integrator3Dpq0GPU(QuadratureProvider &quadrature_provider, const std::array<uint, 4> grid_sizes,
57 const ctype x_extent, const ctype q0_extent, const JSONValue &json)
58 : Integrator3Dpq0GPU(quadrature_provider, grid_sizes, x_extent, q0_extent,
59 json.get_uint("/integration/cudathreadsperblock"))
60 {
61 }
62
63 Integrator3Dpq0GPU(QuadratureProvider &quadrature_provider, std::array<uint, 4> grid_sizes, const ctype x_extent,
64 const ctype q0_extent, const uint max_block_size = 256)
67 pool(rmm::mr::get_current_device_resource(), (device_data_size / 256 + 1) * 256)
68 {
69 if (grid_sizes[2] != grid_sizes[1])
70 throw std::runtime_error("Grid sizes must be currently equal for all angular dimensions!");
71
72 ptr_x_quadrature_p = quadrature_provider.get_device_points<ctype>(grid_sizes[0]);
73 ptr_x_quadrature_w = quadrature_provider.get_device_weights<ctype>(grid_sizes[0]);
74 ptr_ang_quadrature_p = quadrature_provider.get_device_points<ctype>(grid_sizes[1]);
76 ptr_q0_quadrature_p = quadrature_provider.get_device_points<ctype>(grid_sizes[3]);
77 ptr_q0_quadrature_w = quadrature_provider.get_device_weights<ctype>(grid_sizes[3]);
78
79 block_sizes = {max_block_size, max_block_size, max_block_size};
80 // choose block sizes such that the size is both as close to max_block_size as possible and the individual sizes
81 // are as close to each other as possible
82 uint optimize_dim = 2;
83 while (block_sizes[0] * block_sizes[1] * block_sizes[2] > max_block_size || block_sizes[0] > grid_sizes[0] ||
84 block_sizes[1] > grid_sizes[1] || block_sizes[2] > grid_sizes[2] * grid_sizes[3]) {
85 if (block_sizes[optimize_dim] > 1) block_sizes[optimize_dim]--;
86 while (grid_sizes[optimize_dim] % block_sizes[optimize_dim] != 0)
87 block_sizes[optimize_dim]--;
88 optimize_dim = (optimize_dim + 2) % 3;
89 }
90
91 uint blocks1 = grid_sizes[0] / block_sizes[0];
92 uint threads1 = block_sizes[0];
93 uint blocks2 = grid_sizes[1] / block_sizes[1];
94 uint threads2 = block_sizes[1];
95 uint blocks3 = grid_sizes[2] * grid_sizes[3] / block_sizes[2];
96 uint threads3 = block_sizes[2];
97
98 num_blocks = dim3(blocks1, blocks2, blocks3);
99 threads_per_block = dim3(threads1, threads2, threads3);
100 }
101
114
115 template <typename... T> NT get(const ctype k, const T &...t) const
116 {
117 const auto cuda_stream = cuda_stream_pool.get_stream();
118 rmm::device_uvector<NT> device_data(device_data_size, cuda_stream, &pool);
122 check_cuda();
123 return KERNEL::constant(k, t...) + thrust::reduce(thrust::cuda::par.on(cuda_stream.value()), device_data.begin(),
124 device_data.end(), NT(0.), thrust::plus<NT>());
125 }
126
137 template <typename... T> std::future<NT> request(const ctype k, const T &...t) const
138 {
139 const auto cuda_stream = cuda_stream_pool.get_stream();
140 std::shared_ptr<rmm::device_uvector<NT>> device_data =
141 std::make_shared<rmm::device_uvector<NT>>(device_data_size, cuda_stream, &pool);
145 check_cuda();
146 const NT constant = KERNEL::constant(k, t...);
147
148 return std::async(std::launch::deferred, [=, this]() {
149 return constant + thrust::reduce(thrust::cuda::par.on(cuda_stream.value()), (*device_data).begin(),
150 (*device_data).end(), NT(0.), thrust::plus<NT>());
151 });
152 }
153
154 private:
155 const std::array<uint, 4> grid_sizes;
156 std::array<uint, 3> block_sizes;
157
159
166
169
172
173 using PoolMR = rmm::mr::pool_memory_resource<rmm::mr::device_memory_resource>;
174 mutable PoolMR pool;
175 const rmm::cuda_stream_pool cuda_stream_pool;
176 };
177} // namespace DiFfRG
178
179#else
180
181#ifdef USE_CUDA
182
183#else
184
186
187namespace DiFfRG
188{
189 template <typename NT, typename KERNEL> class Integrator3Dpq0GPU : public Integrator3Dpq0TBB<NT, KERNEL>
190 {
191 public:
192 using ctype = typename get_type::ctype<NT>;
193
194 Integrator3Dpq0GPU(QuadratureProvider &quadrature_provider, std::array<uint, 4> grid_sizes, const ctype x_extent,
195 const ctype q0_extent, const uint max_block_size = 256)
196 : Integrator3Dpq0TBB<NT, KERNEL>(quadrature_provider, grid_sizes, x_extent, q0_extent)
197 {
198 }
199
200 Integrator3Dpq0GPU(QuadratureProvider &quadrature_provider, const std::array<uint, 4> grid_sizes,
201 const ctype x_extent, const ctype q0_extent, const JSONValue &)
202 : Integrator3Dpq0TBB<NT, KERNEL>(quadrature_provider, grid_sizes, x_extent, q0_extent)
203 {
204 }
205 };
206} // namespace DiFfRG
207
208#endif
209
210#endif
Definition integrator_3Dpq0_gpu.hh:49
PoolMR pool
Definition integrator_3Dpq0_gpu.hh:174
const std::array< uint, 4 > grid_sizes
Definition integrator_3Dpq0_gpu.hh:155
const rmm::cuda_stream_pool cuda_stream_pool
Definition integrator_3Dpq0_gpu.hh:175
const ctype * ptr_x_quadrature_w
Definition integrator_3Dpq0_gpu.hh:161
rmm::mr::pool_memory_resource< rmm::mr::device_memory_resource > PoolMR
Definition integrator_3Dpq0_gpu.hh:173
const uint device_data_size
Definition integrator_3Dpq0_gpu.hh:158
std::array< uint, 3 > block_sizes
Definition integrator_3Dpq0_gpu.hh:156
Integrator3Dpq0GPU(const Integrator3Dpq0GPU &other)
Definition integrator_3Dpq0_gpu.hh:102
dim3 num_blocks
Definition integrator_3Dpq0_gpu.hh:170
typename get_type::ctype< NT > ctype
Numerical type to be used for integration tasks e.g. the argument or possible jacobians.
Definition integrator_3Dpq0_gpu.hh:54
std::future< NT > request(const ctype k, const T &...t) const
Request a future for the integral of the kernel.
Definition integrator_3Dpq0_gpu.hh:137
Integrator3Dpq0GPU(QuadratureProvider &quadrature_provider, std::array< uint, 4 > grid_sizes, const ctype x_extent, const ctype q0_extent, const uint max_block_size=256)
Definition integrator_3Dpq0_gpu.hh:63
const ctype * ptr_x_quadrature_p
Definition integrator_3Dpq0_gpu.hh:160
const ctype * ptr_ang_quadrature_p
Definition integrator_3Dpq0_gpu.hh:162
NT get(const ctype k, const T &...t) const
Definition integrator_3Dpq0_gpu.hh:115
const ctype x_extent
Definition integrator_3Dpq0_gpu.hh:167
dim3 threads_per_block
Definition integrator_3Dpq0_gpu.hh:171
const ctype * ptr_ang_quadrature_w
Definition integrator_3Dpq0_gpu.hh:163
const ctype q0_extent
Definition integrator_3Dpq0_gpu.hh:168
const ctype * ptr_q0_quadrature_w
Definition integrator_3Dpq0_gpu.hh:165
const ctype * ptr_q0_quadrature_p
Definition integrator_3Dpq0_gpu.hh:164
Integrator3Dpq0GPU(QuadratureProvider &quadrature_provider, const std::array< uint, 4 > grid_sizes, const ctype x_extent, const ctype q0_extent, const JSONValue &json)
Definition integrator_3Dpq0_gpu.hh:56
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_3dpq0(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_3Dpq0_gpu.hh:21