OPALX (Object Oriented Parallel Accelerator Library for Exascale) MINIorX
OPALX
DistributionMoments.h
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#ifndef DISTRIBUTIONMOMENTS_H
19#define DISTRIBUTIONMOMENTS_H
20
21#include "Ippl.h"
22#include <Kokkos_Core.hpp>
24#include "Physics/Physics.h"
25#include "Physics/Units.h"
26
27#include <vector>
28
29template <typename T, unsigned Dim = 3>
31
32typedef typename std::pair<Vector_t<double, 3>, Vector_t<double, 3>> VectorPair_t;
33
34typedef boost::numeric::ublas::matrix<double> matrix_t;
35
36class OpalParticle;
37
39public:
41 void foo();
42 void compute(
43 const std::vector<OpalParticle>::const_iterator&,
44 const std::vector<OpalParticle>::const_iterator&);
48 size_t Np,
49 size_t Nlocal);
53 size_t Np,
54 size_t Nlocal,
55 double density);
56 void computePlasmaParameter(double);
57
69
78
79 double getMeanTime() const;
80 double getStdTime() const;
81 double getMeanGamma() const;
82 double getMeanGammaZ() const;
83 double getMeanKineticEnergy() const;
84 double getTemperature() const;
85 double getDebyeLength() const;
86 double getPlasmaParameter() const;
87 double getStdKineticEnergy() const;
88 double getDx() const;
89 double getDDx() const;
90 double getDy() const;
91 double getDDy() const;
94 matrix_t getMoments6x6() const;
95 double getTotalCharge() const;
96 double getTotalMass() const;
97 double getTotalNumParticles() const;
98
102 size_t Np,
103 size_t Nlocal);
104private:
105 bool isParticleExcluded(const OpalParticle&) const;
106
107 //template <class InputIt>
108 //void computeMeans(const InputIt&, const InputIt&);
109 // template <class InputIt>
110 // void computeStatistics(const InputIt&, const InputIt&);
111 template <class InputIt>
112 void computePercentiles(const InputIt&, const InputIt&);
113 using iterator_t = std::vector<Vector_t<double, 2>>::const_iterator;
114 std::pair<double, iterator_t> determinePercentilesDetail(
115 const iterator_t& begin, const iterator_t& end,
116 const std::vector<int>& globalAccumulatedHistogram,
117 const std::vector<int>& localAccumulatedHistogram, unsigned int dimension,
118 int numRequiredParticles) const;
119 double computeNormalizedEmittance(const iterator_t& begin, const iterator_t& end) const;
120
121 void fillMembers(std::vector<double>&);
122
123 void reset();
124
126
145
147 double stdTime_m;
155
160
164
169};
170
172 return meanR_m;
173}
174
176 return stdR_m;
177}
178
180 return meanP_m;
181}
182
184 return stdP_m;
185}
186
188 return normalizedEps_m;
189}
190
192 return geometricEps_m;
193}
194
196 return stdRP_m;
197}
198
200 return halo_m;
201}
202
204 return minR_m;
205}
206
208 return maxR_m;
209}
210
211inline double DistributionMoments::getMeanTime() const {
212 return meanTime_m;
213}
214
215inline double DistributionMoments::getStdTime() const {
216 return stdTime_m;
217}
218
220 return meanGamma_m;
221}
222
224 return meanGammaZ_m;
225}
226
228 return meanKineticEnergy_m;
229}
230
231// Compute and return the value of temperature in K
234}
236 return debyeLength_m;
237}
239 return plasmaParameter_m;
240}
241
243 return stdKineticEnergy_m;
244}
245
246inline double DistributionMoments::getDx() const {
247 return moments_m(0, 5);
248}
249
250inline double DistributionMoments::getDDx() const {
251 return moments_m(1, 5);
252}
253
254inline double DistributionMoments::getDy() const {
255 return moments_m(2, 5);
256}
257
258inline double DistributionMoments::getDDy() const {
259 return moments_m(3, 5);
260}
261
263 return centroid_m;
264}
265
267 return means_m;
268}
269
271 return moments_m;
272}
273
275 return totalCharge_m;
276}
277
279 return totalMass_m;
280}
281
283 return totalNumParticles_m;
284}
285
288}
289
292}
293
296}
297
300}
301
304}
305
308}
309
312}
313
316}
317
319 Vector_t<double, 3> maxDistance;
320 for (unsigned int i = 0; i < 3; ++i) {
321 maxDistance[i] = std::max(std::abs(maxR_m[i]), std::abs(minR_m[i]));
322 }
323 return maxDistance;
324}
325
326#endif
boost::numeric::ublas::matrix< double > matrix_t
Definition: BoostMatrix.h:23
std::pair< Vector_t< double, 3 >, Vector_t< double, 3 > > VectorPair_t
boost::numeric::ublas::matrix< double > matrix_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 double kB
Boltzman's constant in eV/K.
Definition: Physics.h:60
constexpr double c
The velocity of light in m/s.
Definition: Physics.h:45
constexpr double eV2kg
Definition: Units.h:110
KOKKOS_INLINE_FUNCTION Vector< T, Dim > max(const Vector< T, Dim > &a, const Vector< T, Dim > &b)
Definition: Vector.hpp:226
void fillMembers(std::vector< double > &)
Vector_t< double, 3 > sixtyEightPercentile_m
double getMeanKineticEnergy() const
Vector_t< double, 3 > get68Percentile() const
Vector_t< double, 3 > getMinPosition() const
Vector_t< double, 3 > getMeanMomentum() const
Vector_t< double, 3 > ninetyFivePercentile_m
Vector_t< double, 3 > getStandardDeviationPosition() const
Vector_t< double, 3 > getNormalizedEmittance() const
double getMeanGamma() const
void computeDebyeLength(ippl::ParticleAttrib< Vector_t< double, 3 > >::view_type &Pview, size_t Np, size_t Nlocal, double density)
Vector_t< double, 3 > get99_99Percentile() const
void computePlasmaParameter(double)
std::vector< Vector_t< double, 2 > >::const_iterator iterator_t
double getTotalCharge() const
Vector_t< double, 3 > minR_m
static const double percentileFourSigmasNormalDist_m
double getTemperature() const
matrix_t getMoments6x6() const
double getDebyeLength() const
unsigned int totalNumParticles_m
Vector_t< double, 3 > getHalo() const
Vector_t< double, 3 > ninetyNinePercentile_m
void compute(const std::vector< OpalParticle >::const_iterator &, const std::vector< OpalParticle >::const_iterator &)
double getStdKineticEnergy() const
Vector_t< double, 3 > getNormalizedEmittance99Percentile() const
Vector_t< double, 3 > getStandardDeviationMomentum() const
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 > getMeanPosition() const
Vector_t< double, 3 > getStandardDeviationRP() const
Vector_t< double, 3 > ninetyNine_NinetyNinePercentile_m
double getTotalMass() const
double getPlasmaParameter() const
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
Vector_t< double, 3 > getNormalizedEmittance99_99Percentile() const
double computeNormalizedEmittance(const iterator_t &begin, const iterator_t &end) const
Vector_t< double, 3 > halo_m
double getTotalNumParticles() const
Vector_t< double, 3 > stdP_m
Vector_t< double, 3 > normalizedEps99Percentile_m
Vector_t< double, 3 > getMaxPosition() const
Vector_t< double, 3 > stdR_m
static const double percentileThreeSigmasNormalDist_m
Vector_t< double, 3 > getNormalizedEmittance95Percentile() const
Vector_t< double, 3 > getGeometricEmittance() const
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 > get95Percentile() const
Vector_t< double, 6 > getCentroid() const
Vector_t< double, 3 > getNormalizedEmittance68Percentile() const
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
Computes the percentile and the range of all local particles that are contained therein.
Vector_t< double, 3 > geometricEps_m
static const double percentileOneSigmaNormalDist_m
Vector_t< double, 3 > get99Percentile() const
double getMeanTime() const
Vector_t< double, 6 > getMeans() const
Vector_t< double, 6 > means_m
Vector_t< double, 6 > centroid_m
Vector_t< double, 3 > getMaxR() const
double getMeanGammaZ() const
Vector_t< double, 3 > normalizedEps95Percentile_m
typename detail::ViewType< T, 1, Properties... >::view_type view_type