OPAL (Object Oriented Parallel Accelerator Library)  2021.1.99
OPAL
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 #include "DistributionMoments.h"
19 
20 #include "Utilities/Options.h"
21 #include "Utilities/Util.h"
22 
23 #include "Message/GlobalComm.h"
24 #include "Utility/Inform.h"
25 
26 #include "OpalParticle.h"
27 #include "PartBunchBase.h"
28 
29 extern Inform* gmsg;
30 
32 {
33  reset();
34 }
35 
37 {
38  computeStatistics(bunch.begin(), bunch.end());
39 }
40 
41 void DistributionMoments::compute(const std::vector<OpalParticle>::const_iterator &first,
42  const std::vector<OpalParticle>::const_iterator &last)
43 {
44  computeStatistics(first, last);
45 }
46 
47 template<class InputIt>
48 void DistributionMoments::computeMeans(const InputIt &first, const InputIt &last)
49 {
50  unsigned int localNum = last - first;
51  std::vector<double> localStatistics(10);
52  for (InputIt it = first; it != last; ++ it) {
53  OpalParticle const& particle = *it;
54  if (isParticleExcluded(particle)) {
55  -- localNum;
56  continue;
57  }
58  unsigned int l = 0;
59  for (unsigned int i = 0; i < 6; ++ i, ++ l) {
60  localStatistics[l] += particle[i];
61  }
62 
63  // l = 6
64  localStatistics[l++] += particle.getTime();
65 
66  double gamma = Util::getGamma(particle.getP());
67  double eKin = (gamma - 1.0) * particle.getMass();
68  localStatistics[l++] += eKin;
69  localStatistics[l++] += gamma;
70  }
71  localStatistics.back() = localNum;
72 
73  allreduce(localStatistics.data(), localStatistics.size(), std::plus<double>());
74  totalNumParticles_m = localStatistics.back();
75 
76  double perParticle = 1.0 / totalNumParticles_m;
77  unsigned int l = 0;
78 
79  for (; l < 6; ++ l) {
80  centroid_m[l] = localStatistics[l];
81  }
82  for (unsigned int i = 0; i < 3; ++ i) {
83  meanR_m(i) = centroid_m[2 * i] * perParticle;
84  meanP_m(i) = centroid_m[2 * i + 1] * perParticle;
85  }
86 
87  meanTime_m = localStatistics[l++] * perParticle;
88 
89  meanKineticEnergy_m = localStatistics[l++] * perParticle;
90  meanGamma_m = localStatistics[l++] * perParticle;
91 }
92 
93 /* 2 * Dim centroids + Dim * ( 2 * Dim + 1 ) 2nd moments + 2 * Dim (3rd and 4th order moments)
94  * --> 1st order moments: 0, ..., 2 * Dim - 1
95  * --> 2nd order moments: 2 * Dim, ..., Dim * ( 2 * Dim + 1 )
96  * --> 3rd order moments: Dim * ( 2 * Dim + 1 ) + 1, ..., Dim * ( 2 * Dim + 1 ) + Dim
97  * (only, <x^3>, <y^3> and <z^3>)
98  * --> 4th order moments: Dim * ( 2 * Dim + 1 ) + Dim + 1, ..., Dim * ( 2 * Dim + 1 ) + 2 * Dim
99  *
100  * For a 6x6 matrix we have each 2nd order moment (except diagonal
101  * entries) twice. We only store the upper half of the matrix.
102  */
103 template<class InputIt>
104 void DistributionMoments::computeStatistics(const InputIt &first, const InputIt &last)
105 {
106  reset();
107 
108  computeMeans(first, last);
109 
110  double perParticle = 1.0 / totalNumParticles_m;
111  std::vector<double> localStatistics(37);
112  for (InputIt it = first; it != last; ++ it) {
113  OpalParticle const& particle = *it;
114  if (isParticleExcluded(particle)) {
115  continue;
116  }
117  unsigned int l = 6;
118  for (unsigned int i = 0; i < 6; ++ i) {
119  localStatistics[i] += std::pow(particle[i] - centroid_m[i] * perParticle, 2);
120  for (unsigned int j = 0; j <= i; ++ j, ++ l) {
121  localStatistics[l] += particle[i] * particle[j];
122  }
123  }
124 
125  localStatistics[l++] += std::pow(particle.getTime() - meanTime_m, 2);
126 
127  for (unsigned int i = 0; i < 3; ++ i, l += 2) {
128  double r2 = std::pow(particle[i], 2);
129  localStatistics[l] += r2 * particle[i];
130  localStatistics[l + 1] += r2 * r2;
131  }
132 
133  double eKin = Util::getKineticEnergy(particle.getP(), particle.getMass());
134  localStatistics[l++] += std::pow(eKin - meanKineticEnergy_m, 2);
135  localStatistics[l++] += particle.getCharge();
136  localStatistics[l++] += particle.getMass();
137  }
138 
139  allreduce(localStatistics.data(), localStatistics.size(), std::plus<double>());
140 
141  fillMembers(localStatistics);
142 }
143 
144 void DistributionMoments::fillMembers(std::vector<double> const& localMoments) {
145  Vector_t squaredEps, fac, sumRP;
146  double perParticle = 1.0 / totalNumParticles_m;
147 
148  unsigned int l = 0;
149  for (; l < 6; l += 2) {
150  stdR_m(l / 2) = std::sqrt(localMoments[l] * perParticle);
151  stdP_m(l / 2) = std::sqrt(localMoments[l + 1] * perParticle);
152  }
153 
154  for (unsigned int i = 0; i < 6; ++ i) {
155  for (unsigned int j = 0; j <= i; ++ j, ++ l) {
156  moments_m(i, j) = localMoments[l];
157  moments_m(j, i) = moments_m(i, j);
158  }
159  }
160 
161  stdTime_m = localMoments[l++] * perParticle;
162 
163  for (unsigned int i = 0; i < 3; ++ i, l += 2) {
164  double w1 = centroid_m[2 * i] * perParticle;
165  double w2 = moments_m(2 * i , 2 * i) * perParticle;
166  double w3 = localMoments[l] * perParticle;
167  double w4 = localMoments[l + 1] * perParticle;
168  double tmp = w2 - std::pow(w1, 2);
169 
170  halo_m(i) = (w4 + w1 * (-4 * w3 + 3 * w1 * (tmp + w2))) / tmp;
172  }
173 
174  stdKineticEnergy_m = std::sqrt(localMoments[l++] * perParticle);
175 
176  totalCharge_m = localMoments[l++];
177  totalMass_m = localMoments[l++];
178 
179  for (unsigned int i = 0; i < 3; ++ i) {
180  sumRP(i) = moments_m(2 * i, 2 * i + 1) * perParticle - meanR_m(i) * meanP_m(i);
181  stdRP_m(i) = sumRP(i) / (stdR_m(i) * stdP_m(i));
182  squaredEps(i) = std::pow(stdR_m(i) * stdP_m(i), 2) - std::pow(sumRP(i), 2);
183  normalizedEps_m(i) = std::sqrt(std::max(squaredEps(i), 0.0));
184  }
185 
186  double betaGamma = std::sqrt(std::pow(meanGamma_m, 2) - 1.0);
187  geometricEps_m = normalizedEps_m / Vector_t(betaGamma);
188 }
189 
191 {
192  double data[] = {0.0, 0.0};
193  for (OpalParticle const& particle: bunch) {
194  data[0] += Util::getKineticEnergy(particle.getP(), particle.getMass());
195  }
196  data[1] = bunch.getLocalNum();
197  allreduce(data, 2, std::plus<double>());
198 
199  meanKineticEnergy_m = data[0] / data[1];
200 }
201 
203 {
204  std::fill(std::begin(centroid_m), std::end(centroid_m), 0.0);
205  meanR_m = 0.0;
206  meanP_m = 0.0;
207  stdR_m = 0.0;
208  stdP_m = 0.0;
209  stdRP_m = 0.0;
210  normalizedEps_m = 0.0;
211  geometricEps_m = 0.0;
212  halo_m = 0.0;
213 
214  meanKineticEnergy_m = 0.0;
215  stdKineticEnergy_m = 0.0;
216  Dx_m = 0.0;
217  DDx_m = 0.0;
218  Dy_m = 0.0;
219  DDy_m = 0.0;
221 
222  totalCharge_m = 0.0;
223  totalMass_m = 0.0;
225 }
226 
228 {
229  //FIXME After issue 287 is resolved this shouldn't be necessary anymore
230  return !Options::amr
232  && particle.getId() == 0;
233 }
Tps< T > pow(const Tps< T > &x, int y)
Integer power.
Definition: TpsMath.h:76
Tps< T > sqrt(const Tps< T > &x)
Square root.
Definition: TpsMath.h:91
Inform * gmsg
Definition: Main.cpp:62
PartBunchBase< T, Dim >::ConstIterator end(PartBunchBase< T, Dim > const &bunch)
PartBunchBase< T, Dim >::ConstIterator begin(PartBunchBase< T, Dim > const &bunch)
T::PETE_Expr_t::PETE_Return_t max(const PETE_Expr< T > &expr, NDIndex< D > &loc)
Definition: ReductionLoc.h:84
void allreduce(const T *input, T *output, int count, Op op)
Definition: GlobalComm.hpp:510
double haloShift
The constant parameter C to shift halo, by < w^4 > / < w^2 > ^2 - C (w=x,y,z)
Definition: Options.cpp:107
bool amr
Enable AMR if true.
Definition: Options.cpp:99
double getKineticEnergy(Vector_t p, double mass)
Definition: Util.h:37
double getGamma(Vector_t p)
Definition: Util.h:27
size_t getLocalNum() const
ConstIterator end() const
ConstIterator begin() const
bool isInOPALCyclMode()
Definition: OpalData.cpp:271
static OpalData * getInstance()
Definition: OpalData.cpp:195
FMatrix< double, 6, 6 > moments_m
void computeMeans(const InputIt &, const InputIt &)
void computeMeanKineticEnergy(PartBunchBase< double, 3 > const &)
unsigned int totalNumParticles_m
void compute(const std::vector< OpalParticle >::const_iterator &, const std::vector< OpalParticle >::const_iterator &)
void computeStatistics(const InputIt &, const InputIt &)
bool isParticleExcluded(const OpalParticle &) const
void fillMembers(std::vector< double > const &)
const Vector_t & getP() const
Get momentum.
Definition: OpalParticle.h:231
double getCharge() const
Get charge in Coulomb.
Definition: OpalParticle.h:243
double getMass() const
Get mass in GeV/c^2.
Definition: OpalParticle.h:249
double getTime() const
Get time.
Definition: OpalParticle.h:237
int64_t getId() const
Get the id of the particle.
Definition: OpalParticle.h:176
Definition: Inform.h:42
Vektor< double, 3 > Vector_t
Definition: Vektor.h:6