Geant4  10.01
G4SPSEneDistribution.hh
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
27 //
28 // MODULE: G4SPSEneDistribution.hh
29 //
30 // Version: 1.0
31 // Date: 5/02/04
32 // Author: Fan Lei
33 // Organisation: QinetiQ ltd.
34 // Customer: ESA/ESTEC
35 //
37 //
38 // CHANGE HISTORY
39 // --------------
40 //
41 //
42 // Version 1.0, 05/02/2004, Fan Lei, Created.
43 // Based on the G4GeneralParticleSource class in Geant4 v6.0
44 //
45 //
46 // 26/03/2014, Andrew Green.
47 // Modification to used STL vectors instead of C-style arrays. This should save some space,
48 // particularly when the blackbody function is not used. Also moved to dynamically allocated
49 // memory in the LinearInterpolation, ExpInterpolation and LogInterpolation functions. Again,
50 // this will save space if these functions are unused.
51 //
52 // 06/06/2014 A Dotti
53 // For thread safety: this is a shared object,
54 // mutex has been added to control access to shared resources (data members).
55 // in Getters and Setters, mutex is NOT used in GenerateOne because it is
56 // assumed that properties are not changed during event loop.
57 //
59 //
60 //
61 // Class Description:
62 //
63 // To generate the energy of a primary vertex according to the defined distribution
64 //
66 //
67 // MEMBER FUNCTIONS
68 // ----------------
69 //
70 // G4SPSEneDistribution ()
71 // Constructor: Initializes variables
72 //
73 // ~G4SPSEneDistribution ()
74 // Destructor:
75 //
76 // void SetEnergyDisType(G4String)
77 // Allows the user to choose the energy distribution type. The arguments
78 // are Mono (mono-energetic), Lin (linear), Pow (power-law), Exp
79 // (exponential), Gauss (gaussian), Brem (bremsstrahlung), BBody (black-body), Cdg
80 // (cosmic diffuse gamma-ray), User (user-defined), Arb (arbitrary
81 // point-wise), Epn (energy per nucleon).
82 //
83 // void SetEmin(G4double)
84 // Sets the minimum energy.
85 //
86 // void SetEmax(G4double)
87 // Sets the maximum energy.
88 //
89 // void SetMonoEnergy(G4double)
90 // Sets energy for mono-energetic distribution.
91 //
92 // void SetAlpha(G4double)
93 // Sets alpha for a power-law distribution.
94 //
95 // void SetTemp(G4double)
96 // Sets Temperature for a Brem or BBody distributions.
97 //
98 // void SetEzero(G4double)
99 // Sets Ezero for an exponential distribution.
100 //
101 // void SetGradient(G4double)
102 // Sets gradient for a linear distribution.
103 //
104 // void SetInterCept(G4double)
105 // Sets intercept for a linear distribution.
106 //
107 // void UserEnergyHisto(G4ThreeVector)
108 // Allows user to defined a histogram for the energy distribution.
109 //
110 // void ArbEnergyHisto(G4ThreeVector)
111 // Allows the user to define an Arbitrary set of points for the
112 // energy distribution.
113 //
114 // void EpnEnergyHisto(G4ThreeVector)
115 // Allows the user to define an Energy per nucleon histogram.
116 //
117 // void Calculate()
118 // Controls the calculation of Integral PDF for the Cdg and BBody
119 // distributions.
120 //
121 // void InputEnergySpectra(G4bool)
122 // Allows the user to choose between momentum and energy histograms
123 // for user-defined histograms and arbitrary point-wise spectr.
124 // The default is true (energy).
125 //
126 // void InputDifferentialSpectra(G4bool)
127 // Allows the user to choose between integral and differential
128 // distributions when using the arbitrary point-wise option.
129 //
130 // void ArbInterpolate(G4String)
131 // ArbInterpolate allows the user to specify the type of function to
132 // interpolate the Arbitrary points spectrum with.
133 //
134 // void SetBiasRndm (G4SPSRandomGenerator* a)
135 // Sets the biased random number generator
136 //
137 // G4double GenerateOne(G4ParticleDefinition*);
138 // Generate one random energy for the specified particle
139 //
140 // void ReSetHist(G4String);
141 // Re-sets the histogram for user defined distribution
142 //
143 // void SetVerbosity(G4int)
144 // Sets the verbosity level.
145 //
147 
148 #ifndef G4SPSEneDistribution_h
149 #define G4SPSEneDistribution_h 1
150 
152 #include "G4ParticleMomentum.hh"
153 #include "G4ParticleDefinition.hh"
154 #include "G4DataInterpolation.hh"
155 #include "G4Threading.hh"
156 #include "G4Cache.hh"
157 #include <vector>
158 
159 #include "G4SPSRandomGenerator.hh"
160 
162 public:
165 
166  void SetEnergyDisType(G4String);//
168  void SetEmin(G4double);//
169  G4double GetEmin();//
170  G4double GetArbEmin();//
171  void SetEmax(G4double);//
172  G4double GetEmax();//
173  G4double GetArbEmax();//
174  void SetMonoEnergy(G4double);//
175  void SetAlpha(G4double);//
176  void SetBiasAlpha(G4double);//
177  void SetTemp(G4double);//
179  void SetEzero(G4double);//
180  void SetGradient(G4double);//
181  void SetInterCept(G4double);//
186 
187  void InputEnergySpectra(G4bool);//
189  void ArbInterpolate(G4String);//
190  G4String GetIntType();//
191 
192  void Calculate();//
193 
195  // method to re-set the histograms
196  void ReSetHist(G4String);//
197  // Set the verbosity level.
198  void SetVerbosity(G4int a);//
199 
200  //x
202 
203  G4double GetMonoEnergy(); //Mono-energteic energy
204  G4double GetSE();// Standard deviation for Gaussion distrbution in energy
205  G4double Getalpha(); // alpha (pow)
206  G4double GetEzero(); // E0 (exp)
207  G4double GetTemp(); // Temp (bbody,brem)
208  G4double Getgrad(); // gradient and intercept for linear spectra
209  G4double Getcept(); //
210 
213 
216 
217 
218 private:
219  void LinearInterpolation();//
220  void LogInterpolation();//
221  void ExpInterpolation();//
222  void SplineInterpolation();//
223  void CalculateCdgSpectrum();//
224  void CalculateBbodySpectrum();//
225 
226  // The following methods generate energies according to the spectral
227  // parameters defined above.
228  void GenerateMonoEnergetic();//G4double& outputEne);//
229  void GenerateBiasPowEnergies();//G4double& outputEne,G4double& outputWeight);//
230  void GenerateGaussEnergies();//
231  void GenerateBremEnergies();//
232  void GenerateBbodyEnergies();//
233  void GenerateCdgEnergies();//
234  void GenUserHistEnergies();//
235  void GenEpnHistEnergies();//
236  void GenArbPointEnergies();//<<<<<<<<<<< DOES NOT WORK, REQUIRES UPDATE OF DATA MEMBERS.
237  void GenerateExpEnergies(G4bool);//
239  void GeneratePowEnergies(G4bool);//
240 
241  // converts energy per nucleon to energy.
242  void ConvertEPNToEnergy();
243 
244  void InitHists();//
245 
246 private:
247 
248  G4String EnergyDisType; // energy dis type Variable - Mono,Lin,Exp,etc
249  G4double weight; // particle weight //// NOT INVARIANT
250  G4double MonoEnergy; //Mono-energteic energy
251  G4double SE; // Standard deviation for Gaussion distrbution in energy
252  //Non invariant data members become G4Cache
253  G4double Emin, Emax; // emin and emax ////// NOT INVARIANT
254  G4double alpha, Ezero;// alpha (pow), E0 (exp) ////// NOT INVARIANT
255  G4double Temp; // Temp (bbody,brem)
256  G4double biasalpha; // biased power index
257  G4double grad, cept; // gradient and intercept for linear spectra ////// NOT INVARIANT
258  G4double prob_norm; // normalisation factor use in calculate the probability
259  G4bool Biased; // true - biased to power-law
260  G4bool EnergySpec; // true - energy spectra, false - momentum spectra
261  G4bool DiffSpec; // true - differential spec, false integral spec
262  //G4bool ApplyRig; // false no rigidity cutoff, true then apply one
263  //G4double ERig; // energy of rigidity cutoff
270  G4double CDGhist[3]; // cumulative histo for cdg
271 
272  //AG: Begin edit to use STL vectors.
273 // G4double BBHist[10001], Bbody_x[10001];
274  std::vector<G4double>* BBHist;
275  std::vector<G4double>* Bbody_x;
278 
279  //AG: Edit here to use dynamic memory, will save space inless these functions are used.
280  G4String IntType; // Interpolation type
281 // G4double Arb_grad[1024], Arb_cept[1024]; // grad and cept for 1024 segments AG: Switched to DMA
285 // G4double Arb_alpha[1024], Arb_Const[1024]; // alpha and constants AG: Switched to DMA
289 // G4double Arb_ezero[1024]; // ezero AG: Switched to DMA
292  G4double ArbEmin, ArbEmax; // Emin and Emax for the whole arb distribution used primarily for debug.
293 
295 
297 
298  // Verbosity
300 
302 
303  std::vector<G4DataInterpolation*> SplineInt;//[1024]; // holds Spline stuff required for sampling
304  G4DataInterpolation *Splinetemp; // holds a temp Spline used for calculating area
305 
306  G4Mutex mutex; // protect access to shared resources
307  //Thread local data (non-invariant during event loop).
308  //These are copied from master one at the beginning of generation
309  //of each event
310  struct threadLocal_t {
320  };
322 };
323 
324 #endif
325 
std::vector< G4DataInterpolation * > SplineInt
CLHEP::Hep3Vector G4ThreeVector
void ArbEnergyHisto(G4ThreeVector)
G4PhysicsOrderedFreeVector GetUserDefinedEnergyHisto()
G4double a
Definition: TRTMaterials.hh:39
std::vector< G4double > * Bbody_x
int G4int
Definition: G4Types.hh:78
G4PhysicsOrderedFreeVector GetArbEnergyHisto()
G4PhysicsOrderedFreeVector IPDFEnergyH
bool G4bool
Definition: G4Types.hh:79
G4double GetProbability(G4double)
G4DataInterpolation * Splinetemp
G4double GenerateOne(G4ParticleDefinition *)
void SetBiasRndm(G4SPSRandomGenerator *a)
std::vector< G4double > * BBHist
void EpnEnergyHisto(G4ThreeVector)
G4PhysicsOrderedFreeVector ArbEnergyH
G4int G4Mutex
Definition: G4Threading.hh:161
G4SPSRandomGenerator * eneRndm
G4PhysicsOrderedFreeVector UDefEnergyH
void UserEnergyHisto(G4ThreeVector)
void InputDifferentialSpectra(G4bool)
G4PhysicsOrderedFreeVector EpnEnergyH
void ArbEnergyHistoFile(G4String)
double G4double
Definition: G4Types.hh:76
G4PhysicsOrderedFreeVector IPDFArbEnergyH
G4Cache< threadLocal_t > threadLocalData
G4PhysicsOrderedFreeVector ZeroPhysVector