source: trunk/source/processes/electromagnetic/xrays/include/G4Scintillation.hh @ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 11.4 KB
Line 
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,v 1.21 2010/10/28 23:29:21 gum Exp $
28// GEANT4 tag $Name: xrays-V09-03-05 $
29//
30//
31////////////////////////////////////////////////////////////////////////
32// Scintillation Light Class Definition
33////////////////////////////////////////////////////////////////////////
34//
35// File:        G4Scintillation.hh 
36// Description: Discrete Process - Generation of Scintillation Photons
37// Version:     1.0
38// Created:     1998-11-07
39// Author:      Peter Gumplinger
40// Updated:     2010-10-20 Allow the scintillation yield to be a function
41//                         of energy deposited by particle type
42//                         Thanks to Zach Hartwig (Department of Nuclear
43//                         Science and Engineeering - MIT)
44//              2005-07-28 add G4ProcessType to constructor
45//              2002-11-21 change to user G4Poisson for small MeanNumPotons
46//              2002-11-07 allow for fast and slow scintillation
47//              2002-11-05 make use of constant material properties
48//              2002-05-16 changed to inherit from VRestDiscreteProcess
49//              2002-05-09 changed IsApplicable method
50//              1999-10-29 add method and class descriptors
51//
52// mail:        gum@triumf.ca
53//
54////////////////////////////////////////////////////////////////////////
55
56#ifndef G4Scintillation_h
57#define G4Scintillation_h 1
58
59/////////////
60// Includes
61/////////////
62
63#include "globals.hh"
64#include "templates.hh"
65#include "Randomize.hh"
66#include "G4Poisson.hh"
67#include "G4ThreeVector.hh"
68#include "G4ParticleMomentum.hh"
69#include "G4Step.hh"
70#include "G4VRestDiscreteProcess.hh"
71#include "G4OpticalPhoton.hh"
72#include "G4DynamicParticle.hh"
73#include "G4Material.hh"
74#include "G4PhysicsTable.hh"
75#include "G4MaterialPropertiesTable.hh"
76#include "G4PhysicsOrderedFreeVector.hh"
77
78#include "G4EmSaturation.hh"
79
80// Class Description:
81// RestDiscrete Process - Generation of Scintillation Photons.
82// Class inherits publicly from G4VRestDiscreteProcess.
83// Class Description - End:
84
85/////////////////////
86// Class Definition
87/////////////////////
88
89class G4Scintillation : public G4VRestDiscreteProcess
90{
91
92private:
93
94        //////////////
95        // Operators
96        //////////////
97
98        // G4Scintillation& operator=(const G4Scintillation &right);
99
100public: // Without description
101
102        ////////////////////////////////
103        // Constructors and Destructor
104        ////////////////////////////////
105
106        G4Scintillation(const G4String& processName = "Scintillation",
107                                 G4ProcessType type = fElectromagnetic);
108
109        // G4Scintillation(const G4Scintillation &right);
110
111        ~G4Scintillation();     
112
113        ////////////
114        // Methods
115        ////////////
116
117public: // With description
118
119        // G4Scintillation Process has both PostStepDoIt (for energy
120        // deposition of particles in flight) and AtRestDoIt (for energy
121        // given to the medium by particles at rest)
122
123        G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
124        // Returns true -> 'is applicable', for any particle type except
125        // for an 'opticalphoton' and for short-lived particles
126
127        G4double GetMeanFreePath(const G4Track& aTrack,
128                                       G4double ,
129                                       G4ForceCondition* );
130        // Returns infinity; i. e. the process does not limit the step,
131        // but sets the 'StronglyForced' condition for the DoIt to be
132        // invoked at every step.
133
134        G4double GetMeanLifeTime(const G4Track& aTrack,
135                                 G4ForceCondition* );
136        // Returns infinity; i. e. the process does not limit the time,
137        // but sets the 'StronglyForced' condition for the DoIt to be
138        // invoked at every step.
139
140        G4VParticleChange* PostStepDoIt(const G4Track& aTrack, 
141                                        const G4Step&  aStep);
142        G4VParticleChange* AtRestDoIt (const G4Track& aTrack,
143                                       const G4Step& aStep);
144
145        // These are the methods implementing the scintillation process.
146
147        void SetTrackSecondariesFirst(const G4bool state);
148        // If set, the primary particle tracking is interrupted and any
149        // produced scintillation photons are tracked next. When all
150        // have been tracked, the tracking of the primary resumes.
151
152        void SetFiniteRiseTime(const G4bool state);
153        // If set, the G4Scintillation process expects the user to have
154        // set the constant material property FAST/SLOWSCINTILLATIONRISETIME.
155
156        G4bool GetTrackSecondariesFirst() const;
157        // Returns the boolean flag for tracking secondaries first.
158
159        G4bool GetFiniteRiseTime() const;
160        // Returns the boolean flag for a finite scintillation rise time.
161       
162        void SetScintillationYieldFactor(const G4double yieldfactor);
163        // Called to set the scintillation photon yield factor, needed when
164        // the yield is different for different types of particles. This
165        // scales the yield obtained from the G4MaterialPropertiesTable.
166
167        G4double GetScintillationYieldFactor() const;
168        // Returns the photon yield factor.
169
170        void SetScintillationExcitationRatio(const G4double excitationratio);
171        // Called to set the scintillation exciation ratio, needed when
172        // the scintillation level excitation is different for different
173        // types of particles. This overwrites the YieldRatio obtained
174        // from the G4MaterialPropertiesTable.
175
176        G4double GetScintillationExcitationRatio() const;
177        // Returns the scintillation level excitation ratio.
178
179        G4PhysicsTable* GetFastIntegralTable() const;
180        // Returns the address of the fast scintillation integral table.
181
182        G4PhysicsTable* GetSlowIntegralTable() const;
183        // Returns the address of the slow scintillation integral table.
184
185        void AddSaturation(G4EmSaturation* sat) { emSaturation = sat; }
186        // Adds Birks Saturation to the process.
187
188        void RemoveSaturation() { emSaturation = NULL; }
189        // Removes the Birks Saturation from the process.
190
191        G4EmSaturation* GetSaturation() const { return emSaturation; }
192        // Returns the Birks Saturation.
193
194        void SetScintillationByParticleType(const G4bool );
195        // Called by the user to set the scintillation yield as a function
196        // of energy deposited by particle type
197
198        G4bool GetScintillationByParticleType() const
199        { return scintillationByParticleType; }
200        // Return the boolean that determines the method of scintillation
201        // production
202
203        void DumpPhysicsTable() const;
204        // Prints the fast and slow scintillation integral tables.
205
206protected:
207
208        void BuildThePhysicsTable();
209        // It builds either the fast or slow scintillation integral table;
210        // or both.
211
212        ///////////////////////
213        // Class Data Members
214        ///////////////////////
215
216
217        G4PhysicsTable* theSlowIntegralTable;
218        G4PhysicsTable* theFastIntegralTable;
219
220
221
222        G4bool fTrackSecondariesFirst;
223        G4bool fFiniteRiseTime;
224
225        G4double YieldFactor;
226
227        G4double ExcitationRatio;
228
229        G4bool scintillationByParticleType;
230
231private:
232
233        G4double single_exp(G4double t, G4double tau2);
234        G4double bi_exp(G4double t, G4double tau1, G4double tau2);
235
236        // emission time distribution when there is a finite rise time
237        G4double sample_time(G4double tau1, G4double tau2);
238
239        G4EmSaturation* emSaturation;
240
241};
242
243////////////////////
244// Inline methods
245////////////////////
246
247inline 
248G4bool G4Scintillation::IsApplicable(const G4ParticleDefinition& aParticleType)
249{
250       if (aParticleType.GetParticleName() == "opticalphoton") return false;
251       if (aParticleType.IsShortLived()) return false;
252
253       return true;
254}
255
256inline 
257void G4Scintillation::SetTrackSecondariesFirst(const G4bool state) 
258{ 
259        fTrackSecondariesFirst = state;
260}
261
262inline
263void G4Scintillation::SetFiniteRiseTime(const G4bool state)
264{
265        fFiniteRiseTime = state;
266}
267
268inline
269G4bool G4Scintillation::GetTrackSecondariesFirst() const
270{
271        return fTrackSecondariesFirst;
272}
273
274inline 
275G4bool G4Scintillation::GetFiniteRiseTime() const
276{
277        return fFiniteRiseTime;
278}
279
280inline
281void G4Scintillation::SetScintillationYieldFactor(const G4double yieldfactor)
282{
283        YieldFactor = yieldfactor;
284}
285
286inline
287G4double G4Scintillation::GetScintillationYieldFactor() const
288{
289        return YieldFactor;
290}
291
292inline
293void G4Scintillation::SetScintillationExcitationRatio(const G4double excitationratio)
294{
295        ExcitationRatio = excitationratio;
296}
297
298inline
299G4double G4Scintillation::GetScintillationExcitationRatio() const
300{
301        return ExcitationRatio;
302}
303
304inline
305G4PhysicsTable* G4Scintillation::GetSlowIntegralTable() const
306{
307        return theSlowIntegralTable;
308}
309
310inline
311G4PhysicsTable* G4Scintillation::GetFastIntegralTable() const
312{
313        return theFastIntegralTable;
314}
315
316inline
317void G4Scintillation::DumpPhysicsTable() const
318{
319        if (theFastIntegralTable) {
320           G4int PhysicsTableSize = theFastIntegralTable->entries();
321           G4PhysicsOrderedFreeVector *v;
322
323           for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
324           {
325                v = (G4PhysicsOrderedFreeVector*)(*theFastIntegralTable)[i];
326                v->DumpValues();
327           }
328         }
329
330        if (theSlowIntegralTable) {
331           G4int PhysicsTableSize = theSlowIntegralTable->entries();
332           G4PhysicsOrderedFreeVector *v;
333
334           for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
335           {
336                v = (G4PhysicsOrderedFreeVector*)(*theSlowIntegralTable)[i];
337                v->DumpValues();
338           }
339         }
340}
341
342inline
343G4double G4Scintillation::single_exp(G4double t, G4double tau2)
344{
345         return std::exp(-1.0*t/tau2)/tau2;
346}
347
348inline
349G4double G4Scintillation::bi_exp(G4double t, G4double tau1, G4double tau2)
350{
351         return std::exp(-1.0*t/tau2)*(1-std::exp(-1.0*t/tau1))/tau2/tau2*(tau1+tau2);
352}
353
354#endif /* G4Scintillation_h */
Note: See TracBrowser for help on using the repository browser.