Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4StackManager Class Reference

#include <G4StackManager.hh>

Public Member Functions

 G4StackManager ()
 
 ~G4StackManager ()
 
const G4StackManageroperator= (const G4StackManager &)=delete
 
G4bool operator== (const G4StackManager &) const =delete
 
G4bool operator!= (const G4StackManager &) const =delete
 
G4int PushOneTrack (G4Track *newTrack, G4VTrajectory *newTrajectory=nullptr)
 
G4TrackPopNextTrack (G4VTrajectory **newTrajectory)
 
G4int PrepareNewEvent (G4Event *currentEvent)
 
void ReClassify ()
 
void SetNumberOfAdditionalWaitingStacks (G4int iAdd)
 
void TransferStackedTracks (G4ClassificationOfNewTrack origin, G4ClassificationOfNewTrack destination)
 
void TransferOneStackedTrack (G4ClassificationOfNewTrack origin, G4ClassificationOfNewTrack destination)
 
void RegisterSubEventType (G4int ty, G4int maxEnt)
 
void SetDefaultClassification (G4TrackStatus, G4ClassificationOfNewTrack, G4ExceptionSeverity es=G4ExceptionSeverity::IgnoreTheIssue)
 
void SetDefaultClassification (const G4ParticleDefinition *, G4ClassificationOfNewTrack, G4ExceptionSeverity es=G4ExceptionSeverity::IgnoreTheIssue)
 
G4ClassificationOfNewTrack GetDefaultClassification ()
 
void ReleaseSubEvent (G4int ty)
 
std::size_t GetNSubEventTypes ()
 
G4int GetSubEventType (std::size_t i)
 
void clear ()
 
void ClearUrgentStack ()
 
void ClearWaitingStack (G4int i=0)
 
void ClearPostponeStack ()
 
G4int GetNTotalTrack () const
 
G4int GetNUrgentTrack () const
 
G4int GetNWaitingTrack (G4int i=0) const
 
G4int GetNPostponedTrack () const
 
void SetVerboseLevel (G4int const value)
 
void SetUserStackingAction (G4UserStackingAction *value)
 

Detailed Description

Definition at line 69 of file G4StackManager.hh.

Constructor & Destructor Documentation

◆ G4StackManager()

G4StackManager::G4StackManager ( )

Definition at line 48 of file G4StackManager.cc.

49{
50 theMessenger = new G4StackingMessenger(this);
51#ifdef G4_USESMARTSTACK
52 urgentStack = new G4SmartTrackStack;
53 // G4cout << "+++ G4StackManager uses G4SmartTrackStack. +++" << G4endl;
54#else
55 urgentStack = new G4TrackStack(5000);
56 // G4cout << "+++ G4StackManager uses ordinary G4TrackStack. +++" << G4endl;
57#endif
58 waitingStack = new G4TrackStack(1000);
59 postponeStack = new G4TrackStack(1000);
60}

◆ ~G4StackManager()

G4StackManager::~G4StackManager ( )

Definition at line 62 of file G4StackManager.cc.

63{
64 delete userStackingAction;
65
66#ifdef G4VERBOSE
67 if(verboseLevel>0)
68 {
69 G4cout << "++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++" << G4endl;
70 G4cout << " Maximum number of tracks in the urgent stack : " << urgentStack->GetMaxNTrack() << G4endl;
71 G4cout << "++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++" << G4endl;
72 }
73#endif
74 delete urgentStack;
75 delete waitingStack;
76 delete postponeStack;
77 delete theMessenger;
78 if(numberOfAdditionalWaitingStacks>0)
79 {
80 for(G4int i=0; i<numberOfAdditionalWaitingStacks; ++i)
81 {
82 delete additionalWaitingStacks[i];
83 }
84 }
85}
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
std::size_t GetMaxNTrack() const

Member Function Documentation

◆ clear()

void G4StackManager::clear ( )

Definition at line 617 of file G4StackManager.cc.

618{
621 for(G4int i=1; i<=numberOfAdditionalWaitingStacks; ++i)
622 {
624 }
625}
void ClearWaitingStack(G4int i=0)

Referenced by G4EventManager::AbortCurrentEvent().

◆ ClearPostponeStack()

