#include <G4Scintillation.hh>
Inheritance diagram for G4Scintillation:
Definition at line 88 of file G4Scintillation.hh.
G4Scintillation::G4Scintillation | ( | const G4String & | processName = "Scintillation" , |
|
G4ProcessType | type = fElectromagnetic | |||
) |
Definition at line 96 of file G4Scintillation.cc.
References BuildThePhysicsTable(), ExcitationRatio, fFiniteRiseTime, fScintillation, fTrackSecondariesFirst, G4cout, G4endl, G4VProcess::GetProcessName(), scintillationByParticleType, G4VProcess::SetProcessSubType(), theFastIntegralTable, theSlowIntegralTable, G4VProcess::verboseLevel, and YieldFactor.
00098 : G4VRestDiscreteProcess(processName, type) 00099 { 00100 SetProcessSubType(fScintillation); 00101 00102 fTrackSecondariesFirst = false; 00103 fFiniteRiseTime = false; 00104 00105 YieldFactor = 1.0; 00106 ExcitationRatio = 1.0; 00107 00108 scintillationByParticleType = false; 00109 00110 theFastIntegralTable = NULL; 00111 theSlowIntegralTable = NULL; 00112 00113 if (verboseLevel>0) { 00114 G4cout << GetProcessName() << " is created " << G4endl; 00115 } 00116 00117 BuildThePhysicsTable(); 00118 00119 emSaturation = NULL; 00120 }
G4Scintillation::~G4Scintillation | ( | ) |
Definition at line 126 of file G4Scintillation.cc.
References G4PhysicsTable::clearAndDestroy(), theFastIntegralTable, and theSlowIntegralTable.
00127 { 00128 if (theFastIntegralTable != NULL) { 00129 theFastIntegralTable->clearAndDestroy(); 00130 delete theFastIntegralTable; 00131 } 00132 if (theSlowIntegralTable != NULL) { 00133 theSlowIntegralTable->clearAndDestroy(); 00134 delete theSlowIntegralTable; 00135 } 00136 }
void G4Scintillation::AddSaturation | ( | G4EmSaturation * | sat | ) | [inline] |
Definition at line 183 of file G4Scintillation.hh.
Referenced by G4OpticalPhysics::AddScintillationSaturation(), and G4OpticalPhysics::ConstructProcess().
G4VParticleChange * G4Scintillation::AtRestDoIt | ( | const G4Track & | aTrack, | |
const G4Step & | aStep | |||
) | [virtual] |
Reimplemented from G4VRestDiscreteProcess.
Definition at line 146 of file G4Scintillation.cc.
References PostStepDoIt().
00151 { 00152 return G4Scintillation::PostStepDoIt(aTrack, aStep); 00153 }
void G4Scintillation::BuildThePhysicsTable | ( | ) | [protected] |
Definition at line 537 of file G4Scintillation.cc.
References G4PhysicsVector::Energy(), G4Material::GetMaterialPropertiesTable(), G4Material::GetMaterialTable(), G4Material::GetNumberOfMaterials(), G4MaterialPropertiesTable::GetProperty(), G4PhysicsVector::GetVectorLength(), G4PhysicsTable::insertAt(), theFastIntegralTable, and theSlowIntegralTable.
Referenced by G4Scintillation().
00538 { 00539 if (theFastIntegralTable && theSlowIntegralTable) return; 00540 00541 const G4MaterialTable* theMaterialTable = 00542 G4Material::GetMaterialTable(); 00543 G4int numOfMaterials = G4Material::GetNumberOfMaterials(); 00544 00545 // create new physics table 00546 00547 if(!theFastIntegralTable)theFastIntegralTable = new G4PhysicsTable(numOfMaterials); 00548 if(!theSlowIntegralTable)theSlowIntegralTable = new G4PhysicsTable(numOfMaterials); 00549 00550 // loop for materials 00551 00552 for (G4int i=0 ; i < numOfMaterials; i++) 00553 { 00554 G4PhysicsOrderedFreeVector* aPhysicsOrderedFreeVector = 00555 new G4PhysicsOrderedFreeVector(); 00556 G4PhysicsOrderedFreeVector* bPhysicsOrderedFreeVector = 00557 new G4PhysicsOrderedFreeVector(); 00558 00559 // Retrieve vector of scintillation wavelength intensity for 00560 // the material from the material's optical properties table. 00561 00562 G4Material* aMaterial = (*theMaterialTable)[i]; 00563 00564 G4MaterialPropertiesTable* aMaterialPropertiesTable = 00565 aMaterial->GetMaterialPropertiesTable(); 00566 00567 if (aMaterialPropertiesTable) { 00568 00569 G4MaterialPropertyVector* theFastLightVector = 00570 aMaterialPropertiesTable->GetProperty("FASTCOMPONENT"); 00571 00572 if (theFastLightVector) { 00573 00574 // Retrieve the first intensity point in vector 00575 // of (photon energy, intensity) pairs 00576 00577 G4double currentIN = (*theFastLightVector)[0]; 00578 00579 if (currentIN >= 0.0) { 00580 00581 // Create first (photon energy, Scintillation 00582 // Integral pair 00583 00584 G4double currentPM = theFastLightVector->Energy(0); 00585 00586 G4double currentCII = 0.0; 00587 00588 aPhysicsOrderedFreeVector-> 00589 InsertValues(currentPM , currentCII); 00590 00591 // Set previous values to current ones prior to loop 00592 00593 G4double prevPM = currentPM; 00594 G4double prevCII = currentCII; 00595 G4double prevIN = currentIN; 00596 00597 // loop over all (photon energy, intensity) 00598 // pairs stored for this material 00599 00600 for (size_t ii = 1; 00601 ii < theFastLightVector->GetVectorLength(); 00602 ++ii) 00603 { 00604 currentPM = theFastLightVector->Energy(ii); 00605 currentIN = (*theFastLightVector)[ii]; 00606 00607 currentCII = 0.5 * (prevIN + currentIN); 00608 00609 currentCII = prevCII + 00610 (currentPM - prevPM) * currentCII; 00611 00612 aPhysicsOrderedFreeVector-> 00613 InsertValues(currentPM, currentCII); 00614 00615 prevPM = currentPM; 00616 prevCII = currentCII; 00617 prevIN = currentIN; 00618 } 00619 00620 } 00621 } 00622 00623 G4MaterialPropertyVector* theSlowLightVector = 00624 aMaterialPropertiesTable->GetProperty("SLOWCOMPONENT"); 00625 00626 if (theSlowLightVector) { 00627 00628 // Retrieve the first intensity point in vector 00629 // of (photon energy, intensity) pairs 00630 00631 G4double currentIN = (*theSlowLightVector)[0]; 00632 00633 if (currentIN >= 0.0) { 00634 00635 // Create first (photon energy, Scintillation 00636 // Integral pair 00637 00638 G4double currentPM = theSlowLightVector->Energy(0); 00639 00640 G4double currentCII = 0.0; 00641 00642 bPhysicsOrderedFreeVector-> 00643 InsertValues(currentPM , currentCII); 00644 00645 // Set previous values to current ones prior to loop 00646 00647 G4double prevPM = currentPM; 00648 G4double prevCII = currentCII; 00649 G4double prevIN = currentIN; 00650 00651 // loop over all (photon energy, intensity) 00652 // pairs stored for this material 00653 00654 for (size_t ii = 1; 00655 ii < theSlowLightVector->GetVectorLength(); 00656 ++ii) 00657 { 00658 currentPM = theSlowLightVector->Energy(ii); 00659 currentIN = (*theSlowLightVector)[ii]; 00660 00661 currentCII = 0.5 * (prevIN + currentIN); 00662 00663 currentCII = prevCII + 00664 (currentPM - prevPM) * currentCII; 00665 00666 bPhysicsOrderedFreeVector-> 00667 InsertValues(currentPM, currentCII); 00668 00669 prevPM = currentPM; 00670 prevCII = currentCII; 00671 prevIN = currentIN; 00672 } 00673 00674 } 00675 } 00676 } 00677 00678 // The scintillation integral(s) for a given material 00679 // will be inserted in the table(s) according to the 00680 // position of the material in the material table. 00681 00682 theFastIntegralTable->insertAt(i,aPhysicsOrderedFreeVector); 00683 theSlowIntegralTable->insertAt(i,bPhysicsOrderedFreeVector); 00684 00685 } 00686 }
void G4Scintillation::DumpPhysicsTable | ( | ) | const [inline] |
Definition at line 312 of file G4Scintillation.hh.
References G4PhysicsOrderedFreeVector::DumpValues(), G4PhysicsTable::entries(), theFastIntegralTable, and theSlowIntegralTable.
00313 { 00314 if (theFastIntegralTable) { 00315 G4int PhysicsTableSize = theFastIntegralTable->entries(); 00316 G4PhysicsOrderedFreeVector *v; 00317 00318 for (G4int i = 0 ; i < PhysicsTableSize ; i++ ) 00319 { 00320 v = (G4PhysicsOrderedFreeVector*)(*theFastIntegralTable)[i]; 00321 v->DumpValues(); 00322 } 00323 } 00324 00325 if (theSlowIntegralTable) { 00326 G4int PhysicsTableSize = theSlowIntegralTable->entries(); 00327 G4PhysicsOrderedFreeVector *v; 00328 00329 for (G4int i = 0 ; i < PhysicsTableSize ; i++ ) 00330 { 00331 v = (G4PhysicsOrderedFreeVector*)(*theSlowIntegralTable)[i]; 00332 v->DumpValues(); 00333 } 00334 } 00335 }
G4PhysicsTable * G4Scintillation::GetFastIntegralTable | ( | ) | const [inline] |
Definition at line 306 of file G4Scintillation.hh.
References theFastIntegralTable.
00307 { 00308 return theFastIntegralTable; 00309 }
G4bool G4Scintillation::GetFiniteRiseTime | ( | ) | const [inline] |
Definition at line 270 of file G4Scintillation.hh.
References fFiniteRiseTime.
00271 { 00272 return fFiniteRiseTime; 00273 }
G4double G4Scintillation::GetMeanFreePath | ( | const G4Track & | aTrack, | |
G4double | , | |||
G4ForceCondition * | ||||
) | [virtual] |
Implements G4VRestDiscreteProcess.
Definition at line 705 of file G4Scintillation.cc.
References DBL_MAX, and StronglyForced.
00708 { 00709 *condition = StronglyForced; 00710 00711 return DBL_MAX; 00712 00713 }
G4double G4Scintillation::GetMeanLifeTime | ( | const G4Track & | aTrack, | |
G4ForceCondition * | ||||
) | [virtual] |
Implements G4VRestDiscreteProcess.
Definition at line 719 of file G4Scintillation.cc.
References DBL_MAX, and Forced.
G4EmSaturation* G4Scintillation::GetSaturation | ( | ) | const [inline] |
G4bool G4Scintillation::GetScintillationByParticleType | ( | ) | const [inline] |
Definition at line 196 of file G4Scintillation.hh.
References scintillationByParticleType.
00197 { return scintillationByParticleType; }
G4double G4Scintillation::GetScintillationExcitationRatio | ( | ) | const [inline] |
Definition at line 294 of file G4Scintillation.hh.
References ExcitationRatio.
00295 { 00296 return ExcitationRatio; 00297 }
G4double G4Scintillation::GetScintillationYieldFactor | ( | ) | const [inline] |
Definition at line 282 of file G4Scintillation.hh.
References YieldFactor.
00283 { 00284 return YieldFactor; 00285 }
G4PhysicsTable * G4Scintillation::GetSlowIntegralTable | ( | ) | const [inline] |
Definition at line 300 of file G4Scintillation.hh.
References theSlowIntegralTable.
00301 { 00302 return theSlowIntegralTable; 00303 }
G4bool G4Scintillation::GetTrackSecondariesFirst | ( | ) | const [inline] |
Definition at line 264 of file G4Scintillation.hh.
References fTrackSecondariesFirst.
00265 { 00266 return fTrackSecondariesFirst; 00267 }
G4bool G4Scintillation::IsApplicable | ( | const G4ParticleDefinition & | aParticleType | ) | [inline, virtual] |
Reimplemented from G4VProcess.
Definition at line 243 of file G4Scintillation.hh.
References G4ParticleDefinition::GetParticleName(), and G4ParticleDefinition::IsShortLived().
Referenced by TLBE< T >::ConstructOp(), and G4OpticalPhysics::ConstructProcess().
00244 { 00245 if (aParticleType.GetParticleName() == "opticalphoton") return false; 00246 if (aParticleType.IsShortLived()) return false; 00247 00248 return true; 00249 }
G4VParticleChange * G4Scintillation::PostStepDoIt | ( | const G4Track & | aTrack, | |
const G4Step & | aStep | |||
) | [virtual] |
Reimplemented from G4VRestDiscreteProcess.
Definition at line 159 of file G4Scintillation.cc.
References G4ParticleChange::AddSecondary(), G4Alpha::AlphaDefinition(), G4VProcess::aParticleChange, G4Deuteron::DeuteronDefinition(), G4Electron::ElectronDefinition(), ExcitationRatio, fAlive, FatalException, fFiniteRiseTime, fSuspend, fTrackSecondariesFirst, G4cout, G4endl, G4Exception(), G4Poisson(), G4UniformRand, G4Gamma::GammaDefinition(), G4DynamicParticle::GetDefinition(), G4Step::GetDeltaPosition(), G4Track::GetDynamicParticle(), G4PhysicsOrderedFreeVector::GetEnergy(), G4StepPoint::GetGlobalTime(), G4Material::GetIndex(), G4Track::GetMaterial(), G4Material::GetMaterialPropertiesTable(), G4PhysicsOrderedFreeVector::GetMaxValue(), G4VParticleChange::GetNumberOfSecondaries(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetParticleType(), G4ParticleDefinition::GetPDGCharge(), G4StepPoint::GetPosition(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4MaterialPropertiesTable::GetProperty(), G4Step::GetStepLength(), G4Step::GetTotalEnergyDeposit(), G4StepPoint::GetTouchableHandle(), G4StepPoint::GetVelocity(), G4ParticleChange::Initialize(), G4Neutron::NeutronDefinition(), ns, G4OpticalPhoton::OpticalPhoton(), G4VRestDiscreteProcess::PostStepDoIt(), G4VParticleChange::ProposeTrackStatus(), G4Proton::ProtonDefinition(), scintillationByParticleType, G4DynamicParticle::SetKineticEnergy(), G4VParticleChange::SetNumberOfSecondaries(), G4Track::SetParentID(), G4DynamicParticle::SetPolarization(), G4Track::SetTouchableHandle(), theFastIntegralTable, theSlowIntegralTable, G4Triton::TritonDefinition(), G4PhysicsVector::Value(), G4VProcess::verboseLevel, G4EmSaturation::VisibleEnergyDeposition(), and YieldFactor.
Referenced by AtRestDoIt().
00166 { 00167 aParticleChange.Initialize(aTrack); 00168 00169 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle(); 00170 const G4Material* aMaterial = aTrack.GetMaterial(); 00171 00172 G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint(); 00173 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint(); 00174 00175 G4ThreeVector x0 = pPreStepPoint->GetPosition(); 00176 G4ThreeVector p0 = aStep.GetDeltaPosition().unit(); 00177 G4double t0 = pPreStepPoint->GetGlobalTime(); 00178 00179 G4double TotalEnergyDeposit = aStep.GetTotalEnergyDeposit(); 00180 00181 G4MaterialPropertiesTable* aMaterialPropertiesTable = 00182 aMaterial->GetMaterialPropertiesTable(); 00183 if (!aMaterialPropertiesTable) 00184 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00185 00186 G4MaterialPropertyVector* Fast_Intensity = 00187 aMaterialPropertiesTable->GetProperty("FASTCOMPONENT"); 00188 G4MaterialPropertyVector* Slow_Intensity = 00189 aMaterialPropertiesTable->GetProperty("SLOWCOMPONENT"); 00190 00191 if (!Fast_Intensity && !Slow_Intensity ) 00192 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00193 00194 G4int nscnt = 1; 00195 if (Fast_Intensity && Slow_Intensity) nscnt = 2; 00196 00197 G4double ScintillationYield = 0.; 00198 00199 if (scintillationByParticleType) { 00200 // The scintillation response is a function of the energy 00201 // deposited by particle types. 00202 00203 // Get the definition of the current particle 00204 G4ParticleDefinition *pDef = aParticle->GetDefinition(); 00205 G4MaterialPropertyVector *Scint_Yield_Vector = NULL; 00206 00207 // Obtain the G4MaterialPropertyVectory containing the 00208 // scintillation light yield as a function of the deposited 00209 // energy for the current particle type 00210 00211 // Protons 00212 if(pDef==G4Proton::ProtonDefinition()) 00213 Scint_Yield_Vector = aMaterialPropertiesTable-> 00214 GetProperty("PROTONSCINTILLATIONYIELD"); 00215 00216 // Deuterons 00217 else if(pDef==G4Deuteron::DeuteronDefinition()) 00218 Scint_Yield_Vector = aMaterialPropertiesTable-> 00219 GetProperty("DEUTERONSCINTILLATIONYIELD"); 00220 00221 // Tritons 00222 else if(pDef==G4Triton::TritonDefinition()) 00223 Scint_Yield_Vector = aMaterialPropertiesTable-> 00224 GetProperty("TRITONSCINTILLATIONYIELD"); 00225 00226 // Alphas 00227 else if(pDef==G4Alpha::AlphaDefinition()) 00228 Scint_Yield_Vector = aMaterialPropertiesTable-> 00229 GetProperty("ALPHASCINTILLATIONYIELD"); 00230 00231 // Ions (particles derived from G4VIon and G4Ions) 00232 // and recoil ions below tracking cut from neutrons after hElastic 00233 else if(pDef->GetParticleType()== "nucleus" || 00234 pDef==G4Neutron::NeutronDefinition()) 00235 Scint_Yield_Vector = aMaterialPropertiesTable-> 00236 GetProperty("IONSCINTILLATIONYIELD"); 00237 00238 // Electrons (must also account for shell-binding energy 00239 // attributed to gamma from standard PhotoElectricEffect) 00240 else if(pDef==G4Electron::ElectronDefinition() || 00241 pDef==G4Gamma::GammaDefinition()) 00242 Scint_Yield_Vector = aMaterialPropertiesTable-> 00243 GetProperty("ELECTRONSCINTILLATIONYIELD"); 00244 00245 // Default for particles not enumerated/listed above 00246 else 00247 Scint_Yield_Vector = aMaterialPropertiesTable-> 00248 GetProperty("ELECTRONSCINTILLATIONYIELD"); 00249 00250 // If the user has not specified yields for (p,d,t,a,carbon) 00251 // then these unspecified particles will default to the 00252 // electron's scintillation yield 00253 if(!Scint_Yield_Vector){ 00254 Scint_Yield_Vector = aMaterialPropertiesTable-> 00255 GetProperty("ELECTRONSCINTILLATIONYIELD"); 00256 } 00257 00258 // Throw an exception if no scintillation yield is found 00259 if (!Scint_Yield_Vector) { 00260 G4ExceptionDescription ed; 00261 ed << "\nG4Scintillation::PostStepDoIt(): " 00262 << "Request for scintillation yield for energy deposit and particle type without correct entry in MaterialPropertiesTable\n" 00263 << "ScintillationByParticleType requires at minimum that ELECTRONSCINTILLATIONYIELD is set by the user\n" 00264 << G4endl; 00265 G4String comments = "Missing MaterialPropertiesTable entry - No correct entry in MaterialPropertiesTable"; 00266 G4Exception("G4Scintillation::PostStepDoIt","Scint01", 00267 FatalException,ed,comments); 00268 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00269 } 00270 00271 if (verboseLevel>1) { 00272 G4cout << "\n" 00273 << "Particle = " << pDef->GetParticleName() << "\n" 00274 << "Energy Dep. = " << TotalEnergyDeposit/MeV << "\n" 00275 << "Yield = " 00276 << Scint_Yield_Vector->Value(TotalEnergyDeposit) 00277 << "\n" << G4endl; 00278 } 00279 00280 // Obtain the scintillation yield using the total energy 00281 // deposited by the particle in this step. 00282 00283 // Units: [# scintillation photons] 00284 ScintillationYield = Scint_Yield_Vector-> 00285 Value(TotalEnergyDeposit); 00286 } else { 00287 // The default linear scintillation process 00288 ScintillationYield = aMaterialPropertiesTable-> 00289 GetConstProperty("SCINTILLATIONYIELD"); 00290 00291 // Units: [# scintillation photons / MeV] 00292 ScintillationYield *= YieldFactor; 00293 } 00294 00295 G4double ResolutionScale = aMaterialPropertiesTable-> 00296 GetConstProperty("RESOLUTIONSCALE"); 00297 00298 // Birks law saturation: 00299 00300 //G4double constBirks = 0.0; 00301 00302 //constBirks = aMaterial->GetIonisation()->GetBirksConstant(); 00303 00304 G4double MeanNumberOfPhotons; 00305 00306 // Birk's correction via emSaturation and specifying scintillation by 00307 // by particle type are physically mutually exclusive 00308 00309 if (scintillationByParticleType) 00310 MeanNumberOfPhotons = ScintillationYield; 00311 else if (emSaturation) 00312 MeanNumberOfPhotons = ScintillationYield* 00313 (emSaturation->VisibleEnergyDeposition(&aStep)); 00314 else 00315 MeanNumberOfPhotons = ScintillationYield*TotalEnergyDeposit; 00316 00317 G4int NumPhotons; 00318 00319 if (MeanNumberOfPhotons > 10.) 00320 { 00321 G4double sigma = ResolutionScale * std::sqrt(MeanNumberOfPhotons); 00322 NumPhotons = G4int(G4RandGauss::shoot(MeanNumberOfPhotons,sigma)+0.5); 00323 } 00324 else 00325 { 00326 NumPhotons = G4int(G4Poisson(MeanNumberOfPhotons)); 00327 } 00328 00329 if (NumPhotons <= 0) 00330 { 00331 // return unchanged particle and no secondaries 00332 00333 aParticleChange.SetNumberOfSecondaries(0); 00334 00335 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00336 } 00337 00339 00340 aParticleChange.SetNumberOfSecondaries(NumPhotons); 00341 00342 if (fTrackSecondariesFirst) { 00343 if (aTrack.GetTrackStatus() == fAlive ) 00344 aParticleChange.ProposeTrackStatus(fSuspend); 00345 } 00346 00348 00349 G4int materialIndex = aMaterial->GetIndex(); 00350 00351 // Retrieve the Scintillation Integral for this material 00352 // new G4PhysicsOrderedFreeVector allocated to hold CII's 00353 00354 G4int Num = NumPhotons; 00355 00356 for (G4int scnt = 1; scnt <= nscnt; scnt++) { 00357 00358 G4double ScintillationTime = 0.*ns; 00359 G4double ScintillationRiseTime = 0.*ns; 00360 G4PhysicsOrderedFreeVector* ScintillationIntegral = NULL; 00361 00362 if (scnt == 1) { 00363 if (nscnt == 1) { 00364 if(Fast_Intensity){ 00365 ScintillationTime = aMaterialPropertiesTable-> 00366 GetConstProperty("FASTTIMECONSTANT"); 00367 if (fFiniteRiseTime) { 00368 ScintillationRiseTime = aMaterialPropertiesTable-> 00369 GetConstProperty("FASTSCINTILLATIONRISETIME"); 00370 } 00371 ScintillationIntegral = 00372 (G4PhysicsOrderedFreeVector*)((*theFastIntegralTable)(materialIndex)); 00373 } 00374 if(Slow_Intensity){ 00375 ScintillationTime = aMaterialPropertiesTable-> 00376 GetConstProperty("SLOWTIMECONSTANT"); 00377 if (fFiniteRiseTime) { 00378 ScintillationRiseTime = aMaterialPropertiesTable-> 00379 GetConstProperty("SLOWSCINTILLATIONRISETIME"); 00380 } 00381 ScintillationIntegral = 00382 (G4PhysicsOrderedFreeVector*)((*theSlowIntegralTable)(materialIndex)); 00383 } 00384 } 00385 else { 00386 G4double YieldRatio = aMaterialPropertiesTable-> 00387 GetConstProperty("YIELDRATIO"); 00388 if ( ExcitationRatio == 1.0 ) { 00389 Num = G4int (std::min(YieldRatio,1.0) * NumPhotons); 00390 } 00391 else { 00392 Num = G4int (std::min(ExcitationRatio,1.0) * NumPhotons); 00393 } 00394 ScintillationTime = aMaterialPropertiesTable-> 00395 GetConstProperty("FASTTIMECONSTANT"); 00396 if (fFiniteRiseTime) { 00397 ScintillationRiseTime = aMaterialPropertiesTable-> 00398 GetConstProperty("FASTSCINTILLATIONRISETIME"); 00399 } 00400 ScintillationIntegral = 00401 (G4PhysicsOrderedFreeVector*)((*theFastIntegralTable)(materialIndex)); 00402 } 00403 } 00404 else { 00405 Num = NumPhotons - Num; 00406 ScintillationTime = aMaterialPropertiesTable-> 00407 GetConstProperty("SLOWTIMECONSTANT"); 00408 if (fFiniteRiseTime) { 00409 ScintillationRiseTime = aMaterialPropertiesTable-> 00410 GetConstProperty("SLOWSCINTILLATIONRISETIME"); 00411 } 00412 ScintillationIntegral = 00413 (G4PhysicsOrderedFreeVector*)((*theSlowIntegralTable)(materialIndex)); 00414 } 00415 00416 if (!ScintillationIntegral) continue; 00417 00418 // Max Scintillation Integral 00419 00420 G4double CIImax = ScintillationIntegral->GetMaxValue(); 00421 00422 for (G4int i = 0; i < Num; i++) { 00423 00424 // Determine photon energy 00425 00426 G4double CIIvalue = G4UniformRand()*CIImax; 00427 G4double sampledEnergy = 00428 ScintillationIntegral->GetEnergy(CIIvalue); 00429 00430 if (verboseLevel>1) { 00431 G4cout << "sampledEnergy = " << sampledEnergy << G4endl; 00432 G4cout << "CIIvalue = " << CIIvalue << G4endl; 00433 } 00434 00435 // Generate random photon direction 00436 00437 G4double cost = 1. - 2.*G4UniformRand(); 00438 G4double sint = std::sqrt((1.-cost)*(1.+cost)); 00439 00440 G4double phi = twopi*G4UniformRand(); 00441 G4double sinp = std::sin(phi); 00442 G4double cosp = std::cos(phi); 00443 00444 G4double px = sint*cosp; 00445 G4double py = sint*sinp; 00446 G4double pz = cost; 00447 00448 // Create photon momentum direction vector 00449 00450 G4ParticleMomentum photonMomentum(px, py, pz); 00451 00452 // Determine polarization of new photon 00453 00454 G4double sx = cost*cosp; 00455 G4double sy = cost*sinp; 00456 G4double sz = -sint; 00457 00458 G4ThreeVector photonPolarization(sx, sy, sz); 00459 00460 G4ThreeVector perp = photonMomentum.cross(photonPolarization); 00461 00462 phi = twopi*G4UniformRand(); 00463 sinp = std::sin(phi); 00464 cosp = std::cos(phi); 00465 00466 photonPolarization = cosp * photonPolarization + sinp * perp; 00467 00468 photonPolarization = photonPolarization.unit(); 00469 00470 // Generate a new photon: 00471 00472 G4DynamicParticle* aScintillationPhoton = 00473 new G4DynamicParticle(G4OpticalPhoton::OpticalPhoton(), 00474 photonMomentum); 00475 aScintillationPhoton->SetPolarization 00476 (photonPolarization.x(), 00477 photonPolarization.y(), 00478 photonPolarization.z()); 00479 00480 aScintillationPhoton->SetKineticEnergy(sampledEnergy); 00481 00482 // Generate new G4Track object: 00483 00484 G4double rand; 00485 00486 if (aParticle->GetDefinition()->GetPDGCharge() != 0) { 00487 rand = G4UniformRand(); 00488 } else { 00489 rand = 1.0; 00490 } 00491 00492 G4double delta = rand * aStep.GetStepLength(); 00493 G4double deltaTime = delta / 00494 ((pPreStepPoint->GetVelocity()+ 00495 pPostStepPoint->GetVelocity())/2.); 00496 00497 // emission time distribution 00498 if (ScintillationRiseTime==0.0) { 00499 deltaTime = deltaTime - 00500 ScintillationTime * std::log( G4UniformRand() ); 00501 } else { 00502 deltaTime = deltaTime + 00503 sample_time(ScintillationRiseTime, ScintillationTime); 00504 } 00505 00506 G4double aSecondaryTime = t0 + deltaTime; 00507 00508 G4ThreeVector aSecondaryPosition = 00509 x0 + rand * aStep.GetDeltaPosition(); 00510 00511 G4Track* aSecondaryTrack = 00512 new G4Track(aScintillationPhoton,aSecondaryTime,aSecondaryPosition); 00513 00514 aSecondaryTrack->SetTouchableHandle( 00515 aStep.GetPreStepPoint()->GetTouchableHandle()); 00516 // aSecondaryTrack->SetTouchableHandle((G4VTouchable*)0); 00517 00518 aSecondaryTrack->SetParentID(aTrack.GetTrackID()); 00519 00520 aParticleChange.AddSecondary(aSecondaryTrack); 00521 00522 } 00523 } 00524 00525 if (verboseLevel>0) { 00526 G4cout << "\n Exiting from G4Scintillation::DoIt -- NumberOfSecondaries = " 00527 << aParticleChange.GetNumberOfSecondaries() << G4endl; 00528 } 00529 00530 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep); 00531 }
void G4Scintillation::RemoveSaturation | ( | ) | [inline] |
void G4Scintillation::SetFiniteRiseTime | ( | const G4bool | state | ) | [inline] |
Definition at line 258 of file G4Scintillation.hh.
References fFiniteRiseTime.
Referenced by G4OpticalPhysics::ConstructProcess(), and G4OpticalPhysics::SetFiniteRiseTime().
00259 { 00260 fFiniteRiseTime = state; 00261 }
void G4Scintillation::SetScintillationByParticleType | ( | const | G4bool | ) |
Definition at line 691 of file G4Scintillation.cc.
References G4Exception(), JustWarning, RemoveSaturation(), and scintillationByParticleType.
Referenced by G4OpticalPhysics::ConstructProcess(), and G4OpticalPhysics::SetScintillationByParticleType().
00692 { 00693 if (emSaturation) { 00694 G4Exception("G4Scintillation::SetScintillationByParticleType", "Scint02", 00695 JustWarning, "Redefinition: Birks Saturation is replaced by ScintillationByParticleType!"); 00696 RemoveSaturation(); 00697 } 00698 scintillationByParticleType = scintType; 00699 }
void G4Scintillation::SetScintillationExcitationRatio | ( | const G4double | excitationratio | ) | [inline] |
Definition at line 288 of file G4Scintillation.hh.
References ExcitationRatio.
Referenced by TLBE< T >::ConstructOp(), G4OpticalPhysics::ConstructProcess(), and G4OpticalPhysics::SetScintillationExcitationRatio().
00289 { 00290 ExcitationRatio = excitationratio; 00291 }
void G4Scintillation::SetScintillationYieldFactor | ( | const G4double | yieldfactor | ) | [inline] |
Definition at line 276 of file G4Scintillation.hh.
References YieldFactor.
Referenced by TLBE< T >::ConstructOp(), G4OpticalPhysics::ConstructProcess(), and G4OpticalPhysics::SetScintillationYieldFactor().
00277 { 00278 YieldFactor = yieldfactor; 00279 }
void G4Scintillation::SetTrackSecondariesFirst | ( | const G4bool | state | ) | [inline] |
Definition at line 252 of file G4Scintillation.hh.
References fTrackSecondariesFirst.
Referenced by TLBE< T >::ConstructOp(), and G4OpticalPhysics::SetTrackSecondariesFirst().
00253 { 00254 fTrackSecondariesFirst = state; 00255 }
G4double G4Scintillation::ExcitationRatio [protected] |
Definition at line 222 of file G4Scintillation.hh.
Referenced by G4Scintillation(), GetScintillationExcitationRatio(), PostStepDoIt(), and SetScintillationExcitationRatio().
G4bool G4Scintillation::fFiniteRiseTime [protected] |
Definition at line 218 of file G4Scintillation.hh.
Referenced by G4Scintillation(), GetFiniteRiseTime(), PostStepDoIt(), and SetFiniteRiseTime().
G4bool G4Scintillation::fTrackSecondariesFirst [protected] |
Definition at line 217 of file G4Scintillation.hh.
Referenced by G4Scintillation(), GetTrackSecondariesFirst(), PostStepDoIt(), and SetTrackSecondariesFirst().
G4bool G4Scintillation::scintillationByParticleType [protected] |
Definition at line 224 of file G4Scintillation.hh.
Referenced by G4Scintillation(), GetScintillationByParticleType(), PostStepDoIt(), and SetScintillationByParticleType().
G4PhysicsTable* G4Scintillation::theFastIntegralTable [protected] |
Definition at line 215 of file G4Scintillation.hh.
Referenced by BuildThePhysicsTable(), DumpPhysicsTable(), G4Scintillation(), GetFastIntegralTable(), PostStepDoIt(), and ~G4Scintillation().
G4PhysicsTable* G4Scintillation::theSlowIntegralTable [protected] |
Definition at line 214 of file G4Scintillation.hh.
Referenced by BuildThePhysicsTable(), DumpPhysicsTable(), G4Scintillation(), GetSlowIntegralTable(), PostStepDoIt(), and ~G4Scintillation().
G4double G4Scintillation::YieldFactor [protected] |
Definition at line 220 of file G4Scintillation.hh.
Referenced by G4Scintillation(), GetScintillationYieldFactor(), PostStepDoIt(), and SetScintillationYieldFactor().