Geant4-11
Public Member Functions | Private Member Functions | Private Attributes
G4StackManager Class Reference

#include <G4StackManager.hh>

Public Member Functions

void clear ()
 
void ClearPostponeStack ()
 
void ClearUrgentStack ()
 
void ClearWaitingStack (G4int i=0)
 
 G4StackManager ()
 
G4int GetNPostponedTrack () const
 
G4int GetNTotalTrack () const
 
G4int GetNUrgentTrack () const
 
G4int GetNWaitingTrack (G4int i=0) const
 
G4bool operator!= (const G4StackManager &) const =delete
 
const G4StackManageroperator= (const G4StackManager &)=delete
 
G4bool operator== (const G4StackManager &) const =delete
 
G4TrackPopNextTrack (G4VTrajectory **newTrajectory)
 
G4int PrepareNewEvent ()
 
G4int PushOneTrack (G4Track *newTrack, G4VTrajectory *newTrajectory=nullptr)
 
void ReClassify ()
 
void SetNumberOfAdditionalWaitingStacks (G4int iAdd)
 
void SetUserStackingAction (G4UserStackingAction *value)
 
void SetVerboseLevel (G4int const value)
 
void TransferOneStackedTrack (G4ClassificationOfNewTrack origin, G4ClassificationOfNewTrack destination)
 
void TransferStackedTracks (G4ClassificationOfNewTrack origin, G4ClassificationOfNewTrack destination)
 
 ~G4StackManager ()
 

Private Member Functions

G4ClassificationOfNewTrack DefaultClassification (G4Track *aTrack)
 

Private Attributes

std::vector< G4TrackStack * > additionalWaitingStacks
 
G4int numberOfAdditionalWaitingStacks = 0
 
G4TrackStackpostponeStack = nullptr
 
G4StackingMessengertheMessenger = nullptr
 
G4TrackStackurgentStack = nullptr
 
G4UserStackingActionuserStackingAction = nullptr
 
G4int verboseLevel = 0
 
G4TrackStackwaitingStack = nullptr
 

Detailed Description

Definition at line 61 of file G4StackManager.hh.

Constructor & Destructor Documentation

◆ G4StackManager()

G4StackManager::G4StackManager ( )

Definition at line 45 of file G4StackManager.cc.

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

References postponeStack, theMessenger, urgentStack, and waitingStack.

◆ ~G4StackManager()

G4StackManager::~G4StackManager ( )

Definition at line 59 of file G4StackManager.cc.

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

References additionalWaitingStacks, G4cout, G4endl, G4TrackStack::GetMaxNTrack(), numberOfAdditionalWaitingStacks, postponeStack, theMessenger, urgentStack, userStackingAction, verboseLevel, and waitingStack.

Member Function Documentation

◆ clear()

void G4StackManager::clear ( )

Definition at line 547 of file G4StackManager.cc.

548{
551 for(G4int i=1; i<=numberOfAdditionalWaitingStacks; ++i)
552 {
554 }
555}
void ClearWaitingStack(G4int i=0)

References ClearUrgentStack(), ClearWaitingStack(), and numberOfAdditionalWaitingStacks.

Referenced by G4EventManager::AbortCurrentEvent(), and export_G4StackManager().

◆ ClearPostponeStack()

void G4StackManager::ClearPostponeStack ( )

◆ ClearUrgentStack()

void G4StackManager::ClearUrgentStack ( )

◆ ClearWaitingStack()

void G4StackManager::ClearWaitingStack ( G4int  i = 0)

Definition at line 562 of file G4StackManager.cc.

563{
564 if(i==0)
565 {
567 }
568 else
569 {
571 {
572 additionalWaitingStacks[i-1]->clearAndDestroy();
573 }
574 }
575}

References additionalWaitingStacks, G4TrackStack::clearAndDestroy(), numberOfAdditionalWaitingStacks, and waitingStack.

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

◆ DefaultClassification()

G4ClassificationOfNewTrack G4StackManager::DefaultClassification ( G4Track aTrack)
private

Definition at line 634 of file G4StackManager.cc.

636{
637 G4ClassificationOfNewTrack classification = fUrgent;
638 if( aTrack->GetTrackStatus() == fPostponeToNextEvent )
639 {
640 classification = fPostpone;
641 }
642 return classification;
643}
@ fPostponeToNextEvent
G4TrackStatus GetTrackStatus() const

References fPostpone, fPostponeToNextEvent, fUrgent, and G4Track::GetTrackStatus().