void G4StackManager::ClearPostponeStack ( )

Definition at line 647 of file G4StackManager.cc.

648{
649 postponeStack->clearAndDestroy();
650}
void clearAndDestroy()

Referenced by G4StackingMessenger::SetNewValue().

◆ ClearUrgentStack()

void G4StackManager::ClearUrgentStack ( )

Definition at line 627 of file G4StackManager.cc.

628{
629 urgentStack->clearAndDestroy();
630}

Referenced by clear(), and G4StackingMessenger::SetNewValue().

◆ ClearWaitingStack()

void G4StackManager::ClearWaitingStack ( G4int i = 0)

Definition at line 632 of file G4StackManager.cc.

633{
634 if(i==0)
635 {
636 waitingStack->clearAndDestroy();
637 }
638 else
639 {
640 if(i<=numberOfAdditionalWaitingStacks)
641 {
642 additionalWaitingStacks[i-1]->clearAndDestroy();
643 }
644 }
645}

Referenced by clear(), and G4StackingMessenger::SetNewValue().

◆ GetDefaultClassification()

G4ClassificationOfNewTrack G4StackManager::GetDefaultClassification ( )
inline

Definition at line 133 of file G4StackManager.hh.

134 { return fDefaultClassification; }

◆ GetNPostponedTrack()

G4int G4StackManager::GetNPostponedTrack ( ) const

Definition at line 684 of file G4StackManager.cc.

685{
686 return (G4int)postponeStack->GetNTrack();
687}
std::size_t GetNTrack() const

Referenced by PrepareNewEvent(), and G4StackingMessenger::SetNewValue().

◆ GetNSubEventTypes()

std::size_t G4StackManager::GetNSubEventTypes ( )
inline

Definition at line 138 of file G4StackManager.hh.

139 { return subEvtTypes.size(); }

◆ GetNTotalTrack()

G4int G4StackManager::GetNTotalTrack ( ) const

Definition at line 652 of file G4StackManager.cc.

653{
654 std::size_t n = urgentStack->GetNTrack()
655 + waitingStack->GetNTrack()
656 + postponeStack->GetNTrack();
657 for(G4int i=1; i<=numberOfAdditionalWaitingStacks; ++i)
658 {
659 n += additionalWaitingStacks[i-1]->GetNTrack();
660 }
661 return G4int(n);
662}

◆ GetNUrgentTrack()

G4int G4StackManager::GetNUrgentTrack ( ) const

Definition at line 664 of file G4StackManager.cc.

665{
666 return (G4int)urgentStack->GetNTrack();
667}

Referenced by PopNextTrack(), PushOneTrack(), ReClassify(), and G4StackingMessenger::SetNewValue().

◆ GetNWaitingTrack()

G4int G4StackManager::GetNWaitingTrack ( G4int i = 0) const

Definition at line 669 of file G4StackManager.cc.

670{
671 if(i==0)
672 {
673 return (G4int)waitingStack->GetNTrack();
674 }
675
676 if(i<=numberOfAdditionalWaitingStacks)
677 {
678 return (G4int)additionalWaitingStacks[i-1]->GetNTrack();
679 }
680
681 return 0;
682}

Referenced by PopNextTrack(), and G4StackingMessenger::SetNewValue().

◆ GetSubEventType()

G4int G4StackManager::GetSubEventType ( std::size_t i)
inline

Definition at line 140 of file G4StackManager.hh.

141 { return subEvtTypes[i]; }

◆ operator!=()

G4bool G4StackManager::operator!= ( const G4StackManager & ) const
delete

◆ operator=()

const G4StackManager & G4StackManager::operator= ( const G4StackManager & )
delete

◆ operator==()

G4bool G4StackManager::operator== ( const G4StackManager & ) const
delete

◆ PopNextTrack()

G4Track * G4StackManager::PopNextTrack ( G4VTrajectory ** newTrajectory)

Definition at line 166 of file G4StackManager.cc.

