DiFfRG
Loading...
Searching...
No Matches
integrator_1D_cartesian_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{
20 template <typename NT, typename KERNEL> class Integrator1DCartesianTBB
21 {
22 public:
26 using ctype = typename get_type::ctype<NT>;
27
38 Integrator1DCartesianTBB(QuadratureProvider &quadrature_provider, const std::array<uint, 1> grid_size,
39 const ctype x_extent = 0., const uint max_block_size = 0, const ctype qx_min = -M_PI,
40 const ctype qx_max = M_PI)
41 : Integrator1DCartesianTBB(quadrature_provider, grid_size[0], x_extent, qx_min, qx_max)
42 {
43 }
44
56 Integrator1DCartesianTBB(QuadratureProvider &quadrature_provider, const std::array<uint, 1> grid_size,
57 const JSONValue &json)
58 : Integrator1DCartesianTBB(quadrature_provider, grid_size[0], 0.,
59 json.get_double("/discretization/integration/qx_min", -M_PI),
60 json.get_double("/discretization/integration/qx_max", M_PI))
61 {
62 }
63
73 Integrator1DCartesianTBB(QuadratureProvider &quadrature_provider, const std::array<uint, 1> grid_size,
74 const ctype x_extent, const JSONValue &json)
75 : Integrator1DCartesianTBB(quadrature_provider, grid_size[0], x_extent,
76 json.get_double("/discretization/integration/qx_min", -M_PI),
77 json.get_double("/discretization/integration/qx_max", M_PI))
78 {
79 }
80
90 Integrator1DCartesianTBB(QuadratureProvider &quadrature_provider, const uint grid_size, const ctype x_extent = 0.,
91 const ctype qx_min = -M_PI, const ctype qx_max = M_PI)
92 : grid_size(grid_size), x_quadrature_p(quadrature_provider.get_points<ctype>(grid_size)),
93 x_quadrature_w(quadrature_provider.get_weights<ctype>(grid_size))
94 {
95 this->qx_min = qx_min;
96 this->qx_extent = qx_max - qx_min;
97 }
98
103 {
104 this->qx_extent = this->qx_extent - qx_min + this->qx_min;
105 this->qx_min = qx_min;
106 }
107
111 void set_qx_max(const ctype qx_max) { this->qx_extent = qx_max - qx_min; }
112
121 template <typename... T> NT get(const ctype k, const T &...t) const
122 {
123 constexpr int d = 1;
124 constexpr ctype int_element = powr<-d>(2 * (ctype)M_PI); // fourier factor
125
126 return KERNEL::constant(k, t...) + tbb::parallel_reduce(
127 tbb::blocked_range<uint>(0, grid_size), NT(0),
128 [&](const tbb::blocked_range<uint> &r, NT value) -> NT {
129 for (uint idx_x = r.begin(); idx_x != r.end(); ++idx_x) {
130 const ctype q = qx_min + qx_extent * x_quadrature_p[idx_x];
131 const ctype weight = qx_extent * x_quadrature_w[idx_x];
132 value += int_element * weight * KERNEL::kernel(q, k, t...);
133 }
134 return value;
135 },
136 [&](NT x, NT y) -> NT { return x + y; });
137 }
138
147 template <typename... T> std::future<NT> request(const ctype k, const T &...t) const
148 {
149 return std::async(std::launch::deferred, [=, this]() { return get(k, t...); });
150 }
151
152 private:
154
155 ctype qx_min = -M_PI;
156 ctype qx_extent = 2. * M_PI;
157
158 const std::vector<ctype> &x_quadrature_p;
159 const std::vector<ctype> &x_quadrature_w;
160 };
161} // namespace DiFfRG
Integration of an arbitrary 1D function from qx_min to qx_max using TBB.
Definition integrator_1D_cartesian_cpu.hh:21
void set_qx_max(const ctype qx_max)
Set the maximum value of the qx integration range.
Definition integrator_1D_cartesian_cpu.hh:111
const uint grid_size
Definition integrator_1D_cartesian_cpu.hh:153
Integrator1DCartesianTBB(QuadratureProvider &quadrature_provider, const uint grid_size, const ctype x_extent=0., const ctype qx_min=-M_PI, const ctype qx_max=M_PI)
Construct a new Integrator1DCartesianTBB object.
Definition integrator_1D_cartesian_cpu.hh:90
const std::vector< ctype > & x_quadrature_p
Definition integrator_1D_cartesian_cpu.hh:158
void set_qx_min(const ctype qx_min)
Set the minimum value of the qx integration range.
Definition integrator_1D_cartesian_cpu.hh:102
ctype qx_extent
Definition integrator_1D_cartesian_cpu.hh:156
const std::vector< ctype > & x_quadrature_w
Definition integrator_1D_cartesian_cpu.hh:159
Integrator1DCartesianTBB(QuadratureProvider &quadrature_provider, const std::array< uint, 1 > grid_size, const ctype x_extent, const JSONValue &json)
Construct a new Integrator1DCartesianTBB object.
Definition integrator_1D_cartesian_cpu.hh:73
Integrator1DCartesianTBB(QuadratureProvider &quadrature_provider, const std::array< uint, 1 > grid_size, const ctype x_extent=0., const uint max_block_size=0, const ctype qx_min=-M_PI, const ctype qx_max=M_PI)
Construct a new Integrator1DCartesianTBB object.
Definition integrator_1D_cartesian_cpu.hh:38
ctype qx_min
Definition integrator_1D_cartesian_cpu.hh:155
NT get(const ctype k, const T &...t) const
Get the integration result.
Definition integrator_1D_cartesian_cpu.hh:121
std::future< NT > request(const ctype k, const T &...t) const
Get the integration result asynchronously.
Definition integrator_1D_cartesian_cpu.hh:147
Integrator1DCartesianTBB(QuadratureProvider &quadrature_provider, const std::array< uint, 1 > grid_size, const JSONValue &json)
Construct a new Integrator1DCartesianTBB object.
Definition integrator_1D_cartesian_cpu.hh:56
typename get_type::ctype< NT > ctype
Numerical type to be used for integration tasks e.g. the argument or possible jacobians.
Definition integrator_1D_cartesian_cpu.hh:26
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