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
00029
00030
00031
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
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068 #include "DsG4Scintillation.h"
00069 #include "G4UnitsTable.hh"
00070 #include "G4LossTableManager.hh"
00071 #include "G4MaterialCutsCouple.hh"
00072 #include "G4Gamma.hh"
00073 #include "G4Electron.hh"
00074 #include "globals.hh"
00075
00076
00077
00079
00080 using namespace std;
00081
00083
00085
00087
00089
00090
00091
00092
00093
00095
00097
00098 DsG4Scintillation::DsG4Scintillation(const G4String& processName,
00099 G4ProcessType type)
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 }
00125
00127
00129
00130 DsG4Scintillation::~DsG4Scintillation()
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 }
00145
00147
00149
00150
00151
00152
00153 G4VParticleChange*
00154 DsG4Scintillation::AtRestDoIt(const G4Track& aTrack, const G4Step& aStep)
00155
00156
00157
00158
00159 {
00160 return DsG4Scintillation::PostStepDoIt(aTrack, aStep);
00161 }
00162
00163
00164
00165
00166 G4VParticleChange*
00167 DsG4Scintillation::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
00168
00169
00170
00171
00172
00173
00174 {
00175 aParticleChange.Initialize(aTrack);
00176
00177 if (m_noop) {
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
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
00267
00268
00269
00270
00271
00272
00273
00274
00275
00276
00277
00278
00279
00280 G4double ScintillationYield = 0;
00281 {
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 {
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
00308 }
00309
00310 G4double delta = dE_dx/aMaterial->GetDensity();
00311
00312 G4double birk1 = birksConstant1;
00313 if(abs(aParticle->GetCharge())>1.5)
00314 birk1 = 0.57*birk1;
00315
00316 G4double birk2 = 0;
00317
00318 birk2 = birksConstant2;
00319
00320 G4double QuenchedTotalEnergyDeposit
00321 = TotalEnergyDeposit/(1+birk1*delta+birk2*delta*delta);
00322
00323
00324 if(FastMu300nsTrick) {
00325
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
00340
00341
00342
00343
00344
00345
00346
00347
00348
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
00361 if (NumTracks <= 0) {
00362
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
00386
00387
00388 G4int Num = NumTracks;
00389
00390 G4double fastTimeConstant = 0.0;
00391 {
00392 const G4MaterialPropertyVector* ptable =
00393 aMaterialPropertiesTable->GetProperty("FASTTIMECONSTANT");
00394 if (ptable)
00395 fastTimeConstant = ptable->GetProperty(0);
00396 }
00397
00398 G4double slowTimeConstant = 0.0;
00399 {
00400 const G4MaterialPropertyVector* ptable =
00401 aMaterialPropertiesTable->GetProperty("SLOWTIMECONSTANT");
00402 if (ptable)
00403 slowTimeConstant = ptable->GetProperty(0);
00404 }
00405
00406 G4double YieldRatio = 0.0;
00407 {
00408 const G4MaterialPropertyVector* ptable =
00409 aMaterialPropertiesTable->GetProperty("YIELDRATIO");
00410 if (ptable)
00411 YieldRatio = ptable->GetProperty(0);
00412 }
00413
00414
00415 for (G4int scnt = 1; scnt <= nscnt; scnt++) {
00416
00417 G4double ScintillationTime = 0.*ns;
00418 G4PhysicsOrderedFreeVector* ScintillationIntegral = NULL;
00419
00420 if (scnt == 1) {
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 {
00446 Num = NumTracks - Num;
00447 ScintillationTime = slowTimeConstant;
00448 ScintillationIntegral =
00449 (G4PhysicsOrderedFreeVector*)((*theSlowIntegralTable)(materialIndex));
00450 }
00451
00452 if (!ScintillationIntegral) continue;
00453
00454
00455
00456 for (G4int i = 0; i < Num; i++) {
00457
00458
00459 G4double sampledEnergy;
00460 if ( !flagReemission ) {
00461
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
00475 G4double CIIvalue = G4UniformRand()*
00476 ScintillationIntegral->GetMaxValue();
00477 if (CIIvalue == 0.0) {
00478
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
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
00505
00506 G4ParticleMomentum photonMomentum(px, py, pz);
00507
00508
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
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
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
00583
00584 aSecondaryTrack->SetParentID(aTrack.GetTrackID());
00585
00586
00587 aParticleChange.SetSecondaryWeightByProcess( true );
00588 aParticleChange.AddSecondary(aSecondaryTrack);
00589
00590 aSecondaryTrack->SetWeight( weight );
00591
00592
00593
00594
00595
00596
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 }
00607
00608
00609
00610
00611
00612 void DsG4Scintillation::BuildThePhysicsTable()
00613 {
00614 if (theFastIntegralTable && theSlowIntegralTable && theReemissionIntegralTable) return;
00615
00616 const G4MaterialTable* theMaterialTable =
00617 G4Material::GetMaterialTable();
00618 G4int numOfMaterials = G4Material::GetNumberOfMaterials();
00619
00620
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
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
00638
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
00653
00654
00655 theFastLightVector->ResetIterator();
00656 ++(*theFastLightVector);
00657
00658 G4double currentIN = theFastLightVector->
00659 GetProperty();
00660
00661 if (currentIN >= 0.0) {
00662
00663
00664
00665
00666 G4double currentPM = theFastLightVector->
00667 GetPhotonEnergy();
00668
00669 G4double currentCII = 0.0;
00670
00671 aPhysicsOrderedFreeVector->
00672 InsertValues(currentPM , currentCII);
00673
00674
00675
00676 G4double prevPM = currentPM;
00677 G4double prevCII = currentCII;
00678 G4double prevIN = currentIN;
00679
00680
00681
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
00712
00713
00714 theSlowLightVector->ResetIterator();
00715 ++(*theSlowLightVector);
00716
00717 G4double currentIN = theSlowLightVector->
00718 GetProperty();
00719
00720 if (currentIN >= 0.0) {
00721
00722
00723
00724
00725 G4double currentPM = theSlowLightVector->
00726 GetPhotonEnergy();
00727
00728 G4double currentCII = 0.0;
00729
00730 bPhysicsOrderedFreeVector->
00731 InsertValues(currentPM , currentCII);
00732
00733
00734
00735 G4double prevPM = currentPM;
00736 G4double prevCII = currentCII;
00737 G4double prevIN = currentIN;
00738
00739
00740
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
00771
00772
00773 theReemissionVector->ResetIterator();
00774 ++(*theReemissionVector);
00775
00776 G4double currentIN = theReemissionVector->
00777 GetProperty();
00778
00779 if (currentIN >= 0.0) {
00780
00781
00782
00783
00784 G4double currentPM = theReemissionVector->
00785 GetPhotonEnergy();
00786
00787 G4double currentCII = 0.0;
00788
00789 cPhysicsOrderedFreeVector->
00790 InsertValues(currentPM , currentCII);
00791
00792
00793
00794 G4double prevPM = currentPM;
00795 G4double prevCII = currentCII;
00796 G4double prevIN = currentIN;
00797
00798
00799
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
00827
00828
00829
00830 theFastIntegralTable->insertAt(i,aPhysicsOrderedFreeVector);
00831 theSlowIntegralTable->insertAt(i,bPhysicsOrderedFreeVector);
00832 theReemissionIntegralTable->insertAt(i,cPhysicsOrderedFreeVector);
00833 }
00834 }
00835
00836
00837
00838
00839
00840 G4double DsG4Scintillation::GetMeanFreePath(const G4Track&,
00841 G4double ,
00842 G4ForceCondition* condition)
00843 {
00844 *condition = StronglyForced;
00845
00846 return DBL_MAX;
00847
00848 }
00849
00850
00851
00852
00853
00854 G4double DsG4Scintillation::GetMeanLifeTime(const G4Track&,
00855 G4ForceCondition* condition)
00856 {
00857 *condition = Forced;
00858
00859 return DBL_MAX;
00860
00861 }