167{
168#ifdef G4VERBOSE
169 if( verboseLevel > 1 )
170 {
171 G4cout << "### pop requested out of "
172 << GetNUrgentTrack() << " stacked tracks." << G4endl;
173 }
174#endif
175
176 while( GetNUrgentTrack() == 0 )
177 {
178#ifdef G4VERBOSE
179 if( verboseLevel > 1 )
180 {
181 G4cout << "### " << GetNWaitingTrack()
182 << " waiting tracks are re-classified to" << G4endl;
183 }
184#endif
185 waitingStack->TransferTo(urgentStack);
186 if(numberOfAdditionalWaitingStacks>0)
187 {
188 for(G4int i=0; i<numberOfAdditionalWaitingStacks; ++i)
189 {
190 if(i==0)
191 {
192 additionalWaitingStacks[0]->TransferTo(waitingStack);
193 }
194 else
195 {
196 additionalWaitingStacks[i]->TransferTo(additionalWaitingStacks[i-1]);
197 }
198 }
199 }
200 if(userStackingAction != nullptr)
201 {
202 userStackingAction->NewStage();
203 }
204
205#ifdef G4VERBOSE
206 if( verboseLevel > 1 )
207 G4cout << " " << GetNUrgentTrack()
208 << " urgent tracks and " << GetNWaitingTrack()
209 << " waiting tracks." << G4endl;
210#endif
211 if( ( GetNUrgentTrack()==0 ) && ( GetNWaitingTrack()==0 ) )
212 return nullptr;
213 }
214
215 G4StackedTrack selectedStackedTrack = urgentStack->PopFromStack();
216 G4Track * selectedTrack = selectedStackedTrack.GetTrack();
217 *newTrajectory = selectedStackedTrack.GetTrajectory();
218
219#ifdef G4VERBOSE
220 if( verboseLevel > 2 )
221 {
222 G4cout << "Selected G4StackedTrack : " << &selectedStackedTrack
223 << " with G4Track " << selectedStackedTrack.GetTrack()
224 << " (trackID " << selectedStackedTrack.GetTrack()->GetTrackID()
225 << ", parentID " << selectedStackedTrack.GetTrack()->GetParentID()
226 << ")" << G4endl;
227 }
228#endif
229
230 return selectedTrack;
231}
G4int GetNUrgentTrack() const
G4int GetNWaitingTrack(G4int i=0) const
G4Track * GetTrack() const
G4VTrajectory * GetTrajectory() const
void TransferTo(G4TrackStack *aStack)
G4StackedTrack PopFromStack()
G4int GetTrackID() const
G4int GetParentID() const

◆ PrepareNewEvent()

G4int G4StackManager::PrepareNewEvent ( G4Event * currentEvent)

Definition at line 270 of file G4StackManager.cc.

271{
272 if(userStackingAction != nullptr)
273 {
274 userStackingAction->PrepareNewEvent();
275 }
276
277 // Set the urgentStack in a defined state. Not doing it would
278 // affect reproducibility
279 //
280 urgentStack->clearAndDestroy();
281
282 G4int n_passedFromPrevious = 0;
283
284 if( GetNPostponedTrack() > 0 )
285 {
286#ifdef G4VERBOSE
287 if( verboseLevel > 1 )
288 {
290 << " postponed tracked are now shifted to the stack." << G4endl;
291 }
292#endif
293
294 G4StackedTrack aStackedTrack;
295 G4TrackStack tmpStack;
296
297 postponeStack->TransferTo(&tmpStack);
298
299 while( tmpStack.GetNTrack() > 0 )
300 {
301 aStackedTrack=tmpStack.PopFromStack();
302 G4Track* aTrack = aStackedTrack.GetTrack();
303 DefineDefaultClassification( aTrack );
304 G4ClassificationOfNewTrack classification = fDefaultClassification;
305 if(userStackingAction!=nullptr)
306 {
307 classification = userStackingAction->ClassifyNewTrack( aTrack );
308 if(classification != fDefaultClassification)
309 {
310 if(fExceptionSeverity!=G4ExceptionSeverity::IgnoreTheIssue)
311 {
313 ed << "UserStackingAction has changed the track classification from "
314 << fDefaultClassification << " to " << classification << ". ";
315 G4Exception("G4StackManager::PushOneTrack","Event10052",
316 fExceptionSeverity,ed);
317 }
318 }
319 }
320 if(classification!=fKill)
321 {
322 aTrack->SetParentID(-1);
323 aTrack->SetTrackID(-(++n_passedFromPrevious));
324 }
325 else if(aTrack->GetTrackStatus() == fSuspendAndWait && classification > 0)
326 // to avoid this track sent to Waiting stack again
327 { aTrack->SetTrackStatus( fSuspend ); }
328
329 SortOut(aStackedTrack,classification);
330 }
331 }
332
333 // Reset sub-event stacks for a new event
334 for(auto& ses : subEvtStackMap)
335 { ses.second->PrepareNewEvent(currentEvent); }
336
337 return n_passedFromPrevious;
338}
@ IgnoreTheIssue
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
@ fSuspend
@ fSuspendAndWait
G4int GetNPostponedTrack() const
G4TrackStatus GetTrackStatus() const
void SetTrackStatus(const G4TrackStatus aTrackStatus)
void SetTrackID(const G4int aValue)
void SetParentID(const G4int aValue)
virtual G4ClassificationOfNewTrack ClassifyNewTrack(const G4Track *aTrack)

