source: trunk/source/tracking/src/G4SteppingManager2.cc @ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 22.2 KB
Line 
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,v 1.38 2010/07/19 13:41:21 gcosmo Exp $
28// GEANT4 tag $Name: tracking-V09-03-08 $
29//
30//---------------------------------------------------------------
31//
32// G4SteppingManager2.cc
33//
34// Description:
35//   This class represents the manager who steers to move the give
36//   particle from the TrackingManger by one Step.
37//
38// Contact:
39//   Questions and comments to this code should be sent to
40//     Katsuya Amako  (e-mail: Katsuya.Amako@kek.jp)
41//     Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp)
42//
43//---------------------------------------------------------------
44
45//#define debug
46
47#include "G4UImanager.hh"
48#include "G4ForceCondition.hh"
49#include "G4GPILSelection.hh"
50#include "G4SteppingControl.hh"
51#include "G4TransportationManager.hh"
52#include "G4SteppingManager.hh"
53#include "G4LossTableManager.hh"
54
55/////////////////////////////////////////////////
56void G4SteppingManager::GetProcessNumber()
57/////////////////////////////////////////////////
58{
59#ifdef debug
60  G4cout<<"G4SteppingManager::GetProcessNumber: is called track="
61        <<fTrack<<G4endl;
62#endif
63
64  G4ProcessManager* pm= fTrack->GetDefinition()->GetProcessManager();
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
127/////////////////////////////////////////////////////////
128 void G4SteppingManager::DefinePhysicalStepLength()
129/////////////////////////////////////////////////////////
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             break;
191         case Forced:
192             (*fSelectedPostStepDoItVector)[np] = Forced;
193             break;
194         case StronglyForced:
195             (*fSelectedPostStepDoItVector)[np] = StronglyForced;
196             break;
197         default:
198             (*fSelectedPostStepDoItVector)[np] = InActivated;
199             break;
200     }
201       
202
203
204     if (fCondition==ExclusivelyForced) { 
205         for(size_t nrest=np+1; nrest < MAXofPostStepLoops; nrest++){ 
206             (*fSelectedPostStepDoItVector)[nrest] = InActivated; 
207         } 
208         return;  // Take note the 'return' at here !!!
209     } 
210     else{
211         if(physIntLength < PhysicalStep ){
212             PhysicalStep = physIntLength;
213             fStepStatus = fPostStepDoItProc;
214             fPostStepDoItProcTriggered = G4int(np);
215             fStep->GetPostStepPoint()
216                 ->SetProcessDefinedStep(fCurrentProcess);
217         }
218     }
219     
220
221   }
222
223   if (fPostStepDoItProcTriggered<MAXofPostStepLoops) {
224       if ((*fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] == 
225           InActivated) {
226           (*fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] = 
227               NotForced;
228       }
229   }
230
231// GPIL for AlongStep
232   proposedSafety = DBL_MAX;
233   G4double safetyProposedToAndByProcess = proposedSafety;
234
235   for(size_t kp=0; kp < MAXofAlongStepLoops; kp++){
236     fCurrentProcess = (*fAlongStepGetPhysIntVector)[kp];
237     if (fCurrentProcess== 0) continue;
238         // NULL means the process is inactivated by a user on fly.
239
240     physIntLength = fCurrentProcess->
241                     AlongStepGPIL( *fTrack, fPreviousStepSize,
242                                     PhysicalStep,
243                                     safetyProposedToAndByProcess,
244                                    &fGPILSelection );
245#ifdef G4VERBOSE
246                         // !!!!! Verbose
247           if(verboseLevel>0) fVerbose->DPSLAlongStep();
248#endif
249     if(physIntLength < PhysicalStep){
250       PhysicalStep = physIntLength;
251
252       // Check if the process wants to be the GPIL winner. For example,
253       // multi-scattering proposes Step limit, but won't be the winner.
254       if(fGPILSelection==CandidateForSelection){
255          fStepStatus = fAlongStepDoItProc;
256          fStep->GetPostStepPoint()
257               ->SetProcessDefinedStep(fCurrentProcess);
258       }
259
260          // Transportation is assumed to be the last process in the vector
261       if(kp == MAXofAlongStepLoops-1) {
262           if (fTrack->GetNextVolume() != 0)
263               fStepStatus = fGeomBoundary;
264           else
265               fStepStatus = fWorldBoundary;   
266       }
267     }
268
269     // Make sure to check the safety, even if Step is not limited
270     //  by this process.                      J. Apostolakis, June 20, 1998
271     //
272     if (safetyProposedToAndByProcess < proposedSafety)
273        // proposedSafety keeps the smallest value:
274        proposedSafety               = safetyProposedToAndByProcess;
275     else
276        // safetyProposedToAndByProcess always proposes a valid safety:
277        safetyProposedToAndByProcess = proposedSafety;
278     
279   } 
280} // void G4SteppingManager::DefinePhysicalStepLength() //
281
282
283//////////////////////////////////////////////////////
284void G4SteppingManager::InvokeAtRestDoItProcs()
285//////////////////////////////////////////////////////
286{
287// Select the rest process which has the shortest time before
288// it is invoked. In rest processes, GPIL()
289// returns the time before a process occurs.
290   G4double lifeTime, shortestLifeTime;
291
292   fAtRestDoItProcTriggered = 0;
293   shortestLifeTime = DBL_MAX;
294
295   unsigned int NofInactiveProc=0;
296   for( size_t ri=0 ; ri < MAXofAtRestLoops ; ri++ ){
297     fCurrentProcess = (*fAtRestGetPhysIntVector)[ri];
298     if (fCurrentProcess== 0) {
299       (*fSelectedAtRestDoItVector)[ri] = InActivated;
300       NofInactiveProc++;
301       continue;
302     }   // NULL means the process is inactivated by a user on fly.
303
304     lifeTime =
305       fCurrentProcess->AtRestGPIL( *fTrack, &fCondition );
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          (*fSelectedAtRestDoItVector)[fAtRestDoItProcTriggered] = NotForced;
316       }
317     }
318   }
319
320// at least one process is necessary to destroy the particle 
321// exit with warning
322   if(NofInactiveProc==MAXofAtRestLoops){ 
323     G4cerr << "ERROR - G4SteppingManager::InvokeAtRestDoItProcs()" << G4endl
324            << "        No AtRestDoIt process is active!" << G4endl;
325     // G4Exception("G4SteppingManager::InvokeAtRestDoItProcs", "Tracking0013",
326     //             FatalException, "No AtRestDoIt process is active." );
327   }
328
329   fStep->SetStepLength( 0. );  //the particle has stopped
330   fTrack->SetStepLength( 0. );
331
332// invoke selected process
333   for(size_t np=0; np < MAXofAtRestLoops; np++){
334   //
335   // Note: DoItVector has inverse order against GetPhysIntVector
336   //       and SelectedAtRestDoItVector.
337   //
338     if( (*fSelectedAtRestDoItVector)[MAXofAtRestLoops-np-1] != InActivated){
339
340       fCurrentProcess = (*fAtRestDoItVector)[np];
341       fParticleChange
342         = fCurrentProcess->AtRestDoIt( *fTrack, *fStep);
343                               
344       // Set the current process as a process which defined this Step length
345       fStep->GetPostStepPoint()
346              ->SetProcessDefinedStep(fCurrentProcess);
347
348       // Update Step 
349       fParticleChange->UpdateStepForAtRest(fStep);
350
351       // Now Store the secondaries from ParticleChange to SecondaryList
352       G4Track* tempSecondaryTrack;
353       G4int    num2ndaries;
354
355       num2ndaries = fParticleChange->GetNumberOfSecondaries();
356
357       for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
358         tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
359
360         if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
361         { ApplyProductionCut(tempSecondaryTrack); }
362
363         // Set parentID
364         tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
365
366         // Set the process pointer which created this track
367         tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
368         
369         // If this 2ndry particle has 'zero' kinetic energy, make sure
370         // it invokes a rest process at the beginning of the tracking
371         if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
372           G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
373           if (pm->GetAtRestProcessVector()->entries()>0){
374             tempSecondaryTrack->SetTrackStatus( fStopButAlive );
375             fSecondary->push_back( tempSecondaryTrack );
376             fN2ndariesAtRestDoIt++;
377           } else {
378             delete tempSecondaryTrack;
379           }
380         } else {
381           fSecondary->push_back( tempSecondaryTrack );
382           fN2ndariesAtRestDoIt++;
383         }     
384       } //end of loop on secondary
385
386
387       // clear ParticleChange
388       fParticleChange->Clear();
389
390     } //if(fSelectedAtRestDoItVector[np] != InActivated){
391   } //for(size_t np=0; np < MAXofAtRestLoops; np++){
392   fStep->UpdateTrack();
393
394   fTrack->SetTrackStatus( fStopAndKill );
395}
396
397
398/////////////////////////////////////////////////////////
399void G4SteppingManager::InvokeAlongStepDoItProcs()
400/////////////////////////////////////////////////////////
401{
402
403// If the current Step is defined by a 'ExclusivelyForced'
404// PostStepDoIt, then don't invoke any AlongStepDoIt
405   if(fStepStatus == fExclusivelyForcedProc){
406     return;               // Take note 'return' at here !!!
407   }
408
409// Invoke the all active continuous processes
410   for( size_t ci=0 ; ci<MAXofAlongStepLoops ; ci++ ){
411     fCurrentProcess = (*fAlongStepDoItVector)[ci];
412     if (fCurrentProcess== 0) continue;
413         // NULL means the process is inactivated by a user on fly.
414
415     fParticleChange
416       = fCurrentProcess->AlongStepDoIt( *fTrack, *fStep );
417
418     // Update the PostStepPoint of Step according to ParticleChange
419     fParticleChange->UpdateStepForAlongStep(fStep);
420#ifdef G4VERBOSE
421                         // !!!!! Verbose
422               if(verboseLevel>0) fVerbose->AlongStepDoItOneByOne();
423#endif
424
425     // Now Store the secondaries from ParticleChange to SecondaryList
426     G4Track* tempSecondaryTrack;
427     G4int    num2ndaries;
428
429     num2ndaries = fParticleChange->GetNumberOfSecondaries();
430
431     for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
432         tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
433
434         if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
435         { ApplyProductionCut(tempSecondaryTrack); }
436
437         // Set parentID
438         tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
439
440         // Set the process pointer which created this track
441         tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
442
443         // If this 2ndry particle has 'zero' kinetic energy, make sure
444         // it invokes a rest process at the beginning of the tracking
445         if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
446           G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
447           if (pm->GetAtRestProcessVector()->entries()>0){
448             tempSecondaryTrack->SetTrackStatus( fStopButAlive );
449             fSecondary->push_back( tempSecondaryTrack );
450             fN2ndariesAlongStepDoIt++;
451           } else {
452             delete tempSecondaryTrack;
453           }
454         } else {
455           fSecondary->push_back( tempSecondaryTrack );
456           fN2ndariesAlongStepDoIt++;
457         }
458     } //end of loop on secondary
459     
460     // Set the track status according to what the process defined
461     // if kinetic energy >0, otherwise set  fStopButAlive
462     fTrack->SetTrackStatus( fParticleChange->GetTrackStatus() );
463     
464     // clear ParticleChange
465     fParticleChange->Clear();
466   }
467
468   fStep->UpdateTrack();
469   G4TrackStatus fNewStatus = fTrack->GetTrackStatus();
470
471   if ( fNewStatus == fAlive && fTrack->GetKineticEnergy() <= DBL_MIN ) {
472     if(MAXofAtRestLoops>0) fNewStatus = fStopButAlive;
473     else                   fNewStatus = fStopAndKill;
474     fTrack->SetTrackStatus( fNewStatus );
475   }
476
477}
478
479////////////////////////////////////////////////////////
480void G4SteppingManager::InvokePostStepDoItProcs()
481////////////////////////////////////////////////////////
482{
483
484// Invoke the specified discrete processes
485   for(size_t np=0; np < MAXofPostStepLoops; np++){
486   //
487   // Note: DoItVector has inverse order against GetPhysIntVector
488   //       and SelectedPostStepDoItVector.
489   //
490     G4int Cond = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops-np-1];
491     if(Cond != InActivated){
492       if( ((Cond == NotForced) && (fStepStatus == fPostStepDoItProc)) ||
493           ((Cond == Forced) && (fStepStatus != fExclusivelyForcedProc)) ||
494           ((Cond == Conditionally) && (fStepStatus == fAlongStepDoItProc)) ||
495           ((Cond == ExclusivelyForced) && (fStepStatus == fExclusivelyForcedProc)) || 
496           ((Cond == StronglyForced) ) 
497          ) {
498
499         InvokePSDIP(np);
500       }
501     } //if(*fSelectedPostStepDoItVector(np)........
502
503     // Exit from PostStepLoop if the track has been killed,
504     // but extra treatment for processes with Strongly Forced flag
505     if(fTrack->GetTrackStatus() == fStopAndKill) {
506       for(size_t np1=np+1; np1 < MAXofPostStepLoops; np1++){ 
507         G4int Cond2 = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops-np1-1];
508         if (Cond2 == StronglyForced) {
509           InvokePSDIP(np1);
510         }
511       }
512       break;
513     }
514   } //for(size_t np=0; np < MAXofPostStepLoops; np++){
515}
516
517
518
519void G4SteppingManager::InvokePSDIP(size_t np)
520{
521         fCurrentProcess = (*fPostStepDoItVector)[np];
522         fParticleChange
523            = fCurrentProcess->PostStepDoIt( *fTrack, *fStep);
524
525         // Update PostStepPoint of Step according to ParticleChange
526         fParticleChange->UpdateStepForPostStep(fStep);
527#ifdef G4VERBOSE
528                 // !!!!! Verbose
529           if(verboseLevel>0) fVerbose->PostStepDoItOneByOne();
530#endif
531         // Update G4Track according to ParticleChange after each PostStepDoIt
532         fStep->UpdateTrack();
533
534         // Update safety after each invocation of PostStepDoIts
535         fStep->GetPostStepPoint()->SetSafety( CalculateSafety() );
536
537         // Now Store the secondaries from ParticleChange to SecondaryList
538         G4Track* tempSecondaryTrack;
539         G4int    num2ndaries;
540
541         num2ndaries = fParticleChange->GetNumberOfSecondaries();
542
543         for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
544            tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
545   
546            if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
547            { ApplyProductionCut(tempSecondaryTrack); }
548
549            // Set parentID
550            tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
551           
552            // Set the process pointer which created this track
553            tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
554
555            // If this 2ndry particle has 'zero' kinetic energy, make sure
556            // it invokes a rest process at the beginning of the tracking
557            if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
558              G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
559              if (pm->GetAtRestProcessVector()->entries()>0){
560                tempSecondaryTrack->SetTrackStatus( fStopButAlive );
561                fSecondary->push_back( tempSecondaryTrack );
562                fN2ndariesPostStepDoIt++;
563              } else {
564                delete tempSecondaryTrack;
565              }
566            } else {
567              fSecondary->push_back( tempSecondaryTrack );
568              fN2ndariesPostStepDoIt++;
569            }
570         } //end of loop on secondary
571
572         // Set the track status according to what the process defined
573         fTrack->SetTrackStatus( fParticleChange->GetTrackStatus() );
574
575         // clear ParticleChange
576         fParticleChange->Clear();
577}
578
579#include "G4EnergyLossTables.hh"
580#include "G4ProductionCuts.hh"
581#include "G4ProductionCutsTable.hh"
582
583void G4SteppingManager::ApplyProductionCut(G4Track* aSecondary)
584{
585  G4bool tBelowCutEnergyAndSafety = false;
586  G4int tPtclIdx
587    = G4ProductionCuts::GetIndex(aSecondary->GetDefinition());
588  if (tPtclIdx<0)  { return; }
589  G4ProductionCutsTable* tCutsTbl
590    = G4ProductionCutsTable::GetProductionCutsTable();
591  G4int tCoupleIdx
592    = tCutsTbl->GetCoupleIndex(fPreStepPoint->GetMaterialCutsCouple());
593  G4double tProdThreshold
594    = (*(tCutsTbl->GetEnergyCutsVector(tPtclIdx)))[tCoupleIdx];
595  if( aSecondary->GetKineticEnergy()<tProdThreshold )
596  {
597    tBelowCutEnergyAndSafety = true;
598    if(std::abs(aSecondary->GetDynamicParticle()->GetCharge()) > DBL_MIN)
599    {
600      G4double currentRange
601        = G4LossTableManager::Instance()->GetRange(aSecondary->GetDefinition(),
602                 aSecondary->GetKineticEnergy(),
603                 fPreStepPoint->GetMaterialCutsCouple());
604      tBelowCutEnergyAndSafety = (currentRange < CalculateSafety() );
605    }
606  }
607
608  if( tBelowCutEnergyAndSafety )
609  {
610    if( !(aSecondary->IsGoodForTracking()) )
611    {
612//G4cout << "!! Warning - G4SteppingManager:" << G4endl
613//<< " This physics process generated a secondary"
614//<< " of which energy is below cut but"
615//<< " GoodForTracking is off !!!!!" << G4endl;
616    // Add kinetic energy to the total energy deposit
617      fStep->AddTotalEnergyDeposit(
618      aSecondary->GetKineticEnergy() );
619      aSecondary->SetKineticEnergy(0.0);
620    }
621  }
622}
623
Note: See TracBrowser for help on using the repository browser.