Referenced by PrepareNewEvent(), and PushOneTrack().

◆ GetNPostponedTrack()

G4int G4StackManager::GetNPostponedTrack ( ) const

Definition at line 615 of file G4StackManager.cc.

616{
617 return postponeStack->GetNTrack();
618}
std::size_t GetNTrack() const
Definition: G4TrackStack.hh:68

References G4TrackStack::GetNTrack(), and postponeStack.

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

◆ GetNTotalTrack()

G4int G4StackManager::GetNTotalTrack ( ) const

◆ GetNUrgentTrack()

G4int G4StackManager::GetNUrgentTrack ( ) const

◆ GetNWaitingTrack()

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

Definition at line 599 of file G4StackManager.cc.

600{
601 if(i==0)
602 {
603 return waitingStack->GetNTrack();
604 }
605 else
606 {
608 {
609 return additionalWaitingStacks[i-1]->GetNTrack();
610 }
611 }
612 return 0;
613}

References additionalWaitingStacks, G4TrackStack::GetNTrack(), numberOfAdditionalWaitingStacks, and waitingStack.

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

◆ 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 169 of file G4StackManager.cc.

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

References additionalWaitingStacks, G4cout, G4endl, GetNUrgentTrack(), GetNWaitingTrack(), G4Track::GetParentID(), G4StackedTrack::GetTrack(), G4Track::GetTrackID(), G4StackedTrack::GetTrajectory(), G4UserStackingAction::NewStage(), numberOfAdditionalWaitingStacks, G4TrackStack::PopFromStack(), G4TrackStack::TransferTo(), urgentStack, userStackingAction, verboseLevel, and waitingStack.

Referenced by G4EventManager::DoProcessing().

◆ PrepareNewEvent()

G4int G4StackManager::PrepareNewEvent ( )

Definition at line 283 of file G4StackManager.cc.

284{
286 {
288 }
289
290 // Set the urgentStack in a defined state. Not doing it would
291 // affect reproducibility
292 //
294
295 G4int n_passedFromPrevious = 0;
296
297 if( GetNPostponedTrack() > 0 )
298 {
299#ifdef G4VERBOSE
300 if( verboseLevel > 1 )
301 {
303 << " postponed tracked are now shifted to the stack." << G4endl;
304 }
305#endif
306
307 G4StackedTrack aStackedTrack;
308 G4TrackStack tmpStack;
309
310 postponeStack->TransferTo(&tmpStack);
311
312 while( tmpStack.GetNTrack() > 0 )
313 {
314 aStackedTrack=tmpStack.PopFromStack();
315 G4Track* aTrack = aStackedTrack.GetTrack();
316 aTrack->SetParentID(-1);
317 G4ClassificationOfNewTrack classification;
319 {
320 classification = userStackingAction->ClassifyNewTrack( aTrack );
321 }
322 else
323 {
324 classification = DefaultClassification( aTrack );
325 }
326
327 if(classification==fKill)
328 {
329 delete aTrack;
330 delete aStackedTrack.GetTrajectory();
331 }
332 else
333 {
334 aTrack->SetTrackID(-(++n_passedFromPrevious));
335 switch (classification)
336 {
337 case fUrgent:
338 urgentStack->PushToStack( aStackedTrack );
339 break;
340 case fWaiting:
341 waitingStack->PushToStack( aStackedTrack );
342 break;
343 case fPostpone:
344 postponeStack->PushToStack( aStackedTrack );
345 break;
346 default:
347 G4int i = classification - 10;
349 {
351 ED << "invalid classification " << classification << G4endl;
352 G4Exception("G4StackManager::PrepareNewEvent", "Event0053",
353 FatalException, ED);
354 }
355 else
356 {
357 additionalWaitingStacks[i-1]->PushToStack( aStackedTrack );
358 }
359 break;
360 }
361 }
362 }
363 }
364 return n_passedFromPrevious;
365}
@ 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
G4ClassificationOfNewTrack DefaultClassification(G4Track *aTrack)
G4int GetNPostponedTrack() const
void PushToStack(const G4StackedTrack &aStackedTrack)
Definition: G4TrackStack.hh:59
void SetTrackID(const G4int aValue)
void SetParentID(const G4int aValue)
virtual G4ClassificationOfNewTrack ClassifyNewTrack(const G4Track *aTrack)
virtual void PrepareNewEvent()

