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>
20 template <
typename ctype,
int d,
typename NT,
typename KERNEL,
typename... T>
21 __global__
void gridreduce_angle(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)
25 uint len_x = gridDim.x * blockDim.x;
26 uint idx_x = (blockIdx.x * blockDim.x) + threadIdx.x;
27 uint idx_y = (blockIdx.y * blockDim.y) + threadIdx.y;
28 uint idx = idx_y * len_x + idx_x;
30 const ctype cos = 2 * (ang_quadrature_p[idx_y] - (ctype)0.5);
31 const ctype weight = 2 * ang_quadrature_w[idx_y] * x_quadrature_w[idx_x] * x_extent;
32 const ctype q = k * sqrt(x_quadrature_p[idx_x] * x_extent);
34 const ctype int_element =
S_d
39 dest[idx] = int_element * weight * KERNEL::kernel(q, cos, k, t...);
54 static_assert(d > 1,
"IntegratorAngleGPU: d must be greater than 1, otherwise angles are not needed");
65 json.get_uint(
"/integration/cudathreadsperblock"))
70 const uint max_block_size = 256)
82 uint optimize_dim = 1;
88 optimize_dim = (optimize_dim + 1) % 2;
122 template <
typename... T> NT
get(
const ctype k,
const T &...t)
const
130 return KERNEL::constant(k, t...) + thrust::reduce(thrust::cuda::par.on(cuda_stream.value()), device_data.begin(),
131 device_data.end(), NT(0.), thrust::plus<NT>());
144 template <
typename... T> std::future<NT>
request(
const ctype k,
const T &...t)
const
147 std::shared_ptr<rmm::device_uvector<NT>> device_data =
153 const NT constant = KERNEL::constant(k, t...);
155 return std::async(std::launch::deferred, [=,
this]() {
156 return constant + thrust::reduce(thrust::cuda::par.on(cuda_stream.value()), (*device_data).begin(),
157 (*device_data).end(), NT(0.), thrust::plus<NT>());
177 using PoolMR = rmm::mr::pool_memory_resource<rmm::mr::device_memory_resource>;
189 template <u
int d,
typename NT,
typename KERNEL>
class IntegratorAngleGPU;
198 template <u
int d,
typename NT,
typename KERNEL>
class IntegratorAngleGPU :
public IntegratorAngleTBB<d, NT, KERNEL>
204 const uint max_block_size = 256)
GPU integrator for the integration of a function with one angle with CUDA. Calculates.
Definition integrator_angle_gpu.hh:53
const ctype * ptr_ang_quadrature_w
Definition integrator_angle_gpu.hh:170
std::future< NT > request(const ctype k, const T &...t) const
Request a future for the integral of the kernel.
Definition integrator_angle_gpu.hh:144
rmm::mr::pool_memory_resource< rmm::mr::device_memory_resource > PoolMR
Definition integrator_angle_gpu.hh:177
NT get(const ctype k, const T &...t) const
Get the integral of the kernel.
Definition integrator_angle_gpu.hh:122
const ctype x_extent
Definition integrator_angle_gpu.hh:172
IntegratorAngleGPU(QuadratureProvider &quadrature_provider, std::array< uint, 2 > grid_sizes, const ctype x_extent, const uint max_block_size=256)
Definition integrator_angle_gpu.hh:69
const rmm::cuda_stream_pool cuda_stream_pool
Definition integrator_angle_gpu.hh:179
dim3 threads_per_block
Definition integrator_angle_gpu.hh:175
const uint device_data_size
Definition integrator_angle_gpu.hh:165
PoolMR pool
Definition integrator_angle_gpu.hh:178
const ctype * ptr_x_quadrature_w
Definition integrator_angle_gpu.hh:168
const ctype * ptr_ang_quadrature_p
Definition integrator_angle_gpu.hh:169
dim3 num_blocks
Definition integrator_angle_gpu.hh:174
const std::array< uint, 2 > grid_sizes
Definition integrator_angle_gpu.hh:162
const ctype * ptr_x_quadrature_p
Definition integrator_angle_gpu.hh:167
std::array< uint, 2 > block_sizes
Definition integrator_angle_gpu.hh:163
typename get_type::ctype< NT > ctype
Numerical type to be used for integration tasks e.g. the argument or possible jacobians.
Definition integrator_angle_gpu.hh:60
IntegratorAngleGPU(const IntegratorAngleGPU &other)
Definition integrator_angle_gpu.hh:100
IntegratorAngleGPU(QuadratureProvider &quadrature_provider, const std::array< uint, 2 > grid_sizes, const ctype x_extent, const JSONValue &json)
Definition integrator_angle_gpu.hh:62
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
constexpr __forceinline__ __host__ __device__ double S_d(NT d)
Surface of a d-dimensional sphere.
Definition math.hh:91
consteval NT S_d_prec(uint d)
Surface of a d-dimensional sphere (precompiled)
Definition math.hh:104
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_angle(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_angle_gpu.hh:21