Geant4 10.7.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// G4SteppingManager class implementation - part 2
27//
28// Contact:
29// Questions and comments to this code should be sent to
30// Katsuya Amako (e-mail: [email protected])
31// Takashi Sasaki (e-mail: [email protected])
32// --------------------------------------------------------------------
33
34// #define debug
35
36#include "G4UImanager.hh"
37#include "G4ForceCondition.hh"
38#include "G4GPILSelection.hh"
39#include "G4SteppingControl.hh"
41#include "G4SteppingManager.hh"
42#include "G4LossTableManager.hh"
43#include "G4ParticleTable.hh"
44
45#include "G4EnergyLossTables.hh"
46#include "G4ProductionCuts.hh"
48
49/////////////////////////////////////////////////
51/////////////////////////////////////////////////
52{
53 #ifdef debug
54 G4cout << "G4SteppingManager::GetProcessNumber: is called track="
55 << fTrack << G4endl;
56 #endif
57
59 if(pm == nullptr)
60 {
61 G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
62 << " ProcessManager is NULL for particle = "
63 << fTrack->GetDefinition()->GetParticleName() << ", PDG_code = "
64 << fTrack->GetDefinition()->GetPDGEncoding() << G4endl;
65 G4Exception("G4SteppingManager::GetProcessNumber()", "Tracking0011",
66 FatalException, "Process Manager is not found.");
67 return;
68 }
69
70 // AtRestDoits
71 //
72 MAXofAtRestLoops = pm->GetAtRestProcessVector()->entries();
73 fAtRestDoItVector = pm->GetAtRestProcessVector(typeDoIt);
74 fAtRestGetPhysIntVector = pm->GetAtRestProcessVector(typeGPIL);
75
76 #ifdef debug
77 G4cout << "G4SteppingManager::GetProcessNumber: #ofAtRest="
78 << MAXofAtRestLoops << G4endl;
79 #endif
80
81 // AlongStepDoits
82 //
83 MAXofAlongStepLoops = pm->GetAlongStepProcessVector()->entries();
84 fAlongStepDoItVector = pm->GetAlongStepProcessVector(typeDoIt);
85 fAlongStepGetPhysIntVector = pm->GetAlongStepProcessVector(typeGPIL);
86
87 #ifdef debug
88 G4cout << "G4SteppingManager::GetProcessNumber:#ofAlongStp="
89 << MAXofAlongStepLoops << G4endl;
90 #endif
91
92 // PostStepDoits
93 //
94 MAXofPostStepLoops = pm->GetPostStepProcessVector()->entries();
95 fPostStepDoItVector = pm->GetPostStepProcessVector(typeDoIt);
96 fPostStepGetPhysIntVector = pm->GetPostStepProcessVector(typeGPIL);
97
98 #ifdef debug
99 G4cout << "G4SteppingManager::GetProcessNumber: #ofPostStep="
100 << MAXofPostStepLoops << G4endl;
101 #endif
102
103 if (SizeOfSelectedDoItVector<MAXofAtRestLoops ||
104 SizeOfSelectedDoItVector<MAXofAlongStepLoops ||
105 SizeOfSelectedDoItVector<MAXofPostStepLoops )
106 {
107 G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
108 << " SizeOfSelectedDoItVector= " << SizeOfSelectedDoItVector
109 << " ; is smaller then one of MAXofAtRestLoops= "
110 << MAXofAtRestLoops << G4endl
111 << " or MAXofAlongStepLoops= " << MAXofAlongStepLoops
112 << " or MAXofPostStepLoops= " << MAXofPostStepLoops << G4endl;
113 G4Exception("G4SteppingManager::GetProcessNumber()",
114 "Tracking0012", FatalException,
115 "The array size is smaller than the actual No of processes.");
116 }
117}
118
119// ************************************************************************
120//
121// Private Member Functions
122//
123// ************************************************************************
124
125/////////////////////////////////////////////////////////
126void G4SteppingManager::DefinePhysicalStepLength()
127/////////////////////////////////////////////////////////
128{
129 // ReSet the counter etc.
130 //
131 PhysicalStep = DBL_MAX; // Initialize by a huge number
132 physIntLength = DBL_MAX; // Initialize by a huge number
133
134 #ifdef G4VERBOSE
135 if(verboseLevel>0) fVerbose->DPSLStarted();
136 #endif
137
138 // GPIL for PostStep
139 //
140 fPostStepDoItProcTriggered = MAXofPostStepLoops;
141
142 for(std::size_t np=0; np<MAXofPostStepLoops; ++np)
143 {
144 fCurrentProcess = (*fPostStepGetPhysIntVector)(np);
145 if (fCurrentProcess == nullptr)
146 {
147 (*fSelectedPostStepDoItVector)[np] = InActivated;
148 continue;
149 } // NULL means the process is inactivated by a user on fly
150
151 physIntLength = fCurrentProcess->PostStepGPIL( *fTrack,
152 fPreviousStepSize, &fCondition );
153 #ifdef G4VERBOSE
154 if(verboseLevel>0) fVerbose->DPSLPostStep();
155 #endif
156
157 switch (fCondition)
158 {
160 (*fSelectedPostStepDoItVector)[np] = ExclusivelyForced;
161 fStepStatus = fExclusivelyForcedProc;
162 fStep->GetPostStepPoint()->SetProcessDefinedStep(fCurrentProcess);
163 break;
164 case Conditionally:
165 // (*fSelectedPostStepDoItVector)[np] = Conditionally;
166 G4Exception("G4SteppingManager::DefinePhysicalStepLength()",
167 "Tracking1001", FatalException,
168 "This feature no more supported");
169 break;
170 case Forced:
171 (*fSelectedPostStepDoItVector)[np] = Forced;
172 break;
173 case StronglyForced:
174 (*fSelectedPostStepDoItVector)[np] = StronglyForced;
175 break;
176 default:
177 (*fSelectedPostStepDoItVector)[np] = InActivated;
178 break;
179 }
180
181 if (fCondition==ExclusivelyForced)
182 {
183 for(std::size_t nrest=np+1; nrest<MAXofPostStepLoops; ++nrest)
184 {
185 (*fSelectedPostStepDoItVector)[nrest] = InActivated;
186 }
187 return; // Take note the 'return' at here !!!
188 }
189 else
190 {
191 if(physIntLength < PhysicalStep )
192 {
193 PhysicalStep = physIntLength;
194 fStepStatus = fPostStepDoItProc;
195 fPostStepDoItProcTriggered = G4int(np);
196 fStep->GetPostStepPoint()->SetProcessDefinedStep(fCurrentProcess);
197 }
198 }
199 }
200
201 if (fPostStepDoItProcTriggered<MAXofPostStepLoops)
202 {
203 if ((*fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] == InActivated)
204 {
205 (*fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] = NotForced;
206 }
207 }
208
209 // GPIL for AlongStep
210 //
211 proposedSafety = DBL_MAX;
212 G4double safetyProposedToAndByProcess = proposedSafety;
213
214 for(std::size_t kp=0; kp<MAXofAlongStepLoops; ++kp)
215 {
216 fCurrentProcess = (*fAlongStepGetPhysIntVector)[kp];
217 if (fCurrentProcess == nullptr) continue;
218 // NULL means the process is inactivated by a user on fly
219
220 physIntLength = fCurrentProcess->AlongStepGPIL( *fTrack,
221 fPreviousStepSize, PhysicalStep,
222 safetyProposedToAndByProcess,
223 &fGPILSelection );
224 #ifdef G4VERBOSE
225 if(verboseLevel>0) fVerbose->DPSLAlongStep();
226 #endif
227
228 if(physIntLength < PhysicalStep)
229 {
230 PhysicalStep = physIntLength;
231
232 // Check if the process wants to be the GPIL winner. For example,
233 // multi-scattering proposes Step limit, but won't be the winner
234 //
235 if(fGPILSelection==CandidateForSelection)
236 {
237 fStepStatus = fAlongStepDoItProc;
238 fStep->GetPostStepPoint()->SetProcessDefinedStep(fCurrentProcess);
239 }
240
241 // Transportation is assumed to be the last process in the vector
242 //
243 if(kp == MAXofAlongStepLoops-1) { fStepStatus = fGeomBoundary; }
244 }
245
246 // Make sure to check the safety, even if Step is not limited
247 // by this process
248 //
249 if (safetyProposedToAndByProcess < proposedSafety)
250 {
251 // proposedSafety keeps the smallest value
252 //
253 proposedSafety = safetyProposedToAndByProcess;
254 }
255 else
256 {
257 // safetyProposedToAndByProcess always proposes a valid safety
258 //
259 safetyProposedToAndByProcess = proposedSafety;
260 }
261 }
262}
263
264//////////////////////////////////////////////////////
265void G4SteppingManager::InvokeAtRestDoItProcs()
266//////////////////////////////////////////////////////
267{
268 // Select the rest process which has the shortest time before
269 // it is invoked. In rest processes, GPIL()
270 // returns the time before a process occurs
271
272 G4double lifeTime, shortestLifeTime;
273
274 fAtRestDoItProcTriggered = 0;
275 shortestLifeTime = DBL_MAX;
276
277 unsigned int NofInactiveProc = 0;
278 for( std::size_t ri=0 ; ri < MAXofAtRestLoops ; ++ri )
279 {
280 fCurrentProcess = (*fAtRestGetPhysIntVector)[ri];
281 if (fCurrentProcess == nullptr)
282 {
283 (*fSelectedAtRestDoItVector)[ri] = InActivated;
284 ++NofInactiveProc;
285 continue;
286 } // NULL means the process is inactivated by a user on fly
287
288 lifeTime = fCurrentProcess->AtRestGPIL( *fTrack, &fCondition );
289
290 if(fCondition == Forced)
291 {
292 (*fSelectedAtRestDoItVector)[ri] = Forced;
293 }
294 else
295 {
296 (*fSelectedAtRestDoItVector)[ri] = InActivated;
297 if(lifeTime < shortestLifeTime )
298 {
299 shortestLifeTime = lifeTime;
300 fAtRestDoItProcTriggered = G4int(ri);
301 }
302 }
303 }
304
305 (*fSelectedAtRestDoItVector)[fAtRestDoItProcTriggered] = NotForced;
306
307 // at least one process is necessary to destroy the particle exit with warning
308 //
309 if(NofInactiveProc==MAXofAtRestLoops)
310 {
311 G4Exception("G4SteppingManager::InvokeAtRestDoItProcs()", "Tracking0013",
312 JustWarning, "No AtRestDoIt process is active!" );
313 }
314
315 fStep->SetStepLength( 0. ); // the particle has stopped
316 fTrack->SetStepLength( 0. );
317
318 // invoke selected process
319 //
320 for(std::size_t np=0; np<MAXofAtRestLoops; ++np)
321 {
322 //
323 // Note: DoItVector has inverse order against GetPhysIntVector
324 // and SelectedAtRestDoItVector.
325 //
326 if( (*fSelectedAtRestDoItVector)[MAXofAtRestLoops-np-1] != InActivated)
327 {
328 fCurrentProcess = (*fAtRestDoItVector)[np];
329 fParticleChange = fCurrentProcess->AtRestDoIt(*fTrack, *fStep);
330
331 // Set the current process as a process which defined this Step length
332 //
333 fStep->GetPostStepPoint()->SetProcessDefinedStep(fCurrentProcess);
334
335 // Update Step
336 //
337 fParticleChange->UpdateStepForAtRest(fStep);
338
339 // Now Store the secondaries from ParticleChange to SecondaryList
340
341 G4Track* tempSecondaryTrack;
342 G4int num2ndaries;
343
344 num2ndaries = fParticleChange->GetNumberOfSecondaries();
345
346 for(G4int DSecLoop=0; DSecLoop< num2ndaries; ++DSecLoop)
347 {
348 tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
349
350 if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
351 {
352 ApplyProductionCut(tempSecondaryTrack);
353 }
354
355 // Set parentID
356 tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
357
358 // Set the process pointer which created this track
359 tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
360
361 // If this 2ndry particle has 'zero' kinetic energy, make sure
362 // it invokes a rest process at the beginning of the tracking
363 //
364 if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN)
365 {
366 G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()
368 if(pm == nullptr)
369 {
371 ED << "A track without proper process manager is pushed\n"
372 << "into the track stack.\n"
373 << " Particle name : "
374 << tempSecondaryTrack->GetDefinition()->GetParticleName()
375 << " -- ";
376 if(tempSecondaryTrack->GetParentID()<0)
377 {
378 ED << "created by a primary particle generator.";
379 }
380 else
381 {
382 const G4VProcess* vp = tempSecondaryTrack->GetCreatorProcess();
383 if(vp != nullptr)
384 { ED << "created by " << vp->GetProcessName() << "."; }
385 else
386 { ED << "creaded by unknown process."; }
387 }
388 G4Exception("G4SteppingManager::InvokeAtRestDoItProcs()",
389 "Tracking10051", FatalException, ED);
390 }
391 if (pm->GetAtRestProcessVector()->entries()>0)
392 {
393 tempSecondaryTrack->SetTrackStatus( fStopButAlive );
394 fSecondary->push_back( tempSecondaryTrack );
395 ++fN2ndariesAtRestDoIt;
396 }
397 else
398 {
399 delete tempSecondaryTrack;
400 }
401 }
402 else
403 {
404 fSecondary->push_back( tempSecondaryTrack );
405 ++fN2ndariesAtRestDoIt;
406 }
407 } //end of loop on secondary
408
409 // clear ParticleChange
410 fParticleChange->Clear();
411
412 } // if(fSelectedAtRestDoItVector[np] != InActivated){
413 } // for(std::size_t np=0; np<MAXofAtRestLoops; ++np){
414 fStep->UpdateTrack();
415
416 fTrack->SetTrackStatus( fStopAndKill );
417}
418
419/////////////////////////////////////////////////////////
420void G4SteppingManager::InvokeAlongStepDoItProcs()
421/////////////////////////////////////////////////////////
422{
423 // If the current Step is defined by a 'ExclusivelyForced'
424 // PostStepDoIt, then don't invoke any AlongStepDoIt
425 //
426 if(fStepStatus == fExclusivelyForcedProc)
427 {
428 return; // Take note 'return' is here !!!
429 }
430
431 // Invoke all active continuous processes
432 //
433 for( std::size_t ci=0; ci<MAXofAlongStepLoops; ++ci )
434 {
435 fCurrentProcess = (*fAlongStepDoItVector)[ci];
436 if (fCurrentProcess== 0) continue;
437 // NULL means the process is inactivated by a user on fly.
438
439 fParticleChange = fCurrentProcess->AlongStepDoIt( *fTrack, *fStep );
440
441 // Update the PostStepPoint of Step according to ParticleChange
442 fParticleChange->UpdateStepForAlongStep(fStep);
443
444 #ifdef G4VERBOSE
445 if(verboseLevel>0) fVerbose->AlongStepDoItOneByOne();
446 #endif
447
448 // Now Store the secondaries from ParticleChange to SecondaryList
449
450 G4Track* tempSecondaryTrack;
451 G4int num2ndaries;
452
453 num2ndaries = fParticleChange->GetNumberOfSecondaries();
454
455 for(G4int DSecLoop=0; DSecLoop<num2ndaries; ++DSecLoop)
456 {
457 tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
458
459 if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
460 {
461 ApplyProductionCut(tempSecondaryTrack);
462 }
463
464 // Set parentID
465 tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
466
467 // Set the process pointer which created this track
468 tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
469
470 // If this secondary particle has 'zero' kinetic energy, make sure
471 // it invokes a rest process at the beginning of the tracking
472 //
473 if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN)
474 {
475 G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()
477 if(pm == nullptr)
478 {
480 ED << "A track without proper process manager is pushed\n"
481 << "into the track stack.\n"
482 << " Particle name : "
483 << tempSecondaryTrack->GetDefinition()->GetParticleName()
484 << " -- ";
485 if(tempSecondaryTrack->GetParentID()<0)
486 {
487 ED << "created by a primary particle generator.";
488 }
489 else
490 {
491 const G4VProcess* vp = tempSecondaryTrack->GetCreatorProcess();
492 if(vp != nullptr)
493 { ED << "created by " << vp->GetProcessName() << "."; }
494 else
495 { ED << "creaded by unknown process."; }
496 }
497 G4Exception("G4SteppingManager::InvokeAlongStepDoItProcs()",
498 "Tracking10051", FatalException, ED);
499 }
500 if (pm->GetAtRestProcessVector()->entries()>0)
501 {
502 tempSecondaryTrack->SetTrackStatus( fStopButAlive );
503 fSecondary->push_back( tempSecondaryTrack );
504 ++fN2ndariesAlongStepDoIt;
505 }
506 else
507 {
508 delete tempSecondaryTrack;
509 }
510 }
511 else
512 {
513 fSecondary->push_back( tempSecondaryTrack );
514 ++fN2ndariesAlongStepDoIt;
515 }
516 } //end of loop on secondary
517
518 // Set the track status according to what the process defined
519 // if kinetic energy >0, otherwise set fStopButAlive
520 //
521 fTrack->SetTrackStatus( fParticleChange->GetTrackStatus() );
522
523 // clear ParticleChange
524 fParticleChange->Clear();
525 }
526
527 fStep->UpdateTrack();
528 G4TrackStatus fNewStatus = fTrack->GetTrackStatus();
529
530 if ( fNewStatus == fAlive && fTrack->GetKineticEnergy() <= DBL_MIN )
531 {
532 if(MAXofAtRestLoops>0) fNewStatus = fStopButAlive;
533 else fNewStatus = fStopAndKill;
534 fTrack->SetTrackStatus( fNewStatus );
535 }
536}
537
538////////////////////////////////////////////////////////
539void G4SteppingManager::InvokePostStepDoItProcs()
540////////////////////////////////////////////////////////
541{
542 // Invoke the specified discrete processes
543 //
544 for(std::size_t np=0; np<MAXofPostStepLoops; ++np)
545 {
546 //
547 // Note: DoItVector has inverse order against GetPhysIntVector
548 // and SelectedPostStepDoItVector.
549 //
550 G4int Cond = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops-np-1];
551 if(Cond != InActivated)
552 {
553 if( ((Cond == NotForced) && (fStepStatus == fPostStepDoItProc)) ||
554 ((Cond == Forced) && (fStepStatus != fExclusivelyForcedProc)) ||
555 ((Cond == ExclusivelyForced) && (fStepStatus == fExclusivelyForcedProc)) ||
556 ((Cond == StronglyForced) ) )
557 {
558 InvokePSDIP(np);
559 if ((np==0) && (fTrack->GetNextVolume() == nullptr))
560 {
561 fStepStatus = fWorldBoundary;
562 fStep->GetPostStepPoint()->SetStepStatus( fStepStatus );
563 }
564 }
565 }
566
567 // Exit from PostStepLoop if the track has been killed,
568 // but extra treatment for processes with Strongly Forced flag
569 //
570 if(fTrack->GetTrackStatus() == fStopAndKill)
571 {
572 for(std::size_t np1=np+1; np1<MAXofPostStepLoops; ++np1)
573 {
574 G4int Cond2 = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops-np1-1];
575 if (Cond2 == StronglyForced)
576 {
577 InvokePSDIP(np1);
578 }
579 }
580 break;
581 }
582 }
583}
584
585////////////////////////////////////////////////////////
586void G4SteppingManager::InvokePSDIP(size_t np)
587////////////////////////////////////////////////////////
588{
589 fCurrentProcess = (*fPostStepDoItVector)[np];
590 fParticleChange = fCurrentProcess->PostStepDoIt( *fTrack, *fStep);
591
592 // Update PostStepPoint of Step according to ParticleChange
593 fParticleChange->UpdateStepForPostStep(fStep);
594
595 #ifdef G4VERBOSE
596 if(verboseLevel>0) fVerbose->PostStepDoItOneByOne();
597 #endif
598
599 // Update G4Track according to ParticleChange after each PostStepDoIt
600 fStep->UpdateTrack();
601
602 // Update safety after each invocation of PostStepDoIts
603 fStep->GetPostStepPoint()->SetSafety( CalculateSafety() );
604
605 // Now Store the secondaries from ParticleChange to SecondaryList
606
607 G4Track* tempSecondaryTrack;
608 G4int num2ndaries;
609
610 num2ndaries = fParticleChange->GetNumberOfSecondaries();
611
612 for(G4int DSecLoop=0; DSecLoop<num2ndaries; ++DSecLoop)
613 {
614 tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
615
616 if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
617 {
618 ApplyProductionCut(tempSecondaryTrack);
619 }
620
621 // Set parentID
622 tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
623
624 // Set the process pointer which created this track
625 tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
626
627 // If this secondary particle has 'zero' kinetic energy, make sure
628 // it invokes a rest process at the beginning of the tracking
629 //
630 if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN)
631 {
632 G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()
634 if(pm == nullptr)
635 {
637 ED << "A track without proper process manager is pushed\n"
638 << "into the track stack.\n"
639 << " Particle name : "
640 << tempSecondaryTrack->GetDefinition()->GetParticleName()
641 << " -- ";
642 if(tempSecondaryTrack->GetParentID()<0)
643 {
644 ED << "created by a primary particle generator.";
645 }
646 else
647 {
648 const G4VProcess* vp = tempSecondaryTrack->GetCreatorProcess();
649 if(vp != nullptr)
650 { ED << "created by " << vp->GetProcessName() << "."; }
651 else
652 { ED << "creaded by unknown process."; }
653 }
654 G4Exception("G4SteppingManager::InvokePSDIP()","Tracking10053",
655 FatalException,ED);
656 }
657 if (pm->GetAtRestProcessVector()->entries()>0)
658 {
659 tempSecondaryTrack->SetTrackStatus( fStopButAlive );
660 fSecondary->push_back( tempSecondaryTrack );
661 ++fN2ndariesPostStepDoIt;
662 }
663 else
664 {
665 delete tempSecondaryTrack;
666 }
667 }
668 else
669 {
670 fSecondary->push_back( tempSecondaryTrack );
671 ++fN2ndariesPostStepDoIt;
672 }
673 } //end of loop on secondary
674
675 // Set the track status according to what the process defined
676 fTrack->SetTrackStatus( fParticleChange->GetTrackStatus() );
677
678 // clear ParticleChange
679 fParticleChange->Clear();
680}
681
682////////////////////////////////////////////////////////
683void G4SteppingManager::ApplyProductionCut(G4Track* aSecondary)
684////////////////////////////////////////////////////////
685{
686 G4bool tBelowCutEnergyAndSafety = false;
687 G4int tPtclIdx
689 if (tPtclIdx<0) { return; }
690 G4ProductionCutsTable* tCutsTbl
692 G4int tCoupleIdx
693 = tCutsTbl->GetCoupleIndex(fPreStepPoint->GetMaterialCutsCouple());
694 if (tCoupleIdx<0) { return; }
695 G4double tProdThreshold
696 = (*(tCutsTbl->GetEnergyCutsVector(tPtclIdx)))[tCoupleIdx];
697 if( aSecondary->GetKineticEnergy()<tProdThreshold )
698 {
699 tBelowCutEnergyAndSafety = true;
700 if(std::abs(aSecondary->GetDynamicParticle()->GetCharge()) > DBL_MIN)
701 {
702 G4double currentRange
704 aSecondary->GetKineticEnergy(),
705 fPreStepPoint->GetMaterialCutsCouple());
706 tBelowCutEnergyAndSafety = (currentRange < CalculateSafety() );
707 }
708 }
709
710 if( tBelowCutEnergyAndSafety )
711 {
712 if( !(aSecondary->IsGoodForTracking()) )
713 {
714 // Add kinetic energy to the total energy deposit
715 fStep->AddTotalEnergyDeposit( aSecondary->GetKineticEnergy() );
716 aSecondary->SetKineticEnergy(0.0);
717 }
718 }
719}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
@ InActivated
@ StronglyForced
@ NotForced
@ Conditionally
@ ExclusivelyForced
@ Forced
@ CandidateForSelection
@ typeGPIL
@ typeDoIt
@ fGeomBoundary
Definition: G4StepStatus.hh:43
@ fWorldBoundary
Definition: G4StepStatus.hh:41
@ fPostStepDoItProc
Definition: G4StepStatus.hh:49
@ fAlongStepDoItProc
Definition: G4StepStatus.hh:47
@ fExclusivelyForcedProc
Definition: G4StepStatus.hh:53
G4TrackStatus
@ fAlive
@ fStopAndKill
@ fStopButAlive
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL 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
std::size_t entries() const
const std::vector< G4double > * GetEnergyCutsVector(std::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
const G4VProcess * GetCreatorProcess() const
G4VPhysicalVolume * GetNextVolume() const
G4ParticleDefinition * GetDefinition() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
G4int GetParentID() 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 PostStepGPIL(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
Definition: G4VProcess.hh:479
G4double AtRestGPIL(const G4Track &track, G4ForceCondition *condition)
Definition: G4VProcess.hh:472
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)=0
G4double AlongStepGPIL(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
Definition: G4VProcess.hh:461
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
virtual void AlongStepDoItOneByOne()=0
virtual void DPSLPostStep()=0
virtual void DPSLAlongStep()=0
virtual void PostStepDoItOneByOne()=0
virtual void DPSLStarted()=0
#define DBL_MIN
Definition: templates.hh:54
#define DBL_MAX
Definition: templates.hh:62