References additionalWaitingStacks, G4UserStackingAction::ClassifyNewTrack(), G4TrackStack::clearAndDestroy(), DefaultClassification(), FatalException, fKill, fPostpone, fUrgent, fWaiting, G4cout, G4endl, G4Exception(), GetNPostponedTrack(), G4TrackStack::GetNTrack(), G4StackedTrack::GetTrack(), G4StackedTrack::GetTrajectory(), numberOfAdditionalWaitingStacks, G4TrackStack::PopFromStack(), postponeStack, G4UserStackingAction::PrepareNewEvent(), G4TrackStack::PushToStack(), G4Track::SetParentID(), G4Track::SetTrackID(), G4TrackStack::TransferTo(), urgentStack, userStackingAction, verboseLevel, and waitingStack.

Referenced by G4EventManager::DoProcessing().

◆ PushOneTrack()

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

Definition at line 84 of file G4StackManager.cc.

86{
87 const G4ParticleDefinition* pd = newTrack->GetParticleDefinition();
88 if(pd->GetParticleDefinitionID() < 0)
89 {
91 ED << "A track without proper process manager is pushed \
92 into the track stack.\n"
93 << " Particle name : " << pd->GetParticleName() << " -- ";
94 if(newTrack->GetParentID()<0)
95 {
96 ED << "created by a primary particle generator.";
97 }
98 else
99 {
100 const G4VProcess* vp = newTrack->GetCreatorProcess();
101 if(vp)
102 {
103 ED << "created by " << vp->GetProcessName() << ".";
104 }
105 else
106 {
107 ED << "creaded by unknown process.";
108 }
109 }
110 G4Exception("G4StackManager::PushOneTrack","Event10051",
111 FatalException,ED);
112 delete newTrack;
113 return GetNUrgentTrack();
114 }
115
116 G4ClassificationOfNewTrack classification = DefaultClassification( newTrack );
117 if(userStackingAction != nullptr)
118 {
119 classification = userStackingAction->ClassifyNewTrack( newTrack );
120 }
121
122 if(classification==fKill) // delete newTrack without stacking
123 {
124#ifdef G4VERBOSE
125 if( verboseLevel > 1 )
126 {
127 G4cout << " ---> G4Track " << newTrack << " (trackID "
128 << newTrack->GetTrackID() << ", parentID "
129 << newTrack->GetParentID() << ") is not to be stored." << G4endl;
130 }
131#endif
132 delete newTrack;
133 delete newTrajectory;
134 }
135 else
136 {
137 G4StackedTrack newStackedTrack( newTrack, newTrajectory );
138 switch (classification)
139 {
140 case fUrgent:
141 urgentStack->PushToStack( newStackedTrack );
142 break;
143 case fWaiting:
144 waitingStack->PushToStack( newStackedTrack );
145 break;
146 case fPostpone:
147 postponeStack->PushToStack( newStackedTrack );
148 break;
149 default:
150 G4int i = classification - 10;
152 {
154 ED << "invalid classification " << classification << G4endl;
155 G4Exception("G4StackManager::PushOneTrack", "Event0051",
156 FatalException,ED);
157 }
158 else
159 {
160 additionalWaitingStacks[i-1]->PushToStack( newStackedTrack );
161 }
162 break;
163 }
164 }
165 return GetNUrgentTrack();
166}
G4int GetParticleDefinitionID() const
const G4String & GetParticleName() const
const G4ParticleDefinition * GetParticleDefinition() const
const G4VProcess * GetCreatorProcess() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

References additionalWaitingStacks, G4UserStackingAction::ClassifyNewTrack(), DefaultClassification(), FatalException, fKill, fPostpone, fUrgent, fWaiting, G4cout, G4endl, G4Exception(), G4Track::GetCreatorProcess(), GetNUrgentTrack(), G4Track::GetParentID(), G4Track::GetParticleDefinition(), G4ParticleDefinition::GetParticleDefinitionID(), G4ParticleDefinition::GetParticleName(), G4VProcess::GetProcessName(), G4Track::GetTrackID(), numberOfAdditionalWaitingStacks, postponeStack, G4TrackStack::PushToStack(), urgentStack, userStackingAction, verboseLevel, and waitingStack.

Referenced by G4EventManager::DoProcessing(), and G4EventManager::StackTracks().

◆ ReClassify()

void G4StackManager::ReClassify ( )

Definition at line 236 of file G4StackManager.cc.

