Geant4 9.6.0
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 ()
 
G4int PushOneTrack (G4Track *newTrack, G4VTrajectory *newTrajectory=0)
 
G4TrackPopNextTrack (G4VTrajectory **newTrajectory)
 
G4int PrepareNewEvent ()
 
void ReClassify ()
 
void SetNumberOfAdditionalWaitingStacks (G4int iAdd)
 
void TransferStackedTracks (G4ClassificationOfNewTrack origin, G4ClassificationOfNewTrack destination)
 
void TransferOneStackedTrack (G4ClassificationOfNewTrack origin, G4ClassificationOfNewTrack destination)
 
void clear ()
 
void ClearUrgentStack ()
 
void ClearWaitingStack (int i=0)
 
void ClearPostponeStack ()
 
G4int GetNTotalTrack () const
 
G4int GetNUrgentTrack () const
 
G4int GetNWaitingTrack (int i=0) const
 
G4int GetNPostponedTrack () const
 
void SetVerboseLevel (G4int const value)
 
void SetUserStackingAction (G4UserStackingAction *value)
 

Detailed Description

Definition at line 62 of file G4StackManager.hh.

Constructor & Destructor Documentation

◆ G4StackManager()

G4StackManager::G4StackManager ( )

Definition at line 39 of file G4StackManager.cc.

40:userStackingAction(0),verboseLevel(0),numberOfAdditionalWaitingStacks(0)
41{
42 theMessenger = new G4StackingMessenger(this);
43#ifdef G4_USESMARTSTACK
44 urgentStack = new G4SmartTrackStack;
45 // G4cout<<"+++ G4StackManager uses G4SmartTrackStack. +++"<<G4endl;
46#else
47 urgentStack = new G4TrackStack(5000);
48// G4cout<<"+++ G4StackManager uses ordinary G4TrackStack. +++"<<G4endl;
49#endif
50 waitingStack = new G4TrackStack(1000);
51 postponeStack = new G4TrackStack(1000);
52}

◆ ~G4StackManager()

G4StackManager::~G4StackManager ( )

Definition at line 54 of file G4StackManager.cc.

55{
56 if(userStackingAction) delete userStackingAction;
57
58 if(verboseLevel>0)
59 {
60 G4cout << "++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++" << G4endl;
61 G4cout << " Maximum number of tracks in the urgent stack : " << urgentStack->GetMaxNTrack() << G4endl;
62 G4cout << "++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++" << G4endl;
63 }
64 delete urgentStack;
65 delete waitingStack;
66 delete postponeStack;
67 delete theMessenger;
68 if(numberOfAdditionalWaitingStacks>0) {
69 for(int i=0;i<numberOfAdditionalWaitingStacks;i++) {
70 delete additionalWaitingStacks[i];
71 }
72 }
73}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4int GetMaxNTrack() const
Definition: G4TrackStack.hh:75

Member Function Documentation

◆ clear()

void G4StackManager::clear ( )

Definition at line 454 of file G4StackManager.cc.

455{
458 for(int i=1;i<=numberOfAdditionalWaitingStacks;i++) {ClearWaitingStack(i);}
459}
void ClearWaitingStack(int i=0)

Referenced by G4EventManager::AbortCurrentEvent().

◆ ClearPostponeStack()

void G4StackManager::ClearPostponeStack ( )

Definition at line 475 of file G4StackManager.cc.

476{
477 postponeStack->clearAndDestroy();
478}
void clearAndDestroy()
Definition: G4TrackStack.cc:40

Referenced by G4StackingMessenger::SetNewValue().

◆ ClearUrgentStack()

void G4StackManager::ClearUrgentStack ( )

Definition at line 461 of file G4StackManager.cc.

462{
463 urgentStack->clearAndDestroy();
464}

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

◆ ClearWaitingStack()

void G4StackManager::ClearWaitingStack ( int  i = 0)

