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 #include "G4UImanager.hh"
00047 #include "G4ForceCondition.hh"
00048 #include "G4GPILSelection.hh"
00049 #include "G4SteppingControl.hh"
00050 #include "G4TransportationManager.hh"
00051 #include "G4SteppingManager.hh"
00052 #include "G4LossTableManager.hh"
00053
00055 void G4SteppingManager::GetProcessNumber()
00057 {
00058 #ifdef debug
00059 G4cout<<"G4SteppingManager::GetProcessNumber: is called track="
00060 <<fTrack<<G4endl;
00061 #endif
00062
00063 G4ProcessManager* pm= fTrack->GetDefinition()->GetProcessManager();
00064 if(!pm)
00065 {
00066 G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
00067 << " ProcessManager is NULL for particle = "
00068 << fTrack->GetDefinition()->GetParticleName() << ", PDG_code = "
00069 << fTrack->GetDefinition()->GetPDGEncoding() << G4endl;
00070 G4Exception("G4SteppingManager::GetProcessNumber()", "Tracking0011",
00071 FatalException, "Process Manager is not found.");
00072 return;
00073 }
00074
00075
00076 MAXofAtRestLoops = pm->GetAtRestProcessVector()->entries();
00077 fAtRestDoItVector = pm->GetAtRestProcessVector(typeDoIt);
00078 fAtRestGetPhysIntVector = pm->GetAtRestProcessVector(typeGPIL);
00079 #ifdef debug
00080 G4cout << "G4SteppingManager::GetProcessNumber: #ofAtRest="
00081 << MAXofAtRestLoops << G4endl;
00082 #endif
00083
00084
00085 MAXofAlongStepLoops = pm->GetAlongStepProcessVector()->entries();
00086 fAlongStepDoItVector = pm->GetAlongStepProcessVector(typeDoIt);
00087 fAlongStepGetPhysIntVector = pm->GetAlongStepProcessVector(typeGPIL);
00088 #ifdef debug
00089 G4cout << "G4SteppingManager::GetProcessNumber:#ofAlongStp="
00090 << MAXofAlongStepLoops << G4endl;
00091 #endif
00092
00093
00094 MAXofPostStepLoops = pm->GetPostStepProcessVector()->entries();
00095 fPostStepDoItVector = pm->GetPostStepProcessVector(typeDoIt);
00096 fPostStepGetPhysIntVector = pm->GetPostStepProcessVector(typeGPIL);
00097 #ifdef debug
00098 G4cout << "G4SteppingManager::GetProcessNumber: #ofPostStep="
00099 << MAXofPostStepLoops << G4endl;
00100 #endif
00101
00102 if (SizeOfSelectedDoItVector<MAXofAtRestLoops ||
00103 SizeOfSelectedDoItVector<MAXofAlongStepLoops ||
00104 SizeOfSelectedDoItVector<MAXofPostStepLoops )
00105 {
00106 G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
00107 << " SizeOfSelectedDoItVector= " << SizeOfSelectedDoItVector
00108 << " ; is smaller then one of MAXofAtRestLoops= "
00109 << MAXofAtRestLoops << G4endl
00110 << " or MAXofAlongStepLoops= " << MAXofAlongStepLoops
00111 << " or MAXofPostStepLoops= " << MAXofPostStepLoops << G4endl;
00112 G4Exception("G4SteppingManager::GetProcessNumber()",
00113 "Tracking0012", FatalException,
00114 "The array size is smaller than the actual No of processes.");
00115 }
00116 }
00117
00118
00119
00120
00121
00122
00123
00124
00125
00127 void G4SteppingManager::DefinePhysicalStepLength()
00129 {
00130
00131
00132 PhysicalStep = DBL_MAX;
00133 physIntLength = DBL_MAX;
00134 #ifdef G4VERBOSE
00135
00136 if(verboseLevel>0) fVerbose->DPSLStarted();
00137 #endif
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162 fPostStepDoItProcTriggered = MAXofPostStepLoops;
00163
00164 for(size_t np=0; np < MAXofPostStepLoops; np++){
00165 fCurrentProcess = (*fPostStepGetPhysIntVector)(np);
00166 if (fCurrentProcess== 0) {
00167 (*fSelectedPostStepDoItVector)[np] = InActivated;
00168 continue;
00169 }
00170
00171 physIntLength = fCurrentProcess->
00172 PostStepGPIL( *fTrack,
00173 fPreviousStepSize,
00174 &fCondition );
00175 #ifdef G4VERBOSE
00176
00177 if(verboseLevel>0) fVerbose->DPSLPostStep();
00178 #endif
00179
00180 switch (fCondition) {
00181 case ExclusivelyForced:
00182 (*fSelectedPostStepDoItVector)[np] = ExclusivelyForced;
00183 fStepStatus = fExclusivelyForcedProc;
00184 fStep->GetPostStepPoint()
00185 ->SetProcessDefinedStep(fCurrentProcess);
00186 break;
00187 case Conditionally:
00188
00189 G4Exception("G4SteppingManager::DefinePhysicalStepLength()", "Tracking1001", FatalException, "This feature no more supported");
00190
00191 break;
00192 case Forced:
00193 (*fSelectedPostStepDoItVector)[np] = Forced;
00194 break;
00195 case StronglyForced:
00196 (*fSelectedPostStepDoItVector)[np] = StronglyForced;
00197 break;
00198 default:
00199 (*fSelectedPostStepDoItVector)[np] = InActivated;
00200 break;
00201 }
00202
00203
00204
00205 if (fCondition==ExclusivelyForced) {
00206 for(size_t nrest=np+1; nrest < MAXofPostStepLoops; nrest++){
00207 (*fSelectedPostStepDoItVector)[nrest] = InActivated;
00208 }
00209 return;
00210 }
00211 else{
00212 if(physIntLength < PhysicalStep ){
00213 PhysicalStep = physIntLength;
00214 fStepStatus = fPostStepDoItProc;
00215 fPostStepDoItProcTriggered = G4int(np);
00216 fStep->GetPostStepPoint()
00217 ->SetProcessDefinedStep(fCurrentProcess);
00218 }
00219 }
00220
00221
00222 }
00223
00224 if (fPostStepDoItProcTriggered<MAXofPostStepLoops) {
00225 if ((*fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] ==
00226 InActivated) {
00227 (*fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] =
00228 NotForced;
00229 }
00230 }
00231
00232
00233 proposedSafety = DBL_MAX;
00234 G4double safetyProposedToAndByProcess = proposedSafety;
00235
00236 for(size_t kp=0; kp < MAXofAlongStepLoops; kp++){
00237 fCurrentProcess = (*fAlongStepGetPhysIntVector)[kp];
00238 if (fCurrentProcess== 0) continue;
00239
00240
00241 physIntLength = fCurrentProcess->
00242 AlongStepGPIL( *fTrack, fPreviousStepSize,
00243 PhysicalStep,
00244 safetyProposedToAndByProcess,
00245 &fGPILSelection );
00246 #ifdef G4VERBOSE
00247
00248 if(verboseLevel>0) fVerbose->DPSLAlongStep();
00249 #endif
00250 if(physIntLength < PhysicalStep){
00251 PhysicalStep = physIntLength;
00252
00253
00254
00255 if(fGPILSelection==CandidateForSelection){
00256 fStepStatus = fAlongStepDoItProc;
00257 fStep->GetPostStepPoint()
00258 ->SetProcessDefinedStep(fCurrentProcess);
00259 }
00260
00261
00262 if(kp == MAXofAlongStepLoops-1)
00263 fStepStatus = fGeomBoundary;
00264 }
00265
00266
00267
00268
00269 if (safetyProposedToAndByProcess < proposedSafety)
00270
00271 proposedSafety = safetyProposedToAndByProcess;
00272 else
00273
00274 safetyProposedToAndByProcess = proposedSafety;
00275
00276 }
00277 }
00278
00279
00281 void G4SteppingManager::InvokeAtRestDoItProcs()
00283 {
00284
00285
00286
00287 G4double lifeTime, shortestLifeTime;
00288
00289 fAtRestDoItProcTriggered = 0;
00290 shortestLifeTime = DBL_MAX;
00291
00292 unsigned int NofInactiveProc=0;
00293 for( size_t ri=0 ; ri < MAXofAtRestLoops ; ri++ ){
00294 fCurrentProcess = (*fAtRestGetPhysIntVector)[ri];
00295 if (fCurrentProcess== 0) {
00296 (*fSelectedAtRestDoItVector)[ri] = InActivated;
00297 NofInactiveProc++;
00298 continue;
00299 }
00300
00301 lifeTime =
00302 fCurrentProcess->AtRestGPIL( *fTrack, &fCondition );
00303
00304
00305
00306 if(fCondition==Forced && fCurrentProcess){
00307 (*fSelectedAtRestDoItVector)[ri] = Forced;
00308 }
00309 else{
00310 (*fSelectedAtRestDoItVector)[ri] = InActivated;
00311 if(lifeTime < shortestLifeTime ){
00312 shortestLifeTime = lifeTime;
00313 fAtRestDoItProcTriggered = G4int(ri);
00314 }
00315 }
00316 }
00317
00318 (*fSelectedAtRestDoItVector)[fAtRestDoItProcTriggered] = NotForced;
00319
00320
00321
00322 if(NofInactiveProc==MAXofAtRestLoops){
00323 G4Exception("G4SteppingManager::InvokeAtRestDoItProcs()", "Tracking0013",
00324 JustWarning, "No AtRestDoIt process is active!" );
00325 }
00326
00327 fStep->SetStepLength( 0. );
00328 fTrack->SetStepLength( 0. );
00329
00330
00331 for(size_t np=0; np < MAXofAtRestLoops; np++){
00332
00333
00334
00335
00336 if( (*fSelectedAtRestDoItVector)[MAXofAtRestLoops-np-1] != InActivated){
00337
00338 fCurrentProcess = (*fAtRestDoItVector)[np];
00339 fParticleChange
00340 = fCurrentProcess->AtRestDoIt( *fTrack, *fStep);
00341
00342
00343 fStep->GetPostStepPoint()
00344 ->SetProcessDefinedStep(fCurrentProcess);
00345
00346
00347 fParticleChange->UpdateStepForAtRest(fStep);
00348
00349
00350 G4Track* tempSecondaryTrack;
00351 G4int num2ndaries;
00352
00353 num2ndaries = fParticleChange->GetNumberOfSecondaries();
00354
00355 for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
00356 tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
00357
00358 if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
00359 { ApplyProductionCut(tempSecondaryTrack); }
00360
00361
00362 tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
00363
00364
00365 tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
00366
00367
00368
00369 if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
00370 G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
00371 if (pm->GetAtRestProcessVector()->entries()>0){
00372 tempSecondaryTrack->SetTrackStatus( fStopButAlive );
00373 fSecondary->push_back( tempSecondaryTrack );
00374 fN2ndariesAtRestDoIt++;
00375 } else {
00376 delete tempSecondaryTrack;
00377 }
00378 } else {
00379 fSecondary->push_back( tempSecondaryTrack );
00380 fN2ndariesAtRestDoIt++;
00381 }
00382 }
00383
00384
00385
00386 fParticleChange->Clear();
00387
00388 }
00389 }
00390 fStep->UpdateTrack();
00391
00392 fTrack->SetTrackStatus( fStopAndKill );
00393 }
00394
00395
00397 void G4SteppingManager::InvokeAlongStepDoItProcs()
00399 {
00400
00401
00402
00403 if(fStepStatus == fExclusivelyForcedProc){
00404 return;
00405 }
00406
00407
00408 for( size_t ci=0 ; ci<MAXofAlongStepLoops ; ci++ ){
00409 fCurrentProcess = (*fAlongStepDoItVector)[ci];
00410 if (fCurrentProcess== 0) continue;
00411
00412
00413 fParticleChange
00414 = fCurrentProcess->AlongStepDoIt( *fTrack, *fStep );
00415
00416
00417 fParticleChange->UpdateStepForAlongStep(fStep);
00418 #ifdef G4VERBOSE
00419
00420 if(verboseLevel>0) fVerbose->AlongStepDoItOneByOne();
00421 #endif
00422
00423
00424 G4Track* tempSecondaryTrack;
00425 G4int num2ndaries;
00426
00427 num2ndaries = fParticleChange->GetNumberOfSecondaries();
00428
00429 for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
00430 tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
00431
00432 if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
00433 { ApplyProductionCut(tempSecondaryTrack); }
00434
00435
00436 tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
00437
00438
00439 tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
00440
00441
00442
00443 if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
00444 G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
00445 if (pm->GetAtRestProcessVector()->entries()>0){
00446 tempSecondaryTrack->SetTrackStatus( fStopButAlive );
00447 fSecondary->push_back( tempSecondaryTrack );
00448 fN2ndariesAlongStepDoIt++;
00449 } else {
00450 delete tempSecondaryTrack;
00451 }
00452 } else {
00453 fSecondary->push_back( tempSecondaryTrack );
00454 fN2ndariesAlongStepDoIt++;
00455 }
00456 }
00457
00458
00459
00460 fTrack->SetTrackStatus( fParticleChange->GetTrackStatus() );
00461
00462
00463 fParticleChange->Clear();
00464 }
00465
00466 fStep->UpdateTrack();
00467 G4TrackStatus fNewStatus = fTrack->GetTrackStatus();
00468
00469 if ( fNewStatus == fAlive && fTrack->GetKineticEnergy() <= DBL_MIN ) {
00470 if(MAXofAtRestLoops>0) fNewStatus = fStopButAlive;
00471 else fNewStatus = fStopAndKill;
00472 fTrack->SetTrackStatus( fNewStatus );
00473 }
00474
00475 }
00476
00478 void G4SteppingManager::InvokePostStepDoItProcs()
00480 {
00481
00482
00483 for(size_t np=0; np < MAXofPostStepLoops; np++){
00484
00485
00486
00487
00488 G4int Cond = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops-np-1];
00489 if(Cond != InActivated){
00490 if( ((Cond == NotForced) && (fStepStatus == fPostStepDoItProc)) ||
00491 ((Cond == Forced) && (fStepStatus != fExclusivelyForcedProc)) ||
00492
00493 ((Cond == ExclusivelyForced) && (fStepStatus == fExclusivelyForcedProc)) ||
00494 ((Cond == StronglyForced) )
00495 ) {
00496
00497 InvokePSDIP(np);
00498 if ((np==0) && (fTrack->GetNextVolume() == 0)){
00499 fStepStatus = fWorldBoundary;
00500 fStep->GetPostStepPoint()->SetStepStatus( fStepStatus );
00501 }
00502 }
00503 }
00504
00505
00506
00507 if(fTrack->GetTrackStatus() == fStopAndKill) {
00508 for(size_t np1=np+1; np1 < MAXofPostStepLoops; np1++){
00509 G4int Cond2 = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops-np1-1];
00510 if (Cond2 == StronglyForced) {
00511 InvokePSDIP(np1);
00512 }
00513 }
00514 break;
00515 }
00516 }
00517 }
00518
00519
00520
00521 void G4SteppingManager::InvokePSDIP(size_t np)
00522 {
00523 fCurrentProcess = (*fPostStepDoItVector)[np];
00524 fParticleChange
00525 = fCurrentProcess->PostStepDoIt( *fTrack, *fStep);
00526
00527
00528 fParticleChange->UpdateStepForPostStep(fStep);
00529 #ifdef G4VERBOSE
00530
00531 if(verboseLevel>0) fVerbose->PostStepDoItOneByOne();
00532 #endif
00533
00534 fStep->UpdateTrack();
00535
00536
00537 fStep->GetPostStepPoint()->SetSafety( CalculateSafety() );
00538
00539
00540 G4Track* tempSecondaryTrack;
00541 G4int num2ndaries;
00542
00543 num2ndaries = fParticleChange->GetNumberOfSecondaries();
00544
00545 for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
00546 tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
00547
00548 if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
00549 { ApplyProductionCut(tempSecondaryTrack); }
00550
00551
00552 tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
00553
00554
00555 tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
00556
00557
00558
00559 if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
00560 G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
00561 if (pm->GetAtRestProcessVector()->entries()>0){
00562 tempSecondaryTrack->SetTrackStatus( fStopButAlive );
00563 fSecondary->push_back( tempSecondaryTrack );
00564 fN2ndariesPostStepDoIt++;
00565 } else {
00566 delete tempSecondaryTrack;
00567 }
00568 } else {
00569 fSecondary->push_back( tempSecondaryTrack );
00570 fN2ndariesPostStepDoIt++;
00571 }
00572 }
00573
00574
00575 fTrack->SetTrackStatus( fParticleChange->GetTrackStatus() );
00576
00577
00578 fParticleChange->Clear();
00579 }
00580
00581 #include "G4EnergyLossTables.hh"
00582 #include "G4ProductionCuts.hh"
00583 #include "G4ProductionCutsTable.hh"
00584
00585 void G4SteppingManager::ApplyProductionCut(G4Track* aSecondary)
00586 {
00587 G4bool tBelowCutEnergyAndSafety = false;
00588 G4int tPtclIdx
00589 = G4ProductionCuts::GetIndex(aSecondary->GetDefinition());
00590 if (tPtclIdx<0) { return; }
00591 G4ProductionCutsTable* tCutsTbl
00592 = G4ProductionCutsTable::GetProductionCutsTable();
00593 G4int tCoupleIdx
00594 = tCutsTbl->GetCoupleIndex(fPreStepPoint->GetMaterialCutsCouple());
00595 G4double tProdThreshold
00596 = (*(tCutsTbl->GetEnergyCutsVector(tPtclIdx)))[tCoupleIdx];
00597 if( aSecondary->GetKineticEnergy()<tProdThreshold )
00598 {
00599 tBelowCutEnergyAndSafety = true;
00600 if(std::abs(aSecondary->GetDynamicParticle()->GetCharge()) > DBL_MIN)
00601 {
00602 G4double currentRange
00603 = G4LossTableManager::Instance()->GetRange(aSecondary->GetDefinition(),
00604 aSecondary->GetKineticEnergy(),
00605 fPreStepPoint->GetMaterialCutsCouple());
00606 tBelowCutEnergyAndSafety = (currentRange < CalculateSafety() );
00607 }
00608 }
00609
00610 if( tBelowCutEnergyAndSafety )
00611 {
00612 if( !(aSecondary->IsGoodForTracking()) )
00613 {
00614
00615
00616
00617
00618
00619 fStep->AddTotalEnergyDeposit(
00620 aSecondary->GetKineticEnergy() );
00621 aSecondary->SetKineticEnergy(0.0);
00622 }
00623 }
00624 }
00625