237{
238 G4StackedTrack aStackedTrack;
239 G4TrackStack tmpStack;
240
241 if( userStackingAction == nullptr ) return;
242 if( GetNUrgentTrack() == 0 ) return;
243
244 urgentStack->TransferTo(&tmpStack);
245 while( tmpStack.GetNTrack() > 0 )
246 {
247 aStackedTrack=tmpStack.PopFromStack();
248 G4ClassificationOfNewTrack classification =
249 userStackingAction->ClassifyNewTrack( aStackedTrack.GetTrack() );
250 switch (classification)
251 {
252 case fKill:
253 delete aStackedTrack.GetTrack();
254 delete aStackedTrack.GetTrajectory();
255 break;
256 case fUrgent:
257 urgentStack->PushToStack( aStackedTrack );
258 break;
259 case fWaiting:
260 waitingStack->PushToStack( aStackedTrack );
261 break;
262 case fPostpone:
263 postponeStack->PushToStack( aStackedTrack );
264 break;
265 default:
266 G4int i = classification - 10;
268 {
270 ED << "invalid classification " << classification << G4endl;
271 G4Exception("G4StackManager::ReClassify", "Event0052",
272 FatalException, ED);
273 }
274 else
275 {
276 additionalWaitingStacks[i-1]->PushToStack( aStackedTrack );
277 }
278 break;
279 }
280 }
281}

References additionalWaitingStacks, G4UserStackingAction::ClassifyNewTrack(), FatalException, fKill, fPostpone, fUrgent, fWaiting, G4endl, G4Exception(), G4TrackStack::GetNTrack(), GetNUrgentTrack(), G4StackedTrack::GetTrack(), G4StackedTrack::GetTrajectory(), numberOfAdditionalWaitingStacks, G4TrackStack::PopFromStack(), postponeStack, G4TrackStack::PushToStack(), G4TrackStack::TransferTo(), urgentStack, userStackingAction, and waitingStack.

Referenced by export_G4StackManager(), and G4AdjointStackingAction::NewStage().

◆ SetNumberOfAdditionalWaitingStacks()

void G4StackManager::SetNumberOfAdditionalWaitingStacks ( G4int  iAdd)

Definition at line 367 of file G4StackManager.cc.

368{
370 {
371 for(G4int i=numberOfAdditionalWaitingStacks; i<iAdd; ++i)
372 {
373 G4TrackStack* newStack = new G4TrackStack;
374 additionalWaitingStacks.push_back(newStack);
375 }
377 }
378 else if (iAdd < numberOfAdditionalWaitingStacks)
379 {
380 for(G4int i=numberOfAdditionalWaitingStacks; i>iAdd; --i)
381 {
382 delete additionalWaitingStacks[i];
383 }
384 }
385}

References additionalWaitingStacks, and numberOfAdditionalWaitingStacks.

Referenced by G4EventManager::SetNumberOfAdditionalWaitingStacks().

◆ SetUserStackingAction()

void G4StackManager::SetUserStackingAction ( G4UserStackingAction value)

Definition at line 625 of file G4StackManager.cc.

626{
627 userStackingAction = value;
629 {
631 }
632}
void SetStackManager(G4StackManager *value)

References G4UserStackingAction::SetStackManager(), and userStackingAction.

Referenced by G4EventManager::SetUserAction().

◆ SetVerboseLevel()

void G4StackManager::SetVerboseLevel ( G4int const  value)

Definition at line 620 of file G4StackManager.cc.

621{
622 verboseLevel = value;
623}

References verboseLevel.

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

◆ TransferOneStackedTrack()

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

Definition at line 466 of file G4StackManager.cc.

