Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4SPSEneDistribution.hh
이 파일의 문서화 페이지로 가기
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 //
58 // 24/11/2017 Fan Lei
59 // Added cutoff power-law distribution option. Implementation is similar to that of the BlackBody one.
60 //
62 //
63 //
64 // Class Description:
65 //
66 // To generate the energy of a primary vertex according to the defined distribution
67 //
69 //
70 // MEMBER FUNCTIONS
71 // ----------------
72 //
73 // G4SPSEneDistribution ()
74 // Constructor: Initializes variables
75 //
76 // ~G4SPSEneDistribution ()
77 // Destructor:
78 //
79 // void SetEnergyDisType(G4String)
80 // Allows the user to choose the energy distribution type. The arguments
81 // are Mono (mono-energetic), Lin (linear), Pow (power-law), Exp
82 // (exponential), Gauss (gaussian), Brem (bremsstrahlung), BBody (black-body), Cdg
83 // (cosmic diffuse gamma-ray), User (user-defined), Arb (arbitrary
84 // point-wise), Epn (energy per nucleon).
85 //
86 // void SetEmin(G4double)
87 // Sets the minimum energy.
88 //
89 // void SetEmax(G4double)
90 // Sets the maximum energy.
91 //
92 // void SetMonoEnergy(G4double)
93 // Sets energy for mono-energetic distribution.
94 //
95 // void SetAlpha(G4double)
96 // Sets alpha for a power-law distribution.
97 //
98 // void SetTemp(G4double)
99 // Sets Temperature for a Brem or BBody distributions.
100 //
101 // void SetEzero(G4double)
102 // Sets Ezero for an exponential distribution.
103 //
104 // void SetGradient(G4double)
105 // Sets gradient for a linear distribution.
106 //
107 // void SetInterCept(G4double)
108 // Sets intercept for a linear distribution.
109 //
110 // void UserEnergyHisto(G4ThreeVector)
111 // Allows user to defined a histogram for the energy distribution.
112 //
113 // void ArbEnergyHisto(G4ThreeVector)
114 // Allows the user to define an Arbitrary set of points for the
115 // energy distribution.
116 //
117 // void EpnEnergyHisto(G4ThreeVector)
118 // Allows the user to define an Energy per nucleon histogram.
119 //
120 // void Calculate()
121 // Controls the calculation of Integral PDF for the Cdg and BBody
122 // distributions.
123 //
124 // void InputEnergySpectra(G4bool)
125 // Allows the user to choose between momentum and energy histograms
126 // for user-defined histograms and arbitrary point-wise spectr.
127 // The default is true (energy).
128 //
129 // void InputDifferentialSpectra(G4bool)
130 // Allows the user to choose between integral and differential
131 // distributions when using the arbitrary point-wise option.
132 //
133 // void ArbInterpolate(G4String)
134 // ArbInterpolate allows the user to specify the type of function to
135 // interpolate the Arbitrary points spectrum with.
136 //
137 // void SetBiasRndm (G4SPSRandomGenerator* a)
138 // Sets the biased random number generator
139 //
140 // G4double GenerateOne(G4ParticleDefinition*);
141 // Generate one random energy for the specified particle
142 //
143 // void ReSetHist(G4String);
144 // Re-sets the histogram for user defined distribution
145 //
146 // void SetVerbosity(G4int)
147 // Sets the verbosity level.
148 //
150 
151 #ifndef G4SPSEneDistribution_h
152 #define G4SPSEneDistribution_h 1
153 
155 #include "G4ParticleMomentum.hh"
156 #include "G4ParticleDefinition.hh"
157 #include "G4DataInterpolation.hh"
158 #include "G4Threading.hh"
159 #include "G4Cache.hh"
160 #include <vector>
161 
162 #include "G4SPSRandomGenerator.hh"
163 
175 public:
178 
179  void SetEnergyDisType(G4String);//
181  void SetEmin(G4double);//
182  G4double GetEmin();//
183  G4double GetArbEmin();//
184  void SetEmax(G4double);//
185  G4double GetEmax();//
186  G4double GetArbEmax();//
187  void SetMonoEnergy(G4double);//
188  void SetAlpha(G4double);//
189  void SetBiasAlpha(G4double);//
190  void SetTemp(G4double);//
192  void SetEzero(G4double);//
193  void SetGradient(G4double);//
194  void SetInterCept(G4double);//
199 
200  void InputEnergySpectra(G4bool);//
202  void ArbInterpolate(G4String);//
203  G4String GetIntType();//
204 
205  void Calculate();//
206 
208  // method to re-set the histograms
209  void ReSetHist(G4String);//
210  // Set the verbosity level.
211  void SetVerbosity(G4int a);//
212 
213  //x
215 
216  G4double GetMonoEnergy(); //Mono-energteic energy
217  G4double GetSE();// Standard deviation for Gaussion distrbution in energy
218  G4double Getalpha(); // alpha (pow)
219  G4double GetEzero(); // E0 (exp)
220  G4double GetTemp(); // Temp (bbody,brem)
221  G4double Getgrad(); // gradient and intercept for linear spectra
222  G4double Getcept(); //
223 
226 
229 
230 
231 private:
232  void LinearInterpolation();//
233  void LogInterpolation();//
234  void ExpInterpolation();//
235  void SplineInterpolation();//
236  void CalculateCdgSpectrum();//
237  void CalculateBbodySpectrum();//
238  void CalculateCPowSpectrum();//
239 
240  // The following methods generate energies according to the spectral
241  // parameters defined above.
242  void GenerateMonoEnergetic();//G4double& outputEne);//
243  void GenerateBiasPowEnergies();//G4double& outputEne,G4double& outputWeight);//
244  void GenerateGaussEnergies();//
245  void GenerateBremEnergies();//
246  void GenerateBbodyEnergies();//
247  void GenerateCdgEnergies();//
248  void GenUserHistEnergies();//
249  void GenEpnHistEnergies();//
250  void GenArbPointEnergies();//<<<<<<<<<<< DOES NOT WORK, REQUIRES UPDATE OF DATA MEMBERS.
251  void GenerateExpEnergies(G4bool);//
253  void GeneratePowEnergies(G4bool);//
254  void GenerateCPowEnergies();//
255 
256  // converts energy per nucleon to energy.
257  void ConvertEPNToEnergy();
258 
259  void BBInitHists();//
260  void CPInitHists();//
261 
262 private:
263 
264  G4String EnergyDisType; // energy dis type Variable - Mono,Lin,Exp,etc
265  G4double weight; // particle weight //// NOT INVARIANT
266  G4double MonoEnergy; //Mono-energteic energy
267  G4double SE; // Standard deviation for Gaussion distrbution in energy
268  //Non invariant data members become G4Cache
269  G4double Emin, Emax; // emin and emax ////// NOT INVARIANT
270  G4double alpha, Ezero;// alpha (pow), E0 (exp) ////// NOT INVARIANT
271  G4double Temp; // Temp (bbody,brem)
272  G4double biasalpha; // biased power index
273  G4double grad, cept; // gradient and intercept for linear spectra ////// NOT INVARIANT
274  G4double prob_norm; // normalisation factor use in calculate the probability
275  G4bool Biased; // true - biased to power-law
276  G4bool EnergySpec; // true - energy spectra, false - momentum spectra
277  G4bool DiffSpec; // true - differential spec, false integral spec
278  //G4bool ApplyRig; // false no rigidity cutoff, true then apply one
279  //G4double ERig; // energy of rigidity cutoff
286  G4double CDGhist[3]; // cumulative histo for cdg
287 
288  //AG: Begin edit to use STL vectors.
289 // G4double BBHist[10001], Bbody_x[10001];
290  std::vector<G4double>* BBHist;
291  std::vector<G4double>* Bbody_x;
294 // For cutoff power-law
295  std::vector<G4double>* CPHist;
296  std::vector<G4double>* CP_x;
299 
300 
301  //AG: Edit here to use dynamic memory, will save space inless these functions are used.
302  G4String IntType; // Interpolation type
303 // G4double Arb_grad[1024], Arb_cept[1024]; // grad and cept for 1024 segments AG: Switched to DMA
307 // G4double Arb_alpha[1024], Arb_Const[1024]; // alpha and constants AG: Switched to DMA
311 // G4double Arb_ezero[1024]; // ezero AG: Switched to DMA
314  G4double ArbEmin, ArbEmax; // Emin and Emax for the whole arb distribution used primarily for debug.
315 
317 
319 
320  // Verbosity
322 
324 
325  std::vector<G4DataInterpolation*> SplineInt;//[1024]; // holds Spline stuff required for sampling
326  G4DataInterpolation *Splinetemp; // holds a temp Spline used for calculating area
327 
328  G4Mutex mutex; // protect access to shared resources
329  //Thread local data (non-invariant during event loop).
330  //These are copied from master one at the beginning of generation
331  //of each event
332  struct threadLocal_t {
342  };
344 };
345 
346 #endif
347 
G4double GetProbability(G4double)
std::vector< G4double > * CPHist
G4PhysicsOrderedFreeVector IPDFEnergyH
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
void ArbEnergyHisto(G4ThreeVector)
void EpnEnergyHisto(G4ThreeVector)
void ArbEnergyHistoFile(G4String)
std::vector< G4DataInterpolation * > SplineInt
G4double GenerateOne(G4ParticleDefinition *)
void SetBiasRndm(G4SPSRandomGenerator *a)
std::vector< G4double > * CP_x
std::vector< G4double > * BBHist
G4PhysicsOrderedFreeVector UDefEnergyH
G4PhysicsOrderedFreeVector ZeroPhysVector
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
std::vector< G4double > * Bbody_x
G4PhysicsOrderedFreeVector IPDFArbEnergyH
G4PhysicsOrderedFreeVector EpnEnergyH
void UserEnergyHisto(G4ThreeVector)
int G4int
Definition: G4Types.hh:78
G4PhysicsOrderedFreeVector ArbEnergyH
G4SPSRandomGenerator * eneRndm
G4Cache< threadLocal_t > threadLocalData
G4PhysicsOrderedFreeVector GetUserDefinedEnergyHisto()
G4DataInterpolation * Splinetemp
void InputDifferentialSpectra(G4bool)
G4PhysicsOrderedFreeVector GetArbEnergyHisto()
std::mutex G4Mutex
Definition: G4Threading.hh:84