Definition at line 466 of file G4StackManager.cc.

467{
468 if(i==0) {
469 waitingStack->clearAndDestroy();
470 } else {
471 if(i<=numberOfAdditionalWaitingStacks) additionalWaitingStacks[i-1]->clearAndDestroy();
472 }
473}

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

◆ GetNPostponedTrack()

G4int G4StackManager::GetNPostponedTrack ( ) const

Definition at line 501 of file G4StackManager.cc.

502{
503 return postponeStack->GetNTrack();
504}
G4int GetNTrack() const
Definition: G4TrackStack.hh:74

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

◆ GetNTotalTrack()

G4int G4StackManager::GetNTotalTrack ( ) const

Definition at line 480 of file G4StackManager.cc.

481{
482 int n = urgentStack->GetNTrack() + waitingStack->GetNTrack() + postponeStack->GetNTrack();
483 for(int i=1;i<=numberOfAdditionalWaitingStacks;i++) {n += additionalWaitingStacks[i-1]->GetNTrack();}
484 return n;
485}

◆ GetNUrgentTrack()

G4int G4StackManager::GetNUrgentTrack ( ) const

Definition at line 487 of file G4StackManager.cc.

488{
489 return urgentStack->GetNTrack();
490}

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

◆ GetNWaitingTrack()

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

Definition at line 492 of file G4StackManager.cc.

493{
494 if(i==0) { return waitingStack->GetNTrack(); }
495 else {
496 if(i<=numberOfAdditionalWaitingStacks) { return additionalWaitingStacks[i-1]->GetNTrack();}
497 }
498 return 0;
499}

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

◆ PopNextTrack()

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

Definition at line 135 of file G4StackManager.cc.

136{
137#ifdef G4VERBOSE
138 if( verboseLevel > 1 )
139 {
140 G4cout << "### pop requested out of "
141 << GetNUrgentTrack() << " stacked tracks." << G4endl;
142 }
143#endif
144
145 while( GetNUrgentTrack() == 0 )
146 {
147#ifdef G4VERBOSE
148 if( verboseLevel > 1 ) G4cout << "### " << GetNWaitingTrack()
149 << " waiting tracks are re-classified to" << G4endl;
150#endif
151 waitingStack->TransferTo(urgentStack);
152 if(numberOfAdditionalWaitingStacks>0) {
153 for(int i=0;i<numberOfAdditionalWaitingStacks;i++) {
154 if(i==0) {
155 additionalWaitingStacks[0]->TransferTo(waitingStack);
156 } else {
157 additionalWaitingStacks[i]->TransferTo(additionalWaitingStacks[i-1]);
158 }
159 }
160 }
161 if(userStackingAction) userStackingAction->NewStage();
162#ifdef G4VERBOSE
163 if( verboseLevel > 1 ) G4cout << " " << GetNUrgentTrack()
164 << " urgent tracks and " << GetNWaitingTrack()
165 << " waiting tracks." << G4endl;
166#endif
167 if( ( GetNUrgentTrack()==0 ) && ( GetNWaitingTrack()==0 ) ) return 0;
168 }
169
170 G4StackedTrack selectedStackedTrack = urgentStack->PopFromStack();
171 G4Track * selectedTrack = selectedStackedTrack.GetTrack();
172 *newTrajectory = selectedStackedTrack.GetTrajectory();
173
174#ifdef G4VERBOSE
175 if( verboseLevel > 2 )
176 {
177 G4cout << "Selected G4StackedTrack : " << &selectedStackedTrack
178 << " with G4Track " << selectedStackedTrack.GetTrack()
179 << " (trackID " << selectedStackedTrack.GetTrack()->GetTrackID()
180 << ", parentID " << selectedStackedTrack.GetTrack()->GetParentID()
181 << ")" << G4endl;
182 }
183#endif
184
185 return selectedTrack;
186}
G4int GetNUrgentTrack() const
G4int GetNWaitingTrack(int i=0) const
G4Track * GetTrack() const
G4VTrajectory * GetTrajectory() const
void TransferTo(G4TrackStack *aStack)
Definition: G4TrackStack.cc:49
G4StackedTrack PopFromStack()
Definition: G4TrackStack.hh:63
G4int GetTrackID() const
G4int GetParentID() const