◆ PushOneTrack()

G4int G4StackManager::PushOneTrack ( G4Track * newTrack,
G4VTrajectory * newTrajectory = nullptr )

Definition at line 87 of file G4StackManager.cc.

89{
90 const G4ParticleDefinition* pd = newTrack->GetParticleDefinition();
91 if(pd->GetParticleDefinitionID() < 0)
92 {
94 ED << "A track without proper process manager is pushed \
95 into the track stack.\n"
96 << " Particle name : " << pd->GetParticleName() << " -- ";
97 if(newTrack->GetParentID()==0)
98 {
99 ED << "created by a primary particle generator.";
100 }
101 else
102 {
103 const G4VProcess* vp = newTrack->GetCreatorProcess();
104 if(vp != nullptr)
105 {
106 ED << "created by " << vp->GetProcessName() << ".";
107 }
108 else
109 {
110 ED << "creaded by unknown process.";
111 }
112 }
113 G4Exception("G4StackManager::PushOneTrack","Event10051",
114 FatalException,ED);
115 delete newTrack;
116 return GetNUrgentTrack();
117 }
118
119 DefineDefaultClassification( newTrack );
120 G4ClassificationOfNewTrack classification = fDefaultClassification;
121 if(userStackingAction!=nullptr)
122 {
123 classification = userStackingAction->ClassifyNewTrack( newTrack );
124 if(classification != fDefaultClassification)
125 {
126 if(fExceptionSeverity!=G4ExceptionSeverity::IgnoreTheIssue)
127 {
129 ed << "UserStackingAction has changed the track classification from "
130 << fDefaultClassification << " to " << classification << ". ";
131 G4Exception("G4StackManager::PushOneTrack","Event10052",
132 fExceptionSeverity,ed);
133 }
134 }
135 }
136 if(newTrack->GetTrackStatus() == fSuspendAndWait && classification > 0)
137 // to avoid this track sent to Waiting stack again
138 { newTrack->SetTrackStatus( fSuspend ); }
139
140#ifdef G4VERBOSE
141 if( verboseLevel > 1 )
142 {
143 G4cout << "### Storing a track ("
144 << newTrack->GetParticleDefinition()->GetParticleName()
145 << ",trackID=" << newTrack->GetTrackID()
146 << ",parentID=" << newTrack->GetParentID() << ") ";
147 if(newTrack->GetParentID()==0)
148 { G4cout << "created by a primary particle generator "; }
149 else
150 {
151 const G4VProcess* vp = newTrack->GetCreatorProcess();
152 if(vp != nullptr)
153 { G4cout << "created by " << vp->GetProcessName() << " "; }
154 else
155 { G4cout << "creaded by unknown process "; }
156 }
157 G4cout << "into stack #" << classification << G4endl;
158 }
159#endif
160 G4StackedTrack newStackedTrack( newTrack, newTrajectory );
161 SortOut(newStackedTrack,classification);
162
163 return GetNUrgentTrack();
164}
@ FatalException
G4int GetParticleDefinitionID() const
const G4String & GetParticleName() const
const G4ParticleDefinition * GetParticleDefinition() const
const G4VProcess * GetCreatorProcess() const
const G4String & GetProcessName() const

