Geant4  10.01.p01
G4Scintillation.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 //
27 // $Id: G4Scintillation.hh 85355 2014-10-28 09:58:59Z gcosmo $
28 //
29 //
31 // Scintillation Light Class Definition
33 //
34 // File: G4Scintillation.hh
35 // Description: Discrete Process - Generation of Scintillation Photons
36 // Version: 1.0
37 // Created: 1998-11-07
38 // Author: Peter Gumplinger
39 // Updated: 2010-10-20 Allow the scintillation yield to be a function
40 // of energy deposited by particle type
41 // Thanks to Zach Hartwig (Department of Nuclear
42 // Science and Engineeering - MIT)
43 // 2005-07-28 add G4ProcessType to constructor
44 // 2002-11-21 change to user G4Poisson for small MeanNumPotons
45 // 2002-11-07 allow for fast and slow scintillation
46 // 2002-11-05 make use of constant material properties
47 // 2002-05-16 changed to inherit from VRestDiscreteProcess
48 // 2002-05-09 changed IsApplicable method
49 // 1999-10-29 add method and class descriptors
50 //
51 // mail: gum@triumf.ca
52 //
54 
55 #ifndef G4Scintillation_h
56 #define G4Scintillation_h 1
57 
59 // Includes
61 
62 #include "globals.hh"
63 #include "templates.hh"
64 #include "Randomize.hh"
65 #include "G4Poisson.hh"
66 #include "G4ThreeVector.hh"
67 #include "G4ParticleMomentum.hh"
68 #include "G4Step.hh"
70 #include "G4OpticalPhoton.hh"
71 #include "G4DynamicParticle.hh"
72 #include "G4Material.hh"
73 #include "G4PhysicsTable.hh"
76 
77 #include "G4EmSaturation.hh"
78 
79 // Class Description:
80 // RestDiscrete Process - Generation of Scintillation Photons.
81 // Class inherits publicly from G4VRestDiscreteProcess.
82 // Class Description - End:
83 
85 // Class Definition
87 
89 {
90 
91 public:
92 
94  // Constructors and Destructor
96 
97  G4Scintillation(const G4String& processName = "Scintillation",
100 
101 private:
102 
104 
106  // Operators
108 
110 
111 public:
112 
114  // Methods
116 
117  // G4Scintillation Process has both PostStepDoIt (for energy
118  // deposition of particles in flight) and AtRestDoIt (for energy
119  // given to the medium by particles at rest)
120 
121  G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
122  // Returns true -> 'is applicable', for any particle type except
123  // for an 'opticalphoton' and for short-lived particles
124 
125  void BuildPhysicsTable(const G4ParticleDefinition& aParticleType);
126  // Build table at the right time
127 
128  G4double GetMeanFreePath(const G4Track& aTrack,
129  G4double ,
130  G4ForceCondition* );
131  // Returns infinity; i. e. the process does not limit the step,
132  // but sets the 'StronglyForced' condition for the DoIt to be
133  // invoked at every step.
134 
135  G4double GetMeanLifeTime(const G4Track& aTrack,
136  G4ForceCondition* );
137  // Returns infinity; i. e. the process does not limit the time,
138  // but sets the 'StronglyForced' condition for the DoIt to be
139  // invoked at every step.
140 
141  G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
142  const G4Step& aStep);
143  G4VParticleChange* AtRestDoIt (const G4Track& aTrack,
144  const G4Step& aStep);
145 
147  const G4Step &aStep);
148  // Returns the number of scintillation photons calculated when
149  // scintillation depends on the particle type and energy
150  // deposited (includes nonlinear dependendency)
151 
152  // These are the methods implementing the scintillation process.
153 
154  void SetTrackSecondariesFirst(const G4bool state);
155  // If set, the primary particle tracking is interrupted and any
156  // produced scintillation photons are tracked next. When all
157  // have been tracked, the tracking of the primary resumes.
158 
159  void SetFiniteRiseTime(const G4bool state);
160  // If set, the G4Scintillation process expects the user to have
161  // set the constant material property FAST/SLOWSCINTILLATIONRISETIME.
162 
164  // Returns the boolean flag for tracking secondaries first.
165 
166  G4bool GetFiniteRiseTime() const;
167  // Returns the boolean flag for a finite scintillation rise time.
168 
169  void SetScintillationYieldFactor(const G4double yieldfactor);
170  // Called to set the scintillation photon yield factor, needed when
171  // the yield is different for different types of particles. This
172  // scales the yield obtained from the G4MaterialPropertiesTable.
173 
175  // Returns the photon yield factor.
176 
177  void SetScintillationExcitationRatio(const G4double ratio);
178  // Called to set the scintillation exciation ratio, needed when
179  // the scintillation level excitation is different for different
180  // types of particles. This overwrites the YieldRatio obtained
181  // from the G4MaterialPropertiesTable.
182 
184  // Returns the scintillation level excitation ratio.
185 
187  // Returns the address of the fast scintillation integral table.
188 
190  // Returns the address of the slow scintillation integral table.
191 
193  // Adds Birks Saturation to the process.
194 
195  void RemoveSaturation();
196  // Removes the Birks Saturation from the process.
197 
199  // Returns the Birks Saturation.
200 
202  // Called by the user to set the scintillation yield as a function
203  // of energy deposited by particle type
204 
206  { return fScintillationByParticleType; }
207  // Return the boolean that determines the method of scintillation
208  // production
209 
210  void DumpPhysicsTable() const;
211  // Prints the fast and slow scintillation integral tables.
212 
213 protected:
214 
215  void BuildThePhysicsTable();
216  // It builds either the fast or slow scintillation integral table;
217  // or both.
218 
220  // Class Data Members
222 
225 
226 private:
227 
230 
232 
234 
236 
237 #ifdef G4DEBUG_SCINTILLATION
238  G4double ScintTrackEDep, ScintTrackYield;
239 #endif
240 
242  G4double bi_exp(G4double t, G4double tau1, G4double tau2);
243 
244  // emission time distribution when there is a finite rise time
246 
248 
249 };
250 
252 // Inline methods
254 
255 inline
257 {
258  if (aParticleType.GetParticleName() == "opticalphoton") return false;
259  if (aParticleType.IsShortLived()) return false;
260 
261  return true;
262 }
263 
264 inline
266 {
267  return fTrackSecondariesFirst;
268 }
269 
270 inline
272 {
273  return fFiniteRiseTime;
274 }
275 
276 inline
278 {
279  return fYieldFactor;
280 }
281 
282 inline
284 {
285  return fExcitationRatio;
286 }
287 
288 inline
290 {
291  return fSlowIntegralTable;
292 }
293 
294 inline
296 {
297  return fFastIntegralTable;
298 }
299 
300 inline
302 {
303  if (fFastIntegralTable) {
304  G4int PhysicsTableSize = fFastIntegralTable->entries();
306 
307  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
308  {
310  v->DumpValues();
311  }
312  }
313 
314  if (fSlowIntegralTable) {
315  G4int PhysicsTableSize = fSlowIntegralTable->entries();
317 
318  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
319  {
321  v->DumpValues();
322  }
323  }
324 }
325 
326 inline
328 {
329  return std::exp(-1.0*t/tau2)/tau2;
330 }
331 
332 inline
334 {
335  return std::exp(-1.0*t/tau2)*(1-std::exp(-1.0*t/tau1))/tau2/tau2*(tau1+tau2);
336 }
337 
338 #endif /* G4Scintillation_h */
void SetScintillationByParticleType(const G4bool)
G4Scintillation(const G4String &processName="Scintillation", G4ProcessType type=fElectromagnetic)
G4double bi_exp(G4double t, G4double tau1, G4double tau2)
void SetFiniteRiseTime(const G4bool state)
G4bool fTrackSecondariesFirst
G4bool GetFiniteRiseTime() const
G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *)
G4double single_exp(G4double t, G4double tau2)
G4PhysicsTable * fFastIntegralTable
G4EmSaturation * fEmSaturation
G4double GetScintillationYieldFactor() const
G4bool GetTrackSecondariesFirst() const
G4PhysicsTable * fSlowIntegralTable
G4PhysicsTable * GetSlowIntegralTable() const
G4VParticleChange * AtRestDoIt(const G4Track &aTrack, const G4Step &aStep)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
void SetScintillationYieldFactor(const G4double yieldfactor)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
void AddSaturation(G4EmSaturation *)
bool G4bool
Definition: G4Types.hh:79
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4Step.hh:76
void SetScintillationExcitationRatio(const G4double ratio)
G4double fExcitationRatio
G4Scintillation & operator=(const G4Scintillation &right)
void DumpPhysicsTable() const
void SetTrackSecondariesFirst(const G4bool state)
G4EmSaturation * GetSaturation() const
G4PhysicsTable * GetFastIntegralTable() const
G4double GetScintillationExcitationRatio() const
double G4double
Definition: G4Types.hh:76
G4ForceCondition
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)
G4bool GetScintillationByParticleType() const
G4double sample_time(G4double tau1, G4double tau2)
G4double GetScintillationYieldByParticleType(const G4Track &aTrack, const G4Step &aStep)
size_t entries() const
G4ProcessType
G4bool fScintillationByParticleType