Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4SteppingManager2.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4SteppingManager2.cc 69894 2013-05-17 09:47:23Z gcosmo $
28 //
29 //---------------------------------------------------------------
30 //
31 // G4SteppingManager2.cc
32 //
33 // Description:
34 // This class represents the manager who steers to move the give
35 // particle from the TrackingManger by one Step.
36 //
37 // Contact:
38 // Questions and comments to this code should be sent to
39 // Katsuya Amako (e-mail: Katsuya.Amako@kek.jp)
40 // Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp)
41 //
42 //---------------------------------------------------------------
43 
44 //#define debug
45 
46 #include "G4UImanager.hh"
47 #include "G4ForceCondition.hh"
48 #include "G4GPILSelection.hh"
49 #include "G4SteppingControl.hh"
51 #include "G4SteppingManager.hh"
52 #include "G4LossTableManager.hh"
53 
57 {
58 #ifdef debug
59  G4cout<<"G4SteppingManager::GetProcessNumber: is called track="
60  <<fTrack<<G4endl;
61 #endif
62 
64  if(!pm)
65  {
66  G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
67  << " ProcessManager is NULL for particle = "
68  << fTrack->GetDefinition()->GetParticleName() << ", PDG_code = "
69  << fTrack->GetDefinition()->GetPDGEncoding() << G4endl;
70  G4Exception("G4SteppingManager::GetProcessNumber()", "Tracking0011",
71  FatalException, "Process Manager is not found.");
72  return;
73  }
74 
75 // AtRestDoits
76  MAXofAtRestLoops = pm->GetAtRestProcessVector()->entries();
77  fAtRestDoItVector = pm->GetAtRestProcessVector(typeDoIt);
78  fAtRestGetPhysIntVector = pm->GetAtRestProcessVector(typeGPIL);
79 #ifdef debug
80  G4cout << "G4SteppingManager::GetProcessNumber: #ofAtRest="
81  << MAXofAtRestLoops << G4endl;
82 #endif
83 
84 // AlongStepDoits
85  MAXofAlongStepLoops = pm->GetAlongStepProcessVector()->entries();
86  fAlongStepDoItVector = pm->GetAlongStepProcessVector(typeDoIt);
87  fAlongStepGetPhysIntVector = pm->GetAlongStepProcessVector(typeGPIL);
88 #ifdef debug
89  G4cout << "G4SteppingManager::GetProcessNumber:#ofAlongStp="
90  << MAXofAlongStepLoops << G4endl;
91 #endif
92 
93 // PostStepDoits
94  MAXofPostStepLoops = pm->GetPostStepProcessVector()->entries();
95  fPostStepDoItVector = pm->GetPostStepProcessVector(typeDoIt);
96  fPostStepGetPhysIntVector = pm->GetPostStepProcessVector(typeGPIL);
97 #ifdef debug
98  G4cout << "G4SteppingManager::GetProcessNumber: #ofPostStep="
99  << MAXofPostStepLoops << G4endl;
100 #endif
101 
102  if (SizeOfSelectedDoItVector<MAXofAtRestLoops ||
103  SizeOfSelectedDoItVector<MAXofAlongStepLoops ||
104  SizeOfSelectedDoItVector<MAXofPostStepLoops )
105  {
106  G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
107  << " SizeOfSelectedDoItVector= " << SizeOfSelectedDoItVector
108  << " ; is smaller then one of MAXofAtRestLoops= "
109  << MAXofAtRestLoops << G4endl
110  << " or MAXofAlongStepLoops= " << MAXofAlongStepLoops
111  << " or MAXofPostStepLoops= " << MAXofPostStepLoops << G4endl;
112  G4Exception("G4SteppingManager::GetProcessNumber()",
113  "Tracking0012", FatalException,
114  "The array size is smaller than the actual No of processes.");
115  }
116 }
117 
118 
119 // ************************************************************************
120 //
121 // Private Member Functions
122 //
123 // ************************************************************************
124 
125 
127  void G4SteppingManager::DefinePhysicalStepLength()
129 {
130 
131 // ReSet the counter etc.
132  PhysicalStep = DBL_MAX; // Initialize by a huge number
133  physIntLength = DBL_MAX; // Initialize by a huge number
134 #ifdef G4VERBOSE
135  // !!!!! Verbose
136  if(verboseLevel>0) fVerbose->DPSLStarted();
137 #endif
138 
139 // Obtain the user defined maximum allowed Step in the volume
140 // 1997.12.13 adds argument for GetMaxAllowedStep by K.Kurashige
141 // 2004.01.20 This block will be removed by Geant4 7.0
142 // G4UserLimits* ul= fCurrentVolume->GetLogicalVolume()->GetUserLimits();
143 // if (ul) {
144 // physIntLength = ul->GetMaxAllowedStep(*fTrack);
145 //#ifdef G4VERBOSE
146 // // !!!!! Verbose
147 // if(verboseLevel>0) fVerbose->DPSLUserLimit();
148 //#endif
149 // }
150 //
151 // if(physIntLength < PhysicalStep ){
152 // PhysicalStep = physIntLength;
153 // fStepStatus = fUserDefinedLimit;
154 // fStep->GetPostStepPoint()
155 // ->SetProcessDefinedStep(0);
156 // // Take note that the process pointer is 'NULL' if the Step
157 // // is defined by the user defined limit.
158 // }
159 // 2004.01.20 This block will be removed by Geant4 7.0
160 
161 // GPIL for PostStep
162  fPostStepDoItProcTriggered = MAXofPostStepLoops;
163 
164  for(size_t np=0; np < MAXofPostStepLoops; np++){
165  fCurrentProcess = (*fPostStepGetPhysIntVector)(np);
166  if (fCurrentProcess== 0) {
167  (*fSelectedPostStepDoItVector)[np] = InActivated;
168  continue;
169  } // NULL means the process is inactivated by a user on fly.
170 
171  physIntLength = fCurrentProcess->
172  PostStepGPIL( *fTrack,
173  fPreviousStepSize,
174  &fCondition );
175 #ifdef G4VERBOSE
176  // !!!!! Verbose
177  if(verboseLevel>0) fVerbose->DPSLPostStep();
178 #endif
179 
180  switch (fCondition) {
181  case ExclusivelyForced:
182  (*fSelectedPostStepDoItVector)[np] = ExclusivelyForced;
183  fStepStatus = fExclusivelyForcedProc;
184  fStep->GetPostStepPoint()
185  ->SetProcessDefinedStep(fCurrentProcess);
186  break;
187  case Conditionally:
188  // (*fSelectedPostStepDoItVector)[np] = Conditionally;
189  G4Exception("G4SteppingManager::DefinePhysicalStepLength()", "Tracking1001", FatalException, "This feature no more supported");
190 
191  break;
192  case Forced:
193  (*fSelectedPostStepDoItVector)[np] = Forced;
194  break;
195  case StronglyForced:
196  (*fSelectedPostStepDoItVector)[np] = StronglyForced;
197  break;
198  default:
199  (*fSelectedPostStepDoItVector)[np] = InActivated;
200  break;
201  }
202 
203 
204 
205  if (fCondition==ExclusivelyForced) {
206  for(size_t nrest=np+1; nrest < MAXofPostStepLoops; nrest++){
207  (*fSelectedPostStepDoItVector)[nrest] = InActivated;
208  }
209  return; // Take note the 'return' at here !!!
210  }
211  else{
212  if(physIntLength < PhysicalStep ){
213  PhysicalStep = physIntLength;
214  fStepStatus = fPostStepDoItProc;
215  fPostStepDoItProcTriggered = G4int(np);
216  fStep->GetPostStepPoint()
217  ->SetProcessDefinedStep(fCurrentProcess);
218  }
219  }
220 
221 
222  }
223 
224  if (fPostStepDoItProcTriggered<MAXofPostStepLoops) {
225  if ((*fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] ==
226  InActivated) {
227  (*fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] =
228  NotForced;
229  }
230  }
231 
232 // GPIL for AlongStep
233  proposedSafety = DBL_MAX;
234  G4double safetyProposedToAndByProcess = proposedSafety;
235 
236  for(size_t kp=0; kp < MAXofAlongStepLoops; kp++){
237  fCurrentProcess = (*fAlongStepGetPhysIntVector)[kp];
238  if (fCurrentProcess== 0) continue;
239  // NULL means the process is inactivated by a user on fly.
240 
241  physIntLength = fCurrentProcess->
242  AlongStepGPIL( *fTrack, fPreviousStepSize,
243  PhysicalStep,
244  safetyProposedToAndByProcess,
245  &fGPILSelection );
246 #ifdef G4VERBOSE
247  // !!!!! Verbose
248  if(verboseLevel>0) fVerbose->DPSLAlongStep();
249 #endif
250  if(physIntLength < PhysicalStep){
251  PhysicalStep = physIntLength;
252 
253  // Check if the process wants to be the GPIL winner. For example,
254  // multi-scattering proposes Step limit, but won't be the winner.
255  if(fGPILSelection==CandidateForSelection){
256  fStepStatus = fAlongStepDoItProc;
257  fStep->GetPostStepPoint()
258  ->SetProcessDefinedStep(fCurrentProcess);
259  }
260 
261  // Transportation is assumed to be the last process in the vector
262  if(kp == MAXofAlongStepLoops-1)
263  fStepStatus = fGeomBoundary;
264  }
265 
266  // Make sure to check the safety, even if Step is not limited
267  // by this process. J. Apostolakis, June 20, 1998
268  //
269  if (safetyProposedToAndByProcess < proposedSafety)
270  // proposedSafety keeps the smallest value:
271  proposedSafety = safetyProposedToAndByProcess;
272  else
273  // safetyProposedToAndByProcess always proposes a valid safety:
274  safetyProposedToAndByProcess = proposedSafety;
275 
276  }
277 } // void G4SteppingManager::DefinePhysicalStepLength() //
278 
279 
281 void G4SteppingManager::InvokeAtRestDoItProcs()
283 {
284 // Select the rest process which has the shortest time before
285 // it is invoked. In rest processes, GPIL()
286 // returns the time before a process occurs.
287  G4double lifeTime, shortestLifeTime;
288 
289  fAtRestDoItProcTriggered = 0;
290  shortestLifeTime = DBL_MAX;
291 
292  unsigned int NofInactiveProc=0;
293  for( size_t ri=0 ; ri < MAXofAtRestLoops ; ri++ ){
294  fCurrentProcess = (*fAtRestGetPhysIntVector)[ri];
295  if (fCurrentProcess== 0) {
296  (*fSelectedAtRestDoItVector)[ri] = InActivated;
297  NofInactiveProc++;
298  continue;
299  } // NULL means the process is inactivated by a user on fly.
300 
301  lifeTime =
302  fCurrentProcess->AtRestGPIL( *fTrack, &fCondition );
303 
304 
305 
306  if(fCondition==Forced && fCurrentProcess){
307  (*fSelectedAtRestDoItVector)[ri] = Forced;
308  }
309  else{
310  (*fSelectedAtRestDoItVector)[ri] = InActivated;
311  if(lifeTime < shortestLifeTime ){
312  shortestLifeTime = lifeTime;
313  fAtRestDoItProcTriggered = G4int(ri);
314  }
315  }
316  }
317 
318  (*fSelectedAtRestDoItVector)[fAtRestDoItProcTriggered] = NotForced;
319 
320 // at least one process is necessary to destroy the particle
321 // exit with warning
322  if(NofInactiveProc==MAXofAtRestLoops){
323  G4Exception("G4SteppingManager::InvokeAtRestDoItProcs()", "Tracking0013",
324  JustWarning, "No AtRestDoIt process is active!" );
325  }
326 
327  fStep->SetStepLength( 0. ); //the particle has stopped
328  fTrack->SetStepLength( 0. );
329 
330 // invoke selected process
331  for(size_t np=0; np < MAXofAtRestLoops; np++){
332  //
333  // Note: DoItVector has inverse order against GetPhysIntVector
334  // and SelectedAtRestDoItVector.
335  //
336  if( (*fSelectedAtRestDoItVector)[MAXofAtRestLoops-np-1] != InActivated){
337 
338  fCurrentProcess = (*fAtRestDoItVector)[np];
339  fParticleChange
340  = fCurrentProcess->AtRestDoIt( *fTrack, *fStep);
341 
342  // Set the current process as a process which defined this Step length
343  fStep->GetPostStepPoint()
344  ->SetProcessDefinedStep(fCurrentProcess);
345 
346  // Update Step
347  fParticleChange->UpdateStepForAtRest(fStep);
348 
349  // Now Store the secondaries from ParticleChange to SecondaryList
350  G4Track* tempSecondaryTrack;
351  G4int num2ndaries;
352 
353  num2ndaries = fParticleChange->GetNumberOfSecondaries();
354 
355  for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
356  tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
357 
358  if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
359  { ApplyProductionCut(tempSecondaryTrack); }
360 
361  // Set parentID
362  tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
363 
364  // Set the process pointer which created this track
365  tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
366 
367  // If this 2ndry particle has 'zero' kinetic energy, make sure
368  // it invokes a rest process at the beginning of the tracking
369  if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
370  G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
371  if (pm->GetAtRestProcessVector()->entries()>0){
372  tempSecondaryTrack->SetTrackStatus( fStopButAlive );
373  fSecondary->push_back( tempSecondaryTrack );
374  fN2ndariesAtRestDoIt++;
375  } else {
376  delete tempSecondaryTrack;
377  }
378  } else {
379  fSecondary->push_back( tempSecondaryTrack );
380  fN2ndariesAtRestDoIt++;
381  }
382  } //end of loop on secondary
383 
384 
385  // clear ParticleChange
386  fParticleChange->Clear();
387 
388  } //if(fSelectedAtRestDoItVector[np] != InActivated){
389  } //for(size_t np=0; np < MAXofAtRestLoops; np++){
390  fStep->UpdateTrack();
391 
392  fTrack->SetTrackStatus( fStopAndKill );
393 }
394 
395 
397 void G4SteppingManager::InvokeAlongStepDoItProcs()
399 {
400 
401 // If the current Step is defined by a 'ExclusivelyForced'
402 // PostStepDoIt, then don't invoke any AlongStepDoIt
403  if(fStepStatus == fExclusivelyForcedProc){
404  return; // Take note 'return' at here !!!
405  }
406 
407 // Invoke the all active continuous processes
408  for( size_t ci=0 ; ci<MAXofAlongStepLoops ; ci++ ){
409  fCurrentProcess = (*fAlongStepDoItVector)[ci];
410  if (fCurrentProcess== 0) continue;
411  // NULL means the process is inactivated by a user on fly.
412 
413  fParticleChange
414  = fCurrentProcess->AlongStepDoIt( *fTrack, *fStep );
415 
416  // Update the PostStepPoint of Step according to ParticleChange
417  fParticleChange->UpdateStepForAlongStep(fStep);
418 #ifdef G4VERBOSE
419  // !!!!! Verbose
420  if(verboseLevel>0) fVerbose->AlongStepDoItOneByOne();
421 #endif
422 
423  // Now Store the secondaries from ParticleChange to SecondaryList
424  G4Track* tempSecondaryTrack;
425  G4int num2ndaries;
426 
427  num2ndaries = fParticleChange->GetNumberOfSecondaries();
428 
429  for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
430  tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
431 
432  if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
433  { ApplyProductionCut(tempSecondaryTrack); }
434 
435  // Set parentID
436  tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
437 
438  // Set the process pointer which created this track
439  tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
440 
441  // If this 2ndry particle has 'zero' kinetic energy, make sure
442  // it invokes a rest process at the beginning of the tracking
443  if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
444  G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
445  if (pm->GetAtRestProcessVector()->entries()>0){
446  tempSecondaryTrack->SetTrackStatus( fStopButAlive );
447  fSecondary->push_back( tempSecondaryTrack );
448  fN2ndariesAlongStepDoIt++;
449  } else {
450  delete tempSecondaryTrack;
451  }
452  } else {
453  fSecondary->push_back( tempSecondaryTrack );
454  fN2ndariesAlongStepDoIt++;
455  }
456  } //end of loop on secondary
457 
458  // Set the track status according to what the process defined
459  // if kinetic energy >0, otherwise set fStopButAlive
460  fTrack->SetTrackStatus( fParticleChange->GetTrackStatus() );
461 
462  // clear ParticleChange
463  fParticleChange->Clear();
464  }
465 
466  fStep->UpdateTrack();
467  G4TrackStatus fNewStatus = fTrack->GetTrackStatus();
468 
469  if ( fNewStatus == fAlive && fTrack->GetKineticEnergy() <= DBL_MIN ) {
470  if(MAXofAtRestLoops>0) fNewStatus = fStopButAlive;
471  else fNewStatus = fStopAndKill;
472  fTrack->SetTrackStatus( fNewStatus );
473  }
474 
475 }
476 
478 void G4SteppingManager::InvokePostStepDoItProcs()
480 {
481 
482 // Invoke the specified discrete processes
483  for(size_t np=0; np < MAXofPostStepLoops; np++){
484  //
485  // Note: DoItVector has inverse order against GetPhysIntVector
486  // and SelectedPostStepDoItVector.
487  //
488  G4int Cond = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops-np-1];
489  if(Cond != InActivated){
490  if( ((Cond == NotForced) && (fStepStatus == fPostStepDoItProc)) ||
491  ((Cond == Forced) && (fStepStatus != fExclusivelyForcedProc)) ||
492  // ((Cond == Conditionally) && (fStepStatus == fAlongStepDoItProc)) ||
493  ((Cond == ExclusivelyForced) && (fStepStatus == fExclusivelyForcedProc)) ||
494  ((Cond == StronglyForced) )
495  ) {
496 
497  InvokePSDIP(np);
498  if ((np==0) && (fTrack->GetNextVolume() == 0)){
499  fStepStatus = fWorldBoundary;
500  fStep->GetPostStepPoint()->SetStepStatus( fStepStatus );
501  }
502  }
503  } //if(*fSelectedPostStepDoItVector(np)........
504 
505  // Exit from PostStepLoop if the track has been killed,
506  // but extra treatment for processes with Strongly Forced flag
507  if(fTrack->GetTrackStatus() == fStopAndKill) {
508  for(size_t np1=np+1; np1 < MAXofPostStepLoops; np1++){
509  G4int Cond2 = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops-np1-1];
510  if (Cond2 == StronglyForced) {
511  InvokePSDIP(np1);
512  }
513  }
514  break;
515  }
516  } //for(size_t np=0; np < MAXofPostStepLoops; np++){
517 }
518 
519 
520 
521 void G4SteppingManager::InvokePSDIP(size_t np)
522 {
523  fCurrentProcess = (*fPostStepDoItVector)[np];
524  fParticleChange
525  = fCurrentProcess->PostStepDoIt( *fTrack, *fStep);
526 
527  // Update PostStepPoint of Step according to ParticleChange
528  fParticleChange->UpdateStepForPostStep(fStep);
529 #ifdef G4VERBOSE
530  // !!!!! Verbose
531  if(verboseLevel>0) fVerbose->PostStepDoItOneByOne();
532 #endif
533  // Update G4Track according to ParticleChange after each PostStepDoIt
534  fStep->UpdateTrack();
535 
536  // Update safety after each invocation of PostStepDoIts
537  fStep->GetPostStepPoint()->SetSafety( CalculateSafety() );
538 
539  // Now Store the secondaries from ParticleChange to SecondaryList
540  G4Track* tempSecondaryTrack;
541  G4int num2ndaries;
542 
543  num2ndaries = fParticleChange->GetNumberOfSecondaries();
544 
545  for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
546  tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
547 
548  if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
549  { ApplyProductionCut(tempSecondaryTrack); }
550 
551  // Set parentID
552  tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
553 
554  // Set the process pointer which created this track
555  tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
556 
557  // If this 2ndry particle has 'zero' kinetic energy, make sure
558  // it invokes a rest process at the beginning of the tracking
559  if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
560  G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
561  if (pm->GetAtRestProcessVector()->entries()>0){
562  tempSecondaryTrack->SetTrackStatus( fStopButAlive );
563  fSecondary->push_back( tempSecondaryTrack );
564  fN2ndariesPostStepDoIt++;
565  } else {
566  delete tempSecondaryTrack;
567  }
568  } else {
569  fSecondary->push_back( tempSecondaryTrack );
570  fN2ndariesPostStepDoIt++;
571  }
572  } //end of loop on secondary
573 
574  // Set the track status according to what the process defined
575  fTrack->SetTrackStatus( fParticleChange->GetTrackStatus() );
576 
577  // clear ParticleChange
578  fParticleChange->Clear();
579 }
580 
581 #include "G4EnergyLossTables.hh"
582 #include "G4ProductionCuts.hh"
583 #include "G4ProductionCutsTable.hh"
584 
585 void G4SteppingManager::ApplyProductionCut(G4Track* aSecondary)
586 {
587  G4bool tBelowCutEnergyAndSafety = false;
588  G4int tPtclIdx
589  = G4ProductionCuts::GetIndex(aSecondary->GetDefinition());
590  if (tPtclIdx<0) { return; }
591  G4ProductionCutsTable* tCutsTbl
593  G4int tCoupleIdx
594  = tCutsTbl->GetCoupleIndex(fPreStepPoint->GetMaterialCutsCouple());
595  G4double tProdThreshold
596  = (*(tCutsTbl->GetEnergyCutsVector(tPtclIdx)))[tCoupleIdx];
597  if( aSecondary->GetKineticEnergy()<tProdThreshold )
598  {
599  tBelowCutEnergyAndSafety = true;
600  if(std::abs(aSecondary->GetDynamicParticle()->GetCharge()) > DBL_MIN)
601  {
602  G4double currentRange
604  aSecondary->GetKineticEnergy(),
605  fPreStepPoint->GetMaterialCutsCouple());
606  tBelowCutEnergyAndSafety = (currentRange < CalculateSafety() );
607  }
608  }
609 
610  if( tBelowCutEnergyAndSafety )
611  {
612  if( !(aSecondary->IsGoodForTracking()) )
613  {
614 //G4cout << "!! Warning - G4SteppingManager:" << G4endl
615 //<< " This physics process generated a secondary"
616 //<< " of which energy is below cut but"
617 //<< " GoodForTracking is off !!!!!" << G4endl;
618  // Add kinetic energy to the total energy deposit
619  fStep->AddTotalEnergyDeposit(
620  aSecondary->GetKineticEnergy() );
621  aSecondary->SetKineticEnergy(0.0);
622  }
623  }
624 }
625