Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
28//
29//
30////////////////////////////////////////////////////////////////////////
31// Scintillation Light Class Definition
32////////////////////////////////////////////////////////////////////////
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: [email protected]
52//
53////////////////////////////////////////////////////////////////////////
54
55#ifndef G4Scintillation_h
56#define G4Scintillation_h 1
57
58/////////////
59// Includes
60/////////////
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
84/////////////////////
85// Class Definition
86/////////////////////
87
89{
90
91public:
92
93 ////////////////////////////////
94 // Constructors and Destructor
95 ////////////////////////////////
96
97 G4Scintillation(const G4String& processName = "Scintillation",
100
101private:
102
103 G4Scintillation(const G4Scintillation &right);
104
105 //////////////
106 // Operators
107 //////////////
108
109 G4Scintillation& operator=(const G4Scintillation &right);
110
111public:
112
113 ////////////
114 // Methods
115 ////////////
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 G4double GetMeanFreePath(const G4Track& aTrack,
126 G4double ,
128 // Returns infinity; i. e. the process does not limit the step,
129 // but sets the 'StronglyForced' condition for the DoIt to be
130 // invoked at every step.
131
132 G4double GetMeanLifeTime(const G4Track& aTrack,
134 // Returns infinity; i. e. the process does not limit the time,
135 // but sets the 'StronglyForced' condition for the DoIt to be
136 // invoked at every step.
137
138 G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
139 const G4Step& aStep);
140 G4VParticleChange* AtRestDoIt (const G4Track& aTrack,
141 const G4Step& aStep);
142
143 // These are the methods implementing the scintillation process.
144
145 void SetTrackSecondariesFirst(const G4bool state);
146 // If set, the primary particle tracking is interrupted and any
147 // produced scintillation photons are tracked next. When all
148 // have been tracked, the tracking of the primary resumes.
149
150 void SetFiniteRiseTime(const G4bool state);
151 // If set, the G4Scintillation process expects the user to have
152 // set the constant material property FAST/SLOWSCINTILLATIONRISETIME.
153
155 // Returns the boolean flag for tracking secondaries first.
156
158 // Returns the boolean flag for a finite scintillation rise time.
159
160 void SetScintillationYieldFactor(const G4double yieldfactor);
161 // Called to set the scintillation photon yield factor, needed when
162 // the yield is different for different types of particles. This
163 // scales the yield obtained from the G4MaterialPropertiesTable.
164
166 // Returns the photon yield factor.
167
168 void SetScintillationExcitationRatio(const G4double excitationratio);
169 // Called to set the scintillation exciation ratio, needed when
170 // the scintillation level excitation is different for different
171 // types of particles. This overwrites the YieldRatio obtained
172 // from the G4MaterialPropertiesTable.
173
175 // Returns the scintillation level excitation ratio.
176
178 // Returns the address of the fast scintillation integral table.
179
181 // Returns the address of the slow scintillation integral table.
182
183 void AddSaturation(G4EmSaturation* sat) { emSaturation = sat; }
184 // Adds Birks Saturation to the process.
185
186 void RemoveSaturation() { emSaturation = NULL; }
187 // Removes the Birks Saturation from the process.
188
189 G4EmSaturation* GetSaturation() const { return emSaturation; }
190 // Returns the Birks Saturation.
191
193 // Called by the user to set the scintillation yield as a function
194 // of energy deposited by particle type
195
198 // Return the boolean that determines the method of scintillation
199 // production
200
201 void DumpPhysicsTable() const;
202 // Prints the fast and slow scintillation integral tables.
203
204protected:
205
207 // It builds either the fast or slow scintillation integral table;
208 // or both.
209
210 ///////////////////////
211 // Class Data Members
212 ///////////////////////
213
216
219
221
223
225
226private:
227
228 G4double single_exp(G4double t, G4double tau2);
229 G4double bi_exp(G4double t, G4double tau1, G4double tau2);
230
231 // emission time distribution when there is a finite rise time
232 G4double sample_time(G4double tau1, G4double tau2);
233
234 G4EmSaturation* emSaturation;
235
236};
237
238////////////////////
239// Inline methods
240////////////////////
241
242inline
244{
245 if (aParticleType.GetParticleName() == "opticalphoton") return false;
246 if (aParticleType.IsShortLived()) return false;
247
248 return true;
249}
250
251inline
253{
255}
256
257inline
259{
260 fFiniteRiseTime = state;
261}
262
263inline
265{
267}
268
269inline
271{
272 return fFiniteRiseTime;
273}
274
275inline
277{
278 YieldFactor = yieldfactor;
279}
280
281inline
283{
284 return YieldFactor;
285}
286
287inline
289{
290 ExcitationRatio = excitationratio;
291}
292
293inline
295{
296 return ExcitationRatio;
297}
298
299inline
301{
303}
304
305inline
307{
309}
310
311inline
313{
315 G4int PhysicsTableSize = theFastIntegralTable->entries();
317
318 for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
319 {
321 v->DumpValues();
322 }
323 }
324
326 G4int PhysicsTableSize = theSlowIntegralTable->entries();
328
329 for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
330 {
332 v->DumpValues();
333 }
334 }
335}
336
337inline
338G4double G4Scintillation::single_exp(G4double t, G4double tau2)
339{
340 return std::exp(-1.0*t/tau2)/tau2;
341}
342
343inline
344G4double G4Scintillation::bi_exp(G4double t, G4double tau1, G4double tau2)
345{
346 return std::exp(-1.0*t/tau2)*(1-std::exp(-1.0*t/tau1))/tau2/tau2*(tau1+tau2);
347}
348
349#endif /* G4Scintillation_h */
G4ForceCondition
G4ProcessType
@ fElectromagnetic
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
const G4String & GetParticleName() const
size_t entries() const
G4EmSaturation * GetSaturation() const
void SetScintillationYieldFactor(const G4double yieldfactor)
void SetTrackSecondariesFirst(const G4bool state)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4bool scintillationByParticleType
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
G4bool GetFiniteRiseTime() const
G4VParticleChange * AtRestDoIt(const G4Track &aTrack, const G4Step &aStep)
void SetScintillationExcitationRatio(const G4double excitationratio)
void AddSaturation(G4EmSaturation *sat)
G4double GetScintillationYieldFactor() const
void DumpPhysicsTable() const
G4double ExcitationRatio
G4bool fTrackSecondariesFirst
void SetFiniteRiseTime(const G4bool state)
G4double GetScintillationExcitationRatio() const
G4PhysicsTable * GetFastIntegralTable() const
G4bool GetTrackSecondariesFirst() const
G4bool GetScintillationByParticleType() const
G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *)
G4PhysicsTable * GetSlowIntegralTable() const
void SetScintillationByParticleType(const G4bool)
G4PhysicsTable * theSlowIntegralTable
G4PhysicsTable * theFastIntegralTable
Definition: G4Step.hh:78