OPAL (Object Oriented Parallel Accelerator Library)  2021.1.99
OPAL
Distribution.h
Go to the documentation of this file.
1 //
2 // Class Distribution
3 // This class defines the initial beam that is injected or emitted into the simulation.
4 //
5 // Copyright (c) 2008 - 2020, Paul Scherrer Institut, Villigen PSI, Switzerland
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 OPAL_Distribution_HH
19 #define OPAL_Distribution_HH
20 
21 #include <fstream>
22 #include <string>
23 #include <vector>
24 
26 #include "Algorithms/PartData.h"
27 #include "Algorithms/Vektor.h"
28 #include "AppTypes/SymTenzor.h"
29 #include "Attributes/Attributes.h"
30 
32 
33 #include <gsl/gsl_histogram.h>
34 #include <gsl/gsl_qrng.h>
35 #include <gsl/gsl_rng.h>
36 
37 #ifdef WITH_UNIT_TESTS
38 #include <gtest/gtest_prod.h>
39 #endif
40 
41 class Beam;
42 class Beamline;
43 
44 template <class T, unsigned Dim>
45 class PartBunchBase;
46 
47 class PartBins;
48 class LaserProfile;
49 class H5PartWrapper;
50 
51 namespace DistrTypeT
52 {
62  };
63 }
64 
65 namespace EmissionModelT
66 {
69  NONEQUIL
70  };
71 }
72 
74 {
76  EVOVERC
77  };
78 }
79 
80 namespace Attrib
81 {
82  namespace Distribution
83  {
84  enum AttributesT {
95  W,
96  FE,
125  MX,
126  MY,
127  MZ,
128  MT,
140  R, // the correlation matrix (a la transport)
157  EX, // below is for the matched distribution
158  EY,
159  ET,
160  SECTOR, // Matched-Gauss: single sector or full machine
161  NSECTORS, // Matched-Gauss: number of sectors to average field
162  NSTEPS, // Matched-Gauss: number of steps for closed orbit finder
163  RGUESS, // Matched-Gauss: guess for closed orbit finder
164  DENERGY, // Matched-Gauss: energy step for closed orbit finder
170  // E2,
171  ID1, // special particle that the user can set
172  ID2, // special particle that the user can set
174  SIZE
175  };
176  }
177 
178  namespace Legacy
179  {
180  namespace Distribution
181  {
183  // DESCRIPTION OF THE DISTRIBUTION:
185  // DEBIN,
189  T,
190  PT,
191  // ALPHAX,
192  // ALPHAY,
193  // BETAX,
194  // BETAY,
195  // DX,
196  // DDX,
197  // DY,
198  // DDY,
199  SIZE
200  };
201  }
202  }
203 }
204 
205 
206 class Distribution: public Definition {
207 
208 public:
209 
210  Distribution();
211  virtual ~Distribution();
212 
213  virtual bool canReplaceBy(Object *object);
214  virtual Distribution *clone(const std::string &name);
215  virtual void execute();
216  virtual void update();
217  size_t getNumOfLocalParticlesToCreate(size_t n);
219  size_t numberOfParticles,
220  double current, const Beamline &bl);
222  std::vector<Distribution *> addedDistributions,
223  size_t &numberOfParticles);
224  void createOpalT(PartBunchBase<double, 3> *beam, size_t &numberOfParticles);
225  void doRestartOpalT(PartBunchBase<double, 3> *p, size_t Np, int restartStep, H5PartWrapper *h5wrapper);
226  void doRestartOpalCycl(PartBunchBase<double, 3> *p, size_t Np, int restartStep,
227  const int specifiedNumBunch, H5PartWrapper *h5wrapper);
228  size_t emitParticles(PartBunchBase<double, 3> *beam, double eZ);
229  double getPercentageEmitted() const;
230  static Distribution *find(const std::string &name);
231 
232  bool getIfDistEmitting();
234  double getMaxTOrZ();
235  double getMinTOrZ();
236  size_t getNumberOfEmissionSteps();
237  int getNumberOfEnergyBins();
238  double getEmissionDeltaT();
239  double getEnergyBinDeltaT();
240  double getWeight();
241 
242  double getTEmission();
243 
244  Vector_t get_pmean() const;
245 
246  std::string getTypeofDistribution();
247 
248  Inform &printInfo(Inform &os) const;
249 
250  bool Rebin();
251  void setDistToEmitted(bool emitted);
252  void setDistType();
253  void setSigmaR_m();
254  void setSigmaP_m(double massIneV);
255  void shiftBeam(double &maxTOrZ, double &minTOrZ);
256  double getEmissionTimeShift() const;
257 
259 
261 private:
262 #ifdef WITH_UNIT_TESTS
263  FRIEND_TEST(GaussTest, FullSigmaTest1);
264  FRIEND_TEST(GaussTest, FullSigmaTest2);
265  FRIEND_TEST(BinomialTest, FullSigmaTest1);
266  FRIEND_TEST(BinomialTest, FullSigmaTest2);
267 #endif
268 
269  Distribution(const std::string &name, Distribution *parent);
270 
271  // Not implemented.
272  Distribution(const Distribution &) = delete;
273  void operator=(const Distribution &) = delete;
274 
275  std::vector<double>& getXDist();
276  std::vector<double>& getBGxDist();
277  std::vector<double>& getYDist();
278  std::vector<double>& getBGyDist();
279  std::vector<double>& getTOrZDist();
280  std::vector<double>& getBGzDist();
281  void eraseXDist();
282  void eraseBGxDist();
283  void eraseYDist();
284  void eraseBGyDist();
285  void eraseTOrZDist();
286  void eraseBGzDist();
287 
288  void addDistributions();
289  void applyEmissionModel(double lowEnergyLimit, double &px, double &py, double &pz, std::vector<double> &additionalRNs);
290  void applyEmissModelAstra(double &px, double &py, double &pz, std::vector<double> &additionalRNs);
291  void applyEmissModelNone(double &pz);
292  void applyEmissModelNonEquil(double eZ, double &px, double &py, double &pz, std::vector<double> &additionalRNs);
299  void create(size_t &numberOfParticles, double massIneV, double charge);
300  void calcPartPerDist(size_t numberOfParticles);
302  void checkIfEmitted();
303  void checkParticleNumber(size_t &numberOfParticles);
305  size_t getNumberOfParticlesInFile(std::ifstream &inputFile);
306 
308  public:
310  { }
311 
312  virtual double get(double rand) = 0;
313  };
314 
316  public:
318  ami_m(rhs.ami_m)
319  {}
320 
322  { ami_m = 1.0 / a; }
323 
324  virtual double get(double rand);
325  private:
326  double ami_m;
327  };
328 
330  public:
331  virtual double get(double rand);
332  };
333 
334  void createDistributionBinomial(size_t numberOfParticles, double massIneV);
335  void createDistributionFlattop(size_t numberOfParticles, double massIneV);
336  void createDistributionMultiGauss(size_t numberOfParticles, double massIneV);
337  void createDistributionFromFile(size_t numberOfParticles, double massIneV);
338  void createDistributionGauss(size_t numberOfParticles, double massIneV);
339  void createMatchedGaussDistribution(size_t numberOfParticles,
340  double massIneV, double charge);
341  void sampleUniformDisk(gsl_qrng* quasiRandGen2D, double& x1, double& x2);
342  void fillEBinHistogram();
343  void fillParticleBins();
344  size_t findEBin(double tOrZ);
345  void generateAstraFlattopT(size_t numberOfParticles);
346  void generateBinomial(size_t numberOfParticles);
347  void generateFlattopLaserProfile(size_t numberOfParticles);
348  void generateFlattopT(size_t numberOfParticles);
349  void generateFlattopZ(size_t numberOfParticles);
350  void generateGaussZ(size_t numberOfParticles);
352  size_t numberOfParticles, double massIneV);
353  void generateLongFlattopT(size_t numberOfParticles);
354  void generateTransverseGauss(size_t numberOfParticles);
357  void printDist(Inform &os, size_t numberOfParticles) const;
358  void printDistBinomial(Inform &os) const;
359  void printDistFlattop(Inform &os) const;
360  void printDistMultiGauss(Inform &os) const;
361  void printDistFromFile(Inform &os) const;
362  void printDistGauss(Inform &os) const;
363  void printDistMatchedGauss(Inform &os) const;
364  void printEmissionModel(Inform &os) const;
365  void printEmissionModelAstra(Inform &os) const;
366  void printEmissionModelNone(Inform &os) const;
367  void printEmissionModelNonEquil(Inform &os) const;
368  void printEnergyBins(Inform &os) const;
369  void adjustPhaseSpace(double massIneV);
370  void reflectDistribution(size_t &numberOfParticles);
371  void scaleDistCoordinates();
373  gsl_qrng* selectRandomGenerator(std::string, unsigned int dimension);
374  void setAttributes();
375  void setDistParametersBinomial(double massIneV);
376  void setDistParametersFlattop(double massIneV);
377  void setDistParametersMultiGauss(double massIneV);
378  void setDistParametersGauss(double massIneV);
379  void setEmissionTime(double &maxT, double &minT);
384  void setupEnergyBins(double maxTOrZ, double minTOrZ);
385  void setupParticleBins(double massIneV, PartBunchBase<double, 3> *beam);
386  void shiftDistCoordinates(double massIneV);
387  void writeOutFileHeader();
388  void writeOutFileEmission();
389  void writeOutFileInjection();
390 
391  std::string distT_m;
393 
395 
396  bool emitting_m;
398 
401 
403  std::vector<Distribution *> addedDistributions_m;
404  std::vector<size_t> particlesPerDist_m;
405 
408 
410  double tEmission_m;
411  static const double percentTEmission_m;
415  double tBin_m;
424  gsl_histogram *energyBinHist_m;
426 
427  gsl_rng *randGen_m;
428 
429  // ASTRA and NONE photo emission model.
430  double pTotThermal_m;
432 
433  // NONEQUIL photo emission model.
435  double laserEnergy_m;
437  double cathodeTemp_m;
439 
440  std::vector<std::vector<double> > additionalRNs_m;
441 
444 
445  // Beam coordinate containers.
446  std::vector<double> xDist_m;
447  std::vector<double> pxDist_m;
448  std::vector<double> yDist_m;
449  std::vector<double> pyDist_m;
450  std::vector<double> tOrZDist_m;
451  std::vector<double> pzDist_m;
452 
453  // Initial coordinates for file write.
454  std::vector<double> xWrite_m;
455  std::vector<double> pxWrite_m;
456  std::vector<double> yWrite_m;
457  std::vector<double> pyWrite_m;
458  std::vector<double> tOrZWrite_m;
459  std::vector<double> pzWrite_m;
460  std::vector<size_t> binWrite_m;
461 
462  // for compatibility reasons
463  double avrgpz_m;
464 
465  //Distribution parameters.
467  double sigmaTRise_m;
468  double sigmaTFall_m;
476 
477  // Parameters multiGauss distribution.
478  double sepPeaks_m;
479  unsigned nPeaks_m;
480 
481  // Laser profile.
483  std::string laserImageName_m;
486 
487  // AAA This is for the matched distribution
488  double I_m;
489  double E_m;
490 
492  double tRise_m;
493  double tFall_m;
494  double sigmaRise_m;
495  double sigmaFall_m;
496  double cutoff_m;
497 };
498 
499 inline Inform &operator<<(Inform &os, const Distribution &d) {
500  return d.printInfo(os);
501 }
502 
503 inline
505  return pmean_m;
506 }
507 
508 inline
510  return distrTypeT_m;
511 }
512 
513 inline
515  return (double)totalNumberEmittedParticles_m / (double)totalNumberParticles_m;
516 }
517 
518 inline
521 }
522 
523 #endif // OPAL_Distribution_HH
Inform & operator<<(Inform &os, const Distribution &d)
Definition: Distribution.h:499
std::complex< double > a
const std::string name
std::string getString(const Attribute &attr)
Get string value.
Definition: Attributes.cpp:343
@ GUNGAUSSFLATTOPTH
Definition: Distribution.h:59
The base class for all OPAL definitions.
Definition: Definition.h:30
The base class for all OPAL objects.
Definition: Object.h:48
std::vector< Attribute > itsAttr
The object attributes.
Definition: Object.h:216
Particle reference data.
Definition: PartData.h:35
An abstract sequence of beam line components.
Definition: Beamline.h:34
void shiftDistCoordinates(double massIneV)
static Distribution * find(const std::string &name)
void setDistParametersBinomial(double massIneV)
std::vector< double > & getBGzDist()
std::vector< double > & getTOrZDist()
void fillEBinHistogram()
PartData particleRefData_m
Definition: Distribution.h:399
std::vector< double > xWrite_m
Definition: Distribution.h:454
Vector_t pmean_m
Total thermal momentum.
Definition: Distribution.h:431
int currentEnergyBin_m
Definition: Distribution.h:417
void generateFlattopT(size_t numberOfParticles)
double cathodeTemp_m
Cathode material Fermi energy (eV).
Definition: Distribution.h:437
void createOpalCycl(PartBunchBase< double, 3 > *beam, size_t numberOfParticles, double current, const Beamline &bl)
void adjustPhaseSpace(double massIneV)
std::vector< double > pxWrite_m
Definition: Distribution.h:455
void writeOutFileHeader()
double getMaxTOrZ()
void checkEmissionParameters()
double getMinTOrZ()
double sepPeaks_m
Definition: Distribution.h:478
std::vector< double > tOrZDist_m
Definition: Distribution.h:450
Vector_t get_pmean() const
Definition: Distribution.h:504
double sigmaTFall_m
Definition: Distribution.h:468
double tPulseLengthFWHM_m
Definition: Distribution.h:469
double getWeight()
EmissionModelT::EmissionModelT emissionModel_m
Emission Model.
Definition: Distribution.h:407
int getNumberOfEnergyBins()
virtual void execute()
Execute the command.
gsl_histogram * energyBinHist_m
Distribution energy bins.
Definition: Distribution.h:424
void setupEmissionModelAstra(PartBunchBase< double, 3 > *beam)
void printEnergyBins(Inform &os) const
std::vector< std::vector< double > > additionalRNs_m
Upper limit on emission energy distribution (eV).
Definition: Distribution.h:440
std::vector< double > pxDist_m
Definition: Distribution.h:447
std::vector< double > yWrite_m
Definition: Distribution.h:456
double getEmissionDeltaT()
void printDistBinomial(Inform &os) const
double cathodeFermiEnergy_m
Laser photon energy (eV).
Definition: Distribution.h:436
Distribution(const Distribution &)=delete
double emitEnergyUpperLimit_m
Cathode temperature (K).
Definition: Distribution.h:438
void calcPartPerDist(size_t numberOfParticles)
unsigned int numberOfDistributions_m
and list type for switch statements.
Definition: Distribution.h:394
void addDistributions()
void initializeBeam(PartBunchBase< double, 3 > *beam)
void doRestartOpalT(PartBunchBase< double, 3 > *p, size_t Np, int restartStep, H5PartWrapper *h5wrapper)
void checkIfEmitted()
void applyEmissModelNonEquil(double eZ, double &px, double &py, double &pz, std::vector< double > &additionalRNs)
std::vector< double > & getBGxDist()
void chooseInputMomentumUnits(InputMomentumUnitsT::InputMomentumUnitsT inputMoUnits)
Vector_t cutoffR_m
Definition: Distribution.h:472
void createOpalT(PartBunchBase< double, 3 > *beam, std::vector< Distribution * > addedDistributions, size_t &numberOfParticles)
void generateAstraFlattopT(size_t numberOfParticles)
void fillParticleBins()
size_t getNumberOfParticlesInFile(std::ifstream &inputFile)
DistrTypeT::DistrTypeT getType() const
Definition: Distribution.h:509
DistrTypeT::DistrTypeT distrTypeT_m
Distribution type. Declared as string.
Definition: Distribution.h:392
void createDistributionBinomial(size_t numberOfParticles, double massIneV)
void setupEmissionModelNone(PartBunchBase< double, 3 > *beam)
virtual bool canReplaceBy(Object *object)
Distribution can only be replaced by another distribution.
std::vector< Distribution * > addedDistributions_m
Vector of distributions to be added to this distribution.
Definition: Distribution.h:403
void doRestartOpalCycl(PartBunchBase< double, 3 > *p, size_t Np, int restartStep, const int specifiedNumBunch, H5PartWrapper *h5wrapper)
std::vector< double > tOrZWrite_m
Definition: Distribution.h:458
std::vector< double > & getXDist()
size_t getNumberOfEmissionSteps()
void injectBeam(PartBunchBase< double, 3 > *beam)
double sigmaRise_m
Definition: Distribution.h:494
std::string laserImageName_m
Definition: Distribution.h:483
Vector_t sigmaP_m
Definition: Distribution.h:471
void operator=(const Distribution &)=delete
void create(size_t &numberOfParticles, double massIneV, double charge)
virtual ~Distribution()
std::vector< double > pyDist_m
Definition: Distribution.h:449
int getLastEmittedEnergyBin()
std::string getTypeofDistribution()
Definition: Distribution.h:519
void createDistributionMultiGauss(size_t numberOfParticles, double massIneV)
double laserEnergy_m
Cathode material work function (eV).
Definition: Distribution.h:435
void setEmissionTime(double &maxT, double &minT)
Vector_t mBinomial_m
Definition: Distribution.h:474
gsl_rng * randGen_m
Definition: Distribution.h:427
void setupEnergyBins(double maxTOrZ, double minTOrZ)
void printEmissionModelNone(Inform &os) const
LaserProfile * laserProfile_m
Definition: Distribution.h:485
std::vector< double > yDist_m
Definition: Distribution.h:448
void generateFlattopZ(size_t numberOfParticles)
double cutoff_m
Definition: Distribution.h:496
std::vector< size_t > particlesPerDist_m
Definition: Distribution.h:404
void setNumberOfDistributions(unsigned int n)
Definition: Distribution.h:258
void generateMatchedGauss(const SigmaGenerator::matrix_t &, size_t numberOfParticles, double massIneV)
double getPercentageEmitted() const
Definition: Distribution.h:514
double tRise_m
time binned distribution with thermal energy
Definition: Distribution.h:492
void setSigmaP_m(double massIneV)
std::string laserProfileFileName_m
Definition: Distribution.h:482
std::vector< double > & getYDist()
void printDistFromFile(Inform &os) const
virtual Distribution * clone(const std::string &name)
Return a clone.
double getEmissionTimeShift() const
void printDistMultiGauss(Inform &os) const
void reflectDistribution(size_t &numberOfParticles)
std::vector< double > pzWrite_m
Definition: Distribution.h:459
virtual void update()
Update this object.
double sigmaFall_m
Definition: Distribution.h:495
void eraseTOrZDist()
void setAttributes()
void generateLongFlattopT(size_t numberOfParticles)
double cathodeWorkFunc_m
Definition: Distribution.h:434
void setDistParametersFlattop(double massIneV)
double sigmaTRise_m
Definition: Distribution.h:467
void printDistFlattop(Inform &os) const
void printDistMatchedGauss(Inform &os) const
static const double percentTEmission_m
Definition: Distribution.h:411
void writeOutFileEmission()
PartBins * energyBins_m
Definition: Distribution.h:423
void printEmissionModelAstra(Inform &os) const
void printEmissionModelNonEquil(Inform &os) const
void generateBinomial(size_t numberOfParticles)
std::vector< double > & getBGyDist()
void setupParticleBins(double massIneV, PartBunchBase< double, 3 > *beam)
void sampleUniformDisk(gsl_qrng *quasiRandGen2D, double &x1, double &x2)
double currentEmissionTime_m
Definition: Distribution.h:416
std::vector< size_t > binWrite_m
Definition: Distribution.h:460
double getEnergyBinDeltaT()
void createDistributionFromFile(size_t numberOfParticles, double massIneV)
size_t totalNumberParticles_m
Definition: Distribution.h:442
double pTotThermal_m
Random number generator.
Definition: Distribution.h:430
int numberOfSampleBins_m
Definition: Distribution.h:421
void setDistToEmitted(bool emitted)
void setDistParametersMultiGauss(double massIneV)
void createDistributionFlattop(size_t numberOfParticles, double massIneV)
size_t findEBin(double tOrZ)
void scaleDistCoordinates()
std::vector< double > pyWrite_m
Definition: Distribution.h:457
SymTenzor< double, 6 > correlationMatrix_m
Definition: Distribution.h:475
size_t getNumOfLocalParticlesToCreate(size_t n)
std::vector< double > pzDist_m
Definition: Distribution.h:451
double getTEmission()
double tEmission_m
Emission parameters.
Definition: Distribution.h:410
void generateFlattopLaserProfile(size_t numberOfParticles)
gsl_qrng * selectRandomGenerator(std::string, unsigned int dimension)
Select and allocate gsl random number generator.
void generateGaussZ(size_t numberOfParticles)
void printDist(Inform &os, size_t numberOfParticles) const
bool getIfDistEmitting()
void applyEmissModelNone(double &pz)
void writeOutFileInjection()
std::string distT_m
Definition: Distribution.h:391
size_t emitParticles(PartBunchBase< double, 3 > *beam, double eZ)
void applyEmissionModel(double lowEnergyLimit, double &px, double &py, double &pz, std::vector< double > &additionalRNs)
void applyEmissModelAstra(double &px, double &py, double &pz, std::vector< double > &additionalRNs)
double tFall_m
Definition: Distribution.h:493
void setupEmissionModelNonEquil()
double avrgpz_m
Definition: Distribution.h:463
void setupEmissionModel(PartBunchBase< double, 3 > *beam)
void printDistGauss(Inform &os) const
int currentSampleBin_m
Definition: Distribution.h:418
unsigned nPeaks_m
Definition: Distribution.h:479
InputMomentumUnitsT::InputMomentumUnitsT inputMoUnits_m
Definition: Distribution.h:466
void printEmissionModel(Inform &os) const
void setDistParametersGauss(double massIneV)
Vector_t cutoffP_m
Definition: Distribution.h:473
Inform & printInfo(Inform &os) const
void checkParticleNumber(size_t &numberOfParticles)
void createMatchedGaussDistribution(size_t numberOfParticles, double massIneV, double charge)
void createDistributionGauss(size_t numberOfParticles, double massIneV)
void shiftBeam(double &maxTOrZ, double &minTOrZ)
Vector_t sigmaR_m
Definition: Distribution.h:470
size_t totalNumberEmittedParticles_m
Definition: Distribution.h:443
std::vector< double > xDist_m
Definition: Distribution.h:446
double laserIntensityCut_m
Definition: Distribution.h:484
int numberOfEnergyBins_m
Definition: Distribution.h:419
void generateTransverseGauss(size_t numberOfParticles)
virtual double get(double rand)=0
MDependentBehavior(const MDependentBehavior &rhs)
Definition: Distribution.h:317
virtual double get(double rand)
virtual double get(double rand)
RealDiracMatrix::matrix_t matrix_t
Type for storing maps.
Definition: Beam.h:32
Definition: Inform.h:42