#include <DsG4Scintillation.h>
Public Member Functions | |
DsG4Scintillation (const G4String &processName="Scintillation", G4ProcessType type=fElectromagnetic) | |
~DsG4Scintillation () | |
G4bool | IsApplicable (const G4ParticleDefinition &aParticleType) |
G4double | GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) |
G4double | GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *) |
G4VParticleChange * | PostStepDoIt (const G4Track &aTrack, const G4Step &aStep) |
G4VParticleChange * | AtRestDoIt (const G4Track &aTrack, const G4Step &aStep) |
void | SetTrackSecondariesFirst (const G4bool state) |
G4bool | GetTrackSecondariesFirst () const |
void | SetScintillationYieldFactor (const G4double yieldfactor) |
G4double | GetScintillationYieldFactor () const |
void | SetUseFastMu300nsTrick (const G4bool fastMu300nsTrick) |
G4bool | GetUseFastMu300nsTrick () const |
void | SetScintillationExcitationRatio (const G4double excitationratio) |
G4double | GetScintillationExcitationRatio () const |
G4PhysicsTable * | GetFastIntegralTable () const |
G4PhysicsTable * | GetSlowIntegralTable () const |
G4PhysicsTable * | GetReemissionIntegralTable () const |
void | DumpPhysicsTable () const |
G4double | GetPhotonWeight () const |
void | SetPhotonWeight (G4double weight) |
void | SetDoReemission (bool tf=true) |
bool | GetDoReemission () |
void | SetDoBothProcess (bool tf=true) |
bool | GetDoBothProcess () |
G4bool | GetApplyPreQE () const |
void | SetApplyPreQE (G4bool a) |
G4double | GetPreQE () const |
void | SetPreQE (G4double a) |
void | SetBirksConstant1 (double c1) |
double | GetBirksConstant1 () |
void | SetBirksConstant2 (double c2) |
double | GetBirksConstant2 () |
void | SetNoOp (bool tf=true) |
Protected Attributes | |
G4PhysicsTable * | theSlowIntegralTable |
G4PhysicsTable * | theFastIntegralTable |
G4PhysicsTable * | theReemissionIntegralTable |
G4bool | doReemission |
G4bool | doBothProcess |
double | birksConstant1 |
double | birksConstant2 |
Private Member Functions | |
void | BuildThePhysicsTable () |
Private Attributes | |
G4bool | fTrackSecondariesFirst |
G4double | YieldFactor |
G4bool | FastMu300nsTrick |
G4double | ExcitationRatio |
G4double | fPhotonWeight |
G4bool | fApplyPreQE |
G4double | fPreQE |
bool | m_noop |
Definition at line 101 of file DsG4Scintillation.h.
DsG4Scintillation::DsG4Scintillation | ( | const G4String & | processName = "Scintillation" , |
|
G4ProcessType | type = fElectromagnetic | |||
) |
Definition at line 98 of file DsG4Scintillation.cc.
00100 : G4VRestDiscreteProcess(processName, type) 00101 , doReemission(true) 00102 , doBothProcess(true) 00103 , fPhotonWeight(1.0) 00104 , fApplyPreQE(false) 00105 , fPreQE(1.) 00106 , m_noop(false) 00107 { 00108 SetProcessSubType(fScintillation); 00109 fTrackSecondariesFirst = false; 00110 00111 YieldFactor = 1.0; 00112 ExcitationRatio = 1.0; 00113 00114 theFastIntegralTable = NULL; 00115 theSlowIntegralTable = NULL; 00116 theReemissionIntegralTable = NULL; 00117 00118 if (verboseLevel>0) { 00119 G4cout << GetProcessName() << " is created " << G4endl; 00120 } 00121 00122 BuildThePhysicsTable(); 00123 00124 }
DsG4Scintillation::~DsG4Scintillation | ( | ) |
Definition at line 130 of file DsG4Scintillation.cc.
00131 { 00132 if (theFastIntegralTable != NULL) { 00133 theFastIntegralTable->clearAndDestroy(); 00134 delete theFastIntegralTable; 00135 } 00136 if (theSlowIntegralTable != NULL) { 00137 theSlowIntegralTable->clearAndDestroy(); 00138 delete theSlowIntegralTable; 00139 } 00140 if (theReemissionIntegralTable != NULL) { 00141 theReemissionIntegralTable->clearAndDestroy(); 00142 delete theReemissionIntegralTable; 00143 } 00144 }
G4bool DsG4Scintillation::IsApplicable | ( | const G4ParticleDefinition & | aParticleType | ) | [inline] |
Definition at line 268 of file DsG4Scintillation.h.
00269 { 00270 if (aParticleType.GetParticleName() == "opticalphoton"){ 00271 return true; 00272 } else { 00273 return true; 00274 } 00275 }
G4double DsG4Scintillation::GetMeanFreePath | ( | const G4Track & | aTrack, | |
G4double | , | |||
G4ForceCondition * | ||||
) |
G4double DsG4Scintillation::GetMeanLifeTime | ( | const G4Track & | aTrack, | |
G4ForceCondition * | ||||
) |
G4VParticleChange * DsG4Scintillation::PostStepDoIt | ( | const G4Track & | aTrack, | |
const G4Step & | aStep | |||
) |
Definition at line 167 of file DsG4Scintillation.cc.
00174 { 00175 aParticleChange.Initialize(aTrack); 00176 00177 if (m_noop) { // do nothing, bail 00178 aParticleChange.SetNumberOfSecondaries(0); 00179 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00180 } 00181 00182 00183 G4String pname=""; 00184 G4ThreeVector vertpos; 00185 G4double vertenergy=0.0; 00186 G4double reem_d=0.0; 00187 G4bool flagReemission= false; 00188 if (aTrack.GetDefinition() == G4OpticalPhoton::OpticalPhoton()) { 00189 const G4VProcess* process = aStep.GetTrack()->GetCreatorProcess(); 00190 if(process) pname = process->GetProcessName(); 00191 00192 if (verboseLevel>0) { 00193 G4cout<<"Optical photon!!!!!!!!!!!"<<G4endl; 00194 } 00195 if(doBothProcess) { 00196 flagReemission= doReemission 00197 && aTrack.GetTrackStatus() == fStopAndKill 00198 && aStep.GetPostStepPoint()->GetStepStatus() != fGeomBoundary; 00199 } 00200 else{ 00201 flagReemission= doReemission 00202 && aTrack.GetTrackStatus() == fStopAndKill 00203 && aStep.GetPostStepPoint()->GetStepStatus() != fGeomBoundary 00204 && pname=="Cerenkov"; 00205 } 00206 if(verboseLevel>0) { 00207 G4cout<<"flag of Reemission is "<<flagReemission<<"!!"<<G4endl; 00208 } 00209 if (!flagReemission) { 00210 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00211 } 00212 } 00213 00214 G4double TotalEnergyDeposit = aStep.GetTotalEnergyDeposit(); 00215 00216 if (TotalEnergyDeposit <= 0.0 && !flagReemission) { 00217 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00218 } 00219 00220 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle(); 00221 const G4Material* aMaterial = aTrack.GetMaterial(); 00222 00223 G4MaterialPropertiesTable* aMaterialPropertiesTable = 00224 aMaterial->GetMaterialPropertiesTable(); 00225 if (!aMaterialPropertiesTable) 00226 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00227 00228 const G4MaterialPropertyVector* Fast_Intensity = 00229 aMaterialPropertiesTable->GetProperty("FASTCOMPONENT"); 00230 const G4MaterialPropertyVector* Slow_Intensity = 00231 aMaterialPropertiesTable->GetProperty("SLOWCOMPONENT"); 00232 const G4MaterialPropertyVector* Reemission_Prob = 00233 aMaterialPropertiesTable->GetProperty("REEMISSIONPROB"); 00234 00235 if (!Fast_Intensity && !Slow_Intensity ) 00236 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00237 00238 G4int nscnt = 1; 00239 if (Fast_Intensity && Slow_Intensity) nscnt = 2; 00240 00241 G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint(); 00242 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint(); 00243 00244 G4ThreeVector x0 = pPreStepPoint->GetPosition(); 00245 G4ThreeVector p0 = aStep.GetDeltaPosition().unit(); 00246 G4double t0 = pPreStepPoint->GetGlobalTime(); 00247 00248 //Replace NumPhotons by NumTracks 00249 G4int NumTracks=0; 00250 G4double weight=1.0; 00251 if (flagReemission) { 00252 if(verboseLevel>0){ 00253 G4cout<<"the process name is "<<pname<<"!!"<<G4endl;} 00254 00255 if ( Reemission_Prob == 0) 00256 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00257 G4double p_reemission= 00258 Reemission_Prob->GetProperty(aTrack.GetKineticEnergy()); 00259 if (G4UniformRand() >= p_reemission) 00260 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00261 NumTracks= 1; 00262 weight= aTrack.GetWeight(); 00263 } 00264 else { 00266 // J.B.Birks. The theory and practice of Scintillation Counting. 00267 // Pergamon Press, 1964. 00268 // For particles with energy much smaller than minimum ionization 00269 // energy, the scintillation response is non-linear because of quenching 00270 // effect. The light output is reduced by a parametric factor: 00271 // 1/(1 + birk1*delta + birk2* delta^2). 00272 // Delta is the energy loss per unit mass thickness. birk1 and birk2 00273 // were measured for several organic scintillators. 00274 // Here we use birk1 = 0.0125*g/cm2/MeV and ignore birk2. 00275 // R.L.Craun and D.L.Smith. Nucl. Inst. and Meth., 80:239-244, 1970. 00276 // Liang Zhan 01/27/2006 00277 // ///////////////////////////////////////////////////////////////////// 00278 00279 00280 G4double ScintillationYield = 0; 00281 {// Yield. Material must have this or we lack raisins dayetras 00282 const G4MaterialPropertyVector* ptable = 00283 aMaterialPropertiesTable->GetProperty("SCINTILLATIONYIELD"); 00284 if (!ptable) { 00285 G4cout << "ConstProperty: failed to get SCINTILLATIONYIELD" 00286 << G4endl; 00287 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00288 } 00289 ScintillationYield = ptable->GetProperty(0); 00290 } 00291 00292 G4double ResolutionScale = 1; 00293 {// Resolution Scale 00294 const G4MaterialPropertyVector* ptable = 00295 aMaterialPropertiesTable->GetProperty("RESOLUTIONSCALE"); 00296 if (ptable) 00297 ResolutionScale = ptable->GetProperty(0); 00298 } 00299 00300 G4double dE = TotalEnergyDeposit; 00301 G4double dx = aStep.GetStepLength(); 00302 G4double dE_dx = dE/dx; 00303 if(aTrack.GetDefinition() == G4Gamma::Gamma() && dE > 0) 00304 { 00305 G4LossTableManager* manager = G4LossTableManager::Instance(); 00306 dE_dx = dE/manager->GetRange(G4Electron::Electron(), dE, aTrack.GetMaterialCutsCouple()); 00307 //G4cout<<"gamma dE_dx = "<<dE_dx/(MeV/mm)<<"MeV/mm"<<G4endl; 00308 } 00309 00310 G4double delta = dE_dx/aMaterial->GetDensity();//get scintillator density 00311 //G4double birk1 = 0.0125*g/cm2/MeV; 00312 G4double birk1 = birksConstant1; 00313 if(abs(aParticle->GetCharge())>1.5)//for particle charge greater than 1. 00314 birk1 = 0.57*birk1; 00315 00316 G4double birk2 = 0; 00317 //birk2 = (0.0031*g/MeV/cm2)*(0.0031*g/MeV/cm2); 00318 birk2 = birksConstant2; 00319 00320 G4double QuenchedTotalEnergyDeposit 00321 = TotalEnergyDeposit/(1+birk1*delta+birk2*delta*delta); 00322 00323 //Add 300ns trick for muon simuation, by Haoqi Jan 27, 2011 00324 if(FastMu300nsTrick) { 00325 // cout<<"GlobalTime ="<<aStep.GetTrack()->GetGlobalTime()/ns<<endl; 00326 if(aStep.GetTrack()->GetGlobalTime()/ns>300) { 00327 ScintillationYield = YieldFactor * ScintillationYield; 00328 } 00329 else{ 00330 ScintillationYield=0.; 00331 } 00332 } 00333 else { 00334 ScintillationYield = YieldFactor * ScintillationYield; 00335 } 00336 00337 G4double MeanNumberOfPhotons= ScintillationYield * QuenchedTotalEnergyDeposit; 00338 00339 // Implemented the fast simulation method from GLG4Scint 00340 // Jianglai 09-05-2006 00341 00342 // randomize number of TRACKS (not photons) 00343 // this gets statistics right for number of PE after applying 00344 // boolean random choice to final absorbed track (change from 00345 // old method of applying binomial random choice to final absorbed 00346 // track, which did want poissonian number of photons divided 00347 // as evenly as possible into tracks) 00348 // Note for weight=1, there's no difference between tracks and photons. 00349 G4double MeanNumberOfTracks= MeanNumberOfPhotons/fPhotonWeight; 00350 if ( fApplyPreQE ) MeanNumberOfTracks*=fPreQE; 00351 if (MeanNumberOfTracks > 10.) { 00352 G4double sigma = ResolutionScale * sqrt(MeanNumberOfTracks); 00353 NumTracks = G4int(G4RandGauss::shoot(MeanNumberOfTracks,sigma)+0.5); 00354 } 00355 else { 00356 NumTracks = G4int(G4Poisson(MeanNumberOfTracks)); 00357 } 00358 } 00359 weight*=fPhotonWeight; 00360 // G4cerr<<"Scint weight is "<<weight<<G4endl; 00361 if (NumTracks <= 0) { 00362 // return unchanged particle and no secondaries 00363 aParticleChange.SetNumberOfSecondaries(0); 00364 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00365 } 00366 00368 00369 aParticleChange.SetNumberOfSecondaries(NumTracks); 00370 00371 if (fTrackSecondariesFirst) { 00372 if (!flagReemission) 00373 if (aTrack.GetTrackStatus() == fAlive ) 00374 aParticleChange.ProposeTrackStatus(fSuspend); 00375 } 00376 00378 00379 G4int materialIndex = aMaterial->GetIndex(); 00380 00381 G4PhysicsOrderedFreeVector* ReemissionIntegral = NULL; 00382 ReemissionIntegral = 00383 (G4PhysicsOrderedFreeVector*)((*theReemissionIntegralTable)(materialIndex)); 00384 00385 // Retrieve the Scintillation Integral for this material 00386 // new G4PhysicsOrderedFreeVector allocated to hold CII's 00387 00388 G4int Num = NumTracks; //# tracks is now the loop control 00389 00390 G4double fastTimeConstant = 0.0; 00391 { // Fast Time Constant 00392 const G4MaterialPropertyVector* ptable = 00393 aMaterialPropertiesTable->GetProperty("FASTTIMECONSTANT"); 00394 if (ptable) 00395 fastTimeConstant = ptable->GetProperty(0); 00396 } 00397 00398 G4double slowTimeConstant = 0.0; 00399 { // Slow Time Constant 00400 const G4MaterialPropertyVector* ptable = 00401 aMaterialPropertiesTable->GetProperty("SLOWTIMECONSTANT"); 00402 if (ptable) 00403 slowTimeConstant = ptable->GetProperty(0); 00404 } 00405 00406 G4double YieldRatio = 0.0; 00407 { // Slow Time Constant 00408 const G4MaterialPropertyVector* ptable = 00409 aMaterialPropertiesTable->GetProperty("YIELDRATIO"); 00410 if (ptable) 00411 YieldRatio = ptable->GetProperty(0); 00412 } 00413 00414 //loop over fast/slow scintillations 00415 for (G4int scnt = 1; scnt <= nscnt; scnt++) { 00416 00417 G4double ScintillationTime = 0.*ns; 00418 G4PhysicsOrderedFreeVector* ScintillationIntegral = NULL; 00419 00420 if (scnt == 1) {//fast 00421 if (nscnt == 1) { 00422 if(Fast_Intensity){ 00423 ScintillationTime = fastTimeConstant; 00424 ScintillationIntegral = 00425 (G4PhysicsOrderedFreeVector*)((*theFastIntegralTable)(materialIndex)); 00426 } 00427 if(Slow_Intensity){ 00428 ScintillationTime = slowTimeConstant; 00429 ScintillationIntegral = 00430 (G4PhysicsOrderedFreeVector*)((*theSlowIntegralTable)(materialIndex)); 00431 } 00432 } 00433 else { 00434 if ( ExcitationRatio == 1.0 ) { 00435 Num = G4int (min(YieldRatio,1.0) * NumTracks); 00436 } 00437 else { 00438 Num = G4int (min(ExcitationRatio,1.0) * NumTracks); 00439 } 00440 ScintillationTime = fastTimeConstant; 00441 ScintillationIntegral = 00442 (G4PhysicsOrderedFreeVector*)((*theFastIntegralTable)(materialIndex)); 00443 } 00444 } 00445 else {//slow 00446 Num = NumTracks - Num; 00447 ScintillationTime = slowTimeConstant; 00448 ScintillationIntegral = 00449 (G4PhysicsOrderedFreeVector*)((*theSlowIntegralTable)(materialIndex)); 00450 } 00451 00452 if (!ScintillationIntegral) continue; 00453 00454 // Max Scintillation Integral 00455 00456 for (G4int i = 0; i < Num; i++) { //Num is # of 2ndary tracks now 00457 // Determine photon energy 00458 00459 G4double sampledEnergy; 00460 if ( !flagReemission ) { 00461 // normal scintillation 00462 G4double CIIvalue = G4UniformRand()* 00463 ScintillationIntegral->GetMaxValue(); 00464 sampledEnergy= 00465 ScintillationIntegral->GetEnergy(CIIvalue); 00466 00467 if (verboseLevel>1) 00468 { 00469 G4cout << "sampledEnergy = " << sampledEnergy << G4endl; 00470 G4cout << "CIIvalue = " << CIIvalue << G4endl; 00471 } 00472 } 00473 else { 00474 // reemission, the sample method need modification 00475 G4double CIIvalue = G4UniformRand()* 00476 ScintillationIntegral->GetMaxValue(); 00477 if (CIIvalue == 0.0) { 00478 // return unchanged particle and no secondaries 00479 aParticleChange.SetNumberOfSecondaries(0); 00480 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00481 } 00482 sampledEnergy= 00483 ScintillationIntegral->GetEnergy(CIIvalue); 00484 if (verboseLevel>1) { 00485 G4cout << "oldEnergy = " <<aTrack.GetKineticEnergy() << G4endl; 00486 G4cout << "reemittedSampledEnergy = " << sampledEnergy 00487 << "\nreemittedCIIvalue = " << CIIvalue << G4endl; 00488 } 00489 } 00490 00491 // Generate random photon direction 00492 00493 G4double cost = 1. - 2.*G4UniformRand(); 00494 G4double sint = sqrt((1.-cost)*(1.+cost)); 00495 00496 G4double phi = twopi*G4UniformRand(); 00497 G4double sinp = sin(phi); 00498 G4double cosp = cos(phi); 00499 00500 G4double px = sint*cosp; 00501 G4double py = sint*sinp; 00502 G4double pz = cost; 00503 00504 // Create photon momentum direction vector 00505 00506 G4ParticleMomentum photonMomentum(px, py, pz); 00507 00508 // Determine polarization of new photon 00509 00510 G4double sx = cost*cosp; 00511 G4double sy = cost*sinp; 00512 G4double sz = -sint; 00513 00514 G4ThreeVector photonPolarization(sx, sy, sz); 00515 00516 G4ThreeVector perp = photonMomentum.cross(photonPolarization); 00517 00518 phi = twopi*G4UniformRand(); 00519 sinp = sin(phi); 00520 cosp = cos(phi); 00521 00522 photonPolarization = cosp * photonPolarization + sinp * perp; 00523 00524 photonPolarization = photonPolarization.unit(); 00525 00526 // Generate a new photon: 00527 00528 G4DynamicParticle* aScintillationPhoton = 00529 new G4DynamicParticle(G4OpticalPhoton::OpticalPhoton(), 00530 photonMomentum); 00531 aScintillationPhoton->SetPolarization 00532 (photonPolarization.x(), 00533 photonPolarization.y(), 00534 photonPolarization.z()); 00535 00536 aScintillationPhoton->SetKineticEnergy(sampledEnergy); 00537 00538 // Generate new G4Track object: 00539 00540 G4double rand=0; 00541 G4ThreeVector aSecondaryPosition; 00542 G4double deltaTime; 00543 if (flagReemission) { 00544 deltaTime= pPostStepPoint->GetGlobalTime() - t0 00545 -ScintillationTime * log( G4UniformRand() ); 00546 aSecondaryPosition= pPostStepPoint->GetPosition(); 00547 vertpos = aTrack.GetVertexPosition(); 00548 vertenergy = aTrack.GetKineticEnergy(); 00549 reem_d = 00550 sqrt( pow( aSecondaryPosition.x()-vertpos.x(), 2) 00551 + pow( aSecondaryPosition.y()-vertpos.y(), 2) 00552 + pow( aSecondaryPosition.z()-vertpos.z(), 2) ); 00553 } 00554 else { 00555 if (aParticle->GetDefinition()->GetPDGCharge() != 0) 00556 { 00557 rand = G4UniformRand(); 00558 } 00559 else 00560 { 00561 rand = 1.0; 00562 } 00563 00564 G4double delta = rand * aStep.GetStepLength(); 00565 deltaTime = delta / 00566 ((pPreStepPoint->GetVelocity()+ 00567 pPostStepPoint->GetVelocity())/2.); 00568 00569 deltaTime = deltaTime - 00570 ScintillationTime * log( G4UniformRand() ); 00571 00572 aSecondaryPosition = 00573 x0 + rand * aStep.GetDeltaPosition(); 00574 } 00575 G4double aSecondaryTime = t0 + deltaTime; 00576 00577 G4Track* aSecondaryTrack = 00578 new G4Track(aScintillationPhoton,aSecondaryTime,aSecondaryPosition); 00579 00580 aSecondaryTrack->SetWeight( weight ); 00581 aSecondaryTrack->SetTouchableHandle(aStep.GetPreStepPoint()->GetTouchableHandle()); 00582 // aSecondaryTrack->SetTouchableHandle((G4VTouchable*)0);//this is wrong 00583 00584 aSecondaryTrack->SetParentID(aTrack.GetTrackID()); 00585 00586 // add the secondary to the ParticleChange object 00587 aParticleChange.SetSecondaryWeightByProcess( true ); // recommended 00588 aParticleChange.AddSecondary(aSecondaryTrack); 00589 00590 aSecondaryTrack->SetWeight( weight ); 00591 // The above line is necessary because AddSecondary() 00592 // overrides our setting of the secondary track weight, 00593 // in Geant4.3.1 & earlier. (and also later, at least 00594 // until Geant4.7 (and beyond?) 00595 // -- maybe not required if SetWeightByProcess(true) called, 00596 // but we do both, just to be sure) 00597 } 00598 } 00599 00600 if (verboseLevel>0) { 00601 G4cout << "\n Exiting from G4Scintillation::DoIt -- NumberOfSecondaries = " 00602 << aParticleChange.GetNumberOfSecondaries() << G4endl; 00603 } 00604 00605 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00606 }
G4VParticleChange * DsG4Scintillation::AtRestDoIt | ( | const G4Track & | aTrack, | |
const G4Step & | aStep | |||
) |
Definition at line 154 of file DsG4Scintillation.cc.
00159 { 00160 return DsG4Scintillation::PostStepDoIt(aTrack, aStep); 00161 }
void DsG4Scintillation::SetTrackSecondariesFirst | ( | const G4bool | state | ) | [inline] |
Definition at line 278 of file DsG4Scintillation.h.
00279 { 00280 fTrackSecondariesFirst = state; 00281 }
G4bool DsG4Scintillation::GetTrackSecondariesFirst | ( | ) | const [inline] |
Definition at line 284 of file DsG4Scintillation.h.
00285 { 00286 return fTrackSecondariesFirst; 00287 }
void DsG4Scintillation::SetScintillationYieldFactor | ( | const G4double | yieldfactor | ) | [inline] |
G4double DsG4Scintillation::GetScintillationYieldFactor | ( | ) | const [inline] |
void DsG4Scintillation::SetUseFastMu300nsTrick | ( | const G4bool | fastMu300nsTrick | ) | [inline] |
Definition at line 303 of file DsG4Scintillation.h.
00304 { 00305 FastMu300nsTrick = fastMu300nsTrick; 00306 }
G4bool DsG4Scintillation::GetUseFastMu300nsTrick | ( | ) | const [inline] |
void DsG4Scintillation::SetScintillationExcitationRatio | ( | const G4double | excitationratio | ) | [inline] |
Definition at line 318 of file DsG4Scintillation.h.
00319 { 00320 ExcitationRatio = excitationratio; 00321 }
G4double DsG4Scintillation::GetScintillationExcitationRatio | ( | ) | const [inline] |
G4PhysicsTable * DsG4Scintillation::GetFastIntegralTable | ( | ) | const [inline] |
Definition at line 336 of file DsG4Scintillation.h.
00337 { 00338 return theFastIntegralTable; 00339 }
G4PhysicsTable * DsG4Scintillation::GetSlowIntegralTable | ( | ) | const [inline] |
Definition at line 330 of file DsG4Scintillation.h.
00331 { 00332 return theSlowIntegralTable; 00333 }
G4PhysicsTable * DsG4Scintillation::GetReemissionIntegralTable | ( | ) | const [inline] |
Definition at line 342 of file DsG4Scintillation.h.
00343 { 00344 return theReemissionIntegralTable; 00345 }
void DsG4Scintillation::DumpPhysicsTable | ( | ) | const [inline] |
Definition at line 348 of file DsG4Scintillation.h.
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 }
G4double DsG4Scintillation::GetPhotonWeight | ( | ) | const [inline] |
void DsG4Scintillation::SetPhotonWeight | ( | G4double | weight | ) | [inline] |
void DsG4Scintillation::SetDoReemission | ( | bool | tf = true |
) | [inline] |
bool DsG4Scintillation::GetDoReemission | ( | ) | [inline] |
void DsG4Scintillation::SetDoBothProcess | ( | bool | tf = true |
) | [inline] |
bool DsG4Scintillation::GetDoBothProcess | ( | ) | [inline] |
G4bool DsG4Scintillation::GetApplyPreQE | ( | ) | const [inline] |
void DsG4Scintillation::SetApplyPreQE | ( | G4bool | a | ) | [inline] |
G4double DsG4Scintillation::GetPreQE | ( | ) | const [inline] |
void DsG4Scintillation::SetPreQE | ( | G4double | a | ) | [inline] |
void DsG4Scintillation::SetBirksConstant1 | ( | double | c1 | ) | [inline] |
double DsG4Scintillation::GetBirksConstant1 | ( | ) | [inline] |
void DsG4Scintillation::SetBirksConstant2 | ( | double | c2 | ) | [inline] |
double DsG4Scintillation::GetBirksConstant2 | ( | ) | [inline] |
void DsG4Scintillation::SetNoOp | ( | bool | tf = true |
) | [inline] |
void DsG4Scintillation::BuildThePhysicsTable | ( | ) | [private] |
Definition at line 612 of file DsG4Scintillation.cc.
00613 { 00614 if (theFastIntegralTable && theSlowIntegralTable && theReemissionIntegralTable) return; 00615 00616 const G4MaterialTable* theMaterialTable = 00617 G4Material::GetMaterialTable(); 00618 G4int numOfMaterials = G4Material::GetNumberOfMaterials(); 00619 00620 // create new physics table 00621 00622 if(!theFastIntegralTable)theFastIntegralTable = new G4PhysicsTable(numOfMaterials); 00623 if(!theSlowIntegralTable)theSlowIntegralTable = new G4PhysicsTable(numOfMaterials); 00624 if(!theReemissionIntegralTable)theReemissionIntegralTable 00625 = new G4PhysicsTable(numOfMaterials); 00626 00627 // loop for materials 00628 00629 for (G4int i=0 ; i < numOfMaterials; i++) { 00630 G4PhysicsOrderedFreeVector* aPhysicsOrderedFreeVector = 00631 new G4PhysicsOrderedFreeVector(); 00632 G4PhysicsOrderedFreeVector* bPhysicsOrderedFreeVector = 00633 new G4PhysicsOrderedFreeVector(); 00634 G4PhysicsOrderedFreeVector* cPhysicsOrderedFreeVector = 00635 new G4PhysicsOrderedFreeVector(); 00636 00637 // Retrieve vector of scintillation wavelength intensity for 00638 // the material from the material's optical properties table. 00639 00640 G4Material* aMaterial = (*theMaterialTable)[i]; 00641 00642 G4MaterialPropertiesTable* aMaterialPropertiesTable = 00643 aMaterial->GetMaterialPropertiesTable(); 00644 00645 if (aMaterialPropertiesTable) { 00646 00647 G4MaterialPropertyVector* theFastLightVector = 00648 aMaterialPropertiesTable->GetProperty("FASTCOMPONENT"); 00649 00650 if (theFastLightVector) { 00651 00652 // Retrieve the first intensity point in vector 00653 // of (photon energy, intensity) pairs 00654 00655 theFastLightVector->ResetIterator(); 00656 ++(*theFastLightVector); // advance to 1st entry 00657 00658 G4double currentIN = theFastLightVector-> 00659 GetProperty(); 00660 00661 if (currentIN >= 0.0) { 00662 00663 // Create first (photon energy, Scintillation 00664 // Integral pair 00665 00666 G4double currentPM = theFastLightVector-> 00667 GetPhotonEnergy(); 00668 00669 G4double currentCII = 0.0; 00670 00671 aPhysicsOrderedFreeVector-> 00672 InsertValues(currentPM , currentCII); 00673 00674 // Set previous values to current ones prior to loop 00675 00676 G4double prevPM = currentPM; 00677 G4double prevCII = currentCII; 00678 G4double prevIN = currentIN; 00679 00680 // loop over all (photon energy, intensity) 00681 // pairs stored for this material 00682 00683 while(++(*theFastLightVector)) { 00684 currentPM = theFastLightVector-> 00685 GetPhotonEnergy(); 00686 00687 currentIN=theFastLightVector-> 00688 GetProperty(); 00689 00690 currentCII = 0.5 * (prevIN + currentIN); 00691 00692 currentCII = prevCII + 00693 (currentPM - prevPM) * currentCII; 00694 00695 aPhysicsOrderedFreeVector-> 00696 InsertValues(currentPM, currentCII); 00697 00698 prevPM = currentPM; 00699 prevCII = currentCII; 00700 prevIN = currentIN; 00701 } 00702 00703 } 00704 } 00705 00706 G4MaterialPropertyVector* theSlowLightVector = 00707 aMaterialPropertiesTable->GetProperty("SLOWCOMPONENT"); 00708 00709 if (theSlowLightVector) { 00710 00711 // Retrieve the first intensity point in vector 00712 // of (photon energy, intensity) pairs 00713 00714 theSlowLightVector->ResetIterator(); 00715 ++(*theSlowLightVector); // advance to 1st entry 00716 00717 G4double currentIN = theSlowLightVector-> 00718 GetProperty(); 00719 00720 if (currentIN >= 0.0) { 00721 00722 // Create first (photon energy, Scintillation 00723 // Integral pair 00724 00725 G4double currentPM = theSlowLightVector-> 00726 GetPhotonEnergy(); 00727 00728 G4double currentCII = 0.0; 00729 00730 bPhysicsOrderedFreeVector-> 00731 InsertValues(currentPM , currentCII); 00732 00733 // Set previous values to current ones prior to loop 00734 00735 G4double prevPM = currentPM; 00736 G4double prevCII = currentCII; 00737 G4double prevIN = currentIN; 00738 00739 // loop over all (photon energy, intensity) 00740 // pairs stored for this material 00741 00742 while(++(*theSlowLightVector)) { 00743 currentPM = theSlowLightVector-> 00744 GetPhotonEnergy(); 00745 00746 currentIN=theSlowLightVector-> 00747 GetProperty(); 00748 00749 currentCII = 0.5 * (prevIN + currentIN); 00750 00751 currentCII = prevCII + 00752 (currentPM - prevPM) * currentCII; 00753 00754 bPhysicsOrderedFreeVector-> 00755 InsertValues(currentPM, currentCII); 00756 00757 prevPM = currentPM; 00758 prevCII = currentCII; 00759 prevIN = currentIN; 00760 } 00761 00762 } 00763 } 00764 00765 G4MaterialPropertyVector* theReemissionVector = 00766 aMaterialPropertiesTable->GetProperty("REEMISSIONPROB"); 00767 00768 if (theReemissionVector) { 00769 00770 // Retrieve the first intensity point in vector 00771 // of (photon energy, intensity) pairs 00772 00773 theReemissionVector->ResetIterator(); 00774 ++(*theReemissionVector); // advance to 1st entry 00775 00776 G4double currentIN = theReemissionVector-> 00777 GetProperty(); 00778 00779 if (currentIN >= 0.0) { 00780 00781 // Create first (photon energy, Scintillation 00782 // Integral pair 00783 00784 G4double currentPM = theReemissionVector-> 00785 GetPhotonEnergy(); 00786 00787 G4double currentCII = 0.0; 00788 00789 cPhysicsOrderedFreeVector-> 00790 InsertValues(currentPM , currentCII); 00791 00792 // Set previous values to current ones prior to loop 00793 00794 G4double prevPM = currentPM; 00795 G4double prevCII = currentCII; 00796 G4double prevIN = currentIN; 00797 00798 // loop over all (photon energy, intensity) 00799 // pairs stored for this material 00800 00801 while(++(*theReemissionVector)) { 00802 currentPM = theReemissionVector-> 00803 GetPhotonEnergy(); 00804 00805 currentIN=theReemissionVector-> 00806 GetProperty(); 00807 00808 currentCII = 0.5 * (prevIN + currentIN); 00809 00810 currentCII = prevCII + 00811 (currentPM - prevPM) * currentCII; 00812 00813 cPhysicsOrderedFreeVector-> 00814 InsertValues(currentPM, currentCII); 00815 00816 prevPM = currentPM; 00817 prevCII = currentCII; 00818 prevIN = currentIN; 00819 } 00820 00821 } 00822 } 00823 00824 } 00825 00826 // The scintillation integral(s) for a given material 00827 // will be inserted in the table(s) according to the 00828 // position of the material in the material table. 00829 00830 theFastIntegralTable->insertAt(i,aPhysicsOrderedFreeVector); 00831 theSlowIntegralTable->insertAt(i,bPhysicsOrderedFreeVector); 00832 theReemissionIntegralTable->insertAt(i,cPhysicsOrderedFreeVector); 00833 } 00834 }
G4PhysicsTable* DsG4Scintillation::theSlowIntegralTable [protected] |
Definition at line 235 of file DsG4Scintillation.h.
G4PhysicsTable* DsG4Scintillation::theFastIntegralTable [protected] |
Definition at line 236 of file DsG4Scintillation.h.
G4PhysicsTable* DsG4Scintillation::theReemissionIntegralTable [protected] |
Definition at line 237 of file DsG4Scintillation.h.
G4bool DsG4Scintillation::doReemission [protected] |
Definition at line 240 of file DsG4Scintillation.h.
G4bool DsG4Scintillation::doBothProcess [protected] |
Definition at line 242 of file DsG4Scintillation.h.
double DsG4Scintillation::birksConstant1 [protected] |
Definition at line 245 of file DsG4Scintillation.h.
double DsG4Scintillation::birksConstant2 [protected] |
Definition at line 246 of file DsG4Scintillation.h.
G4bool DsG4Scintillation::fTrackSecondariesFirst [private] |
Definition at line 250 of file DsG4Scintillation.h.
G4double DsG4Scintillation::YieldFactor [private] |
Definition at line 252 of file DsG4Scintillation.h.
G4bool DsG4Scintillation::FastMu300nsTrick [private] |
Definition at line 253 of file DsG4Scintillation.h.
G4double DsG4Scintillation::ExcitationRatio [private] |
Definition at line 254 of file DsG4Scintillation.h.
G4double DsG4Scintillation::fPhotonWeight [private] |
Definition at line 257 of file DsG4Scintillation.h.
G4bool DsG4Scintillation::fApplyPreQE [private] |
Definition at line 258 of file DsG4Scintillation.h.
G4double DsG4Scintillation::fPreQE [private] |
Definition at line 259 of file DsG4Scintillation.h.
bool DsG4Scintillation::m_noop [private] |
Definition at line 260 of file DsG4Scintillation.h.