469{
470 if(origin==destination) return;
471 if(origin==fKill) return;
472 G4TrackStack* originStack = nullptr;
473 switch(origin)
474 {
475 case fUrgent:
476 originStack = nullptr;
477 break;
478 case fWaiting:
479 originStack = waitingStack;
480 break;
481 case fPostpone:
482 originStack = postponeStack;
483 break;
484 default:
485 G4int i = origin - 10;
487 {
488 originStack = additionalWaitingStacks[i-1];
489 }
490 break;
491 }
492
493 G4StackedTrack aStackedTrack;
494 if(destination==fKill)
495 {
496 if( originStack != nullptr && originStack->GetNTrack() )
497 {
498 aStackedTrack = originStack->PopFromStack();
499 delete aStackedTrack.GetTrack();
500 delete aStackedTrack.GetTrajectory();
501 }
502 else if (urgentStack->GetNTrack() )
503 {
504 aStackedTrack = urgentStack->PopFromStack();
505 delete aStackedTrack.GetTrack();
506 delete aStackedTrack.GetTrajectory();
507 }
508 }
509 else
510 {
511 G4TrackStack* targetStack = nullptr;
512 switch(destination)
513 {
514 case fUrgent:
515 targetStack = nullptr;
516 break;
517 case fWaiting:
518 targetStack = waitingStack;
519 break;
520 case fPostpone:
521 targetStack = postponeStack;
522 break;
523 default:
524 G4int i = destination - 10;
526 {
527 targetStack = additionalWaitingStacks[i-1];
528 }
529 break;
530 }
531 if(originStack && originStack->GetNTrack())
532 {
533 aStackedTrack = originStack->PopFromStack();
534 if(targetStack) { targetStack->PushToStack(aStackedTrack); }
535 else { urgentStack->PushToStack(aStackedTrack); }
536 }
537 else if(urgentStack->GetNTrack())
538 {
539 aStackedTrack = urgentStack->PopFromStack();
540 if(targetStack) { targetStack->PushToStack(aStackedTrack); }
541 else { urgentStack->PushToStack(aStackedTrack); }
542 }
543 }
544 return;
545}

References additionalWaitingStacks, fKill, fPostpone, fUrgent, fWaiting, G4TrackStack::GetNTrack(), G4StackedTrack::GetTrack(), G4StackedTrack::GetTrajectory(), numberOfAdditionalWaitingStacks, G4TrackStack::PopFromStack(), postponeStack, G4TrackStack::PushToStack(), urgentStack, and waitingStack.

◆ TransferStackedTracks()

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

Definition at line 387 of file G4StackManager.cc.

390{
391 if(origin==destination) return;
392 if(origin==fKill) return;
393 G4TrackStack* originStack = nullptr;
394 switch(origin)
395 {
396 case fUrgent:
397 originStack = nullptr;
398 break;
399 case fWaiting:
400 originStack = waitingStack;
401 break;
402 case fPostpone:
403 originStack = postponeStack;
404 break;
405 default:
406 G4int i = origin - 10;
408 {
409 originStack = additionalWaitingStacks[i-1];
410 }
411 break;
412 }
413
414 if(destination==fKill)
415 {
416 if(originStack != nullptr)
417 {
418 originStack->clearAndDestroy();
419 }
420 else
421 {
423 }
424 }
425 else
426 {
427 G4TrackStack* targetStack = nullptr;
428 switch(destination)
429 {
430 case fUrgent:
431 targetStack = nullptr;
432 break;
433 case fWaiting:
434 targetStack = waitingStack;
435 break;
436 case fPostpone:
437 targetStack = postponeStack;
438 break;
439 default:
440 G4int i = destination - 10;
442 {
443 targetStack = additionalWaitingStacks[i-1];
444 }
445 break;
446 }
447 if(originStack != nullptr)
448 {
449 if(targetStack != nullptr)
450 {
451 originStack->TransferTo(targetStack);
452 }
453 else
454 {
455 originStack->TransferTo(urgentStack);
456 }
457 }
458 else
459 {
460 urgentStack->TransferTo(targetStack);
461 }
462 }
463 return;
464}

References additionalWaitingStacks, G4TrackStack::clearAndDestroy(), fKill, fPostpone, fUrgent, fWaiting, numberOfAdditionalWaitingStacks, postponeStack, G4TrackStack::TransferTo(), urgentStack, and waitingStack.

Field Documentation

◆ additionalWaitingStacks

std::vector<G4TrackStack*> G4StackManager::additionalWaitingStacks
private

◆ numberOfAdditionalWaitingStacks

G4int G4StackManager::numberOfAdditionalWaitingStacks = 0
private

◆ postponeStack

G4TrackStack* G4StackManager::postponeStack = nullptr
private

◆ theMessenger

G4StackingMessenger* G4StackManager::theMessenger = nullptr
private

Definition at line 136 of file G4StackManager.hh.

Referenced by G4StackManager(), and ~G4StackManager().

◆ urgentStack

G4TrackStack* G4StackManager::urgentStack = nullptr
private

◆ userStackingAction

G4UserStackingAction* G4StackManager::userStackingAction = nullptr
private

◆ verboseLevel

G4int G4StackManager::verboseLevel = 0
private

◆ waitingStack

G4TrackStack* G4StackManager::waitingStack = nullptr
private

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