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

Last change on this file since 1187 was 1011, checked in by garnier, 15 years ago

update

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.35 2008/12/05 22:15:04 asaim Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
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
279// Select the rest process which has the shortest time before
280// it is invoked. In rest processes, GPIL()
281// returns the time before a process occurs.
282   G4double lifeTime, shortestLifeTime;
283
284   fAtRestDoItProcTriggered = 0;
285   shortestLifeTime = DBL_MAX;
286
287   unsigned int NofInactiveProc=0;
288   for( size_t ri=0 ; ri < MAXofAtRestLoops ; ri++ ){
289     fCurrentProcess = (*fAtRestGetPhysIntVector)[ri];
290     if (fCurrentProcess== NULL) {
291       (*fSelectedAtRestDoItVector)[ri] = InActivated;
292       NofInactiveProc++;
293       continue;
294     }   // NULL means the process is inactivated by a user on fly.
295
296     lifeTime = 
297       fCurrentProcess->AtRestGPIL( 
298                                                     *fTrack,
299                                                &fCondition );
300     if(fCondition==Forced){
301       (*fSelectedAtRestDoItVector)[ri] = Forced;
302     }
303     else{
304       (*fSelectedAtRestDoItVector)[ri] = InActivated;
305      if(lifeTime < shortestLifeTime ){
306         shortestLifeTime = lifeTime;
307         fAtRestDoItProcTriggered =  G4int(int(ri)); 
308       }
309     }
310   }
311
312// at least one process is necessary to destory the particle 
313// exit with warning
314   if(NofInactiveProc==MAXofAtRestLoops){ 
315     //     G4Exception("G4SteppingManager::InvokeAtRestDoItProcs: No AtRestDoIt process is active. " );
316     G4cerr << "G4SteppingManager::InvokeAtRestDoItProcs: No AtRestDoIt process is active. " << G4endl;
317   } else {
318        (*fSelectedAtRestDoItVector)[fAtRestDoItProcTriggered] = NotForced;
319   }
320
321   fStep->SetStepLength( 0. );  //the particle has stopped
322   fTrack->SetStepLength( 0. );
323
324// invoke selected process
325   for(size_t np=0; np < MAXofAtRestLoops; np++){
326   //
327   // Note: DoItVector has inverse order against GetPhysIntVector
328   //       and SelectedAtRestDoItVector.
329   //
330     if( (*fSelectedAtRestDoItVector)[MAXofAtRestLoops-np-1] != InActivated){
331
332
333       fCurrentProcess = (*fAtRestDoItVector)[np];
334       fParticleChange
335         = fCurrentProcess->AtRestDoIt( *fTrack, *fStep);
336                               
337       // Set the current process as a process which defined this Step length
338       fStep->GetPostStepPoint()
339              ->SetProcessDefinedStep(fCurrentProcess);
340
341       // Update Step 
342       fParticleChange->UpdateStepForAtRest(fStep);
343
344       // Now Store the secondaries from ParticleChange to SecondaryList
345       G4Track* tempSecondaryTrack;
346       G4int    num2ndaries;
347
348       num2ndaries = fParticleChange->GetNumberOfSecondaries();
349
350       for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
351         tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
352
353         if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
354         { ApplyProductionCut(tempSecondaryTrack); }
355
356         // Set parentID
357         tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
358
359         // Set the process pointer which created this track
360         tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
361         
362         // If this 2ndry particle has 'zero' kinetic energy, make sure
363         // it invokes a rest process at the beginning of the tracking
364         if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
365           G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
366           if (pm->GetAtRestProcessVector()->entries()>0){
367             tempSecondaryTrack->SetTrackStatus( fStopButAlive );
368             fSecondary->push_back( tempSecondaryTrack );
369             fN2ndariesAtRestDoIt++;
370           } else {
371             delete tempSecondaryTrack;
372           }
373         } else {
374           fSecondary->push_back( tempSecondaryTrack );
375           fN2ndariesAtRestDoIt++;
376         }     
377       } //end of loop on secondary
378
379
380       // clear ParticleChange
381       fParticleChange->Clear();
382
383     } //if(fSelectedAtRestDoItVector[np] != InActivated){
384   } //for(size_t np=0; np < MAXofAtRestLoops; np++){
385
386   fStep->UpdateTrack();
387
388   fTrack->SetTrackStatus( fStopAndKill );
389}
390
391
392/////////////////////////////////////////////////////////
393void G4SteppingManager::InvokeAlongStepDoItProcs()
394/////////////////////////////////////////////////////////
395{
396
397// If the current Step is defined by a 'ExclusivelyForced'
398// PostStepDoIt, then don't invoke any AlongStepDoIt
399   if(fStepStatus == fExclusivelyForcedProc){
400     return;               // Take note 'return' at here !!!
401   }
402
403// Invoke the all active continuous processes
404   for( size_t ci=0 ; ci<MAXofAlongStepLoops ; ci++ ){
405     fCurrentProcess = (*fAlongStepDoItVector)[ci];
406     if (fCurrentProcess== NULL) continue;
407         // NULL means the process is inactivated by a user on fly.
408
409     fParticleChange
410       = fCurrentProcess->AlongStepDoIt( *fTrack, *fStep );
411
412     // Update the PostStepPoint of Step according to ParticleChange
413     fParticleChange->UpdateStepForAlongStep(fStep);
414#ifdef G4VERBOSE
415                         // !!!!! Verbose
416               if(verboseLevel>0) fVerbose->AlongStepDoItOneByOne();
417#endif
418
419     // Now Store the secondaries from ParticleChange to SecondaryList
420     G4Track* tempSecondaryTrack;
421     G4int    num2ndaries;
422
423     num2ndaries = fParticleChange->GetNumberOfSecondaries();
424
425     for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
426         tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
427
428         if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
429         { ApplyProductionCut(tempSecondaryTrack); }
430
431         // Set parentID
432         tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
433
434         // Set the process pointer which created this track
435         tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
436
437         // If this 2ndry particle has 'zero' kinetic energy, make sure
438         // it invokes a rest process at the beginning of the tracking
439         if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
440           G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
441           if (pm->GetAtRestProcessVector()->entries()>0){
442             tempSecondaryTrack->SetTrackStatus( fStopButAlive );
443             fSecondary->push_back( tempSecondaryTrack );
444             fN2ndariesAlongStepDoIt++;
445           } else {
446             delete tempSecondaryTrack;
447           }
448         } else {
449           fSecondary->push_back( tempSecondaryTrack );
450           fN2ndariesAlongStepDoIt++;
451         }
452     } //end of loop on secondary
453     
454     // Set the track status according to what the process defined
455     // if kinetic energy >0, otherwise set  fStopButAlive
456     fTrack->SetTrackStatus( fParticleChange->GetTrackStatus() );
457     
458     // clear ParticleChange
459     fParticleChange->Clear();
460   }
461
462   fStep->UpdateTrack();
463   G4TrackStatus fNewStatus = fTrack->GetTrackStatus();
464
465   if ( fNewStatus == fAlive && fTrack->GetKineticEnergy() <= DBL_MIN ) {
466     if(MAXofAtRestLoops>0) fNewStatus = fStopButAlive;
467     else                   fNewStatus = fStopAndKill;
468     fTrack->SetTrackStatus( fNewStatus );
469   }
470
471}
472
473////////////////////////////////////////////////////////
474void G4SteppingManager::InvokePostStepDoItProcs()
475////////////////////////////////////////////////////////
476{
477
478// Invoke the specified discrete processes
479   for(size_t np=0; np < MAXofPostStepLoops; np++){
480   //
481   // Note: DoItVector has inverse order against GetPhysIntVector
482   //       and SelectedPostStepDoItVector.
483   //
484     G4int Cond = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops-np-1];
485     if(Cond != InActivated){
486       if( ((Cond == NotForced) && (fStepStatus == fPostStepDoItProc)) ||
487           ((Cond == Forced) && (fStepStatus != fExclusivelyForcedProc)) ||
488           ((Cond == Conditionally) && (fStepStatus == fAlongStepDoItProc)) ||
489           ((Cond == ExclusivelyForced) && (fStepStatus == fExclusivelyForcedProc)) || 
490           ((Cond == StronglyForced) ) 
491          ) {
492
493         InvokePSDIP(np);
494       }
495     } //if(*fSelectedPostStepDoItVector(np)........
496
497     // Exit from PostStepLoop if the track has been killed,
498     // but extra treatment for processes with Strongly Forced flag
499     if(fTrack->GetTrackStatus() == fStopAndKill) {
500       for(size_t np1=np+1; np1 < MAXofPostStepLoops; np1++){ 
501         G4int Cond2 = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops-np1-1];
502         if (Cond2 == StronglyForced) {
503           InvokePSDIP(np1);
504         }
505       }
506       break;
507     }
508   } //for(size_t np=0; np < MAXofPostStepLoops; np++){
509}
510
511
512
513void G4SteppingManager::InvokePSDIP(size_t np)
514{
515         fCurrentProcess = (*fPostStepDoItVector)[np];
516         fParticleChange
517            = fCurrentProcess->PostStepDoIt( *fTrack, *fStep);
518
519         // Update PostStepPoint of Step according to ParticleChange
520         fParticleChange->UpdateStepForPostStep(fStep);
521#ifdef G4VERBOSE
522                 // !!!!! Verbose
523           if(verboseLevel>0) fVerbose->PostStepDoItOneByOne();
524#endif
525         // Update G4Track according to ParticleChange after each PostStepDoIt
526         fStep->UpdateTrack();
527
528         // Update safety after each invocation of PostStepDoIts
529         fStep->GetPostStepPoint()->SetSafety( CalculateSafety() );
530
531         // Now Store the secondaries from ParticleChange to SecondaryList
532         G4Track* tempSecondaryTrack;
533         G4int    num2ndaries;
534
535         num2ndaries = fParticleChange->GetNumberOfSecondaries();
536
537         for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
538            tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
539   
540            if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
541            { ApplyProductionCut(tempSecondaryTrack); }
542
543            // Set parentID
544            tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
545           
546            // Set the process pointer which created this track
547            tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
548
549            // If this 2ndry particle has 'zero' kinetic energy, make sure
550            // it invokes a rest process at the beginning of the tracking
551            if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
552              G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
553              if (pm->GetAtRestProcessVector()->entries()>0){
554                tempSecondaryTrack->SetTrackStatus( fStopButAlive );
555                fSecondary->push_back( tempSecondaryTrack );
556                fN2ndariesPostStepDoIt++;
557              } else {
558                delete tempSecondaryTrack;
559              }
560            } else {
561              fSecondary->push_back( tempSecondaryTrack );
562              fN2ndariesPostStepDoIt++;
563            }
564         } //end of loop on secondary
565
566         // Set the track status according to what the process defined
567         fTrack->SetTrackStatus( fParticleChange->GetTrackStatus() );
568
569         // clear ParticleChange
570         fParticleChange->Clear();
571}
572
573#include "G4EnergyLossTables.hh"
574#include "G4ProductionCuts.hh"
575#include "G4ProductionCutsTable.hh"
576
577void G4SteppingManager::ApplyProductionCut(G4Track* aSecondary)
578{
579  G4bool tBelowCutEnergyAndSafety = false;
580  G4int tPtclIdx
581    = G4ProductionCuts::GetIndex(aSecondary->GetDefinition());
582  G4ProductionCutsTable* tCutsTbl
583    = G4ProductionCutsTable::GetProductionCutsTable();
584  G4int tCoupleIdx
585    = tCutsTbl->GetCoupleIndex(fPreStepPoint->GetMaterialCutsCouple());
586  G4double tProdThreshold
587    = (*(tCutsTbl->GetEnergyCutsVector(tPtclIdx)))[tCoupleIdx];
588  if( aSecondary->GetKineticEnergy()<tProdThreshold )
589  {
590    tBelowCutEnergyAndSafety = true;
591    if(std::abs(aSecondary->GetDynamicParticle()->GetCharge()) > DBL_MIN)
592    {
593      G4double currentRange
594        = G4LossTableManager::Instance()->GetRange(aSecondary->GetDefinition(),
595                 aSecondary->GetKineticEnergy(),
596                 fPreStepPoint->GetMaterialCutsCouple());
597      tBelowCutEnergyAndSafety = (currentRange < CalculateSafety() );
598    }
599  }
600
601  if( tBelowCutEnergyAndSafety )
602  {
603    if( !(aSecondary->IsGoodForTracking()) )
604    {
605//G4cout << "!! Warning - G4SteppingManager:" << G4endl
606//<< " This physics process generated a secondary"
607//<< " of which energy is below cut but"
608//<< " GoodForTracking is off !!!!!" << G4endl;
609    // Add kinetic energy to the total energy deposit
610      fStep->AddTotalEnergyDeposit(
611      aSecondary->GetKineticEnergy() );
612      aSecondary->SetKineticEnergy(0.0);
613    }
614  }
615}
616
617
Note: See TracBrowser for help on using the repository browser.