DiFfRG
Loading...
Searching...
No Matches
integrator_angle_finiteTx0_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 <int d, typename NT, typename KERNEL> class IntegratorAngleFiniteTx0TBB
15 {
16 static_assert(d >= 2, "dimension must be at least 2");
17
18 public:
19 using ctype = typename get_type::ctype<NT>;
20
22 const ctype x_extent, const ctype x0_extent, const uint x0_summands,
23 const JSONValue &json)
25 json.get_double("/physical/T"), json.get_uint("/integration/cudathreadsperblock"))
26 {
27 }
28
30 const ctype x_extent, const ctype x0_extent, const uint _x0_summands, const ctype T,
31 const uint max_block_size = 0)
33 grid_sizes{{_grid_sizes[0], _grid_sizes[1], _grid_sizes[2] + _x0_summands}}, x_extent(x_extent),
34 x0_extent(x0_extent), original_x0_summands(_x0_summands)
35 {
36 ptr_x_quadrature_p = quadrature_provider.get_points<ctype>(grid_sizes[0]).data();
37 ptr_x_quadrature_w = quadrature_provider.get_weights<ctype>(grid_sizes[0]).data();
38 ptr_ang_quadrature_p = quadrature_provider.get_points<ctype>(grid_sizes[1]).data();
39 ptr_ang_quadrature_w = quadrature_provider.get_weights<ctype>(grid_sizes[1]).data();
40
41 set_T(T);
42
43 (void)max_block_size;
44 }
45
46 void set_T(const ctype T)
47 {
48 m_T = T;
49 if (is_close(T, 0.))
50 x0_summands = 0;
51 else
52 x0_summands = original_x0_summands;
53 ptr_x0_quadrature_p = quadrature_provider.get_points<ctype>(grid_sizes[2] - x0_summands).data();
54 ptr_x0_quadrature_w = quadrature_provider.get_weights<ctype>(grid_sizes[2] - x0_summands).data();
55 }
56
57 void set_x0_extent(const ctype val) { x0_extent = val; }
58
60 : quadrature_provider(other.quadrature_provider), grid_sizes(other.grid_sizes),
61 ptr_x_quadrature_p(other.ptr_x_quadrature_p), ptr_x_quadrature_w(other.ptr_x_quadrature_w),
62 ptr_ang_quadrature_p(other.ptr_ang_quadrature_p), ptr_ang_quadrature_w(other.ptr_ang_quadrature_w),
63 ptr_x0_quadrature_p(other.ptr_x0_quadrature_p), ptr_x0_quadrature_w(other.ptr_x0_quadrature_w),
64 x_extent(other.x_extent), x0_extent(other.x0_extent), original_x0_summands(other.original_x0_summands),
65 x0_summands(other.x0_summands), m_T(other.m_T)
66 {
67 }
68
69 template <typename... T> NT get(const ctype k, const T &...t) const
70 {
71 constexpr ctype S_dm1 = S_d_prec<ctype>(d - 1);
72 using std::sqrt, std::exp, std::log;
73
74 const ctype integral_start = (2 * x0_summands * (ctype)M_PI * m_T) / k;
75 const ctype log_start = log(integral_start + (m_T == 0) * ctype(1e-3));
76 const ctype log_ext = log(x0_extent / (integral_start + (m_T == 0) * ctype(1e-3)));
77
78 const auto constant = KERNEL::constant(k, t...);
79 return constant +
80 tbb::parallel_reduce(
81 tbb::blocked_range3d<uint, uint, uint>(0, grid_sizes[0], 0, grid_sizes[1], 0, grid_sizes[2]), NT(0),
82 [&](const tbb::blocked_range3d<uint, uint, uint> &r, NT value) -> NT {
83 for (uint idx_x = r.pages().begin(); idx_x != r.pages().end(); ++idx_x) {
84 const ctype q = k * sqrt(ptr_x_quadrature_p[idx_x] * x_extent);
85 for (uint idx_y = r.rows().begin(); idx_y != r.rows().end(); ++idx_y) {
86 const ctype cos = 2 * (ptr_ang_quadrature_p[idx_y] - (ctype)0.5);
87 for (uint idx_z = r.cols().begin(); idx_z != r.cols().end(); ++idx_z) {
88 if (idx_z >= x0_summands) {
89 const ctype q0 = k * (exp(log_start + log_ext * ptr_x0_quadrature_p[idx_z - x0_summands]) -
90 (m_T == 0) * ctype(1e-3));
91
92 const ctype int_element = S_dm1 // solid nd angle
93 *
94 (powr<d - 3>(q) / (ctype)2 * powr<2>(k)) // x = p^2 / k^2 integral
95 * (k) // x0 = q0 / k integral
96 * (1 / (ctype)2) // divide the cos integral out
97 / powr<d>(2 * (ctype)M_PI); // fourier factor
98 const ctype weight = 2 * ptr_ang_quadrature_w[idx_y] * ptr_x_quadrature_w[idx_x] * x_extent *
99 (ptr_x0_quadrature_w[idx_z - x0_summands] * log_ext * q0 / k);
100
101 value += int_element * weight *
102 (KERNEL::kernel(q, cos, q0, k, t...) + KERNEL::kernel(q, cos, -q0, k, t...));
103 } else {
104 const ctype q0 = 2 * (ctype)M_PI * m_T * idx_z;
105 const ctype int_element = m_T * S_dm1 // solid nd angle
106 *
107 (powr<d - 3>(q) / (ctype)2 * powr<2>(k)) // x = p^2 / k^2 integral
108 * (1 / (ctype)2) // divide the cos integral out
109 / powr<d - 1>(2 * (ctype)M_PI); // fourier factor
110 const ctype weight = 2 * ptr_ang_quadrature_w[idx_y] * ptr_x_quadrature_w[idx_x] * x_extent;
111 value += int_element * weight *
112 (idx_z == 0
113 ? KERNEL::kernel(q, cos, (ctype)0, k, t...)
114 : KERNEL::kernel(q, cos, q0, k, t...) + KERNEL::kernel(q, cos, -q0, k, t...));
115 }
116 }
117 }
118 }
119 return value;
120 },
121 [&](NT x, NT y) -> NT { return x + y; });
122 }
123
124 template <typename... T> std::future<NT> request(const ctype k, const T &...t) const
125 {
126 return std::async(std::launch::deferred, [=, this]() { return get(k, t...); });
127 }
128
129 private:
131
132 const std::array<uint, 3> grid_sizes;
133
139
146 };
147} // namespace DiFfRG
Definition integrator_angle_finiteTx0_cpu.hh:15
const ctype * ptr_x_quadrature_w
Definition integrator_angle_finiteTx0_cpu.hh:141
void set_T(const ctype T)
Definition integrator_angle_finiteTx0_cpu.hh:46
const ctype * ptr_x0_quadrature_w
Definition integrator_angle_finiteTx0_cpu.hh:145
std::future< NT > request(const ctype k, const T &...t) const
Definition integrator_angle_finiteTx0_cpu.hh:124
QuadratureProvider & quadrature_provider
Definition integrator_angle_finiteTx0_cpu.hh:130
const ctype * ptr_ang_quadrature_w
Definition integrator_angle_finiteTx0_cpu.hh:143
const uint original_x0_summands
Definition integrator_angle_finiteTx0_cpu.hh:136
typename get_type::ctype< NT > ctype
Definition integrator_angle_finiteTx0_cpu.hh:19
const ctype * ptr_x_quadrature_p
Definition integrator_angle_finiteTx0_cpu.hh:140
void set_x0_extent(const ctype val)
Definition integrator_angle_finiteTx0_cpu.hh:57
uint x0_summands
Definition integrator_angle_finiteTx0_cpu.hh:137
IntegratorAngleFiniteTx0TBB(QuadratureProvider &quadrature_provider, const std::array< uint, 3 > grid_sizes, const ctype x_extent, const ctype x0_extent, const uint x0_summands, const JSONValue &json)
Definition integrator_angle_finiteTx0_cpu.hh:21
const ctype * ptr_ang_quadrature_p
Definition integrator_angle_finiteTx0_cpu.hh:142
const ctype * ptr_x0_quadrature_p
Definition integrator_angle_finiteTx0_cpu.hh:144
ctype m_T
Definition integrator_angle_finiteTx0_cpu.hh:138
ctype x0_extent
Definition integrator_angle_finiteTx0_cpu.hh:135
const std::array< uint, 3 > grid_sizes
Definition integrator_angle_finiteTx0_cpu.hh:132
IntegratorAngleFiniteTx0TBB(const IntegratorAngleFiniteTx0TBB &other)
Definition integrator_angle_finiteTx0_cpu.hh:59
NT get(const ctype k, const T &...t) const
Definition integrator_angle_finiteTx0_cpu.hh:69
const ctype x_extent
Definition integrator_angle_finiteTx0_cpu.hh:134
IntegratorAngleFiniteTx0TBB(QuadratureProvider &quadrature_provider, std::array< uint, 3 > _grid_sizes, const ctype x_extent, const ctype x0_extent, const uint _x0_summands, const ctype T, const uint max_block_size=0)
Definition integrator_angle_finiteTx0_cpu.hh:29
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 std::vector< NT > & get_points(const size_t order, const QuadratureType type=QuadratureType::legendre)
Get the quadrature points for a quadrature of size quadrature_size.
Definition quadrature_provider.hh:151
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 auto & get(named_tuple< tuple_type, strs... > &ob)
get a reference to the element with the given name
Definition tuples.hh:82
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
consteval NT S_d_prec(uint d)
Surface of a d-dimensional sphere (precompiled)
Definition math.hh:104
unsigned int uint
Definition utils.hh:22