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

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

tag geant4.9.4 beta 1 + modifs locales

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