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

Last change on this file since 1058 was 1007, checked in by garnier, 15 years ago

update to geant4.9.2

File size: 9.5 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.15 2008/06/13 01:04:49 gum Exp $
28// GEANT4 tag $Name: geant4-09-02 $
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:     2005-07-28 add G4ProcessType to constructor
41//              2002-11-21 change to user G4Poisson for small MeanNumPotons
42//              2002-11-07 allow for fast and slow scintillation
43//              2002-11-05 make use of constant material properties
44//              2002-05-16 changed to inherit from VRestDiscreteProcess
45//              2002-05-09 changed IsApplicable method
46//              1999-10-29 add method and class descriptors
47//
48// mail:        gum@triumf.ca
49//
50////////////////////////////////////////////////////////////////////////
51
52#ifndef G4Scintillation_h
53#define G4Scintillation_h 1
54
55/////////////
56// Includes
57/////////////
58
59#include "globals.hh"
60#include "templates.hh"
61#include "Randomize.hh"
62#include "G4Poisson.hh"
63#include "G4ThreeVector.hh"
64#include "G4ParticleMomentum.hh"
65#include "G4Step.hh"
66#include "G4VRestDiscreteProcess.hh"
67#include "G4OpticalPhoton.hh"
68#include "G4DynamicParticle.hh"
69#include "G4Material.hh"
70#include "G4PhysicsTable.hh"
71#include "G4MaterialPropertiesTable.hh"
72#include "G4PhysicsOrderedFreeVector.hh"
73
74#include "G4EmSaturation.hh"
75
76// Class Description:
77// RestDiscrete Process - Generation of Scintillation Photons.
78// Class inherits publicly from G4VRestDiscreteProcess.
79// Class Description - End:
80
81/////////////////////
82// Class Definition
83/////////////////////
84
85class G4Scintillation : public G4VRestDiscreteProcess
86{
87
88private:
89
90        //////////////
91        // Operators
92        //////////////
93
94        // G4Scintillation& operator=(const G4Scintillation &right);
95
96public: // Without description
97
98        ////////////////////////////////
99        // Constructors and Destructor
100        ////////////////////////////////
101
102        G4Scintillation(const G4String& processName = "Scintillation",
103                                 G4ProcessType type = fElectromagnetic);
104
105        // G4Scintillation(const G4Scintillation &right);
106
107        ~G4Scintillation();     
108
109        ////////////
110        // Methods
111        ////////////
112
113public: // With description
114
115        // G4Scintillation Process has both PostStepDoIt (for energy
116        // deposition of particles in flight) and AtRestDoIt (for energy
117        // given to the medium by particles at rest)
118
119        G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
120        // Returns true -> 'is applicable', for any particle type except
121        // for an 'opticalphoton'
122
123        G4double GetMeanFreePath(const G4Track& aTrack,
124                                       G4double ,
125                                       G4ForceCondition* );
126        // Returns infinity; i. e. the process does not limit the step,
127        // but sets the 'StronglyForced' condition for the DoIt to be
128        // invoked at every step.
129
130        G4double GetMeanLifeTime(const G4Track& aTrack,
131                                 G4ForceCondition* );
132        // Returns infinity; i. e. the process does not limit the time,
133        // but sets the 'StronglyForced' condition for the DoIt to be
134        // invoked at every step.
135
136        G4VParticleChange* PostStepDoIt(const G4Track& aTrack, 
137                                        const G4Step&  aStep);
138        G4VParticleChange* AtRestDoIt (const G4Track& aTrack,
139                                       const G4Step& aStep);
140
141        // These are the methods implementing the scintillation process.
142
143        void SetTrackSecondariesFirst(const G4bool state);
144        // If set, the primary particle tracking is interrupted and any
145        // produced scintillation photons are tracked next. When all
146        // have been tracked, the tracking of the primary resumes.
147
148        G4bool GetTrackSecondariesFirst() const;
149        // Returns the boolean flag for tracking secondaries first.
150       
151        void SetScintillationYieldFactor(const G4double yieldfactor);
152        // Called to set the scintillation photon yield factor, needed when
153        // the yield is different for different types of particles. This
154        // scales the yield obtained from the G4MaterialPropertiesTable.
155
156        G4double GetScintillationYieldFactor() const;
157        // Returns the photon yield factor.
158
159        void SetScintillationExcitationRatio(const G4double excitationratio);
160        // Called to set the scintillation exciation ratio, needed when
161        // the scintillation level excitation is different for different
162        // types of particles. This overwrites the YieldRatio obtained
163        // from the G4MaterialPropertiesTable.
164
165        G4double GetScintillationExcitationRatio() const;
166        // Returns the scintillation level excitation ratio.
167
168        G4PhysicsTable* GetFastIntegralTable() const;
169        // Returns the address of the fast scintillation integral table.
170
171        G4PhysicsTable* GetSlowIntegralTable() const;
172        // Returns the address of the slow scintillation integral table.
173
174        void AddSaturation(G4EmSaturation* sat) { emSaturation = sat; }
175        // Adds Birks Saturation to the process.
176
177        G4EmSaturation* GetSaturation() const { return emSaturation; }
178        // Returns the Birks Saturation.
179
180        void DumpPhysicsTable() const;
181        // Prints the fast and slow scintillation integral tables.
182
183protected:
184
185        void BuildThePhysicsTable();
186        // It builds either the fast or slow scintillation integral table;
187        // or both.
188
189        ///////////////////////
190        // Class Data Members
191        ///////////////////////
192
193
194        G4PhysicsTable* theSlowIntegralTable;
195        G4PhysicsTable* theFastIntegralTable;
196
197
198
199        G4bool fTrackSecondariesFirst;
200
201        G4double YieldFactor;
202
203        G4double ExcitationRatio;
204
205private:
206
207        G4EmSaturation* emSaturation;
208
209};
210
211////////////////////
212// Inline methods
213////////////////////
214
215inline 
216G4bool G4Scintillation::IsApplicable(const G4ParticleDefinition& aParticleType)
217{
218        if (aParticleType.GetParticleName() == "opticalphoton"){
219           return false;
220        } else {
221           return true;
222        }
223}
224
225inline 
226void G4Scintillation::SetTrackSecondariesFirst(const G4bool state) 
227{ 
228        fTrackSecondariesFirst = state;
229}
230
231inline
232G4bool G4Scintillation::GetTrackSecondariesFirst() const
233{
234        return fTrackSecondariesFirst;
235}
236
237inline
238void G4Scintillation::SetScintillationYieldFactor(const G4double yieldfactor)
239{
240        YieldFactor = yieldfactor;
241}
242
243inline
244G4double G4Scintillation::GetScintillationYieldFactor() const
245{
246        return YieldFactor;
247}
248
249inline
250void G4Scintillation::SetScintillationExcitationRatio(const G4double excitationratio)
251{
252        ExcitationRatio = excitationratio;
253}
254
255inline
256G4double G4Scintillation::GetScintillationExcitationRatio() const
257{
258        return ExcitationRatio;
259}
260
261inline
262G4PhysicsTable* G4Scintillation::GetSlowIntegralTable() const
263{
264        return theSlowIntegralTable;
265}
266
267inline
268G4PhysicsTable* G4Scintillation::GetFastIntegralTable() const
269{
270        return theFastIntegralTable;
271}
272
273inline
274void G4Scintillation::DumpPhysicsTable() const
275{
276        if (theFastIntegralTable) {
277           G4int PhysicsTableSize = theFastIntegralTable->entries();
278           G4PhysicsOrderedFreeVector *v;
279
280           for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
281           {
282                v = (G4PhysicsOrderedFreeVector*)(*theFastIntegralTable)[i];
283                v->DumpValues();
284           }
285         }
286
287        if (theSlowIntegralTable) {
288           G4int PhysicsTableSize = theSlowIntegralTable->entries();
289           G4PhysicsOrderedFreeVector *v;
290
291           for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
292           {
293                v = (G4PhysicsOrderedFreeVector*)(*theSlowIntegralTable)[i];
294                v->DumpValues();
295           }
296         }
297}
298
299#endif /* G4Scintillation_h */
Note: See TracBrowser for help on using the repository browser.