00001 // 00002 // ******************************************************************** 00003 // * DISCLAIMER * 00004 // * * 00005 // * The following disclaimer summarizes all the specific disclaimers * 00006 // * of contributors to this software. The specific disclaimers,which * 00007 // * govern, are listed with their locations in: * 00008 // * http://cern.ch/geant4/license * 00009 // * * 00010 // * Neither the authors of this software system, nor their employing * 00011 // * institutes,nor the agencies providing financial support for this * 00012 // * work make any representation or warranty, express or implied, * 00013 // * regarding this software system or assume any liability for its * 00014 // * use. * 00015 // * * 00016 // * This code implementation is the intellectual property of the * 00017 // * GEANT4 collaboration. * 00018 // * By copying, distributing or modifying the Program (or any work * 00019 // * based on the Program) you indicate your acceptance of this * 00020 // * statement, and all its terms. * 00021 // ******************************************************************** 00022 // 00023 // 00024 // 00025 // ////////////////////////////////////////////////////////////////////// 00026 // Scintillation Light Class Implementation 00027 // ////////////////////////////////////////////////////////////////////// 00028 // 00029 // File: G4Scintillation.cc 00030 // Description: RestDiscrete Process - Generation of Scintillation Photons 00031 // Version: 1.0 00032 // Created: 1998-11-07 00033 // Author: Peter Gumplinger 00034 // Updated: 2005-08-17 by Peter Gumplinger 00035 // > change variable name MeanNumPhotons -> MeanNumberOfPhotons 00036 // 2005-07-28 by Peter Gumplinger 00037 // > add G4ProcessType to constructor 00038 // 2004-08-05 by Peter Gumplinger 00039 // > changed StronglyForced back to Forced in GetMeanLifeTime 00040 // 2002-11-21 by Peter Gumplinger 00041 // > change to use G4Poisson for small MeanNumberOfPhotons 00042 // 2002-11-07 by Peter Gumplinger 00043 // > now allow for fast and slow scintillation component 00044 // 2002-11-05 by Peter Gumplinger 00045 // > now use scintillation constants from G4Material 00046 // 2002-05-09 by Peter Gumplinger 00047 // > use only the PostStepPoint location for the origin of 00048 // scintillation photons when energy is lost to the medium 00049 // by a neutral particle 00050 // 2000-09-18 by Peter Gumplinger 00051 // > change: aSecondaryPosition=x0+rand*aStep.GetDeltaPosition(); 00052 // aSecondaryTrack->SetTouchable(0); 00053 // 2001-09-17, migration of Materials to pure STL (mma) 00054 // 2003-06-03, V.Ivanchenko fix compilation warnings 00055 // 00056 //mail: gum@triumf.ca 00057 // 00059 00060 //----------------------------------------------------------------- 00061 // DsG4Scintillation 00062 // 00063 // A class modified from G4Scintillation. 00064 // Birks' law is used to calculate the scintillation photon number 00065 // Author: Liang Zhan, 2006/01/27 00066 // Modified: bv@bnl.gov, 2008/4/16 for DetSim 00067 //-------------------------------------------------------------- 00068 00069 #ifndef DsG4Scintillation_h 00070 #define DsG4Scintillation_h 1 00071 00072 #include "globals.hh" 00073 #include "templates.hh" 00074 #include "Randomize.hh" 00075 #include "G4Poisson.hh" 00076 #include "G4ThreeVector.hh" 00077 #include "G4ParticleMomentum.hh" 00078 #include "G4Step.hh" 00079 #include "G4VRestDiscreteProcess.hh" 00080 #include "G4OpticalPhoton.hh" 00081 #include "G4DynamicParticle.hh" 00082 #include "G4Material.hh" 00083 #include "G4PhysicsTable.hh" 00084 #include "G4MaterialPropertiesTable.hh" 00085 #include "G4PhysicsOrderedFreeVector.hh" 00086 #include "G4UImessenger.hh" 00087 #include "DsPhysConsOptical.h" 00088 class G4UIcommand; 00089 class G4UIdirectory; 00090 00092 // Class Description: 00093 // RestDiscrete Process - Generation of Scintillation Photons. 00094 // Class inherits publicly from G4VRestDiscreteProcess. 00095 // Class Description - End: 00096 00098 // Class Definition 00100 00101 class DsG4Scintillation : public G4VRestDiscreteProcess, public G4UImessenger 00102 { //too lazy to create another UImessenger class 00103 00104 private: 00105 00107 // Operators 00109 00110 // DsG4Scintillation& operator=(const DsG4Scintillation &right); 00111 00112 public: // Without description 00113 00115 // Constructors and Destructor 00117 00118 DsG4Scintillation(const G4String& processName = "Scintillation", 00119 G4ProcessType type = fElectromagnetic); 00120 00121 // DsG4Scintillation(const DsG4Scintillation &right); 00122 00123 ~DsG4Scintillation(); 00124 00126 // Methods 00128 00129 public: // With description 00130 00131 // DsG4Scintillation Process has both PostStepDoIt (for energy 00132 // deposition of particles in flight) and AtRestDoIt (for energy 00133 // given to the medium by particles at rest) 00134 00135 G4bool IsApplicable(const G4ParticleDefinition& aParticleType); 00136 // Returns true -> 'is applicable', for any particle type except 00137 // for an 'opticalphoton' 00138 00139 G4double GetMeanFreePath(const G4Track& aTrack, 00140 G4double , 00141 G4ForceCondition* ); 00142 // Returns infinity; i. e. the process does not limit the step, 00143 // but sets the 'StronglyForced' condition for the DoIt to be 00144 // invoked at every step. 00145 00146 G4double GetMeanLifeTime(const G4Track& aTrack, 00147 G4ForceCondition* ); 00148 // Returns infinity; i. e. the process does not limit the time, 00149 // but sets the 'StronglyForced' condition for the DoIt to be 00150 // invoked at every step. 00151 00152 G4VParticleChange* PostStepDoIt(const G4Track& aTrack, 00153 const G4Step& aStep); 00154 G4VParticleChange* AtRestDoIt (const G4Track& aTrack, 00155 const G4Step& aStep); 00156 00157 // These are the methods implementing the scintillation process. 00158 00159 void SetTrackSecondariesFirst(const G4bool state); 00160 // If set, the primary particle tracking is interrupted and any 00161 // produced scintillation photons are tracked next. When all 00162 // have been tracked, the tracking of the primary resumes. 00163 00164 G4bool GetTrackSecondariesFirst() const; 00165 // Returns the boolean flag for tracking secondaries first. 00166 00167 void SetScintillationYieldFactor(const G4double yieldfactor); 00168 // Called to set the scintillation photon yield factor, needed when 00169 // the yield is different for different types of particles. This 00170 // scales the yield obtained from the G4MaterialPropertiesTable. 00171 G4double GetScintillationYieldFactor() const; 00172 // Returns the photon yield factor. 00173 00174 void SetUseFastMu300nsTrick(const G4bool fastMu300nsTrick); 00175 G4bool GetUseFastMu300nsTrick() const; 00176 00177 00178 void SetScintillationExcitationRatio(const G4double excitationratio); 00179 // Called to set the scintillation exciation ratio, needed when 00180 // the scintillation level excitation is different for different 00181 // types of particles. This overwrites the YieldRatio obtained 00182 // from the G4MaterialPropertiesTable. 00183 00184 G4double GetScintillationExcitationRatio() const; 00185 // Returns the scintillation level excitation ratio. 00186 00187 G4PhysicsTable* GetFastIntegralTable() const; 00188 // Returns the address of the fast scintillation integral table. 00189 00190 G4PhysicsTable* GetSlowIntegralTable() const; 00191 // Returns the address of the slow scintillation integral table. 00192 00193 G4PhysicsTable* GetReemissionIntegralTable() const; 00194 // Returns the address of the reemission integral table. 00195 00196 void DumpPhysicsTable() const; 00197 // Prints the fast and slow scintillation integral tables. 00198 00199 // Configuration 00200 G4double GetPhotonWeight() const { return fPhotonWeight; } 00201 void SetPhotonWeight(G4double weight) { fPhotonWeight = weight; } 00202 void SetDoReemission(bool tf = true) { doReemission = tf; } 00203 bool GetDoReemission() { return doReemission; } 00204 void SetDoBothProcess(bool tf = true) { doBothProcess = tf; } 00205 bool GetDoBothProcess() { return doBothProcess; } 00206 00207 G4bool GetApplyPreQE() const { return fApplyPreQE; } 00208 void SetApplyPreQE(G4bool a) { fApplyPreQE = a; } 00209 00210 G4double GetPreQE() const { return fPreQE; } 00211 void SetPreQE(G4double a) { fPreQE = a; } 00212 00213 void SetBirksConstant1(double c1) {birksConstant1 = c1;} 00214 double GetBirksConstant1() {return birksConstant1;} 00215 void SetBirksConstant2(double c2) {birksConstant2 = c2;} 00216 double GetBirksConstant2() {return birksConstant2;} 00217 00218 // Don't actually do anything. This is needed, as apposed to 00219 // simply not including the scintilation process, because 00220 // w/out this process no optical photons make it into the 00221 // photocathode (???) 00222 void SetNoOp(bool tf = true) { m_noop = tf; } 00223 private: 00224 00225 void BuildThePhysicsTable(); 00226 // It builds either the fast or slow scintillation integral table; 00227 // or both. 00228 00230 // Class Data Members 00232 00233 protected: 00234 00235 G4PhysicsTable* theSlowIntegralTable; 00236 G4PhysicsTable* theFastIntegralTable; 00237 G4PhysicsTable* theReemissionIntegralTable; 00238 00239 // on/off flag for absorbed opticalphoton reemission 00240 G4bool doReemission; 00241 // choose only reemission of Cerenkov or both of Cerenkov and Scintillaton; 00242 G4bool doBothProcess; 00243 00244 // Birks constant C1 and C2 00245 double birksConstant1; 00246 double birksConstant2; 00247 00248 private: 00249 00250 G4bool fTrackSecondariesFirst; 00251 00252 G4double YieldFactor; 00253 G4bool FastMu300nsTrick; 00254 G4double ExcitationRatio; 00255 00256 //mean number of true photons per secondary track in GLG4Scint 00257 G4double fPhotonWeight; 00258 G4bool fApplyPreQE; 00259 G4double fPreQE; 00260 bool m_noop; 00261 }; 00262 00264 // Inline methods 00266 00267 inline 00268 G4bool DsG4Scintillation::IsApplicable(const G4ParticleDefinition& aParticleType) 00269 { 00270 if (aParticleType.GetParticleName() == "opticalphoton"){ 00271 return true; 00272 } else { 00273 return true; 00274 } 00275 } 00276 00277 inline 00278 void DsG4Scintillation::SetTrackSecondariesFirst(const G4bool state) 00279 { 00280 fTrackSecondariesFirst = state; 00281 } 00282 00283 inline 00284 G4bool DsG4Scintillation::GetTrackSecondariesFirst() const 00285 { 00286 return fTrackSecondariesFirst; 00287 } 00288 00289 inline 00290 void DsG4Scintillation::SetScintillationYieldFactor(const G4double yieldfactor) 00291 { 00292 YieldFactor = yieldfactor; 00293 } 00294 00295 00296 inline 00297 G4double DsG4Scintillation::GetScintillationYieldFactor() const 00298 { 00299 return YieldFactor; 00300 } 00301 00302 inline 00303 void DsG4Scintillation::SetUseFastMu300nsTrick(const G4bool fastMu300nsTrick) 00304 { 00305 FastMu300nsTrick = fastMu300nsTrick; 00306 } 00307 00308 inline 00309 G4bool DsG4Scintillation::GetUseFastMu300nsTrick() const 00310 { 00311 return FastMu300nsTrick; 00312 } 00313 00314 00315 00316 00317 inline 00318 void DsG4Scintillation::SetScintillationExcitationRatio(const G4double excitationratio) 00319 { 00320 ExcitationRatio = excitationratio; 00321 } 00322 00323 inline 00324 G4double DsG4Scintillation::GetScintillationExcitationRatio() const 00325 { 00326 return ExcitationRatio; 00327 } 00328 00329 inline 00330 G4PhysicsTable* DsG4Scintillation::GetSlowIntegralTable() const 00331 { 00332 return theSlowIntegralTable; 00333 } 00334 00335 inline 00336 G4PhysicsTable* DsG4Scintillation::GetFastIntegralTable() const 00337 { 00338 return theFastIntegralTable; 00339 } 00340 00341 inline 00342 G4PhysicsTable* DsG4Scintillation::GetReemissionIntegralTable() const 00343 { 00344 return theReemissionIntegralTable; 00345 } 00346 00347 inline 00348 void DsG4Scintillation::DumpPhysicsTable() const 00349 { 00350 if (theFastIntegralTable) { 00351 G4int PhysicsTableSize = theFastIntegralTable->entries(); 00352 G4PhysicsOrderedFreeVector *v; 00353 00354 for (G4int i = 0 ; i < PhysicsTableSize ; i++ ) 00355 { 00356 v = (G4PhysicsOrderedFreeVector*)(*theFastIntegralTable)[i]; 00357 v->DumpValues(); 00358 } 00359 } 00360 00361 if (theSlowIntegralTable) { 00362 G4int PhysicsTableSize = theSlowIntegralTable->entries(); 00363 G4PhysicsOrderedFreeVector *v; 00364 00365 for (G4int i = 0 ; i < PhysicsTableSize ; i++ ) 00366 { 00367 v = (G4PhysicsOrderedFreeVector*)(*theSlowIntegralTable)[i]; 00368 v->DumpValues(); 00369 } 00370 } 00371 00372 if (theReemissionIntegralTable) { 00373 G4int PhysicsTableSize = theReemissionIntegralTable->entries(); 00374 G4PhysicsOrderedFreeVector *v; 00375 00376 for (G4int i = 0 ; i < PhysicsTableSize ; i++ ) 00377 { 00378 v = (G4PhysicsOrderedFreeVector*)(*theReemissionIntegralTable)[i]; 00379 v->DumpValues(); 00380 } 00381 } 00382 } 00383 00384 #endif /* DsG4Scintillation_h */