Geant4  10.03
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 98002 2016-06-30 13:03:36Z 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  explicit G4Scintillation(const G4String& processName = "Scintillation",
100 
101 private:
102 
103  G4Scintillation(const G4Scintillation &right) = delete;
104 
106  // Operators
108 
109  G4Scintillation& operator=(const G4Scintillation &right) = delete;
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 
122  const G4ParticleDefinition& aParticleType) override;
123  // Returns true -> 'is applicable', for any particle type except
124  // for an 'opticalphoton' and for short-lived particles
125 
126  void BuildPhysicsTable(
127  const G4ParticleDefinition& aParticleType) override;
128  // Build table at the right time
129 
130  G4double GetMeanFreePath(const G4Track& aTrack,
131  G4double ,
132  G4ForceCondition* ) override;
133  // Returns infinity; i. e. the process does not limit the step,
134  // but sets the 'StronglyForced' condition for the DoIt to be
135  // invoked at every step.
136 
137  G4double GetMeanLifeTime(const G4Track& aTrack,
138  G4ForceCondition* ) override;
139  // Returns infinity; i. e. the process does not limit the time,
140  // but sets the 'StronglyForced' condition for the DoIt to be
141  // invoked at every step.
142 
143  G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
144  const G4Step& aStep) override;
145  G4VParticleChange* AtRestDoIt (const G4Track& aTrack,
146  const G4Step& aStep) override;
147 
149  const G4Step &aStep);
150  // Returns the number of scintillation photons calculated when
151  // scintillation depends on the particle type and energy
152  // deposited (includes nonlinear dependendency)
153 
154  // These are the methods implementing the scintillation process.
155 
156  void SetTrackSecondariesFirst(const G4bool state);
157  // If set, the primary particle tracking is interrupted and any
158  // produced scintillation photons are tracked next. When all
159  // have been tracked, the tracking of the primary resumes.
160 
162  // Returns the boolean flag for tracking secondaries first.
163 
164  void SetFiniteRiseTime(const G4bool state);
165  // If set, the G4Scintillation process expects the user to have
166  // set the constant material property FAST/SLOWSCINTILLATIONRISETIME.
167 
168  G4bool GetFiniteRiseTime() const;
169  // Returns the boolean flag for a finite scintillation rise time.
170 
171  void SetScintillationYieldFactor(const G4double yieldfactor);
172  // Called to set the scintillation photon yield factor, needed when
173  // the yield is different for different types of particles. This
174  // scales the yield obtained from the G4MaterialPropertiesTable.
175 
177  // Returns the photon yield factor.
178 
179  void SetScintillationExcitationRatio(const G4double ratio);
180  // Called to set the scintillation exciation ratio, needed when
181  // the scintillation level excitation is different for different
182  // types of particles. This overwrites the YieldRatio obtained
183  // from the G4MaterialPropertiesTable.
184 
186  // Returns the scintillation level excitation ratio.
187 
189  // Returns the address of the fast scintillation integral table.
190 
192  // Returns the address of the slow scintillation integral table.
193 
194  void AddSaturation(G4EmSaturation* sat);
195  // Adds Birks Saturation to the process.
196 
197  void RemoveSaturation();
198  // Removes the Birks Saturation from the process.
199 
200  G4EmSaturation* GetSaturation() const;
201  // Returns the Birks Saturation.
202 
204  // Called by the user to set the scintillation yield as a function
205  // of energy deposited by particle type
206 
208  // Return the boolean that determines the method of scintillation
209  // production
210 
211  void SetScintillationTrackInfo(const G4bool trackType);
212  // Call by the user to set the G4ScintillationTrackInformation
213  // to scintillation photon track
214 
216  // Return the boolean for whether or not the
217  // G4ScintillationTrackInformation is set to the scint. photon track
218 
219  void SetStackPhotons(const G4bool );
220  // Call by the user to set the flag for stacking the scint. photons
221 
222  G4bool GetStackPhotons() const;
223  // Return the boolean for whether or not the scint. photons are stacked
224 
225  G4int GetNumPhotons() const;
226  // Returns the current number of scint. photons (after PostStepDoIt)
227 
228  void DumpPhysicsTable() const;
229  // Prints the fast and slow scintillation integral tables.
230 
231 protected:
232 
233  void BuildThePhysicsTable();
234  // It builds either the fast or slow scintillation integral table;
235  // or both.
236 
238  // Class Data Members
240 
243 
244 private:
245 
248 
250 
252 
254 
256 
258 
260 
261 #ifdef G4DEBUG_SCINTILLATION
262  G4double ScintTrackEDep, ScintTrackYield;
263 #endif
264 
266  G4double bi_exp(G4double t, G4double tau1, G4double tau2);
267 
268  // emission time distribution when there is a finite rise time
270 
272 
273 };
274 
276 // Inline methods
278 
279 inline
281 {
282  if (aParticleType.GetParticleName() == "opticalphoton") return false;
283  if (aParticleType.IsShortLived()) return false;
284 
285  return true;
286 }
287 
288 inline
290 {
291  fTrackSecondariesFirst = state;
292 }
293 
294 inline
296 {
297  return fTrackSecondariesFirst;
298 }
299 
300 inline
302 {
303  fFiniteRiseTime = state;
304 }
305 
306 inline
308 {
309  return fFiniteRiseTime;
310 }
311 
312 inline
314 {
315  fYieldFactor = yieldfactor;
316 }
317 
318 inline
320 {
321  return fYieldFactor;
322 }
323 
324 inline
326 {
327  fExcitationRatio = ratio;
328 }
329 
330 inline
332 {
333  return fExcitationRatio;
334 }
335 
336 inline
338 {
339  return fSlowIntegralTable;
340 }
341 
342 inline
344 {
345  return fFastIntegralTable;
346 }
347 
348 inline
350 {
351  fEmSaturation = sat;
352 }
353 
354 inline
356 {
357  fEmSaturation = nullptr;
358 }
359 
360 inline
362 {
363  return fEmSaturation;
364 }
365 
366 inline
368 {
370 }
371 
372 inline
374 {
375  fScintillationTrackInfo = trackType;
376 }
377 
378 inline
380 {
382 }
383 
384 inline
385 void G4Scintillation::SetStackPhotons(const G4bool stackingFlag)
386 {
387  fStackingFlag = stackingFlag;
388 }
389 
390 inline
392 {
393  return fStackingFlag;
394 }
395 
396 inline
398 {
399  return fNumPhotons;
400 }
401 
402 inline
404 {
405  if (fFastIntegralTable) {
406  G4int PhysicsTableSize = fFastIntegralTable->entries();
408 
409  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
410  {
412  v->DumpValues();
413  }
414  }
415 
416  if (fSlowIntegralTable) {
417  G4int PhysicsTableSize = fSlowIntegralTable->entries();
419 
420  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
421  {
423  v->DumpValues();
424  }
425  }
426 }
427 
428 inline
430 {
431  return std::exp(-1.0*t/tau2)/tau2;
432 }
433 
434 inline
436 {
437  return std::exp(-1.0*t/tau2)*(1-std::exp(-1.0*t/tau1))/tau2/tau2*(tau1+tau2);
438 }
439 
440 #endif /* G4Scintillation_h */
void SetScintillationByParticleType(const G4bool)
G4Scintillation(const G4String &processName="Scintillation", G4ProcessType type=fElectromagnetic)
G4bool GetScintillationTrackInfo() const
G4double bi_exp(G4double t, G4double tau1, G4double tau2)
void SetScintillationTrackInfo(const G4bool trackType)
void SetFiniteRiseTime(const G4bool state)
G4bool fTrackSecondariesFirst
void DumpValues(G4double unitE=1.0, G4double unitV=1.0) const
G4bool GetFiniteRiseTime() const
G4double single_exp(G4double t, G4double tau2)
G4PhysicsTable * fFastIntegralTable
G4EmSaturation * fEmSaturation
G4double GetScintillationYieldFactor() const
G4bool GetTrackSecondariesFirst() const
G4PhysicsTable * fSlowIntegralTable
G4PhysicsTable * GetSlowIntegralTable() const
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType) override
int G4int
Definition: G4Types.hh:78
G4bool GetStackPhotons() const
const G4String & GetParticleName() const
void SetScintillationYieldFactor(const G4double yieldfactor)
G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *) override
void SetStackPhotons(const G4bool)
bool G4bool
Definition: G4Types.hh:79
Definition: G4Step.hh:76
void SetScintillationExcitationRatio(const G4double ratio)
G4double fExcitationRatio
void DumpPhysicsTable() const
G4bool fScintillationTrackInfo
void SetTrackSecondariesFirst(const G4bool state)
G4Scintillation & operator=(const G4Scintillation &right)=delete
G4PhysicsTable * GetFastIntegralTable() const
G4double GetScintillationExcitationRatio() const
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
double G4double
Definition: G4Types.hh:76
G4EmSaturation * GetSaturation() const
G4ForceCondition
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
G4bool GetScintillationByParticleType() const
G4double sample_time(G4double tau1, G4double tau2)
G4double GetScintillationYieldByParticleType(const G4Track &aTrack, const G4Step &aStep)
void AddSaturation(G4EmSaturation *sat)
size_t entries() const
G4int GetNumPhotons() const
G4VParticleChange * AtRestDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4ProcessType
G4bool fScintillationByParticleType