◆ PrepareNewEvent()

G4int G4StackManager::PrepareNewEvent ( )

Definition at line 232 of file G4StackManager.cc.

233{
234 if(userStackingAction) userStackingAction->PrepareNewEvent();
235
236 urgentStack->clearAndDestroy(); // Set the urgentStack in a defined state. Not doing it would affect reproducibility.
237
238 G4int n_passedFromPrevious = 0;
239
240 if( GetNPostponedTrack() > 0 )
241 {
242#ifdef G4VERBOSE
243 if( verboseLevel > 1 )
244 {
246 << " postponed tracked are now shifted to the stack." << G4endl;
247 }
248#endif
249
250 G4StackedTrack aStackedTrack;
251 G4TrackStack tmpStack;
252
253 postponeStack->TransferTo(&tmpStack);
254
255 while( tmpStack.GetNTrack() > 0 )
256 {
257 aStackedTrack=tmpStack.PopFromStack();
258 G4Track* aTrack = aStackedTrack.GetTrack();
259 aTrack->SetParentID(-1);
260 G4ClassificationOfNewTrack classification;
261 if(userStackingAction)
262 { classification = userStackingAction->ClassifyNewTrack( aTrack ); }
263 else
264 { classification = DefaultClassification( aTrack ); }
265
266 if(classification==fKill)
267 {
268 delete aTrack;
269 delete aStackedTrack.GetTrajectory();
270 }
271 else
272 {
273 aTrack->SetTrackID(-(++n_passedFromPrevious));
274 switch (classification)
275 {
276 case fUrgent:
277 urgentStack->PushToStack( aStackedTrack );
278 break;
279 case fWaiting:
280 waitingStack->PushToStack( aStackedTrack );
281 break;
282 case fPostpone:
283 postponeStack->PushToStack( aStackedTrack );
284 break;
285 default:
286 G4int i = classification - 10;
287 if(i<1||i>numberOfAdditionalWaitingStacks) {
289 ED << "invalid classification " << classification << G4endl;
290 G4Exception("G4StackManager::PrepareNewEvent","Event0053",
291 FatalException,ED);
292 } else {
293 additionalWaitingStacks[i-1]->PushToStack( aStackedTrack );
294 }
295 break;
296 }
297 }
298 }
299 }
300
301 return n_passedFromPrevious;
302}
@ FatalException
int G4int
Definition: G4Types.hh:66
G4int GetNPostponedTrack() const
void PushToStack(const G4StackedTrack &aStackedTrack)
Definition: G4TrackStack.hh:62
void SetTrackID(const G4int aValue)
void SetParentID(const G4int aValue)
virtual G4ClassificationOfNewTrack ClassifyNewTrack(const G4Track *aTrack)
virtual void PrepareNewEvent()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

◆ PushOneTrack()

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

Definition at line 82 of file G4StackManager.cc.