Referenced by G4EventManager::StackTracks().

◆ ReClassify()

void G4StackManager::ReClassify ( )

Definition at line 233 of file G4StackManager.cc.

234{
235 G4StackedTrack aStackedTrack;
236 G4TrackStack tmpStack;
237
238 if( userStackingAction == nullptr ) return;
239 if( GetNUrgentTrack() == 0 ) return;
240
241 urgentStack->TransferTo(&tmpStack);
242 while( tmpStack.GetNTrack() > 0 )
243 {
244 aStackedTrack=tmpStack.PopFromStack();
245 DefineDefaultClassification( aStackedTrack.GetTrack() );
246 G4ClassificationOfNewTrack classification = fDefaultClassification;
247 if(userStackingAction!=nullptr)
248 {
249 classification = userStackingAction->ClassifyNewTrack( aStackedTrack.GetTrack() );
250 if(classification != fDefaultClassification)
251 {
252 if(fExceptionSeverity!=G4ExceptionSeverity::IgnoreTheIssue)
253 {
255 ed << "UserStackingAction has changed the track classification from "
256 << fDefaultClassification << " to " << classification << ". ";
257 G4Exception("G4StackManager::PushOneTrack","Event10052",
258 fExceptionSeverity,ed);
259 }
260 }
261 }
262 if(aStackedTrack.GetTrack()->GetTrackStatus() == fSuspendAndWait && classification > 0)
263 // to avoid this track sent to Waiting stack again
264 { aStackedTrack.GetTrack()->SetTrackStatus( fSuspend ); }
265
266 SortOut(aStackedTrack,classification);
267 }
268}

Referenced by G4AdjointStackingAction::NewStage().

◆ RegisterSubEventType()

void G4StackManager::RegisterSubEventType ( G4int ty,
G4int maxEnt )

Definition at line 780 of file G4StackManager.cc.

781{
782 if(subEvtStackMap.find(ty)==subEvtStackMap.end())
783 {
784 subEvtStackMap[ty] = new G4SubEventTrackStack(ty,maxEnt);
785 subEvtTypes.push_back(ty);
786#ifdef G4VERBOSE
787 subEvtStackMap[ty]->SetVerboseLevel(verboseLevel);
788 if( verboseLevel > 0 )
789 {
790 G4cout << " ---> New sub-event stack for sub-event type "
791 << ty << " is created. Classification id for this stack is "
792 << subEvtTypes.size() + 99 << "." << G4endl;
793 }
794 }
795 else
796 {
797 if( verboseLevel > 1 )
798 {
799 G4cout << " ---> Sub-event stack for sub-event type "
800 << ty << " already registered." << G4endl;
801 }
802#endif
803 }
804}

◆ ReleaseSubEvent()

void G4StackManager::ReleaseSubEvent ( G4int ty)

Definition at line 806 of file G4StackManager.cc.

807{
808 auto ses = subEvtStackMap.find(ty);
809 if(ses==subEvtStackMap.end())
810 {
812 ED << "Un-registered sub-event type " << ty << " requested.";
813 G4Exception("G4StackManager::PopSubEvent", "SubEvt8001",
814 FatalException,ED);
815 return; // NOLINT: Required to silence Coverity, which does not recognize G4Exception as exit point
816 }
817
818 ses->second->ReleaseSubEvent();
819}

◆ SetDefaultClassification() [1/2]

void G4StackManager::SetDefaultClassification ( const G4ParticleDefinition * pd,
G4ClassificationOfNewTrack val,
G4ExceptionSeverity es = G4ExceptionSeverity::IgnoreTheIssue )

Definition at line 757 of file G4StackManager.cc.

759{
760 auto pdm = defClassPartDef.find(pd);
761 if(pdm==defClassPartDef.end())
762 { defClassPartDef[pd] = std::pair(val,es); }
763 else
764 {
765 if(pdm->second.first!=val)
766 { // alternating default classification
768 ed << "Default classification for particle " << pd->GetParticleName()
769 << " is changed from " << pdm->second.first << " to "
770 << val << ".";
771 G4Exception("G4StackManager::SetDefaultClassification",
772 "Event11052", JustWarning,ed);
773 pdm->second.first = val;
774 }
775 // Change severity if needed.
776 if(pdm->second.second>es) pdm->second.second = es;
777 }
778}
@ JustWarning

