00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00030
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00070
00071 #include "G4ios.hh"
00072 #include "globals.hh"
00073 #include "G4PhysicalConstants.hh"
00074 #include "G4SystemOfUnits.hh"
00075 #include "G4ParticleTypes.hh"
00076 #include "G4EmProcessSubType.hh"
00077
00078 #include "G4Scintillation.hh"
00079
00081
00083
00085
00087
00088
00089
00090
00091
00093
00095
00096 G4Scintillation::G4Scintillation(const G4String& processName,
00097 G4ProcessType type)
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 }
00121
00123
00125
00126 G4Scintillation::~G4Scintillation()
00127 {
00128 if (theFastIntegralTable != NULL) {
00129 theFastIntegralTable->clearAndDestroy();
00130 delete theFastIntegralTable;
00131 }
00132 if (theSlowIntegralTable != NULL) {
00133 theSlowIntegralTable->clearAndDestroy();
00134 delete theSlowIntegralTable;
00135 }
00136 }
00137
00139
00141
00142
00143
00144
00145 G4VParticleChange*
00146 G4Scintillation::AtRestDoIt(const G4Track& aTrack, const G4Step& aStep)
00147
00148
00149
00150
00151 {
00152 return G4Scintillation::PostStepDoIt(aTrack, aStep);
00153 }
00154
00155
00156
00157
00158 G4VParticleChange*
00159 G4Scintillation::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
00160
00161
00162
00163
00164
00165
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
00201
00202
00203
00204 G4ParticleDefinition *pDef = aParticle->GetDefinition();
00205 G4MaterialPropertyVector *Scint_Yield_Vector = NULL;
00206
00207
00208
00209
00210
00211
00212 if(pDef==G4Proton::ProtonDefinition())
00213 Scint_Yield_Vector = aMaterialPropertiesTable->
00214 GetProperty("PROTONSCINTILLATIONYIELD");
00215
00216
00217 else if(pDef==G4Deuteron::DeuteronDefinition())
00218 Scint_Yield_Vector = aMaterialPropertiesTable->
00219 GetProperty("DEUTERONSCINTILLATIONYIELD");
00220
00221
00222 else if(pDef==G4Triton::TritonDefinition())
00223 Scint_Yield_Vector = aMaterialPropertiesTable->
00224 GetProperty("TRITONSCINTILLATIONYIELD");
00225
00226
00227 else if(pDef==G4Alpha::AlphaDefinition())
00228 Scint_Yield_Vector = aMaterialPropertiesTable->
00229 GetProperty("ALPHASCINTILLATIONYIELD");
00230
00231
00232
00233 else if(pDef->GetParticleType()== "nucleus" ||
00234 pDef==G4Neutron::NeutronDefinition())
00235 Scint_Yield_Vector = aMaterialPropertiesTable->
00236 GetProperty("IONSCINTILLATIONYIELD");
00237
00238
00239
00240 else if(pDef==G4Electron::ElectronDefinition() ||
00241 pDef==G4Gamma::GammaDefinition())
00242 Scint_Yield_Vector = aMaterialPropertiesTable->
00243 GetProperty("ELECTRONSCINTILLATIONYIELD");
00244
00245
00246 else
00247 Scint_Yield_Vector = aMaterialPropertiesTable->
00248 GetProperty("ELECTRONSCINTILLATIONYIELD");
00249
00250
00251
00252
00253 if(!Scint_Yield_Vector){
00254 Scint_Yield_Vector = aMaterialPropertiesTable->
00255 GetProperty("ELECTRONSCINTILLATIONYIELD");
00256 }
00257
00258
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
00281
00282
00283
00284 ScintillationYield = Scint_Yield_Vector->
00285 Value(TotalEnergyDeposit);
00286 } else {
00287
00288 ScintillationYield = aMaterialPropertiesTable->
00289 GetConstProperty("SCINTILLATIONYIELD");
00290
00291
00292 ScintillationYield *= YieldFactor;
00293 }
00294
00295 G4double ResolutionScale = aMaterialPropertiesTable->
00296 GetConstProperty("RESOLUTIONSCALE");
00297
00298
00299
00300
00301
00302
00303
00304 G4double MeanNumberOfPhotons;
00305
00306
00307
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
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
00352
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
00419
00420 G4double CIImax = ScintillationIntegral->GetMaxValue();
00421
00422 for (G4int i = 0; i < Num; i++) {
00423
00424
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
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
00449
00450 G4ParticleMomentum photonMomentum(px, py, pz);
00451
00452
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
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
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
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
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 }
00532
00533
00534
00535
00536
00537 void G4Scintillation::BuildThePhysicsTable()
00538 {
00539 if (theFastIntegralTable && theSlowIntegralTable) return;
00540
00541 const G4MaterialTable* theMaterialTable =
00542 G4Material::GetMaterialTable();
00543 G4int numOfMaterials = G4Material::GetNumberOfMaterials();
00544
00545
00546
00547 if(!theFastIntegralTable)theFastIntegralTable = new G4PhysicsTable(numOfMaterials);
00548 if(!theSlowIntegralTable)theSlowIntegralTable = new G4PhysicsTable(numOfMaterials);
00549
00550
00551
00552 for (G4int i=0 ; i < numOfMaterials; i++)
00553 {
00554 G4PhysicsOrderedFreeVector* aPhysicsOrderedFreeVector =
00555 new G4PhysicsOrderedFreeVector();
00556 G4PhysicsOrderedFreeVector* bPhysicsOrderedFreeVector =
00557 new G4PhysicsOrderedFreeVector();
00558
00559
00560
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
00575
00576
00577 G4double currentIN = (*theFastLightVector)[0];
00578
00579 if (currentIN >= 0.0) {
00580
00581
00582
00583
00584 G4double currentPM = theFastLightVector->Energy(0);
00585
00586 G4double currentCII = 0.0;
00587
00588 aPhysicsOrderedFreeVector->
00589 InsertValues(currentPM , currentCII);
00590
00591
00592
00593 G4double prevPM = currentPM;
00594 G4double prevCII = currentCII;
00595 G4double prevIN = currentIN;
00596
00597
00598
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
00629
00630
00631 G4double currentIN = (*theSlowLightVector)[0];
00632
00633 if (currentIN >= 0.0) {
00634
00635
00636
00637
00638 G4double currentPM = theSlowLightVector->Energy(0);
00639
00640 G4double currentCII = 0.0;
00641
00642 bPhysicsOrderedFreeVector->
00643 InsertValues(currentPM , currentCII);
00644
00645
00646
00647 G4double prevPM = currentPM;
00648 G4double prevCII = currentCII;
00649 G4double prevIN = currentIN;
00650
00651
00652
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
00679
00680
00681
00682 theFastIntegralTable->insertAt(i,aPhysicsOrderedFreeVector);
00683 theSlowIntegralTable->insertAt(i,bPhysicsOrderedFreeVector);
00684
00685 }
00686 }
00687
00688
00689
00690
00691 void G4Scintillation::SetScintillationByParticleType(const G4bool scintType)
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 }
00700
00701
00702
00703
00704
00705 G4double G4Scintillation::GetMeanFreePath(const G4Track&,
00706 G4double ,
00707 G4ForceCondition* condition)
00708 {
00709 *condition = StronglyForced;
00710
00711 return DBL_MAX;
00712
00713 }
00714
00715
00716
00717
00718
00719 G4double G4Scintillation::GetMeanLifeTime(const G4Track&,
00720 G4ForceCondition* condition)
00721 {
00722 *condition = Forced;
00723
00724 return DBL_MAX;
00725
00726 }
00727
00728 G4double G4Scintillation::sample_time(G4double tau1, G4double tau2)
00729 {
00730
00731
00732 while(1) {
00733
00734 G4double ran1 = G4UniformRand();
00735 G4double ran2 = G4UniformRand();
00736
00737
00738
00739 G4double d = (tau1+tau2)/tau2;
00740
00741
00742 G4double t = -1.0*tau2*std::log(1-ran1);
00743 G4double gg = d*single_exp(t,tau2);
00744 if (ran2 <= bi_exp(t,tau1,tau2)/gg) return t;
00745 }
00746 return -1.0;
00747 }