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,
typename NT,
typename KERNEL,
typename... T>
22 const ctype *ang_quadrature_p,
const ctype *ang_quadrature_w,
23 const ctype *x0_quadrature_p,
const ctype *x0_quadrature_w,
24 const ctype x_extent,
const ctype x0_extent,
const uint x0_quadrature_size,
25 const uint x0_summands,
const ctype m_T,
const ctype k, T... t)
29 uint len_x = gridDim.x * blockDim.x;
30 uint len_y = gridDim.y * blockDim.y;
31 uint idx_x = (blockIdx.x * blockDim.x) + threadIdx.x;
32 uint idx_y = (blockIdx.y * blockDim.y) + threadIdx.y;
33 uint idx_z = (blockIdx.z * blockDim.z) + threadIdx.z;
34 uint idx = idx_z * len_x * len_y + idx_y * len_x + idx_x;
36 const ctype q = k * sqrt(x_quadrature_p[idx_x] * x_extent);
37 const ctype cos = 2 * (ang_quadrature_p[idx_y] - (ctype)0.5);
38 const ctype phi = 2 * (ctype)M_PI * ang_quadrature_p[idx_z];
40 2 * (ctype)M_PI * ang_quadrature_w[idx_z] * 2 * ang_quadrature_w[idx_y] * x_quadrature_w[idx_x] * x_extent;
45 const ctype int_element_int = (
powr<d - 3>(q) / (ctype)2 *
powr<2>(k))
49 const ctype integral_start = (2 * x0_summands * (ctype)M_PI * m_T) / k;
50 const ctype log_start = log(integral_start + (m_T == 0) * ctype(1e-3));
51 const ctype log_ext = log(x0_extent / (integral_start + (m_T == 0) * ctype(1e-3)));
53 for (
uint idx_0 = 0; idx_0 < x0_quadrature_size; ++idx_0) {
54 const ctype q0 = k * (exp(log_start + log_ext * x0_quadrature_p[idx_0]) - (m_T == 0) * ctype(1e-3));
56 const ctype m_weight = weight * (x0_quadrature_w[idx_0] * log_ext * q0 / k);
58 res += int_element_int * m_weight *
59 (KERNEL::kernel(q, cos, phi, q0, k, t...) + KERNEL::kernel(q, cos, phi, -q0, k, t...));
63 const ctype int_element_sum = m_T
66 for (
uint idx_0 = 0; idx_0 < x0_summands; ++idx_0) {
67 const ctype q0 = 2 * (ctype)M_PI * m_T * idx_0;
68 res += int_element_sum * weight *
69 (idx_0 == 0 ? KERNEL::kernel(q, cos, phi, (ctype)0, k, t...)
70 : KERNEL::kernel(q, cos, phi, q0, k, t...) + KERNEL::kernel(q, cos, phi, -q0, k, t...));
84 json.get_double(
"/physical/T"), json.get_uint(
"/integration/cudathreadsperblock"))
90 const uint max_block_size = 256)
92 grid_sizes({_grid_sizes[0], _grid_sizes[1], _grid_sizes[2], _grid_sizes[3] + _x0_summands}),
98 throw std::runtime_error(
"Grid sizes must be currently equal for all angular dimensions!");
107 block_sizes = {max_block_size, max_block_size, max_block_size};
110 uint optimize_dim = 2;
116 optimize_dim = (optimize_dim + 2) % 3;
157 template <
typename... T> NT
get(
const ctype k,
const T &...t)
const
166 return KERNEL::constant(k, t...) + thrust::reduce(thrust::cuda::par.on(cuda_stream.value()), device_data.begin(),
167 device_data.end(), NT(0.), thrust::plus<NT>());
170 template <
typename... T> std::future<NT>
request(
const ctype k,
const T &...t)
const
173 std::shared_ptr<rmm::device_uvector<NT>> device_data =
180 const NT constant = KERNEL::constant(k, t...);
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>());
212 using PoolMR = rmm::mr::pool_memory_resource<rmm::mr::device_memory_resource>;
224 template <
typename NT,
typename KERNEL>
class Integrator4DFiniteTx0GPU;
233 template <
typename NT,
typename KERNEL>
class Integrator4DFiniteTx0GPU :
public Integrator4DFiniteTx0TBB<NT, KERNEL>
240 const uint max_block_size = 256)
Definition integrator_4D_finiteTx0_gpu.hh:77
const ctype * ptr_x_quadrature_w
Definition integrator_4D_finiteTx0_gpu.hh:197
dim3 num_blocks
Definition integrator_4D_finiteTx0_gpu.hh:209
uint device_data_size
Definition integrator_4D_finiteTx0_gpu.hh:194
const rmm::cuda_stream_pool cuda_stream_pool
Definition integrator_4D_finiteTx0_gpu.hh:214
NT get(const ctype k, const T &...t) const
Definition integrator_4D_finiteTx0_gpu.hh:157
Integrator4DFiniteTx0GPU(const Integrator4DFiniteTx0GPU &other)
Definition integrator_4D_finiteTx0_gpu.hh:143
const ctype x_extent
Definition integrator_4D_finiteTx0_gpu.hh:203
ctype x0_extent
Definition integrator_4D_finiteTx0_gpu.hh:204
const ctype * ptr_x_quadrature_p
Definition integrator_4D_finiteTx0_gpu.hh:196
rmm::mr::pool_memory_resource< rmm::mr::device_memory_resource > PoolMR
Definition integrator_4D_finiteTx0_gpu.hh:212
Integrator4DFiniteTx0GPU(QuadratureProvider &quadrature_provider, const std::array< uint, 4 > _grid_sizes, const ctype x_extent, const ctype x0_extent, const uint _x0_summands, const ctype T, const uint max_block_size=256)
Definition integrator_4D_finiteTx0_gpu.hh:88
dim3 threads_per_block
Definition integrator_4D_finiteTx0_gpu.hh:210
typename get_type::ctype< NT > ctype
Definition integrator_4D_finiteTx0_gpu.hh:79
PoolMR pool
Definition integrator_4D_finiteTx0_gpu.hh:213
const ctype * ptr_ang_quadrature_p
Definition integrator_4D_finiteTx0_gpu.hh:198
std::future< NT > request(const ctype k, const T &...t) const
Definition integrator_4D_finiteTx0_gpu.hh:170
const uint original_x0_summands
Definition integrator_4D_finiteTx0_gpu.hh:205
uint x0_summands
Definition integrator_4D_finiteTx0_gpu.hh:206
const ctype * ptr_ang_quadrature_w
Definition integrator_4D_finiteTx0_gpu.hh:199
const ctype * ptr_x0_quadrature_w
Definition integrator_4D_finiteTx0_gpu.hh:201
std::array< uint, 3 > block_sizes
Definition integrator_4D_finiteTx0_gpu.hh:192
void set_x0_extent(const ctype val)
Definition integrator_4D_finiteTx0_gpu.hh:141
const std::array< uint, 4 > grid_sizes
Definition integrator_4D_finiteTx0_gpu.hh:191
void set_T(const ctype T)
Definition integrator_4D_finiteTx0_gpu.hh:130
QuadratureProvider & quadrature_provider
Definition integrator_4D_finiteTx0_gpu.hh:189
ctype m_T
Definition integrator_4D_finiteTx0_gpu.hh:207
Integrator4DFiniteTx0GPU(QuadratureProvider &quadrature_provider, const std::array< uint, 4 > grid_sizes, const ctype x_extent, const ctype x0_extent, const uint x0_summands, const JSONValue &json)
Definition integrator_4D_finiteTx0_gpu.hh:81
const ctype * ptr_x0_quadrature_p
Definition integrator_4D_finiteTx0_gpu.hh:200
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_4d_finiteTx0(NT *dest, const ctype *x_quadrature_p, const ctype *x_quadrature_w, const ctype *ang_quadrature_p, const ctype *ang_quadrature_w, const ctype *x0_quadrature_p, const ctype *x0_quadrature_w, const ctype x_extent, const ctype x0_extent, const uint x0_quadrature_size, const uint x0_summands, const ctype m_T, const ctype k, T... t)
Definition integrator_4D_finiteTx0_gpu.hh:21
bool __forceinline__ __host__ __device__ is_close(T1 a, T2 b, T3 eps_)
Function to evaluate whether two floats are equal to numerical precision. Tests for both relative and...
Definition math.hh:160
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