DiFfRG
Loading...
Searching...
No Matches
integrator_4D_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{
25 template <typename NT, typename KERNEL> class Integrator4DTBB
26 {
27 public:
31 using ctype = typename get_type::ctype<NT>;
32
33 Integrator4DTBB(QuadratureProvider &quadrature_provider, const std::array<uint, 4> grid_sizes, const ctype x_extent,
34 const JSONValue &)
35 : Integrator4DTBB(quadrature_provider, grid_sizes, x_extent)
36 {
37 }
38
39 Integrator4DTBB(QuadratureProvider &quadrature_provider, std::array<uint, 4> grid_sizes, const ctype x_extent,
40 const uint max_block_size = 0)
42 x_quadrature_p(quadrature_provider.get_points<ctype>(grid_sizes[0])),
43 x_quadrature_w(quadrature_provider.get_weights<ctype>(grid_sizes[0])),
44 ang_quadrature_p1(quadrature_provider.get_points<ctype>(grid_sizes[1])),
45 ang_quadrature_w1(quadrature_provider.get_weights<ctype>(grid_sizes[1])),
46 ang_quadrature_p2(quadrature_provider.get_points<ctype>(grid_sizes[2])),
47 ang_quadrature_w2(quadrature_provider.get_weights<ctype>(grid_sizes[2])),
48 ang_quadrature_p3(quadrature_provider.get_points<ctype>(grid_sizes[3])),
49 ang_quadrature_w3(quadrature_provider.get_weights<ctype>(grid_sizes[3]))
50 {
51 (void)max_block_size;
52 }
53
64 template <typename... T> NT get(const ctype k, const T &...t) const
65 {
66 using std::sqrt;
67
68 const auto constant = KERNEL::constant(k, t...);
69 return constant + tbb::parallel_reduce(
70 tbb::blocked_range3d<uint, uint, uint>(0, grid_sizes[0], 0, grid_sizes[1], 0,
71 grid_sizes[2] * grid_sizes[3]),
72 NT(0),
73 [&](const tbb::blocked_range3d<uint, uint, uint> &r, NT value) -> NT {
74 for (uint idx_x = r.pages().begin(); idx_x != r.pages().end(); ++idx_x) {
75 const ctype q = k * sqrt(x_quadrature_p[idx_x] * x_extent);
76 for (uint idx_y = r.rows().begin(); idx_y != r.rows().end(); ++idx_y) {
77 const ctype cos1 = 2 * (ang_quadrature_p1[idx_y] - (ctype)0.5);
78 const ctype int_element =
79 (powr<2>(q) * (ctype)0.5 * powr<2>(k)) // x = p^2 / k^2 integral
80 * sqrt(1. - powr<2>(cos1)) // cos1 integral jacobian
81 / powr<4>(2 * (ctype)M_PI); // fourier factor
82 for (uint idx_z = r.cols().begin(); idx_z != r.cols().end(); ++idx_z) {
83 const uint cos2_idx = idx_z / grid_sizes[2];
84 const uint phi_idx = idx_z % grid_sizes[3];
85 const ctype cos2 = 2 * (ang_quadrature_p2[cos2_idx] - (ctype)0.5);
86 const ctype phi = 2 * (ctype)M_PI * ang_quadrature_p3[phi_idx];
87 const ctype weight = 2 * (ctype)M_PI * ang_quadrature_w3[phi_idx] // phi weight
88 * 2 * ang_quadrature_w2[cos2_idx] // cos2 weight
89 * 2 * ang_quadrature_w1[idx_y] // cos1 weight
90 * x_quadrature_w[idx_x] * x_extent; // x weight
91 value += int_element * weight * KERNEL::kernel(q, cos1, cos2, phi, k, t...);
92 }
93 }
94 }
95 return value;
96 },
97 [&](NT x, NT y) -> NT { return x + y; });
98 }
99
110 template <typename... T> std::future<NT> request(const ctype k, const T &...t) const
111 {
112 return std::async(std::launch::deferred, [=, this]() { return get(k, t...); });
113 }
114
115 private:
116 const std::array<uint, 4> grid_sizes;
117
119
120 const std::vector<ctype> &x_quadrature_p;
121 const std::vector<ctype> &x_quadrature_w;
122 const std::vector<ctype> &ang_quadrature_p1;
123 const std::vector<ctype> &ang_quadrature_w1;
124 const std::vector<ctype> &ang_quadrature_p2;
125 const std::vector<ctype> &ang_quadrature_w2;
126 const std::vector<ctype> &ang_quadrature_p3;
127 const std::vector<ctype> &ang_quadrature_w3;
128 };
129} // namespace DiFfRG
Integrator for the integration of a 4D function with three angles with TBB. Calculates.
Definition integrator_4D_cpu.hh:26
const std::vector< ctype > & ang_quadrature_w2
Definition integrator_4D_cpu.hh:125
const std::vector< ctype > & ang_quadrature_p2
Definition integrator_4D_cpu.hh:124
typename get_type::ctype< NT > ctype
Numerical type to be used for integration tasks e.g. the argument or possible jacobians.
Definition integrator_4D_cpu.hh:31
NT get(const ctype k, const T &...t) const
Get the integral of the kernel.
Definition integrator_4D_cpu.hh:64
const std::vector< ctype > & x_quadrature_w
Definition integrator_4D_cpu.hh:121
const std::vector< ctype > & ang_quadrature_w1
Definition integrator_4D_cpu.hh:123
const std::vector< ctype > & x_quadrature_p
Definition integrator_4D_cpu.hh:120
std::future< NT > request(const ctype k, const T &...t) const
Request a future for the integral of the kernel.
Definition integrator_4D_cpu.hh:110
const std::vector< ctype > & ang_quadrature_w3
Definition integrator_4D_cpu.hh:127
Integrator4DTBB(QuadratureProvider &quadrature_provider, const std::array< uint, 4 > grid_sizes, const ctype x_extent, const JSONValue &)
Definition integrator_4D_cpu.hh:33
const std::vector< ctype > & ang_quadrature_p1
Definition integrator_4D_cpu.hh:122
const std::vector< ctype > & ang_quadrature_p3
Definition integrator_4D_cpu.hh:126
Integrator4DTBB(QuadratureProvider &quadrature_provider, std::array< uint, 4 > grid_sizes, const ctype x_extent, const uint max_block_size=0)
Definition integrator_4D_cpu.hh:39
const ctype x_extent
Definition integrator_4D_cpu.hh:118
const std::array< uint, 4 > grid_sizes
Definition integrator_4D_cpu.hh:116
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