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

Last change on this file since 842 was 827, checked in by garnier, 16 years ago

import all except CVS

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