Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
26// G4SPSEneDistribution
27//
28// Class Description:
29//
30// To generate the energy of a primary vertex according to the
31// defined distribution. This is a shared class between threads.
32// Only one thread should use the set-methods here.
33// Note that this is exactly what is achieved using UI commands.
34// If you use the set methods to set defaults in your application take
35// care that only one thread is executing them.
36// In addition take care of calling these methods before the run is
37// started. Do not use the setters during the event loop
38
39// Author: Fan Lei, QinetiQ ltd.
40// Customer: ESA/ESTEC
41// History:
42// - 05/02/2004, Fan Lei - Created.
43// Based on the G4GeneralParticleSource class.
44// - 26/03/2014, Andrew Green.
45// Modification to use STL vectors instead of C-style arrays.
46// Also moved to dynamically allocated memory in the LinearInterpolation(),
47// ExpInterpolation() and LogInterpolation() functions.
48// - 06/06/2014, Andrea Dotti.
49// For thread safety: this is a shared object.
50// Added mutex to control access to shared resources (data members).
51// in Getters and Setters, mutex is NOT used in GenerateOne() because it
52// is assumed that properties are not changed during event loop.
53// - 24/11/2017, Fan Lei
54// Added cutoff power-law distribution option. Implementation is similar
55// to that of the BlackBody one.
56// --------------------------------------------------------------------
57#ifndef G4SPSEneDistribution_hh
58#define G4SPSEneDistribution_hh 1
59
61#include "G4ParticleMomentum.hh"
64#include "G4Threading.hh"
65#include "G4Cache.hh"
66#include <vector>
67
69
71{
72 public:
73
75 // Constructor: initializes variables
77 // Destructor
78
79 void SetEnergyDisType(const G4String&);
80 // Allows the user to choose the energy distribution type.
81 // The arguments are: Mono (mono-energetic), Lin (linear),
82 // Pow (power-law), Exp (exponential), Gauss (gaussian),
83 // Brem (bremsstrahlung), BBody (black-body),
84 // Cdg (cosmic diffuse gamma-ray), User (user-defined),
85 // Arb (arbitrary point-wise), Epn (energy per nucleon)
86
88
89 void SetEmin(G4double);
90 // Sets the minimum energy
91
92 G4double GetEmin() const;
94
95 void SetEmax(G4double);
96 // Sets the maximum energy
97
98 G4double GetEmax() const;
100
102 // Sets energy for mono-energetic distribution
103
104 void SetAlpha(G4double);
105 // Sets alpha for a power-law distribution
106
108
109 void SetTemp(G4double);
110 // Sets Temperature for a Brem or BBody distributions
111
113
114 void SetEzero(G4double);
115 // Sets Ezero for an exponential distribution
116
117 void SetGradient(G4double);
118 // Sets gradient for a linear distribution
119
121 // Sets intercept for a linear distribution
122
123 void UserEnergyHisto(const G4ThreeVector&);
124 // Allows user to defined a histogram for the energy distribution
125
126 void ArbEnergyHisto(const G4ThreeVector&);
127 // Allows the user to define an Arbitrary set of points for the
128 // energy distribution
129
130 void ArbEnergyHistoFile(const G4String&);
131
132 void EpnEnergyHisto(const G4ThreeVector&);
133 // Allows the user to define an Energy per nucleon histogram
134
136 // Allows the user to choose between momentum and energy histograms
137 // for user-defined histograms and arbitrary point-wise spectra.
138 // The default is true (energy)
139
141 // Allows the user to choose between integral and differential
142 // distributions when using the arbitrary point-wise option
143
144 void ArbInterpolate(const G4String&);
145 // Allows the user to specify the type of function to
146 // interpolate the Arbitrary points spectrum with
147
148 const G4String& GetIntType();
149
150 void Calculate();
151 // Controls the calculation of Integral PDF for the Cdg and BBody
152 // distributions
153
155 // Sets the biased random number generator
156
157 void ReSetHist(const G4String&);
158 // Resets the histogram for user defined distribution
159
160 void SetVerbosity(G4int a);
161 // Sets the verbosity level
162
163 G4double GetWeight() const;
164
166 // Mono-energetic energy
167
168 G4double GetSE();
169 // Standard deviation for Gaussian distribution in energy
170
171 G4double Getalpha() const;
172 // Alpha (pow)
173
174 G4double GetEzero() const;
175 // E0 (exp)
176
178 // Temp (bbody,brem)
179
180 G4double Getgrad() const;
181 // Gradient and intercept for linear spectra
182
183 G4double Getcept() const;
184
186
188
190 // Generate one random energy for the specified particle
191
193
195
196 inline void ApplyEnergyWeight(G4bool val) { applyEvergyWeight = val; }
197 inline G4bool IfApplyEnergyWeight() const { return applyEvergyWeight; }
198
199 private:
200
201 void LinearInterpolation();
202 void LogInterpolation();
203 void ExpInterpolation();
204 void SplineInterpolation();
205 void CalculateCdgSpectrum();
206 void CalculateBbodySpectrum();
207 void CalculateCPowSpectrum();
208
209 // The following methods generate energies according
210 // to the spectral parameters defined above
211
212 void GenerateMonoEnergetic();
213 void GenerateBiasPowEnergies();
214 void GenerateGaussEnergies();
215 void GenerateBremEnergies();
216 void GenerateBbodyEnergies();
217 void GenerateCdgEnergies();
218 void GenUserHistEnergies();
219 void GenEpnHistEnergies();
220 void GenArbPointEnergies(); // NOTE: REQUIRES UPDATE OF DATA MEMBERS
221 void GenerateExpEnergies(G4bool);
222 void GenerateLinearEnergies(G4bool);
223 void GeneratePowEnergies(G4bool);
224 void GenerateCPowEnergies();
225
226 void ConvertEPNToEnergy();
227 // Converts energy per nucleon to energy
228
229 void BBInitHists();
230 void CPInitHists();
231
232 private: // Non invariant data members become G4Cache
233
234 G4String EnergyDisType; // energy dis type Variable - Mono,Lin,Exp,etc
235 G4double weight; // particle weight //// NOT INVARIANT
236 G4double MonoEnergy; //Mono-energteic energy
237 G4double SE; // Standard deviation for Gaussian distribution in energy
238
239 G4double Emin, Emax; // emin and emax //// NOT INVARIANT
240 G4double alpha, Ezero;// alpha (pow), E0 (exp) //// NOT INVARIANT
241 G4double Temp; // Temp (bbody,brem)
242 G4double biasalpha; // biased power index
243 G4double grad, cept; // gradient and intercept for linear spectra //// NOT INVARIANT
244 G4double prob_norm; // normalisation factor use in calculate the probability
245 G4bool Biased = false; // biased to power-law
246 G4bool EnergySpec = true; // energy spectra, false - momentum spectra
247 G4bool DiffSpec = true; // differential spec, false integral spec
248
249 G4PhysicsOrderedFreeVector UDefEnergyH; // energy hist data
250 G4PhysicsOrderedFreeVector IPDFEnergyH;
251 G4bool IPDFEnergyExist = false, IPDFArbExist = false, Epnflag = false;
252 G4PhysicsOrderedFreeVector ArbEnergyH; // Arb x,y histogram
253 G4PhysicsOrderedFreeVector IPDFArbEnergyH; // IPDF for Arb
255 G4double CDGhist[3]; // cumulative histo for cdg
256
257 std::vector<G4double>* BBHist = nullptr;
258 std::vector<G4double>* Bbody_x = nullptr;
259 G4bool BBhistInit = false;
260 G4bool BBhistCalcd = false;
261
262 // For cutoff power-law
263 //
264 std::vector<G4double>* CPHist = nullptr;
265 std::vector<G4double>* CP_x = nullptr;
266 G4bool CPhistInit = false;
267 G4bool CPhistCalcd = false;
268
269 G4String IntType; // Interpolation type
270 G4double* Arb_grad = nullptr;
271 G4double* Arb_cept = nullptr;
272 G4bool Arb_grad_cept_flag = false;
273 G4double* Arb_alpha = nullptr;
274 G4double* Arb_Const = nullptr;
275 G4bool Arb_alpha_Const_flag = false;
276 G4double* Arb_ezero = nullptr;
277 G4bool Arb_ezero_flag = false;
278
279 G4bool applyEvergyWeight = false;
280
281 G4double ArbEmin, ArbEmax;
282 // Emin and Emax for the whole arb distribution used primarily for debug.
283
284 G4double particle_energy;
285
286 G4SPSRandomGenerator* eneRndm = nullptr;
287
288 G4int verbosityLevel;
289
290 G4PhysicsOrderedFreeVector ZeroPhysVector; // for re-set only
291
292 std::vector<G4DataInterpolation*> SplineInt;
293 // Holds Spline stuff required for sampling
294 G4DataInterpolation* Splinetemp = nullptr;
295 // Holds a temp Spline used for calculating area
296
297 G4Mutex mutex; // protect access to shared resources
298
299 // Thread local data (non-invariant during event loop).
300 // These are copied from master one at the beginning of
301 // generation of each event
302 //
303 struct threadLocal_t
304 {
305 G4double Emin;
306 G4double Emax;
307 G4double alpha;
309 G4double grad;
310 G4double cept;
311 G4ParticleDefinition* particle_definition;
312 G4double weight;
313 G4double particle_energy;
314 };
315 G4Cache<threadLocal_t> threadLocalData;
316};
317
318#endif
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4double GetProbability(G4double)
G4double GetArbEneWeight(G4double)
G4double GenerateOne(G4ParticleDefinition *)
G4PhysicsOrderedFreeVector GetArbEnergyHisto()
void ArbInterpolate(const G4String &)
const G4String & GetEnergyDisType()
void SetEnergyDisType(const G4String &)
void EpnEnergyHisto(const G4ThreeVector &)
void ArbEnergyHistoFile(const G4String &)
G4bool IfApplyEnergyWeight() const
void ReSetHist(const G4String &)
void InputDifferentialSpectra(G4bool)
void SetBiasRndm(G4SPSRandomGenerator *a)
G4PhysicsOrderedFreeVector GetUserDefinedEnergyHisto()
void ApplyEnergyWeight(G4bool val)
const G4String & GetIntType()
void UserEnergyHisto(const G4ThreeVector &)
void ArbEnergyHisto(const G4ThreeVector &)