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