83{
84 G4ClassificationOfNewTrack classification;
85 if(userStackingAction)
86 { classification = userStackingAction->ClassifyNewTrack( newTrack ); }
87 else
88 { classification = DefaultClassification( newTrack ); }
89
90 if(classification==fKill) // delete newTrack without stacking
91 {
92#ifdef G4VERBOSE
93 if( verboseLevel > 1 )
94 {
95 G4cout << " ---> G4Track " << newTrack << " (trackID "
96 << newTrack->GetTrackID() << ", parentID "
97 << newTrack->GetParentID() << ") is not to be stored." << G4endl;
98 }
99#endif
100 delete newTrack;
101 delete newTrajectory;
102 }
103 else
104 {
105 G4StackedTrack newStackedTrack( newTrack, newTrajectory );
106 switch (classification)
107 {
108 case fUrgent:
109 urgentStack->PushToStack( newStackedTrack );
110 break;
111 case fWaiting:
112 waitingStack->PushToStack( newStackedTrack );
113 break;
114 case fPostpone:
115 postponeStack->PushToStack( newStackedTrack );
116 break;
117 default:
118 G4int i = classification - 10;
119 if(i<1||i>numberOfAdditionalWaitingStacks) {
121 ED << "invalid classification " << classification << G4endl;
122 G4Exception("G4StackManager::PushOneTrack","Event0051",
123 FatalException,ED);
124 } else {
125 additionalWaitingStacks[i-1]->PushToStack( newStackedTrack );
126 }
127 break;
128 }
129 }
130
131 return GetNUrgentTrack();
132}

◆ ReClassify()

void G4StackManager::ReClassify ( )

Definition at line 188 of file G4StackManager.cc.

189{
190 G4StackedTrack aStackedTrack;
191 G4TrackStack tmpStack;
192
193 if( !userStackingAction ) return;
194 if( GetNUrgentTrack() == 0 ) return;
195
196 urgentStack->TransferTo(&tmpStack);
197 while( tmpStack.GetNTrack() > 0 )
198 {
199 aStackedTrack=tmpStack.PopFromStack();
200 G4ClassificationOfNewTrack classification =
201 userStackingAction->ClassifyNewTrack( aStackedTrack.GetTrack() );
202 switch (classification)
203 {
204 case fKill:
205 delete aStackedTrack.GetTrack();
206 delete aStackedTrack.GetTrajectory();
207 break;
208 case fUrgent:
209 urgentStack->PushToStack( aStackedTrack );
210 break;
211 case fWaiting:
212 waitingStack->PushToStack( aStackedTrack );
213 break;
214 case fPostpone:
215 postponeStack->PushToStack( aStackedTrack );
216 break;
217 default:
218 G4int i = classification - 10;
219 if(i<1||i>numberOfAdditionalWaitingStacks) {
221 ED << "invalid classification " << classification << G4endl;
222 G4Exception("G4StackManager::ReClassify","Event0052",
223 FatalException,ED);
224 } else {
225 additionalWaitingStacks[i-1]->PushToStack( aStackedTrack );
226 }
227 break;
228 }
229 }
230}

◆ SetNumberOfAdditionalWaitingStacks()

void G4StackManager::SetNumberOfAdditionalWaitingStacks ( G4int  iAdd)

Definition at line 304 of file G4StackManager.cc.

305{
306 if(iAdd > numberOfAdditionalWaitingStacks)
307 {
308 for(int i=numberOfAdditionalWaitingStacks;i<iAdd;i++)
309 {
310 G4TrackStack* newStack = new G4TrackStack;
311 additionalWaitingStacks.push_back(newStack);
312 }
313 numberOfAdditionalWaitingStacks = iAdd;
314 }
315 else if (iAdd < numberOfAdditionalWaitingStacks)
316 {
317 for(int i=numberOfAdditionalWaitingStacks;i>iAdd;i--)
318 {
319 delete additionalWaitingStacks[i];
320 }
321 }
322}

Referenced by G4EventManager::SetNumberOfAdditionalWaitingStacks().

◆ SetUserStackingAction()

void G4StackManager::SetUserStackingAction ( G4UserStackingAction value)

Definition at line 511 of file G4StackManager.cc.

512{
513 userStackingAction = value;
514 if(userStackingAction) userStackingAction->SetStackManager(this);
515}
void SetStackManager(G4StackManager *value)

Referenced by G4EventManager::SetUserAction().

◆ SetVerboseLevel()

void G4StackManager::SetVerboseLevel ( G4int const  value)

Definition at line 506 of file G4StackManager.cc.

