OPALX (Object Oriented Parallel Accelerator Library for Exascal) MINIorX
OPALX
DistributionMoments.cpp
Go to the documentation of this file.
1
2// Class DistributionMoments
3// Computes the statistics of particle distributions.
4//
5// Copyright (c) 2021, Christof Metzger-Kraus
6// All rights reserved
7//
8// This file is part of OPAL.
9//
10// OPAL is free software: you can redistribute it and/or modify
11// it under the terms of the GNU General Public License as published by
12// the Free Software Foundation, either version 3 of the License, or
13// (at your option) any later version.
14//
15// You should have received a copy of the GNU General Public License
16// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
17//
18
20
21#include "Utilities/Options.h"
22#include "Utilities/Util.h"
23
24#include "Utility/Inform.h"
25
27
28#include <gsl/gsl_histogram.h>
29
30#include <boost/numeric/conversion/cast.hpp>
31
32extern Inform* gmsg;
33
34const double DistributionMoments::percentileOneSigmaNormalDist_m = std::erf(1 / sqrt(2));
35const double DistributionMoments::percentileTwoSigmasNormalDist_m = std::erf(2 / sqrt(2));
36const double DistributionMoments::percentileThreeSigmasNormalDist_m = std::erf(3 / sqrt(2));
37const double DistributionMoments::percentileFourSigmasNormalDist_m = std::erf(4 / sqrt(2));
38
40 reset();
42
43 moments_m.resize(6, 6, false);
44 notCentMoments_m.resize(6, 6, false);
45}
46
48 ippl::ParticleAttrib<Vector_t<double,3>>::view_type& Pview,
49 ippl::ParticleAttrib<double>::view_type& Mview,
50 size_t Np,
51 size_t Nlocal) {
52 /*
53 Np is the total number of particles (reduced over ranks). In this function, it is only used to
54 average over the number of total particles. For an empty simulation, this leads to divisions by
55 0. Since, however, some of the computed moments might be needed regardless, we compute them but
56 need to make sure that we do not divide by 0.
57 Solution: Set Np to 1 if it is 0.
58 */
59 Np = (Np == 0) ? 1 : Np;
60
61 const int Dim = 3;
62 double loc_centroid[2 * Dim] = {};
63 double centroid[2 * Dim] = {};
64 double loc_Ekin, loc_gamma, loc_gammaz, gammaz;
65
66 int rank;
67 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
68
69 Kokkos::parallel_reduce(
70 "calc moments of particle distr.", Nlocal,
71 KOKKOS_LAMBDA(
72 const int k, double& ekin, double& gamma, double& gammaz) {
73 double gamma0 = 0.0;
74 double ekin0 = 0.0;
75
76 for(unsigned j=0; j<Dim; j++){
77 gamma0 += Pview(k)[j]*Pview(k)[j];
78 }
79 gamma0 = Kokkos::sqrt(gamma0+1.0);
80 ekin0 = (gamma0-1.0) * Mview(k);
81 gamma += gamma0;
82 ekin += ekin0;
83
84 gammaz += Pview(k)[2];
85 },
86 Kokkos::Sum<double>(loc_Ekin), Kokkos::Sum<double>(loc_gamma), Kokkos::Sum<double>(loc_gammaz));
87 Kokkos::fence();
88
89 for (unsigned i = 0; i < 2 * Dim; ++i) {
90 Kokkos::parallel_reduce(
91 "calc moments of particle distr.", Nlocal,
92 KOKKOS_LAMBDA(
93 const int k, double& cent) {
94 double part[2 * Dim];
95 part[0] = Rview(k)[0];
96 part[1] = Pview(k)[0];
97 part[2] = Rview(k)[1];
98 part[3] = Pview(k)[1];
99 part[4] = Rview(k)[2];
100 part[5] = Pview(k)[2];
101
102 cent += part[i];
103 },
104 Kokkos::Sum<double>(loc_centroid[i]));
105 Kokkos::fence();
106 }
107 ippl::Comm->barrier();
108
109 MPI_Allreduce(
110 loc_centroid, centroid, 2 * Dim, MPI_DOUBLE, MPI_SUM, ippl::Comm->getCommunicator());
111 MPI_Allreduce(
112 &loc_Ekin, &meanKineticEnergy_m, 1, MPI_DOUBLE, MPI_SUM, ippl::Comm->getCommunicator());
113 MPI_Allreduce(
114 &loc_gamma, &meanGamma_m, 1, MPI_DOUBLE, MPI_SUM, ippl::Comm->getCommunicator());
115 MPI_Allreduce(
116 &loc_gammaz, &gammaz, 1, MPI_DOUBLE, MPI_SUM, ippl::Comm->getCommunicator());
117
118 for (unsigned i = 0; i < 2 * Dim; i++) {
119 centroid_m(i) = centroid[i];
120 means_m(i) = centroid_m[i]/Np;
121 }
122
124 meanGamma_m = meanGamma_m / (1.*Np);
125 gammaz = gammaz / (1.*Np);
126 gammaz *= gammaz;
127 gammaz = std::sqrt(gammaz + 1.0);
128 meanGammaZ_m = gammaz;
129
130 // store mean R, mean P, std R, std P in class member variables
131 for (unsigned i = 0; i < Dim; i++) {
132 meanR_m(i) = means_m[2*i];
133 meanP_m(i) = means_m[2*i+1];
134 }
135}
137 ippl::ParticleAttrib<Vector_t<double,3>>::view_type& Pview,
138 ippl::ParticleAttrib<double>::view_type& Mview,
139 size_t Np,
140 size_t Nlocal) {
141 Np = (Np == 0) ? 1 : Np; // Explanation: see DistributionMoments::computeMeans implementation
142
143 reset();
144 computeMeans(Rview, Pview, Mview, Np, Nlocal);
145
146 double meanR_loc[Dim] = {};
147 double meanP_loc[Dim] = {};
148 double loc_moment[2 * Dim][2 * Dim] = {};
149 double moment[2 * Dim][2 * Dim] = {};
150
151 double loc_moment_ncent[2 * Dim][2 * Dim] = {};
152 double moment_ncent[2 * Dim][2 * Dim] = {};
153
154 // compute central moments
155 for (unsigned i = 0; i < Dim; i++) {
156 meanR_loc[i] = meanR_m[i];
157 meanP_loc[i] = meanP_m[i];
158 }
159
160 for (unsigned i = 0; i < 2 * Dim; ++i) {
161 Kokkos::parallel_reduce(
162 "calc moments of particle distr.", Nlocal,
163 KOKKOS_LAMBDA(
164 const int k, double& mom0, double& mom1, double& mom2,
165 double& mom3, double& mom4, double& mom5) {
166 double part[2 * Dim];
167 part[0] = Rview(k)[0]-meanR_loc[0];
168 part[1] = Pview(k)[0]-meanP_loc[0];
169 part[2] = Rview(k)[1]-meanR_loc[1];
170 part[3] = Pview(k)[1]-meanP_loc[1];
171 part[4] = Rview(k)[2]-meanR_loc[2];
172 part[5] = Pview(k)[2]-meanP_loc[2];
173
174 mom0 += part[i] * part[0];
175 mom1 += part[i] * part[1];
176 mom2 += part[i] * part[2];
177 mom3 += part[i] * part[3];
178 mom4 += part[i] * part[4];
179 mom5 += part[i] * part[5];
180 },
181 Kokkos::Sum<double>(loc_moment[i][0]),
182 Kokkos::Sum<double>(loc_moment[i][1]), Kokkos::Sum<double>(loc_moment[i][2]),
183 Kokkos::Sum<double>(loc_moment[i][3]), Kokkos::Sum<double>(loc_moment[i][4]),
184 Kokkos::Sum<double>(loc_moment[i][5]));
185 Kokkos::fence();
186 }
187
188 MPI_Allreduce(
189 loc_moment, moment, 2 * Dim * 2 * Dim, MPI_DOUBLE, MPI_SUM, ippl::Comm->getCommunicator());
190
191 for (unsigned i = 0; i < 2 * Dim; i++) {
192 for (unsigned j = 0; j < 2 * Dim; j++) {
193 moments_m(i,j) = moment[i][j] / Np;
194 }
195 }
196
197 for (unsigned i = 0; i < Dim; i++) {
198 stdR_m(i) = std::sqrt( moments_m(2*i, 2*i) );
199 stdP_m(i) = std::sqrt( moments_m(2*i+1, 2*i+1) );
200 }
201
202 double mekin = meanKineticEnergy_m;
203 double loc_std_mekin;
204
205 // compute non-central moments
206 for (unsigned i = 0; i < 2 * Dim; ++i) {
207 Kokkos::parallel_reduce(
208 "calc moments of particle distr.", Nlocal,
209 KOKKOS_LAMBDA(
210 const int k, double& mom0, double& mom1, double& mom2,
211 double& mom3, double& mom4, double& mom5) {
212 double part[2 * Dim];
213 part[0] = Rview(k)[0];
214 part[1] = Pview(k)[0];
215 part[2] = Rview(k)[1];
216 part[3] = Pview(k)[1];
217 part[4] = Rview(k)[2];
218 part[5] = Pview(k)[2];
219
220 mom0 += part[i] * part[0];
221 mom1 += part[i] * part[1];
222 mom2 += part[i] * part[2];
223 mom3 += part[i] * part[3];
224 mom4 += part[i] * part[4];
225 mom5 += part[i] * part[5];
226 },
227 Kokkos::Sum<double>(loc_moment_ncent[i][0]),
228 Kokkos::Sum<double>(loc_moment_ncent[i][1]), Kokkos::Sum<double>(loc_moment_ncent[i][2]),
229 Kokkos::Sum<double>(loc_moment_ncent[i][3]), Kokkos::Sum<double>(loc_moment_ncent[i][4]),
230 Kokkos::Sum<double>(loc_moment_ncent[i][5]));
231 Kokkos::fence();
232 }
233 ippl::Comm->barrier();
234
235 Kokkos::parallel_reduce(
236 "calc moments of particle distr.", Nlocal,
237 KOKKOS_LAMBDA(
238 const int k, double& ekin) {
239 double gamma0 = 0;
240 double ekin0 = 0.0;
241
242 for(unsigned j=0; j<Dim; j++){
243 gamma0 += Pview(k)[j]*Pview(k)[j];
244 }
245 gamma0 = Kokkos::sqrt(gamma0+1.0);
246 ekin0 = (gamma0-1.0) * Mview(k);
247
248 ekin += (ekin0-mekin)*(ekin0-mekin);
249 }, Kokkos::Sum<double>(loc_std_mekin) );
250 Kokkos::fence();
251
252 MPI_Allreduce(
253 loc_moment_ncent, moment_ncent, 2 * Dim * 2 * Dim, MPI_DOUBLE, MPI_SUM, ippl::Comm->getCommunicator());
254
255 MPI_Allreduce(
256 &loc_std_mekin, &stdKineticEnergy_m, 1, MPI_DOUBLE, MPI_SUM, ippl::Comm->getCommunicator());
257
259
260 for (unsigned i = 0; i < 2 * Dim; i++) {
261 for (unsigned j = 0; j < 2 * Dim; j++) {
262 notCentMoments_m(i,j) = moment_ncent[i][j];
263 }
264 }
265
266 // compute emmitance, halo, ...
267 double perParticle = 1./(1.*Np);
268 Vector_t<double, 3> squaredEps, fac, sumRP;
269 unsigned int l = 0;
270 for (unsigned int i = 0; i < 3; ++ i, l += 2) {
271 double w1 = centroid_m[2 * i] * perParticle;
272 double w2 = moments_m(2 * i , 2 * i);
273 //not clear which components of non-central moments are computed, needs to be checked
274 double w3 = notCentMoments_m(l, l) * perParticle;
275 double w4 = notCentMoments_m(l + 1, l + 1) * perParticle;
276 double tmp = w2 - std::pow(w1, 2);
277
278 halo_m(i) = (w4 + w1 * (-4 * w3 + 3 * w1 * (tmp + w2))) / tmp;
280 }
281
282 //stdKineticEnergy_m = std::sqrt(localMoments[l++] * perParticle);
283 //totalCharge_m = localMoments[l++];
284 //totalMass_m = localMoments[l++];
285
286 for (unsigned int i = 0; i < 3; ++ i) {
287 sumRP(i) = notCentMoments_m(2 * i, 2 * i + 1) * perParticle - meanR_m(i) * meanP_m(i);
288 stdRP_m(i) = sumRP(i) / (stdR_m(i) * stdP_m(i));
289 squaredEps(i) = std::pow(stdR_m(i) * stdP_m(i), 2) - std::pow(sumRP(i), 2);
290 normalizedEps_m(i) = std::sqrt(std::max(squaredEps(i), 0.0));
291 }
292
293 double betaGamma = std::sqrt(std::pow(meanGamma_m, 2) - 1.0);
295
296}
297
298void DistributionMoments::computeMinMaxPosition(ippl::ParticleAttrib<Vector_t<double,3>>::view_type& Rview, size_t Nlocal)
299{
300 const int Dim = 3;
301
302 double rmax_loc[Dim];
303 double rmin_loc[Dim];
304 double rmax[Dim];
305 double rmin[Dim];
306
307 for(int i=0; i<Dim; i++){
308 rmin_loc[i] = 0.;
309 rmax_loc[i] = 0.;
310 }
311
312 for (unsigned d = 0; d < Dim; ++d) {
313 if (Nlocal > 0) {
314 Kokkos::parallel_reduce(
315 "rel max", Nlocal,
316 KOKKOS_LAMBDA(const int i, double& mm) {
317 double tmp_vel = Rview(i)[d];
318 mm = tmp_vel > mm ? tmp_vel : mm;
319 },
320 Kokkos::Max<double>(rmax_loc[d]));
321
322 Kokkos::parallel_reduce(
323 "rel min", ippl::getRangePolicy(Rview),
324 KOKKOS_LAMBDA(const int i, double& mm) {
325 double tmp_vel = Rview(i)[d];
326 mm = tmp_vel < mm ? tmp_vel : mm;
327 },
328 Kokkos::Min<double>(rmin_loc[d]));
329 }
330 }
331 Kokkos::fence();
332 MPI_Allreduce(rmax_loc, rmax, Dim, MPI_DOUBLE, MPI_MAX, ippl::Comm->getCommunicator());
333 MPI_Allreduce(rmin_loc, rmin, Dim, MPI_DOUBLE, MPI_MIN, ippl::Comm->getCommunicator());
334 ippl::Comm->barrier();
335
336 // store min and max R in class member variables
337 for (unsigned int i=0; i<Dim; i++) {
338 minR_m(i) = rmin[i];
339 maxR_m(i) = rmax[i];
340 }
341}
342
344 const std::vector<OpalParticle>::const_iterator& first,
345 const std::vector<OpalParticle>::const_iterator& last) {
346 // computeStatistics(first, last);
347}
348
349/*
350template <class InputIt>
351void DistributionMoments::computeMeans(const InputIt& first, const InputIt& last) {
352 unsigned int localNum = last - first;
353 std::vector<double> localStatistics(10);
354 std::vector<double> maxima(6, std::numeric_limits<double>::lowest());
355 for (InputIt it = first; it != last; ++it) {
356 OpalParticle const& particle = *it;
357 if (isParticleExcluded(particle)) {
358 --localNum;
359 continue;
360 }
361 for (unsigned int i = 0; i < 3; ++i) {
362 maxima[2 * i] = std::max(maxima[2 * i], particle[2 * i]);
363 maxima[2 * i + 1] =
364 std::max(maxima[2 * i + 1], -particle[2 * i]); // calculates the minimum
365 }
366
367 unsigned int l = 0;
368 for (unsigned int i = 0; i < 6; ++i, ++l) {
369 localStatistics[l] += particle[i];
370 }
371
372 // l = 6
373 localStatistics[l++] += particle.getTime();
374
375 double gamma = Util::getGamma(particle.getP());
376 double eKin = (gamma - 1.0) * particle.getMass();
377 localStatistics[l++] += eKin;
378 localStatistics[l++] += gamma;
379 }
380 localStatistics.back() = localNum;
381
382 ippl::Comm->allreduce(localStatistics.data(), localStatistics.size(), std::plus<double>());
383 ippl::Comm->allreduce(maxima.data(), 6, std::greater<double>());
384
385 totalNumParticles_m = localStatistics.back();
386
387 double perParticle = 1.0 / totalNumParticles_m;
388 unsigned int l = 0;
389
390 for (; l < 6; ++l) {
391 centroid_m[l] = localStatistics[l];
392 }
393 for (unsigned int i = 0; i < 3; ++i) {
394 meanR_m(i) = centroid_m[2 * i] * perParticle;
395 meanP_m(i) = centroid_m[2 * i + 1] * perParticle;
396 maxR_m(i) = maxima[2 * i];
397 minR_m(i) = -maxima[2 * i + 1];
398 }
399
400 meanTime_m = localStatistics[l++] * perParticle;
401
402 meanKineticEnergy_m = localStatistics[l++] * perParticle;
403 meanGamma_m = localStatistics[l++] * perParticle;
404}
405*/
406
407/* 2 * Dim centroids + Dim * ( 2 * Dim + 1 ) 2nd moments + 2 * Dim (3rd and 4th order moments)
408 * --> 1st order moments: 0, ..., 2 * Dim - 1
409 * --> 2nd order moments: 2 * Dim, ..., Dim * ( 2 * Dim + 1 )
410 * --> 3rd order moments: Dim * ( 2 * Dim + 1 ) + 1, ..., Dim * ( 2 * Dim + 1 ) + Dim
411 * (only, <x^3>, <y^3> and <z^3>)
412 * --> 4th order moments: Dim * ( 2 * Dim + 1 ) + Dim + 1, ..., Dim * ( 2 * Dim + 1 ) + 2 * Dim
413 *
414 * For a 6x6 matrix we have each 2nd order moment (except diagonal
415 * entries) twice. We only store the upper half of the matrix.
416
417template <class InputIt>
418void DistributionMoments::computeStatistics(const InputIt& first, const InputIt& last) {
419 reset();
420
421 computeMeans(first, last);
422
423 double perParticle = 1.0 / totalNumParticles_m;
424 std::vector<double> localStatistics(37);
425 for (InputIt it = first; it != last; ++it) {
426 OpalParticle const& particle = *it;
427 if (isParticleExcluded(particle)) {
428 continue;
429 }
430 unsigned int l = 6;
431 for (unsigned int i = 0; i < 6; ++i) {
432 localStatistics[i] += std::pow(particle[i] - centroid_m[i] * perParticle, 2);
433 for (unsigned int j = 0; j <= i; ++j, ++l) {
434 localStatistics[l] += particle[i] * particle[j];
435 }
436 }
437
438 localStatistics[l++] += std::pow(particle.getTime() - meanTime_m, 2);
439
440 for (unsigned int i = 0; i < 3; ++i, l += 2) {
441 double r2 = std::pow(particle[i], 2);
442 localStatistics[l] += r2 * particle[i];
443 localStatistics[l + 1] += r2 * r2;
444 }
445
446 double eKin = Util::getKineticEnergy(particle.getP(), particle.getMass());
447 localStatistics[l++] += std::pow(eKin - meanKineticEnergy_m, 2);
448 localStatistics[l++] += particle.getCharge();
449 localStatistics[l++] += particle.getMass();
450 }
451
452 ippl::Comm->allreduce(localStatistics.data(), localStatistics.size(), std::plus<double>());
453
454 fillMembers(localStatistics);
455
456 computePercentiles(first, last);
457}
458*/
459
460template <class InputIt>
461void DistributionMoments::computePercentiles(const InputIt& first, const InputIt& last) {
463 return;
464 }
465
466 std::vector<gsl_histogram*> histograms(3);
467 // For a normal distribution the number of exchanged data between the cores is minimized
468 // if the number of histogram bins follows the following formula. Since we can't know
469 // how many particles are in each bin for the real distribution we use this formula.
470 unsigned int numBins = 3.5 * std::pow(3, std::log10(totalNumParticles_m));
471
473 for (unsigned int d = 0; d < 3; ++d) {
474 maxR(d) = 1.0000001 * std::max(maxR_m[d] - meanR_m[d], meanR_m[d] - minR_m[d]);
475 histograms[d] = gsl_histogram_alloc(numBins);
476 gsl_histogram_set_ranges_uniform(histograms[d], 0.0, maxR(d));
477 }
478 for (InputIt it = first; it != last; ++it) {
479 OpalParticle const& particle = *it;
480 for (unsigned int d = 0; d < 3; ++d) {
481 gsl_histogram_increment(histograms[d], std::abs(particle[2 * d] - meanR_m[d]));
482 }
483 }
484
485 std::vector<int> localHistogramValues(3 * (numBins + 1)),
486 globalHistogramValues(3 * (numBins + 1));
487 for (unsigned int d = 0; d < 3; ++d) {
488 int j = 0;
489 size_t accumulated = 0;
490 std::generate(
491 localHistogramValues.begin() + d * (numBins + 1) + 1,
492 localHistogramValues.begin() + (d + 1) * (numBins + 1),
493 [&histograms, &d, &j, &accumulated]() {
494 accumulated += gsl_histogram_get(histograms[d], j++);
495 return accumulated;
496 });
497
498 gsl_histogram_free(histograms[d]);
499 }
500
501 ippl::Comm->allreduce(
502 localHistogramValues.data(), globalHistogramValues.data(), 3 * (numBins + 1),
503 std::plus<int>());
504
505 int numParticles68 = boost::numeric_cast<int>(
507 int numParticles95 = boost::numeric_cast<int>(
509 int numParticles99 = boost::numeric_cast<int>(
511 int numParticles99_99 = boost::numeric_cast<int>(
513
514 for (int d = 0; d < 3; ++d) {
515 unsigned int localNum = last - first, current = 0;
516 std::vector<Vector_t<double, 2>> oneDPhaseSpace(localNum);
517 for (InputIt it = first; it != last; ++it, ++current) {
518 OpalParticle const& particle = *it;
519 oneDPhaseSpace[current](0) = particle[2 * d];
520 oneDPhaseSpace[current](1) = particle[2 * d + 1];
521 }
522 std::sort(
523 oneDPhaseSpace.begin(), oneDPhaseSpace.end(),
524 [d, this](Vector_t<double, 2>& left, Vector_t<double, 2>& right) {
525 return std::abs(left[0] - meanR_m[d]) < std::abs(right[0] - meanR_m[d]);
526 });
527
528 iterator_t endSixtyEight, endNinetyFive, endNinetyNine, endNinetyNine_NinetyNine;
529 endSixtyEight = endNinetyFive = endNinetyNine = endNinetyNine_NinetyNine =
530 oneDPhaseSpace.end();
531
532 std::tie(sixtyEightPercentile_m[d], endSixtyEight) = determinePercentilesDetail(
533 oneDPhaseSpace.begin(), oneDPhaseSpace.end(), globalHistogramValues,
534 localHistogramValues, d, numParticles68);
535
536 std::tie(ninetyFivePercentile_m[d], endNinetyFive) = determinePercentilesDetail(
537 oneDPhaseSpace.begin(), oneDPhaseSpace.end(), globalHistogramValues,
538 localHistogramValues, d, numParticles95);
539
540 std::tie(ninetyNinePercentile_m[d], endNinetyNine) = determinePercentilesDetail(
541 oneDPhaseSpace.begin(), oneDPhaseSpace.end(), globalHistogramValues,
542 localHistogramValues, d, numParticles99);
543
544 std::tie(ninetyNine_NinetyNinePercentile_m[d], endNinetyNine_NinetyNine) =
546 oneDPhaseSpace.begin(), oneDPhaseSpace.end(), globalHistogramValues,
547 localHistogramValues, d, numParticles99_99);
548
550 computeNormalizedEmittance(oneDPhaseSpace.begin(), endSixtyEight);
552 computeNormalizedEmittance(oneDPhaseSpace.begin(), endNinetyFive);
554 computeNormalizedEmittance(oneDPhaseSpace.begin(), endNinetyNine);
556 computeNormalizedEmittance(oneDPhaseSpace.begin(), endNinetyNine_NinetyNine);
557 }
558}
559
590std::pair<double, DistributionMoments::iterator_t> DistributionMoments::determinePercentilesDetail(
592 const std::vector<int>& globalAccumulatedHistogram,
593 const std::vector<int>& localAccumulatedHistogram, unsigned int dimension,
594 int numRequiredParticles) const {
595 unsigned int numBins = globalAccumulatedHistogram.size() / 3;
596 double percentile = 0.0;
597 iterator_t endPercentile = end;
598 for (unsigned int i = 1; i < numBins; ++i) {
599 unsigned int idx = dimension * numBins + i;
600 if (globalAccumulatedHistogram[idx] > numRequiredParticles) {
601 iterator_t beginBin = begin + localAccumulatedHistogram[idx - 1];
602 iterator_t endBin = begin + localAccumulatedHistogram[idx];
603 unsigned int numMissingParticles =
604 numRequiredParticles - globalAccumulatedHistogram[idx - 1];
605 unsigned int shift = 2;
606 while (numMissingParticles == 0) {
607 beginBin = begin + localAccumulatedHistogram[idx - shift];
608 numMissingParticles =
609 numRequiredParticles - globalAccumulatedHistogram[idx - shift];
610 ++shift;
611 }
612
613 std::vector<unsigned int> numParticlesInBin(ippl::Comm->size() + 1);
614 numParticlesInBin[ippl::Comm->rank() + 1] = endBin - beginBin;
615
616 ippl::Comm->allreduce(
617 &(numParticlesInBin[1]), ippl::Comm->size(), std::plus<unsigned int>());
618
619 std::partial_sum(
620 numParticlesInBin.begin(), numParticlesInBin.end(), numParticlesInBin.begin());
621
622 std::vector<double> positions(numParticlesInBin.back());
623 std::transform(
624 beginBin, endBin, positions.begin() + numParticlesInBin[ippl::Comm->rank()],
625 [&dimension, this](Vector_t<double, 2> const& particle) {
626 return std::abs(particle[0] - meanR_m[dimension]);
627 });
628 ippl::Comm->allreduce(&(positions[0]), positions.size(), std::plus<double>());
629 std::sort(positions.begin(), positions.end());
630
631 percentile = (*(positions.begin() + numMissingParticles - 1)
632 + *(positions.begin() + numMissingParticles))
633 / 2;
634 for (iterator_t it = beginBin; it != endBin; ++it) {
635 if (std::abs((*it)[0] - meanR_m[dimension]) > percentile) {
636 return std::make_pair(percentile, it);
637 }
638 }
639 endPercentile = endBin;
640 break;
641 }
642 }
643 return std::make_pair(percentile, endPercentile);
644}
645
649 double localStatistics[] = {0.0, 0.0, 0.0, 0.0};
650 localStatistics[0] = end - begin;
651 for (iterator_t it = begin; it < end; ++it) {
652 const Vector_t<double, 2>& rp = *it;
653 localStatistics[1] += rp(0);
654 localStatistics[2] += rp(1);
655 localStatistics[3] += rp(0) * rp(1);
656 }
657 ippl::Comm->allreduce(&(localStatistics[0]), 4, std::plus<double>());
658
659 double numParticles = localStatistics[0];
660 double perParticle = 1 / localStatistics[0];
661 double meanR = localStatistics[1] * perParticle;
662 double meanP = localStatistics[2] * perParticle;
663 double RP = localStatistics[3] * perParticle;
664
665 localStatistics[0] = 0.0;
666 localStatistics[1] = 0.0;
667 for (iterator_t it = begin; it < end; ++it) {
668 const Vector_t<double, 2>& rp = *it;
669 localStatistics[0] += std::pow(rp(0) - meanR, 2);
670 localStatistics[1] += std::pow(rp(1) - meanP, 2);
671 }
672 ippl::Comm->allreduce(&(localStatistics[0]), 2, std::plus<double>());
673
674 double stdR = std::sqrt(localStatistics[0] / numParticles);
675 double stdP = std::sqrt(localStatistics[1] / numParticles);
676 double sumRP = RP - meanR * meanP;
677 double squaredEps = std::pow(stdR * stdP, 2) - std::pow(sumRP, 2);
678 double normalizedEps = std::sqrt(std::max(squaredEps, 0.0));
679
680 return normalizedEps;
681}
682
683void DistributionMoments::fillMembers(std::vector<double>& localMoments) {
684 /*
685 Vector_t<double, 3> squaredEps, fac, sumRP;
686 double perParticle = 1.0 / totalNumParticles_m;
687
688 unsigned int l = 0;
689 for (; l < 6; l += 2) {
690 stdR_m(l / 2) = std::sqrt(localMoments[l] / totalNumParticles_m);
691 stdP_m(l / 2) = std::sqrt(localMoments[l + 1] / totalNumParticles_m);
692 }
693
694 for (unsigned i = 0; i < moments_m.size1(); ++i) {
695 for (unsigned j = 0; j <= i; ++j, ++l) {
696 moments_m(i, j) = localMoments[l] * perParticle;
697 moments_m(j, i) = moments_m(i, j);
698 }
699 }
700
701 stdTime_m = localMoments[l++] * perParticle;
702
703 for (unsigned int i = 0; i < 3; ++i, l += 2) {
704 double w1 = centroid_m[2 * i] * perParticle;
705 double w2 = moments_m(2 * i, 2 * i);
706 double w3 = localMoments[l] * perParticle;
707 double w4 = localMoments[l + 1] * perParticle;
708 double tmp = w2 - std::pow(w1, 2);
709
710 halo_m(i) = (w4 + w1 * (-4 * w3 + 3 * w1 * (tmp + w2))) / tmp;
711 halo_m(i) -= Options::haloShift;
712 }
713
714 stdKineticEnergy_m = std::sqrt(localMoments[l++] * perParticle);
715
716 totalCharge_m = localMoments[l++];
717 totalMass_m = localMoments[l++];
718
719 for (unsigned int i = 0; i < 3; ++i) {
720 sumRP(i) = moments_m(2 * i, 2 * i + 1) - meanR_m(i) * meanP_m(i);
721 stdRP_m(i) = sumRP(i) / (stdR_m(i) * stdP_m(i));
722 squaredEps(i) = std::pow(stdR_m(i) * stdP_m(i), 2) - std::pow(sumRP(i), 2);
723 normalizedEps_m(i) = std::sqrt(std::max(squaredEps(i), 0.0));
724 }
725
726 double betaGamma = std::sqrt(std::pow(meanGamma_m, 2) - 1.0);
727 geometricEps_m = normalizedEps_m / Vector_t<double, 3>(betaGamma);
728 */
729}
730
732/*
733 double data[] = {0.0, 0.0};
734 // ada for (OpalParticle const& particle : bunch) {
735 // data[0] += Util::getKineticEnergy(particle.getP(), particle.getMass());
736 // }
737 data[1] = bunch.getLocalNum();
738 ippl::Comm->allreduce(data, 2, std::plus<double>());
739
740 meanKineticEnergy_m = data[0] / data[1];
741*/
742}
743
745 ippl::ParticleAttrib<Vector_t<double,3>>::view_type& Pview,
746 size_t Np,
747 size_t Nlocal,
748 double density){
749 Np = (Np == 0) ? 1 : Np; // Explanation: see DistributionMoments::computeMeans implementation
750
752 double loc_avgVel[3] = {};
753 double avgVel[3] = {};
754 double c = Physics::c;
755
756 // From P in \beta\gamma to get v in m/s: v = (P*c)/\gamma
757 Kokkos::parallel_reduce(
758 "calc moments of particle distr.", Nlocal,
759 KOKKOS_LAMBDA(
760 const int k, double& mom0, double& mom1, double& mom2){
761
762 double gamma0 = 0.0;
763 for(unsigned j=0; j<3; j++){
764 gamma0 += Pview(k)[j]*Pview(k)[j];
765 }
766 gamma0 = Kokkos::sqrt(gamma0+1.0);
767
768 mom0 += Pview(k)[0] * c / gamma0;
769 mom1 += Pview(k)[1] * c / gamma0;
770 mom2 += Pview(k)[2] * c / gamma0;
771 },
772 Kokkos::Sum<double>(loc_avgVel[0]),
773 Kokkos::Sum<double>(loc_avgVel[1]),
774 Kokkos::Sum<double>(loc_avgVel[2]));
775 Kokkos::fence();
776 ippl::Comm->barrier();
777
778 MPI_Allreduce(
779 loc_avgVel, avgVel, Dim, MPI_DOUBLE, MPI_SUM, ippl::Comm->getCommunicator());
780
781 for (unsigned i = 0; i < 3; i++) {
782 avgVel[i] = avgVel[i] / Np;
783 }
784
785 double tempAvg=0;
786 double loc_tempAvg = 0.0;
788 /*
789 for (OpalParticle const& particle_r : bunch_r) {
790 for (unsigned i = 0; i < 3; i++) {
791 tempAvg += std::pow(
792 (((particle_r.getP()[i] * Physics::c) / (Util::getGamma(particle_r.getP())))
793 - avgVel[i]),
794 2);
795 }
796 }*/
797
798 Kokkos::parallel_reduce(
799 "calc moments of particle distr.", Nlocal,
800 KOKKOS_LAMBDA(
801 const int k, double& mom0){
802
803 double gamma0 = 0.0;
804 for(unsigned j=0; j<3; j++){
805 gamma0 += Pview(k)[j]*Pview(k)[j];
806 }
807 gamma0 = Kokkos::sqrt(gamma0+1.0);
808
809 mom0 += Kokkos::pow( Pview(k)[0] * c / gamma0 - avgVel[0], 2);
810 mom0 += Kokkos::pow( Pview(k)[1] * c / gamma0 - avgVel[1], 2);
811 mom0 += Kokkos::pow( Pview(k)[2] * c / gamma0 - avgVel[2], 2);
812 },
813 Kokkos::Sum<double>(loc_tempAvg));
814 Kokkos::fence();
815 ippl::Comm->barrier();
816
817 MPI_Allreduce(
818 &loc_tempAvg, &tempAvg, 1, MPI_DOUBLE, MPI_SUM, ippl::Comm->getCommunicator());
819
820 // Compute the average temperature k_B T in units of kg m^2/s^2, where k_B is
821 // Boltzmann constant
822 temperature_m = (1.0 / 3.0) * Units::eV2kg * Units::GeV2eV * Physics::m_e * (tempAvg / Np);
823
825 std::sqrt((temperature_m * Physics::epsilon_0) / (density * std::pow(Physics::q_e, 2)));
826
827 computePlasmaParameter(density);
828}
829
831 // Plasma parameter: Average number of particles within the Debye sphere
832 plasmaParameter_m = (4.0 / 3.0) * Physics::pi * std::pow(debyeLength_m, 3) * density;
833}
834
836 std::fill(std::begin(centroid_m), std::end(centroid_m), 0.0);
837 meanR_m = 0.0;
838 meanP_m = 0.0;
839 stdR_m = 0.0;
840 stdP_m = 0.0;
841 stdRP_m = 0.0;
842 normalizedEps_m = 0.0;
843 geometricEps_m = 0.0;
844 halo_m = 0.0;
845
847 stdKineticEnergy_m = 0.0;
848
849 totalCharge_m = 0.0;
850 totalMass_m = 0.0;
852
861}
862
868
871 // FIXME After issue 287 is resolved this shouldn't be necessary anymore
872 return true;
873}
Inform * gmsg
Definition changes.cpp:7
ippl::Vector< T, Dim > Vector_t
PartBunch< T, Dim >::ConstIterator end(PartBunch< T, Dim > const &bunch)
PartBunch< T, Dim >::ConstIterator begin(PartBunch< T, Dim > const &bunch)
typename ippl::detail::ViewType< ippl::Vector< double, Dim >, 1 >::view_type view_type
constexpr unsigned Dim
Definition datatypes.h:6
constexpr double epsilon_0
The permittivity of vacuum in As/Vm.
Definition Physics.h:51
constexpr double q_e
The elementary charge in As.
Definition Physics.h:69
constexpr double m_e
The electron rest mass in GeV.
Definition Physics.h:78
constexpr double c
The velocity of light in m/s.
Definition Physics.h:45
constexpr double pi
The value of.
Definition Physics.h:30
constexpr double GeV2eV
Definition Units.h:68
constexpr double eV2kg
Definition Units.h:110
bool computePercentiles
Definition Options.cpp:105
double haloShift
The constant parameter C to shift halo, by < w^4 > / < w^2 > ^2 - C (w=x,y,z).
Definition Options.cpp:101
void fillMembers(std::vector< double > &)
Vector_t< double, 3 > sixtyEightPercentile_m
std::vector< Vector_t< double, 2 > >::const_iterator iterator_t
Vector_t< double, 3 > ninetyFivePercentile_m
Vector_t< double, 3 > minR_m
static const double percentileFourSigmasNormalDist_m
Vector_t< double, 3 > ninetyNinePercentile_m
void compute(const std::vector< OpalParticle >::const_iterator &, const std::vector< OpalParticle >::const_iterator &)
void computeMoments(ippl::ParticleAttrib< Vector_t< double, 3 > >::view_type &Rview, ippl::ParticleAttrib< Vector_t< double, 3 > >::view_type &Pview, ippl::ParticleAttrib< double >::view_type &Mview, size_t Np, size_t Nlocal)
Vector_t< double, 3 > ninetyNine_NinetyNinePercentile_m
Vector_t< double, 3 > normalizedEps_m
Vector_t< double, 3 > normalizedEps68Percentile_m
void computeMinMaxPosition(ippl::ParticleAttrib< Vector_t< double, 3 > >::view_type &Rview, size_t Nlcoal)
Vector_t< double, 3 > meanP_m
void computeDebyeLength(ippl::ParticleAttrib< Vector_t< double, 3 > >::view_type &Rview, ippl::ParticleAttrib< Vector_t< double, 3 > >::view_type &Pview, size_t Np, size_t Nlocal, double density)
double computeNormalizedEmittance(const iterator_t &begin, const iterator_t &end) const
Vector_t< double, 3 > halo_m
Vector_t< double, 3 > stdP_m
Vector_t< double, 3 > normalizedEps99Percentile_m
Vector_t< double, 3 > stdR_m
static const double percentileThreeSigmasNormalDist_m
Vector_t< double, 3 > meanR_m
void computeMeans(ippl::ParticleAttrib< Vector_t< double, 3 > >::view_type &Rview, ippl::ParticleAttrib< Vector_t< double, 3 > >::view_type &Pview, ippl::ParticleAttrib< double >::view_type &Mview, size_t Np, size_t Nlocal)
void computePercentiles(const InputIt &, const InputIt &)
static const double percentileTwoSigmasNormalDist_m
Vector_t< double, 3 > normalizedEps99_99Percentile_m
bool isParticleExcluded(const OpalParticle &) const
Vector_t< double, 3 > maxR_m
Vector_t< double, 3 > stdRP_m
std::pair< double, iterator_t > determinePercentilesDetail(const iterator_t &begin, const iterator_t &end, const std::vector< int > &globalAccumulatedHistogram, const std::vector< int > &localAccumulatedHistogram, unsigned int dimension, int numRequiredParticles) const
Vector_t< double, 3 > geometricEps_m
static const double percentileOneSigmaNormalDist_m
Vector_t< double, 6 > means_m
Vector_t< double, 6 > centroid_m
Vector_t< double, 3 > normalizedEps95Percentile_m