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

Last change on this file since 846 was 819, checked in by garnier, 16 years ago

import all except CVS

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