507{
508 verboseLevel = value;
509}

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

◆ TransferOneStackedTrack()

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

Definition at line 385 of file G4StackManager.cc.

386{
387 if(origin==destination) return;
388 if(origin==fKill) return;
389 G4TrackStack* originStack = 0;
390 switch(origin)
391 {
392 case fUrgent:
393 originStack = 0;
394 break;
395 case fWaiting:
396 originStack = waitingStack;
397 break;
398 case fPostpone:
399 originStack = postponeStack;
400 break;
401 default:
402 int i = origin - 10;
403 if(i<=numberOfAdditionalWaitingStacks) originStack = additionalWaitingStacks[i-1];
404 break;
405 }
406
407 G4StackedTrack aStackedTrack;
408 if(destination==fKill)
409 {
410 if( originStack && originStack->GetNTrack() ) {
411 aStackedTrack = originStack->PopFromStack();
412 delete aStackedTrack.GetTrack();
413 delete aStackedTrack.GetTrajectory();
414 }
415 else if (urgentStack->GetNTrack() ) {
416 aStackedTrack = urgentStack->PopFromStack();
417 delete aStackedTrack.GetTrack();
418 delete aStackedTrack.GetTrajectory();
419 }
420 }
421 else
422 {
423 G4TrackStack* targetStack = 0;
424 switch(destination)
425 {
426 case fUrgent:
427 targetStack = 0;
428 break;
429 case fWaiting:
430 targetStack = waitingStack;
431 break;
432 case fPostpone:
433 targetStack = postponeStack;
434 break;
435 default:
436 int i = destination - 10;
437 if(i<=numberOfAdditionalWaitingStacks) targetStack = additionalWaitingStacks[i-1];
438 break;
439 }
440 if(originStack && originStack->GetNTrack()) {
441 aStackedTrack = originStack->PopFromStack();
442 if(targetStack) { targetStack->PushToStack(aStackedTrack); }
443 else { urgentStack->PushToStack(aStackedTrack); }
444 }
445 else if(urgentStack->GetNTrack()) {
446 aStackedTrack = urgentStack->PopFromStack();
447 if(targetStack) { targetStack->PushToStack(aStackedTrack); }
448 else { urgentStack->PushToStack(aStackedTrack); }
449 }
450 }
451 return;
452}

◆ TransferStackedTracks()

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

Definition at line 324 of file G4StackManager.cc.

325{
326 if(origin==destination) return;
327 if(origin==fKill) return;
328 G4TrackStack* originStack = 0;
329 switch(origin)
330 {
331 case fUrgent:
332 originStack = 0;
333 break;
334 case fWaiting:
335 originStack = waitingStack;
336 break;
337 case fPostpone:
338 originStack = postponeStack;
339 break;
340 default:
341 int i = origin - 10;
342 if(i<=numberOfAdditionalWaitingStacks) originStack = additionalWaitingStacks[i-1];
343 break;
344 }
345
346 if(destination==fKill)
347 {
348 if(originStack)
349 { originStack->clearAndDestroy(); }
350 else
351 { urgentStack->clearAndDestroy(); }
352 }
353 else
354 {
355 G4TrackStack* targetStack = 0;
356 switch(destination)
357 {
358 case fUrgent:
359 targetStack = 0;
360 break;
361 case fWaiting:
362 targetStack = waitingStack;
363 break;
364 case fPostpone:
365 targetStack = postponeStack;
366 break;
367 default:
368 int i = destination - 10;
369 if(i<=numberOfAdditionalWaitingStacks) targetStack = additionalWaitingStacks[i-1];
370 break;
371 }
372 if(originStack)
373 {
374 if(targetStack)
375 { originStack->TransferTo(targetStack); }
376 else
377 { originStack->TransferTo(urgentStack); }
378 }
379 else
380 { urgentStack->TransferTo(targetStack); }
381 }
382 return;
383}

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