DiFfRG
Loading...
Searching...
No Matches
integrator_3Dpq0_cpu.hh
Go to the documentation of this file.
1#pragma once
2
3// standard library
4#include <future>
5
6// external libraries
7#include <tbb/tbb.h>
8
9// DiFfRG
11
12namespace DiFfRG
13{
14 template <typename NT, typename KERNEL> class Integrator3Dpq0TBB
15 {
16 public:
20 using ctype = typename get_type::ctype<NT>;
21
22 Integrator3Dpq0TBB(QuadratureProvider &quadrature_provider, const std::array<uint, 4> grid_sizes,
23 const ctype x_extent, const ctype q0_extent, const JSONValue &)
24 : Integrator3Dpq0TBB(quadrature_provider, grid_sizes, x_extent, q0_extent)
25 {
26 }
27
28 Integrator3Dpq0TBB(QuadratureProvider &quadrature_provider, std::array<uint, 4> grid_sizes, const ctype x_extent,
29 const ctype q0_extent, const uint max_block_size = 0)
31 x_quadrature_p(quadrature_provider.get_points<ctype>(grid_sizes[0])),
32 x_quadrature_w(quadrature_provider.get_weights<ctype>(grid_sizes[0])),
33 ang_quadrature_p(quadrature_provider.get_points<ctype>(grid_sizes[1])),
34 ang_quadrature_w(quadrature_provider.get_weights<ctype>(grid_sizes[1])),
35 q0_quadrature_p(quadrature_provider.get_points<ctype>(grid_sizes[2])),
36 q0_quadrature_w(quadrature_provider.get_weights<ctype>(grid_sizes[2]))
37 {
38 (void)max_block_size;
39 if (grid_sizes[2] != grid_sizes[1])
40 throw std::runtime_error("Grid sizes must be currently equal for all angular dimensions!");
41 }
42
43 template <typename... T> NT get(const ctype k, const T &...t) const
44 {
45 using std::sqrt;
46
47 const auto constant = KERNEL::constant(k, t...);
48 return constant +
49 tbb::parallel_reduce(
50 tbb::blocked_range3d<uint, uint, uint>(0, grid_sizes[0], 0, grid_sizes[1], 0,
51 grid_sizes[2] * grid_sizes[3]),
52 NT(0.),
53 [&](const tbb::blocked_range3d<uint, uint, uint> &r, NT value) -> NT {
54 for (uint idx_x = r.pages().begin(); idx_x != r.pages().end(); ++idx_x) {
55 const ctype q = k * sqrt(x_quadrature_p[idx_x] * x_extent);
56 for (uint idx_y = r.rows().begin(); idx_y != r.rows().end(); ++idx_y) {
57 const ctype cos = 2 * (ang_quadrature_p[idx_y] - (ctype)0.5);
58 const ctype int_element = (powr<1>(q) / (ctype)2 * powr<2>(k)) // x = p^2 / k^2 integral
59 / powr<4>(2 * (ctype)M_PI); // fourier factor
60 for (uint idx_z = r.cols().begin(); idx_z != r.cols().end(); ++idx_z) {
61 const ctype phi = 2 * (ctype)M_PI * ang_quadrature_p[idx_z / grid_sizes[1]];
62 const ctype q0 = q0_quadrature_p[idx_z % grid_sizes[1]] * q0_extent;
63 const ctype weight = q0_quadrature_w[idx_z % grid_sizes[1]] * q0_extent * 2 * (ctype)M_PI *
64 ang_quadrature_w[idx_z / grid_sizes[1]] * 2 * ang_quadrature_w[idx_y] *
65 x_quadrature_w[idx_x] * x_extent;
66 value +=
67 int_element * weight *
68 (KERNEL::kernel(q, cos, phi, q0, k, t...) + KERNEL::kernel(q, cos, phi, -q0, k, t...));
69 }
70 }
71 }
72 return value;
73 },
74 [&](NT x, NT y) -> NT { return x + y; });
75 }
76
87 template <typename... T> std::future<NT> request(const ctype k, const T &...t) const
88 {
89 return std::async(std::launch::deferred, [=, this]() { return get(k, t...); });
90 }
91
92 private:
93 const std::array<uint, 4> grid_sizes;
94
97
98 const std::vector<ctype> &x_quadrature_p;
99 const std::vector<ctype> &x_quadrature_w;
100 const std::vector<ctype> &ang_quadrature_p;
101 const std::vector<ctype> &ang_quadrature_w;
102 const std::vector<ctype> &q0_quadrature_p;
103 const std::vector<ctype> &q0_quadrature_w;
104 };
105} // namespace DiFfRG
Definition integrator_3Dpq0_cpu.hh:15
NT get(const ctype k, const T &...t) const
Definition integrator_3Dpq0_cpu.hh:43
Integrator3Dpq0TBB(QuadratureProvider &quadrature_provider, std::array< uint, 4 > grid_sizes, const ctype x_extent, const ctype q0_extent, const uint max_block_size=0)
Definition integrator_3Dpq0_cpu.hh:28
typename get_type::ctype< NT > ctype
Numerical type to be used for integration tasks e.g. the argument or possible jacobians.
Definition integrator_3Dpq0_cpu.hh:20
const std::vector< ctype > & x_quadrature_w
Definition integrator_3Dpq0_cpu.hh:99
const ctype q0_extent
Definition integrator_3Dpq0_cpu.hh:96
const ctype x_extent
Definition integrator_3Dpq0_cpu.hh:95
std::future< NT > request(const ctype k, const T &...t) const
Request a future for the integral of the kernel.
Definition integrator_3Dpq0_cpu.hh:87
const std::vector< ctype > & x_quadrature_p
Definition integrator_3Dpq0_cpu.hh:98
Integrator3Dpq0TBB(QuadratureProvider &quadrature_provider, const std::array< uint, 4 > grid_sizes, const ctype x_extent, const ctype q0_extent, const JSONValue &)
Definition integrator_3Dpq0_cpu.hh:22
const std::vector< ctype > & ang_quadrature_p
Definition integrator_3Dpq0_cpu.hh:100
const std::vector< ctype > & q0_quadrature_p
Definition integrator_3Dpq0_cpu.hh:102
const std::array< uint, 4 > grid_sizes
Definition integrator_3Dpq0_cpu.hh:93
const std::vector< ctype > & q0_quadrature_w
Definition integrator_3Dpq0_cpu.hh:103
const std::vector< ctype > & ang_quadrature_w
Definition integrator_3Dpq0_cpu.hh:101
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
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
unsigned int uint
Definition utils.hh:22