◆ SetDefaultClassification() [2/2]

void G4StackManager::SetDefaultClassification ( G4TrackStatus ts,
G4ClassificationOfNewTrack val,
G4ExceptionSeverity es = G4ExceptionSeverity::IgnoreTheIssue )

Definition at line 734 of file G4StackManager.cc.

736{
737 auto tsm = defClassTrackStatus.find(ts);
738 if(tsm==defClassTrackStatus.end())
739 { defClassTrackStatus[ts] = std::pair(val,es); }
740 else
741 {
742 if(tsm->second.first!=val)
743 { // alternating default classification
745 ed << "Default classification for track status " << ts
746 << " is changed from " << tsm->second.first << " to "
747 << val << ".";
748 G4Exception("G4StackManager::SetDefaultClassification",
749 "Event11051",JustWarning,ed);
750 tsm->second.first = val;
751 }
752 // Change severity if needed.
753 if(tsm->second.second>es) tsm->second.second = es;
754 }
755}

◆ SetNumberOfAdditionalWaitingStacks()

void G4StackManager::SetNumberOfAdditionalWaitingStacks ( G4int iAdd)

Definition at line 409 of file G4StackManager.cc.

410{
411 if(iAdd > numberOfAdditionalWaitingStacks)
412 {
413 for(G4int i=numberOfAdditionalWaitingStacks; i<iAdd; ++i)
414 {
415 auto* newStack = new G4TrackStack;
416 additionalWaitingStacks.push_back(newStack);
417 }
418 numberOfAdditionalWaitingStacks = iAdd;
419 }
420 else if (iAdd < numberOfAdditionalWaitingStacks)
421 {
422 for(G4int i=numberOfAdditionalWaitingStacks; i>iAdd; --i)
423 {
424 delete additionalWaitingStacks[i];
425 }
426 }
427}

Referenced by G4RunManager::SetNumberOfAdditionalWaitingStacks().

◆ SetUserStackingAction()

void G4StackManager::SetUserStackingAction ( G4UserStackingAction * value)

Definition at line 696 of file G4StackManager.cc.

697{
698 userStackingAction = value;
699 if(userStackingAction != nullptr)
700 {
701 userStackingAction->SetStackManager(this);
702 }
703}
void SetStackManager(G4StackManager *value)

Referenced by G4EventManager::SetUserAction().

◆ SetVerboseLevel()

void G4StackManager::SetVerboseLevel ( G4int const value)

Definition at line 689 of file G4StackManager.cc.

690{
691 verboseLevel = value;
692 for(auto& sets : subEvtStackMap)
693 { sets.second->SetVerboseLevel(value); }
694}

Referenced by G4StackingMessenger::SetNewValue(), and G4EventManager::SetVerboseLevel().

◆ TransferOneStackedTrack()

void G4StackManager::TransferOneStackedTrack ( G4ClassificationOfNewTrack origin,
G4ClassificationOfNewTrack destination )

Definition at line 522 of file G4StackManager.cc.

