DiFfRG
Loading...
Searching...
No Matches
integrator_finiteTq0_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 IntegratorFiniteTq0TBB
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 JSONValue &json)
23 : IntegratorFiniteTq0TBB(quadrature_provider, grid_sizes, x_extent, json.get_double("/physical/T"))
24 {
25 }
26
27 IntegratorFiniteTq0TBB(QuadratureProvider &quadrature_provider, const std::array<uint, 1> _grid_sizes,
28 const ctype x_extent, const ctype T, const uint max_block_size = 0)
29 : quadrature_provider(quadrature_provider), grid_sizes{{_grid_sizes[0], 0}}, x_extent(x_extent), m_T(T),
30 manual_E(false)
31 {
32 ptr_x_quadrature_p = quadrature_provider.get_points<ctype>(grid_sizes[0]).data();
33 ptr_x_quadrature_w = quadrature_provider.get_weights<ctype>(grid_sizes[0]).data();
34
35 set_T(T);
36
37 (void)max_block_size;
38 }
39
47 void set_T(const ctype T, const ctype E = 0)
48 {
49 if (is_close(T, m_T) && (std::abs(E - m_E) / std::max(E, m_E) < 2.5e-2)) return;
50
51 m_T = T;
52 // the default typical energy scale will default the matsubara size to 11.
53 m_E = is_close(E, 0.) ? 10 * m_T : E;
54 manual_E = !is_close(E, 0.);
55
56 grid_sizes[1] = quadrature_provider.get_matsubara_points<ctype>(m_T, m_E).size();
57 ptr_matsubara_quadrature_p = quadrature_provider.get_matsubara_points<ctype>(m_T, m_E).data();
58 ptr_matsubara_quadrature_w = quadrature_provider.get_matsubara_weights<ctype>(m_T, m_E).data();
59 }
60
66 void set_E(const ctype E) { set_T(m_T, E); }
67
69 : quadrature_provider(other.quadrature_provider), grid_sizes(other.grid_sizes),
70 ptr_x_quadrature_p(other.ptr_x_quadrature_p), ptr_x_quadrature_w(other.ptr_x_quadrature_w),
71 ptr_matsubara_quadrature_p(other.ptr_matsubara_quadrature_p),
72 ptr_matsubara_quadrature_w(other.ptr_matsubara_quadrature_w), x_extent(other.x_extent), m_T(other.m_T),
73 m_E(other.m_E), manual_E(other.manual_E)
74 {
75 }
76
77 template <typename... T> NT get(const ctype k, const T &...t)
78 {
79 if (!manual_E && (std::abs(k - m_E) / std::max(k, m_E) > 2.5e-2)) {
80 set_T(m_T, k);
81 manual_E = false;
82 }
83
84 constexpr ctype S_dm1 = S_d_prec<ctype>(d - 1);
85 using std::sqrt, std::exp, std::log;
86
87 const auto constant = KERNEL::constant(k, t...);
88 return constant +
89 tbb::parallel_reduce(
90 tbb::blocked_range2d<uint, uint>(0, grid_sizes[0], 0, grid_sizes[1]), NT(0),
91 [&](const tbb::blocked_range2d<uint, uint> &r, NT value) -> NT {
92 for (uint idx_x = r.rows().begin(); idx_x != r.rows().end(); ++idx_x) {
93 const ctype q = k * sqrt(ptr_x_quadrature_p[idx_x] * x_extent);
94
95 const ctype x_weight = ptr_x_quadrature_w[idx_x] * x_extent;
96 const ctype int_element = S_dm1 // solid nd angle
97 * (powr<d - 3>(q) / (ctype)2 * powr<2>(k)) // x = p^2 / k^2 integral
98 / powr<d - 1>(2 * (ctype)M_PI); // fourier factor
99
100 for (uint idx_y = r.cols().begin(); idx_y != r.cols().end(); ++idx_y) {
101 // integral part
102 const ctype q0 = ptr_matsubara_quadrature_p[idx_y];
103 const ctype weight = x_weight * ptr_matsubara_quadrature_w[idx_y];
104
105 value +=
106 int_element * weight * (KERNEL::kernel(q, q0, k, t...) + KERNEL::kernel(q, -q0, k, t...));
107 if (m_T > 0 && idx_y == 0)
108 value += int_element * x_weight * m_T * KERNEL::kernel(q, (ctype)0, k, t...);
109 }
110 }
111 return value;
112 },
113 [&](NT x, NT y) -> NT { return x + y; });
114 }
115
116 template <typename... T> std::future<NT> request(const ctype k, const T &...t)
117 {
118 return std::async(std::launch::deferred, [=, this]() { return get(k, t...); });
119 }
120
121 private:
123
124 std::array<uint, 2> grid_sizes;
125
127 ctype m_T, m_E;
129
134 };
135} // namespace DiFfRG
Definition integrator_finiteTq0_cpu.hh:15
IntegratorFiniteTq0TBB(QuadratureProvider &quadrature_provider, const std::array< uint, 1 > _grid_sizes, const ctype x_extent, const ctype T, const uint max_block_size=0)
Definition integrator_finiteTq0_cpu.hh:27
const ctype * ptr_matsubara_quadrature_p
Definition integrator_finiteTq0_cpu.hh:132
std::array< uint, 2 > grid_sizes
Definition integrator_finiteTq0_cpu.hh:124
const ctype * ptr_x_quadrature_w
Definition integrator_finiteTq0_cpu.hh:131
QuadratureProvider & quadrature_provider
Definition integrator_finiteTq0_cpu.hh:122
std::future< NT > request(const ctype k, const T &...t)
Definition integrator_finiteTq0_cpu.hh:116
const ctype * ptr_matsubara_quadrature_w
Definition integrator_finiteTq0_cpu.hh:133
void set_T(const ctype T, const ctype E=0)
Set the temperature and typical energy scale of the integrator and recompute the Matsubara quadrature...
Definition integrator_finiteTq0_cpu.hh:47
typename get_type::ctype< NT > ctype
Definition integrator_finiteTq0_cpu.hh:19
IntegratorFiniteTq0TBB(QuadratureProvider &quadrature_provider, const std::array< uint, 1 > grid_sizes, const ctype x_extent, const JSONValue &json)
Definition integrator_finiteTq0_cpu.hh:21
const ctype * ptr_x_quadrature_p
Definition integrator_finiteTq0_cpu.hh:130
IntegratorFiniteTq0TBB(const IntegratorFiniteTq0TBB &other)
Definition integrator_finiteTq0_cpu.hh:68
const ctype x_extent
Definition integrator_finiteTq0_cpu.hh:126
ctype m_E
Definition integrator_finiteTq0_cpu.hh:127
bool manual_E
Definition integrator_finiteTq0_cpu.hh:128
void set_E(const ctype E)
Set the typical energy scale of the integrator and recompute the Matsubara quadrature rule.
Definition integrator_finiteTq0_cpu.hh:66
NT get(const ctype k, const T &...t)
Definition integrator_finiteTq0_cpu.hh:77
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