Geant4  10.03.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 73605 2013-09-02 09:34:55Z 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 #include "G4ParticleTable.hh"
54 
58 {
59 #ifdef debug
60  G4cout<<"G4SteppingManager::GetProcessNumber: is called track="
61  <<fTrack<<G4endl;
62 #endif
63 
65  if(!pm)
66  {
67  G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
68  << " ProcessManager is NULL for particle = "
69  << fTrack->GetDefinition()->GetParticleName() << ", PDG_code = "
70  << fTrack->GetDefinition()->GetPDGEncoding() << G4endl;
71  G4Exception("G4SteppingManager::GetProcessNumber()", "Tracking0011",
72  FatalException, "Process Manager is not found.");
73  return;
74  }
75 
76 // AtRestDoits
77  MAXofAtRestLoops = pm->GetAtRestProcessVector()->entries();
78  fAtRestDoItVector = pm->GetAtRestProcessVector(typeDoIt);
79  fAtRestGetPhysIntVector = pm->GetAtRestProcessVector(typeGPIL);
80 #ifdef debug
81  G4cout << "G4SteppingManager::GetProcessNumber: #ofAtRest="
82  << MAXofAtRestLoops << G4endl;
83 #endif
84 
85 // AlongStepDoits
86  MAXofAlongStepLoops = pm->GetAlongStepProcessVector()->entries();
87  fAlongStepDoItVector = pm->GetAlongStepProcessVector(typeDoIt);
88  fAlongStepGetPhysIntVector = pm->GetAlongStepProcessVector(typeGPIL);
89 #ifdef debug
90  G4cout << "G4SteppingManager::GetProcessNumber:#ofAlongStp="
91  << MAXofAlongStepLoops << G4endl;
92 #endif
93 
94 // PostStepDoits
95  MAXofPostStepLoops = pm->GetPostStepProcessVector()->entries();
96  fPostStepDoItVector = pm->GetPostStepProcessVector(typeDoIt);
97  fPostStepGetPhysIntVector = pm->GetPostStepProcessVector(typeGPIL);
98 #ifdef debug
99  G4cout << "G4SteppingManager::GetProcessNumber: #ofPostStep="
100  << MAXofPostStepLoops << G4endl;
101 #endif
102 
103  if (SizeOfSelectedDoItVector<MAXofAtRestLoops ||
104  SizeOfSelectedDoItVector<MAXofAlongStepLoops ||
105  SizeOfSelectedDoItVector<MAXofPostStepLoops )
106  {
107  G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
108  << " SizeOfSelectedDoItVector= " << SizeOfSelectedDoItVector
109  << " ; is smaller then one of MAXofAtRestLoops= "
110  << MAXofAtRestLoops << G4endl
111  << " or MAXofAlongStepLoops= " << MAXofAlongStepLoops
112  << " or MAXofPostStepLoops= " << MAXofPostStepLoops << G4endl;
113  G4Exception("G4SteppingManager::GetProcessNumber()",
114  "Tracking0012", FatalException,
115  "The array size is smaller than the actual No of processes.");
116  }
117 }
118 
119 
120 // ************************************************************************
121 //
122 // Private Member Functions
123 //
124 // ************************************************************************
125 
126 
128  void G4SteppingManager::DefinePhysicalStepLength()
130 {
131 
132 // ReSet the counter etc.
133  PhysicalStep = DBL_MAX; // Initialize by a huge number
134  physIntLength = DBL_MAX; // Initialize by a huge number
135 #ifdef G4VERBOSE
136  // !!!!! Verbose
137  if(verboseLevel>0) fVerbose->DPSLStarted();
138 #endif
139 
140 // Obtain the user defined maximum allowed Step in the volume
141 // 1997.12.13 adds argument for GetMaxAllowedStep by K.Kurashige
142 // 2004.01.20 This block will be removed by Geant4 7.0
143 // G4UserLimits* ul= fCurrentVolume->GetLogicalVolume()->GetUserLimits();
144 // if (ul) {
145 // physIntLength = ul->GetMaxAllowedStep(*fTrack);
146 //#ifdef G4VERBOSE
147 // // !!!!! Verbose
148 // if(verboseLevel>0) fVerbose->DPSLUserLimit();
149 //#endif
150 // }
151 //
152 // if(physIntLength < PhysicalStep ){
153 // PhysicalStep = physIntLength;
154 // fStepStatus = fUserDefinedLimit;
155 // fStep->GetPostStepPoint()
156 // ->SetProcessDefinedStep(0);
157 // // Take note that the process pointer is 'NULL' if the Step
158 // // is defined by the user defined limit.
159 // }
160 // 2004.01.20 This block will be removed by Geant4 7.0
161 
162 // GPIL for PostStep
163  fPostStepDoItProcTriggered = MAXofPostStepLoops;
164 
165  for(size_t np=0; np < MAXofPostStepLoops; np++){
166  fCurrentProcess = (*fPostStepGetPhysIntVector)(np);
167  if (fCurrentProcess== 0) {
168  (*fSelectedPostStepDoItVector)[np] = InActivated;
169  continue;
170  } // NULL means the process is inactivated by a user on fly.
171 
172  physIntLength = fCurrentProcess->
173  PostStepGPIL( *fTrack,
174  fPreviousStepSize,
175  &fCondition );
176 #ifdef G4VERBOSE
177  // !!!!! Verbose
178  if(verboseLevel>0) fVerbose->DPSLPostStep();
179 #endif
180 
181  switch (fCondition) {
182  case ExclusivelyForced:
183  (*fSelectedPostStepDoItVector)[np] = ExclusivelyForced;
184  fStepStatus = fExclusivelyForcedProc;
185  fStep->GetPostStepPoint()
186  ->SetProcessDefinedStep(fCurrentProcess);
187  break;
188  case Conditionally:
189  // (*fSelectedPostStepDoItVector)[np] = Conditionally;
190  G4Exception("G4SteppingManager::DefinePhysicalStepLength()", "Tracking1001", FatalException, "This feature no more supported");
191 
192  break;
193  case Forced:
194  (*fSelectedPostStepDoItVector)[np] = Forced;
195  break;
196  case StronglyForced:
197  (*fSelectedPostStepDoItVector)[np] = StronglyForced;
198  break;
199  default:
200  (*fSelectedPostStepDoItVector)[np] = InActivated;
201  break;
202  }
203 
204 
205 
206  if (fCondition==ExclusivelyForced) {
207  for(size_t nrest=np+1; nrest < MAXofPostStepLoops; nrest++){
208  (*fSelectedPostStepDoItVector)[nrest] = InActivated;
209  }
210  return; // Take note the 'return' at here !!!
211  }
212  else{
213  if(physIntLength < PhysicalStep ){
214  PhysicalStep = physIntLength;
215  fStepStatus = fPostStepDoItProc;
216  fPostStepDoItProcTriggered = G4int(np);
217  fStep->GetPostStepPoint()
218  ->SetProcessDefinedStep(fCurrentProcess);
219  }
220  }
221 
222 
223  }
224 
225  if (fPostStepDoItProcTriggered<MAXofPostStepLoops) {
226  if ((*fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] ==
227  InActivated) {
228  (*fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] =
229  NotForced;
230  }
231  }
232 
233 // GPIL for AlongStep
234  proposedSafety = DBL_MAX;
235  G4double safetyProposedToAndByProcess = proposedSafety;
236 
237  for(size_t kp=0; kp < MAXofAlongStepLoops; kp++){
238  fCurrentProcess = (*fAlongStepGetPhysIntVector)[kp];
239  if (fCurrentProcess== 0) continue;
240  // NULL means the process is inactivated by a user on fly.
241 
242  physIntLength = fCurrentProcess->
243  AlongStepGPIL( *fTrack, fPreviousStepSize,
244  PhysicalStep,
245  safetyProposedToAndByProcess,
246  &fGPILSelection );
247 #ifdef G4VERBOSE
248  // !!!!! Verbose
249  if(verboseLevel>0) fVerbose->DPSLAlongStep();
250 #endif
251  if(physIntLength < PhysicalStep){
252  PhysicalStep = physIntLength;
253 
254  // Check if the process wants to be the GPIL winner. For example,
255  // multi-scattering proposes Step limit, but won't be the winner.
256  if(fGPILSelection==CandidateForSelection){
257  fStepStatus = fAlongStepDoItProc;
258  fStep->GetPostStepPoint()
259  ->SetProcessDefinedStep(fCurrentProcess);
260  }
261 
262  // Transportation is assumed to be the last process in the vector
263  if(kp == MAXofAlongStepLoops-1)
264  fStepStatus = fGeomBoundary;
265  }
266 
267  // Make sure to check the safety, even if Step is not limited
268  // by this process. J. Apostolakis, June 20, 1998
269  //
270  if (safetyProposedToAndByProcess < proposedSafety)
271  // proposedSafety keeps the smallest value:
272  proposedSafety = safetyProposedToAndByProcess;
273  else
274  // safetyProposedToAndByProcess always proposes a valid safety:
275  safetyProposedToAndByProcess = proposedSafety;
276 
277  }
278 } // void G4SteppingManager::DefinePhysicalStepLength() //
279 
280 
282 void G4SteppingManager::InvokeAtRestDoItProcs()
284 {
285 // Select the rest process which has the shortest time before
286 // it is invoked. In rest processes, GPIL()
287 // returns the time before a process occurs.
288  G4double lifeTime, shortestLifeTime;
289 
290  fAtRestDoItProcTriggered = 0;
291  shortestLifeTime = DBL_MAX;
292 
293  unsigned int NofInactiveProc=0;
294  for( size_t ri=0 ; ri < MAXofAtRestLoops ; ri++ ){
295  fCurrentProcess = (*fAtRestGetPhysIntVector)[ri];
296  if (fCurrentProcess== 0) {
297  (*fSelectedAtRestDoItVector)[ri] = InActivated;
298  NofInactiveProc++;
299  continue;
300  } // NULL means the process is inactivated by a user on fly.
301 
302  lifeTime =
303  fCurrentProcess->AtRestGPIL( *fTrack, &fCondition );
304 
305 
306 
307  if(fCondition==Forced && fCurrentProcess){
308  (*fSelectedAtRestDoItVector)[ri] = Forced;
309  }
310  else{
311  (*fSelectedAtRestDoItVector)[ri] = InActivated;
312  if(lifeTime < shortestLifeTime ){
313  shortestLifeTime = lifeTime;
314  fAtRestDoItProcTriggered = G4int(ri);
315  }
316  }
317  }
318 
319  (*fSelectedAtRestDoItVector)[fAtRestDoItProcTriggered] = NotForced;
320 
321 // at least one process is necessary to destroy the particle
322 // exit with warning
323  if(NofInactiveProc==MAXofAtRestLoops){
324  G4Exception("G4SteppingManager::InvokeAtRestDoItProcs()", "Tracking0013",
325  JustWarning, "No AtRestDoIt process is active!" );
326  }
327 
328  fStep->SetStepLength( 0. ); //the particle has stopped
329  fTrack->SetStepLength( 0. );
330 
331 // invoke selected process
332  for(size_t np=0; np < MAXofAtRestLoops; np++){
333  //
334  // Note: DoItVector has inverse order against GetPhysIntVector
335  // and SelectedAtRestDoItVector.
336  //
337  if( (*fSelectedAtRestDoItVector)[MAXofAtRestLoops-np-1] != InActivated){
338 
339  fCurrentProcess = (*fAtRestDoItVector)[np];
340  fParticleChange
341  = fCurrentProcess->AtRestDoIt( *fTrack, *fStep);
342 
343  // Set the current process as a process which defined this Step length
344  fStep->GetPostStepPoint()
345  ->SetProcessDefinedStep(fCurrentProcess);
346 
347  // Update Step
348  fParticleChange->UpdateStepForAtRest(fStep);
349 
350  // Now Store the secondaries from ParticleChange to SecondaryList
351  G4Track* tempSecondaryTrack;
352  G4int num2ndaries;
353 
354  num2ndaries = fParticleChange->GetNumberOfSecondaries();
355 
356  for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
357  tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
358 
359  if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
360  { ApplyProductionCut(tempSecondaryTrack); }
361 
362  // Set parentID
363  tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
364 
365  // Set the process pointer which created this track
366  tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
367 
368  // If this 2ndry particle has 'zero' kinetic energy, make sure
369  // it invokes a rest process at the beginning of the tracking
370  if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
371  G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
372  if(!pm && tempSecondaryTrack->GetDefinition()->IsGeneralIon())
374  if (pm->GetAtRestProcessVector()->entries()>0){
375  tempSecondaryTrack->SetTrackStatus( fStopButAlive );
376  fSecondary->push_back( tempSecondaryTrack );
377  fN2ndariesAtRestDoIt++;
378  } else {
379  delete tempSecondaryTrack;
380  }
381  } else {
382  fSecondary->push_back( tempSecondaryTrack );
383  fN2ndariesAtRestDoIt++;
384  }
385  } //end of loop on secondary
386 
387 
388  // clear ParticleChange
389  fParticleChange->Clear();
390 
391  } //if(fSelectedAtRestDoItVector[np] != InActivated){
392  } //for(size_t np=0; np < MAXofAtRestLoops; np++){
393  fStep->UpdateTrack();
394 
395  fTrack->SetTrackStatus( fStopAndKill );
396 }
397 
398 
400 void G4SteppingManager::InvokeAlongStepDoItProcs()
402 {
403 
404 // If the current Step is defined by a 'ExclusivelyForced'
405 // PostStepDoIt, then don't invoke any AlongStepDoIt
406  if(fStepStatus == fExclusivelyForcedProc){
407  return; // Take note 'return' at here !!!
408  }
409 
410 // Invoke the all active continuous processes
411  for( size_t ci=0 ; ci<MAXofAlongStepLoops ; ci++ ){
412  fCurrentProcess = (*fAlongStepDoItVector)[ci];
413  if (fCurrentProcess== 0) continue;
414  // NULL means the process is inactivated by a user on fly.
415 
416  fParticleChange
417  = fCurrentProcess->AlongStepDoIt( *fTrack, *fStep );
418 
419  // Update the PostStepPoint of Step according to ParticleChange
420  fParticleChange->UpdateStepForAlongStep(fStep);
421 #ifdef G4VERBOSE
422  // !!!!! Verbose
423  if(verboseLevel>0) fVerbose->AlongStepDoItOneByOne();
424 #endif
425 
426  // Now Store the secondaries from ParticleChange to SecondaryList
427  G4Track* tempSecondaryTrack;
428  G4int num2ndaries;
429 
430  num2ndaries = fParticleChange->GetNumberOfSecondaries();
431 
432  for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
433  tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
434 
435  if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
436  { ApplyProductionCut(tempSecondaryTrack); }
437 
438  // Set parentID
439  tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
440 
441  // Set the process pointer which created this track
442  tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
443 
444  // If this 2ndry particle has 'zero' kinetic energy, make sure
445  // it invokes a rest process at the beginning of the tracking
446  if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
447  G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
448  if(!pm && tempSecondaryTrack->GetDefinition()->IsGeneralIon())
450  if (pm->GetAtRestProcessVector()->entries()>0){
451  tempSecondaryTrack->SetTrackStatus( fStopButAlive );
452  fSecondary->push_back( tempSecondaryTrack );
453  fN2ndariesAlongStepDoIt++;
454  } else {
455  delete tempSecondaryTrack;
456  }
457  } else {
458  fSecondary->push_back( tempSecondaryTrack );
459  fN2ndariesAlongStepDoIt++;
460  }
461  } //end of loop on secondary
462 
463  // Set the track status according to what the process defined
464  // if kinetic energy >0, otherwise set fStopButAlive
465  fTrack->SetTrackStatus( fParticleChange->GetTrackStatus() );
466 
467  // clear ParticleChange
468  fParticleChange->Clear();
469  }
470 
471  fStep->UpdateTrack();
472  G4TrackStatus fNewStatus = fTrack->GetTrackStatus();
473 
474  if ( fNewStatus == fAlive && fTrack->GetKineticEnergy() <= DBL_MIN ) {
475  if(MAXofAtRestLoops>0) fNewStatus = fStopButAlive;
476  else fNewStatus = fStopAndKill;
477  fTrack->SetTrackStatus( fNewStatus );
478  }
479 
480 }
481 
483 void G4SteppingManager::InvokePostStepDoItProcs()
485 {
486 
487 // Invoke the specified discrete processes
488  for(size_t np=0; np < MAXofPostStepLoops; np++){
489  //
490  // Note: DoItVector has inverse order against GetPhysIntVector
491  // and SelectedPostStepDoItVector.
492  //
493  G4int Cond = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops-np-1];
494  if(Cond != InActivated){
495  if( ((Cond == NotForced) && (fStepStatus == fPostStepDoItProc)) ||
496  ((Cond == Forced) && (fStepStatus != fExclusivelyForcedProc)) ||
497  // ((Cond == Conditionally) && (fStepStatus == fAlongStepDoItProc)) ||
498  ((Cond == ExclusivelyForced) && (fStepStatus == fExclusivelyForcedProc)) ||
499  ((Cond == StronglyForced) )
500  ) {
501 
502  InvokePSDIP(np);
503  if ((np==0) && (fTrack->GetNextVolume() == 0)){
504  fStepStatus = fWorldBoundary;
505  fStep->GetPostStepPoint()->SetStepStatus( fStepStatus );
506  }
507  }
508  } //if(*fSelectedPostStepDoItVector(np)........
509 
510  // Exit from PostStepLoop if the track has been killed,
511  // but extra treatment for processes with Strongly Forced flag
512  if(fTrack->GetTrackStatus() == fStopAndKill) {
513  for(size_t np1=np+1; np1 < MAXofPostStepLoops; np1++){
514  G4int Cond2 = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops-np1-1];
515  if (Cond2 == StronglyForced) {
516  InvokePSDIP(np1);
517  }
518  }
519  break;
520  }
521  } //for(size_t np=0; np < MAXofPostStepLoops; np++){
522 }
523 
524 
525 
526 void G4SteppingManager::InvokePSDIP(size_t np)
527 {
528  fCurrentProcess = (*fPostStepDoItVector)[np];
529  fParticleChange
530  = fCurrentProcess->PostStepDoIt( *fTrack, *fStep);
531 
532  // Update PostStepPoint of Step according to ParticleChange
533  fParticleChange->UpdateStepForPostStep(fStep);
534 #ifdef G4VERBOSE
535  // !!!!! Verbose
536  if(verboseLevel>0) fVerbose->PostStepDoItOneByOne();
537 #endif
538  // Update G4Track according to ParticleChange after each PostStepDoIt
539  fStep->UpdateTrack();
540 
541  // Update safety after each invocation of PostStepDoIts
542  fStep->GetPostStepPoint()->SetSafety( CalculateSafety() );
543 
544  // Now Store the secondaries from ParticleChange to SecondaryList
545  G4Track* tempSecondaryTrack;
546  G4int num2ndaries;
547 
548  num2ndaries = fParticleChange->GetNumberOfSecondaries();
549 
550  for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
551  tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
552 
553  if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
554  { ApplyProductionCut(tempSecondaryTrack); }
555 
556  // Set parentID
557  tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
558 
559  // Set the process pointer which created this track
560  tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
561 
562  // If this 2ndry particle has 'zero' kinetic energy, make sure
563  // it invokes a rest process at the beginning of the tracking
564  if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
565  G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
566  if(!pm && tempSecondaryTrack->GetDefinition()->IsGeneralIon())
568  if (pm->GetAtRestProcessVector()->entries()>0){
569  tempSecondaryTrack->SetTrackStatus( fStopButAlive );
570  fSecondary->push_back( tempSecondaryTrack );
571  fN2ndariesPostStepDoIt++;
572  } else {
573  delete tempSecondaryTrack;
574  }
575  } else {
576  fSecondary->push_back( tempSecondaryTrack );
577  fN2ndariesPostStepDoIt++;
578  }
579  } //end of loop on secondary
580 
581  // Set the track status according to what the process defined
582  fTrack->SetTrackStatus( fParticleChange->GetTrackStatus() );
583 
584  // clear ParticleChange
585  fParticleChange->Clear();
586 }
587 
588 #include "G4EnergyLossTables.hh"
589 #include "G4ProductionCuts.hh"
590 #include "G4ProductionCutsTable.hh"
591 
592 void G4SteppingManager::ApplyProductionCut(G4Track* aSecondary)
593 {
594  G4bool tBelowCutEnergyAndSafety = false;
595  G4int tPtclIdx
596  = G4ProductionCuts::GetIndex(aSecondary->GetDefinition());
597  if (tPtclIdx<0) { return; }
598  G4ProductionCutsTable* tCutsTbl
600  G4int tCoupleIdx
601  = tCutsTbl->GetCoupleIndex(fPreStepPoint->GetMaterialCutsCouple());
602  G4double tProdThreshold
603  = (*(tCutsTbl->GetEnergyCutsVector(tPtclIdx)))[tCoupleIdx];
604  if( aSecondary->GetKineticEnergy()<tProdThreshold )
605  {
606  tBelowCutEnergyAndSafety = true;
607  if(std::abs(aSecondary->GetDynamicParticle()->GetCharge()) > DBL_MIN)
608  {
609  G4double currentRange
611  aSecondary->GetKineticEnergy(),
612  fPreStepPoint->GetMaterialCutsCouple());
613  tBelowCutEnergyAndSafety = (currentRange < CalculateSafety() );
614  }
615  }
616 
617  if( tBelowCutEnergyAndSafety )
618  {
619  if( !(aSecondary->IsGoodForTracking()) )
620  {
621 //G4cout << "!! Warning - G4SteppingManager:" << G4endl
622 //<< " This physics process generated a secondary"
623 //<< " of which energy is below cut but"
624 //<< " GoodForTracking is off !!!!!" << G4endl;
625  // Add kinetic energy to the total energy deposit
626  fStep->AddTotalEnergyDeposit(
627  aSecondary->GetKineticEnergy() );
628  aSecondary->SetKineticEnergy(0.0);
629  }
630  }
631 }
632 
void SetTrackStatus(const G4TrackStatus aTrackStatus)
G4ParticleDefinition * GetDefinition() const
void SetStepLength(G4double value)
G4int GetNumberOfSecondaries() const
virtual G4Step * UpdateStepForAlongStep(G4Step *Step)
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
static G4LossTableManager * Instance()
G4Track * GetSecondary(G4int anIndex) const
static G4int GetIndex(const G4String &name)
void SetProcessDefinedStep(const G4VProcess *aValue)
virtual void PostStepDoItOneByOne()=0
const G4DynamicParticle * GetDynamicParticle() const
G4TrackStatus GetTrackStatus() const
G4ParticleDefinition * GetGenericIon() const
int G4int
Definition: G4Types.hh:78
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void SetStepStatus(const G4StepStatus aValue)
const G4String & GetParticleName() const
G4bool IsGeneralIon() const
void SetCreatorProcess(const G4VProcess *aValue)
G4VPhysicalVolume * GetNextVolume() const
G4double AtRestGPIL(const G4Track &track, G4ForceCondition *condition)
Definition: G4VProcess.hh:490
void UpdateTrack()
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4double GetKineticEnergy() const
G4int entries() const
virtual void DPSLAlongStep()=0
G4GLOB_DLL std::ostream G4cout
virtual void DPSLPostStep()=0
bool G4bool
Definition: G4Types.hh:79
void SetStepLength(G4double value)
G4int GetCoupleIndex(const G4MaterialCutsCouple *aCouple) const
G4double GetCharge() const
G4int GetTrackID() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4ProductionCutsTable * GetProductionCutsTable()
static G4ParticleTable * GetParticleTable()
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)=0
G4ProcessManager * GetProcessManager() const
virtual void AlongStepDoItOneByOne()=0
void SetParentID(const G4int aValue)
G4StepPoint * GetPostStepPoint() const
void AddTotalEnergyDeposit(G4double value)
#define DBL_MIN
Definition: templates.hh:75
G4bool IsGoodForTracking() const
#define G4endl
Definition: G4ios.hh:61
G4bool GetApplyCutsFlag() const
G4TrackStatus GetTrackStatus() const
void SetKineticEnergy(const G4double aValue)
double G4double
Definition: G4Types.hh:76
virtual G4Step * UpdateStepForPostStep(G4Step *Step)
void SetSafety(const G4double aValue)
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0
#define DBL_MAX
Definition: templates.hh:83
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4TrackStatus
virtual G4Step * UpdateStepForAtRest(G4Step *Step)
virtual void DPSLStarted()=0
static const size_t SizeOfSelectedDoItVector
G4GLOB_DLL std::ostream G4cerr
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)=0
G4ProcessVector * GetPostStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const