525{
526 if(origin==destination) return;
527 if(origin==fKill) return;
528 G4TrackStack* originStack = nullptr;
529 switch(origin)
530 {
531 case fUrgent:
532 originStack = nullptr;
533 break;
534 case fWaiting:
535 originStack = waitingStack;
536 break;
537 case fPostpone:
538 originStack = postponeStack;
539 break;
540 default:
541 G4int i = origin - 10;
542 if(i<=numberOfAdditionalWaitingStacks)
543 {
544 originStack = additionalWaitingStacks[i-1];
545 }
546 else
547 {
549 ED << "Invalid origin stack ID " << origin;
550 G4Exception("G4StackManager::TransferStackedTracks","Stack0911",
551 FatalException, ED);
552 }
553 break;
554 }
555
556 G4StackedTrack aStackedTrack;
557 if(destination==fKill)
558 {
559 if( originStack != nullptr && (originStack->GetNTrack() != 0u) )
560 {
561 aStackedTrack = originStack->PopFromStack();
562 delete aStackedTrack.GetTrack();
563 delete aStackedTrack.GetTrajectory();
564 }
565 else if (urgentStack->GetNTrack() != 0u )
566 {
567 aStackedTrack = urgentStack->PopFromStack();
568 delete aStackedTrack.GetTrack();
569 delete aStackedTrack.GetTrajectory();
570 }
571 }
572 else
573 {
574 G4TrackStack* targetStack = nullptr;
575 switch(destination)
576 {
577 case fUrgent:
578 targetStack = nullptr;
579 break;
580 case fWaiting:
581 targetStack = waitingStack;
582 break;
583 case fPostpone:
584 targetStack = postponeStack;
585 break;
586 default:
587 G4int i = destination - 10;
588 if(i<=numberOfAdditionalWaitingStacks)
589 {
590 targetStack = additionalWaitingStacks[i-1];
591 }
592 else
593 {
595 ED << "Invalid origin stack ID " << origin;
596 G4Exception("G4StackManager::TransferStackedTracks","Stack0911",
597 FatalException, ED);
598 }
599 break;
600 }
601 if((originStack != nullptr) && (originStack->GetNTrack() != 0u))
602 {
603 aStackedTrack = originStack->PopFromStack();
604 if(targetStack != nullptr) { targetStack->PushToStack(aStackedTrack); }
605 else { urgentStack->PushToStack(aStackedTrack); }
606 }
607 else if(urgentStack->GetNTrack() != 0u)
608 {
609 aStackedTrack = urgentStack->PopFromStack();
610 if(targetStack != nullptr) { targetStack->PushToStack(aStackedTrack); }
611 else { urgentStack->PushToStack(aStackedTrack); }
612 }
613 }
614 return;
615}
void PushToStack(const G4StackedTrack &aStackedTrack)

◆ TransferStackedTracks()

void G4StackManager::TransferStackedTracks ( G4ClassificationOfNewTrack origin,
G4ClassificationOfNewTrack destination )

Definition at line 429 of file G4StackManager.cc.

432{
433 if(origin==destination) return;
434 if(origin==fKill) return;
435 G4TrackStack* originStack = nullptr;
436 switch(origin)
437 {
438 case fUrgent:
439 originStack = nullptr;
440 break;
441 case fWaiting:
442 originStack = waitingStack;
443 break;
444 case fPostpone:
445 originStack = postponeStack;
446 break;
447 default:
448 G4int i = origin - 10;
449 if(i<=numberOfAdditionalWaitingStacks)
450 {
451 originStack = additionalWaitingStacks[i-1];
452 }
453 else
454 {
456 ED << "Invalid origin stack ID " << origin;
457 G4Exception("G4StackManager::TransferStackedTracks","Stack0911",
458 FatalException, ED);
459 }
460 break;
461 }
462
463 if(destination==fKill)
464 {
465 if(originStack != nullptr)
466 {
467 originStack->clearAndDestroy();
468 }
469 else
470 {
471 urgentStack->clearAndDestroy();
472 }
473 }
474 else
475 {
476 G4TrackStack* targetStack = nullptr;
477 switch(destination)
478 {
479 case fUrgent:
480 targetStack = nullptr;
481 break;
482 case fWaiting:
483 targetStack = waitingStack;
484 break;
485 case fPostpone:
486 targetStack = postponeStack;
487 break;
488 default:
489 G4int i = destination - 10;
490 if(i<=numberOfAdditionalWaitingStacks)
491 {
492 targetStack = additionalWaitingStacks[i-1];
493 }
494 else
495 {
497 ED << "Invalid origin stack ID " << origin;
498 G4Exception("G4StackManager::TransferStackedTracks","Stack0911",
499 FatalException, ED);
500 }
501 break;
502 }
503 if(originStack != nullptr)
504 {
505 if(targetStack != nullptr)
506 {
507 originStack->TransferTo(targetStack);
508 }
509 else
510 {
511 originStack->TransferTo(urgentStack);
512 }
513 }
514 else
515 {
516 urgentStack->TransferTo(targetStack);
517 }
518 }
519 return;
520}

The documentation for this class was generated from the following files: