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>
27 template <
typename ctype,
typename NT,
typename KERNEL,
typename... T>
29 const ctype *y_quadrature_p,
const ctype *y_quadrature_w,
const ctype qx_min,
30 const ctype qy_min,
const ctype qx_extent,
const ctype qy_extent,
31 const ctype k, T... t)
33 uint len_x = gridDim.x * blockDim.x;
34 uint idx_x = (blockIdx.x * blockDim.x) + threadIdx.x;
35 uint idx_y = (blockIdx.y * blockDim.y) + threadIdx.y;
36 uint idx = idx_y * len_x + idx_x;
39 constexpr ctype int_element =
powr<-d>(2 * (ctype)M_PI);
41 const ctype qx = qx_min + qx_extent * x_quadrature_p[idx_x];
42 const ctype qy = qy_min + qy_extent * y_quadrature_p[idx_y];
43 const ctype weight = qy_extent * y_quadrature_w[idx_y] * qx_extent * x_quadrature_w[idx_x];
45 dest[idx] = int_element * weight * KERNEL::kernel(qx, qy, k, t...);
73 json.get_uint(
"/integration/cudathreadsperblock"))
102 uint optimize_dim = 1;
108 optimize_dim = (optimize_dim + 1) % 2;
182 template <
typename... T> NT
get(
const ctype k,
const T &...t)
const
190 return KERNEL::constant(k, t...) + thrust::reduce(thrust::cuda::par.on(cuda_stream.value()), device_data.begin(),
191 device_data.end(), NT(0.), thrust::plus<NT>());
201 template <
typename... T> std::future<NT>
request(
const ctype k,
const T &...t)
const
204 std::shared_ptr<rmm::device_uvector<NT>> device_data =
210 const NT constant = KERNEL::constant(k, t...);
212 return std::async(std::launch::deferred, [=,
this]() {
213 return constant + thrust::reduce(thrust::cuda::par.on(cuda_stream.value()), (*device_data).begin(),
214 (*device_data).end(), NT(0.), thrust::plus<NT>());
237 using PoolMR = rmm::mr::pool_memory_resource<rmm::mr::device_memory_resource>;
249 template <
typename NT,
typename KERNEL>
class Integrator2DCartesianGPU;
258 template <u
typename NT,
typename KERNEL>
class Integrator2DCartesianGPU :
public Integrator2DCartesianTBB<NT, KERNEL>
266 : Integrator2DCartesianTBB<NT, KERNEL>(quadrature_provider,
grid_sizes, x_extent, max_block_size,
qx_min,
272 const ctype x_extent,
const JSONValue &json)
273 : Integrator2DCartesianTBB<NT, KERNEL>(quadrature_provider,
grid_sizes, x_extent, json)
Integration of an arbitrary 2D function from (qx_min, qy_min) to (qx_max, qy_max) using TBB.
Definition integrator_2D_cartesian_gpu.hh:55
const std::array< uint, 2 > grid_sizes
Definition integrator_2D_cartesian_gpu.hh:219
void set_qy_max(const ctype qy_max)
Set the maximum value of the qy integration range.
Definition integrator_2D_cartesian_gpu.hh:173
const ctype * ptr_y_quadrature_w
Definition integrator_2D_cartesian_gpu.hh:232
ctype qx_min
Definition integrator_2D_cartesian_gpu.hh:224
dim3 num_blocks
Definition integrator_2D_cartesian_gpu.hh:234
const uint device_data_size
Definition integrator_2D_cartesian_gpu.hh:222
void set_qx_min(const ctype qx_min)
Set the minimum value of the qx integration range.
Definition integrator_2D_cartesian_gpu.hh:150
const ctype * ptr_x_quadrature_p
Definition integrator_2D_cartesian_gpu.hh:229
rmm::mr::pool_memory_resource< rmm::mr::device_memory_resource > PoolMR
Definition integrator_2D_cartesian_gpu.hh:237
Integrator2DCartesianGPU(QuadratureProvider &quadrature_provider, const std::array< uint, 2 > grid_sizes, const ctype x_extent, const JSONValue &json)
Construct a new Integrator2DCartesianGPU object.
Definition integrator_2D_cartesian_gpu.hh:70
const ctype * ptr_x_quadrature_w
Definition integrator_2D_cartesian_gpu.hh:230
Integrator2DCartesianGPU(QuadratureProvider &quadrature_provider, std::array< uint, 2 > grid_sizes, const ctype x_extent=0., const uint max_block_size=256, const ctype qx_min=-M_PI, const ctype qy_min=-M_PI, const ctype qx_max=M_PI, const ctype qy_max=M_PI)
Construct a new Integrator2DCartesianGPU object.
Definition integrator_2D_cartesian_gpu.hh:88
const rmm::cuda_stream_pool cuda_stream_pool
Definition integrator_2D_cartesian_gpu.hh:239
typename get_type::ctype< NT > ctype
Numerical type to be used for integration tasks e.g. the argument or possible jacobians.
Definition integrator_2D_cartesian_gpu.hh:60
Integrator2DCartesianGPU(const Integrator2DCartesianGPU &other)
Copy a Integrator2DCartesianGPU object.
Definition integrator_2D_cartesian_gpu.hh:130
dim3 threads_per_block
Definition integrator_2D_cartesian_gpu.hh:235
ctype qy_min
Definition integrator_2D_cartesian_gpu.hh:225
const ctype * ptr_y_quadrature_p
Definition integrator_2D_cartesian_gpu.hh:231
std::future< NT > request(const ctype k, const T &...t) const
Get the result of the integration asynchronously.
Definition integrator_2D_cartesian_gpu.hh:201
NT get(const ctype k, const T &...t) const
Get the result of the integration.
Definition integrator_2D_cartesian_gpu.hh:182
ctype qy_extent
Definition integrator_2D_cartesian_gpu.hh:227
void set_qx_max(const ctype qx_max)
Set the maximum value of the qx integration range.
Definition integrator_2D_cartesian_gpu.hh:168
PoolMR pool
Definition integrator_2D_cartesian_gpu.hh:238
std::array< uint, 2 > block_sizes
Definition integrator_2D_cartesian_gpu.hh:220
void set_qy_min(const ctype qy_min)
Set the minimum value of the qy integration range.
Definition integrator_2D_cartesian_gpu.hh:159
ctype qx_extent
Definition integrator_2D_cartesian_gpu.hh:226
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.
__global__ void gridreduce_2d_cartesian(NT *dest, const ctype *x_quadrature_p, const ctype *x_quadrature_w, const ctype *y_quadrature_p, const ctype *y_quadrature_w, const ctype qx_min, const ctype qy_min, const ctype qx_extent, const ctype qy_extent, const ctype k, T... t)
GPU kernel for the integration of an arbitrary 2D function from qx_min to qx_max and qy_min to qy_max...
Definition integrator_2D_cartesian_gpu.hh:28
unsigned int uint
Definition utils.hh:22