IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
InverseTransformSampling.h
Go to the documentation of this file.
1// Class InverseTransformSampling
2// This class can be used for generating samples of a given distribution function class using
3// inverse of its cumulative distribution on host or device.
4//
5
6#ifndef IPPL_INVERSE_TRANSFORM_SAMPLING_H
7#define IPPL_INVERSE_TRANSFORM_SAMPLING_H
8
9#include "Random/Utility.h"
10
11namespace ippl {
12 namespace random {
25 template <typename T, unsigned Dim, class DeviceType, class Distribution>
27 public:
30
34
44 template <class RegionLayout>
46 Vector<T, Dim>& rmin_r, RegionLayout& rlayout_r,
47 size_type& ntotal_r)
48 : dist_m(dist_r)
49 , ntotal_m(ntotal_r) {
50 const typename RegionLayout::host_mirror_type regions =
51 rlayout_r.gethLocalRegions();
52 int rank = ippl::Comm->rank();
53
54 Vector<T, Dim> locrmax, locrmin;
55 for (unsigned d = 0; d < Dim; ++d) {
56 locrmax[d] = regions(rank)[d].max();
57 locrmin[d] = regions(rank)[d].min();
58 }
59
60 updateBounds(rmax_r, rmin_r, locrmax, locrmin);
61 }
62
75 Vector<T, Dim>& rmin_r, Vector<T, Dim>& locrmax_r,
76 Vector<T, Dim>& locrmin_r, size_type& ntotal_r)
77 : dist_m(dist_r)
78 , ntotal_m(ntotal_r) {
79 updateBounds(rmax_r, rmin_r, locrmax_r, locrmin_r);
80 }
81
92 Vector<T, Dim>& rmin_r, size_type& ntotal_r)
93 : dist_m(dist_r)
94 , ntotal_m(ntotal_r) {
95 updateBounds(rmax_r, rmin_r);
96
98 }
99
118 Vector<T, Dim>& locrmin) {
119 int rank = ippl::Comm->rank();
120 Vector<T, Dim> nr_m, dr_m;
121 for (unsigned d = 0; d < Dim; ++d) {
122 nr_m[d] = dist_m.getCdf(locrmax[d], d) - dist_m.getCdf(locrmin[d], d);
123 dr_m[d] = dist_m.getCdf(rmax[d], d) - dist_m.getCdf(rmin[d], d);
124 umin_m[d] = dist_m.getCdf(locrmin[d], d);
125 umax_m[d] = dist_m.getCdf(locrmax[d], d);
126 }
127 T pnr = std::accumulate(nr_m.begin(), nr_m.end(), 1.0, std::multiplies<T>());
128 T pdr = std::accumulate(dr_m.begin(), dr_m.end(), 1.0, std::multiplies<T>());
129
130 T factor = pnr / pdr;
131 nlocal_m = (size_type)(factor * ntotal_m);
132
133 size_type nglobal = 0;
134 ippl::Comm->allreduce(&nlocal_m, &nglobal, 1, std::plus<size_type>());
135
136 int rest = (int)(ntotal_m - nglobal);
137 if (rank < rest) {
138 ++nlocal_m;
139 }
140 }
141
156 Vector<T, Dim> nr_m, dr_m;
157 for (unsigned d = 0; d < Dim; ++d) {
158 umin_m[d] = dist_m.getCdf(rmin[d], d);
159 umax_m[d] = dist_m.getCdf(rmax[d], d);
160 }
161 }
162
167
171 template <class GeneratorPool>
172 struct fill_random {
173 using value_type = T;
176 GeneratorPool pool_m;
179 unsigned int dim_m;
180
190 KOKKOS_FUNCTION
191 fill_random(Distribution& dist_r, view_type& x_r, GeneratorPool& rand_pool_r,
192 Vector<T, Dim>& umin_r, Vector<T, Dim>& umax_r, unsigned int& d_r)
193 : targetdist_m(dist_r)
194 , sample_m(x_r)
195 , pool_m(rand_pool_r)
196 , minbound_m(umin_r)
197 , maxbound_m(umax_r)
198 , dim_m(d_r) {}
199
205 KOKKOS_INLINE_FUNCTION void operator()(const size_t i) const {
206 typename GeneratorPool::generator_type rand_gen = pool_m.get_state();
207
208 value_type u = 0.0;
209
210 u = rand_gen.drand(minbound_m[dim_m], maxbound_m[dim_m]);
211
212 // first guess for Newton-Raphson
213 sample_m(i)[dim_m] = targetdist_m.getEstimate(u, dim_m);
214
215 // solve
218
219 solver.solve(dim_m, sample_m(i)[dim_m], u);
220
221 pool_m.free_state(rand_gen);
222 }
223 };
224
230 KOKKOS_INLINE_FUNCTION size_type getLocalSamplesNum() const { return nlocal_m; }
231
237 KOKKOS_INLINE_FUNCTION void setLocalSamplesNum(size_type nlocal) { nlocal_m = nlocal; }
238
245 void generate(view_type view, Kokkos::Random_XorShift64_Pool<> rand_pool64) {
246 Vector<T, Dim> minbound_m = umin_m;
247 Vector<T, Dim> maxbound_m = umax_m;
248 Distribution targetdist_m = dist_m;
249 size_type numlocal_m = nlocal_m;
250 for (unsigned d = 0; d < Dim; ++d) {
251 Kokkos::parallel_for(numlocal_m, fill_random<Kokkos::Random_XorShift64_Pool<>>(
252 targetdist_m, view, rand_pool64,
253 minbound_m, maxbound_m, d));
254 Kokkos::fence();
255 }
256 }
257
267 void generate(view_type view, size_type startIndex, size_type endIndex,
268 Kokkos::Random_XorShift64_Pool<> rand_pool64) {
269 Vector<T, Dim> minbound_m = umin_m;
270 Vector<T, Dim> maxbound_m = umax_m;
271 Distribution targetdist_m = dist_m;
272 for (unsigned d = 0; d < Dim; ++d) {
273 Kokkos::parallel_for(
274 Kokkos::RangePolicy<>(startIndex, endIndex),
275 fill_random<Kokkos::Random_XorShift64_Pool<>>(
276 targetdist_m, view, rand_pool64, minbound_m, maxbound_m, d));
277 Kokkos::fence();
278 }
279 }
280
281 private:
283 };
284
285 } // namespace random
286} // namespace ippl
287
288#endif
constexpr unsigned Dim
Definition Archive.h:20
std::unique_ptr< mpi::Communicator > Comm
Definition Ippl.h:22
std::size_t size_type
Definition IpplTypes.h:13
The class that represents a distribution.
InverseTransformSampling(Distribution &dist_r, Vector< T, Dim > &rmax_r, Vector< T, Dim > &rmin_r, RegionLayout &rlayout_r, size_type &ntotal_r)
Constructor for InverseTransformSampling class with domain decomposition.
KOKKOS_INLINE_FUNCTION size_type getLocalSamplesNum() const
Get the local number of samples.
KOKKOS_INLINE_FUNCTION void setLocalSamplesNum(size_type nlocal)
Set the local number of particles.
typename ippl::detail::ViewType< Vector< T, Dim >, 1 >::view_type view_type
InverseTransformSampling(Distribution &dist_r, Vector< T, Dim > &rmax_r, Vector< T, Dim > &rmin_r, size_type &ntotal_r)
Constructor for InverseTransformSampling class. In this method, we do not consider any domain decompo...
void updateBounds(Vector< T, Dim > &rmax, Vector< T, Dim > &rmin)
Updates the sampling bounds using the CDF without any domain decomposition.
InverseTransformSampling(Distribution &dist_r, Vector< T, Dim > &rmax_r, Vector< T, Dim > &rmin_r, Vector< T, Dim > &locrmax_r, Vector< T, Dim > &locrmin_r, size_type &ntotal_r)
Constructor for InverseTransformSampling class without applying domain decomposition....
void generate(view_type view, Kokkos::Random_XorShift64_Pool<> rand_pool64)
Generate random samples using inverse transform sampling.
void generate(view_type view, size_type startIndex, size_type endIndex, Kokkos::Random_XorShift64_Pool<> rand_pool64)
Generate random samples using inverse transform sampling for a specific range of particles.
void updateBounds(Vector< T, Dim > &rmax, Vector< T, Dim > &rmin, Vector< T, Dim > &locrmax, Vector< T, Dim > &locrmin)
Updates the sampling bounds and reinitializes internal variables.
~InverseTransformSampling()
Deconstructor for InverseTransformSampling class.
Functor that is used for generating samples.
KOKKOS_INLINE_FUNCTION void operator()(const size_t i) const
Operator to fill random values.
KOKKOS_FUNCTION fill_random(Distribution &dist_r, view_type &x_r, GeneratorPool &rand_pool_r, Vector< T, Dim > &umin_r, Vector< T, Dim > &umax_r, unsigned int &d_r)
Constructor for the fill_random functor.
Functor for solving equations using the Newton-Raphson method.
Definition Utility.h:28
KOKKOS_INLINE_FUNCTION void solve(unsigned int d, T &x, T &u)
Solve an equation using the Newton-Raphson method.
Definition Utility.h:52
KOKKOS_INLINE_FUNCTION constexpr iterator begin()
Definition Vector.hpp:160
KOKKOS_INLINE_FUNCTION constexpr iterator end()
Definition Vector.hpp:165