Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4StackManager.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// G4StackManager class implementation
27//
28// Author: Makoto Asai, 1996
29// --------------------------------------------------------------------
30
31#include "G4StackManager.hh"
33#include "G4VTrajectory.hh"
34#include "G4ios.hh"
35
37#include "G4VProcess.hh"
38
39// Needed for temporal service
40//
41#include "G4ParticleTable.hh"
42#include "G4ProcessManager.hh"
43#include "G4ProcessVector.hh"
44
46{
47 theMessenger = new G4StackingMessenger(this);
48#ifdef G4_USESMARTSTACK
49 urgentStack = new G4SmartTrackStack;
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}
58
60{
61 delete userStackingAction;
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;
75 if(numberOfAdditionalWaitingStacks>0)
76 {
77 for(G4int i=0; i<numberOfAdditionalWaitingStacks; ++i)
78 {
79 delete additionalWaitingStacks[i];
80 }
81 }
82}
83
85PushOneTrack(G4Track* newTrack, G4VTrajectory* newTrajectory)
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 != nullptr)
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;
151 if(i<1 || i>numberOfAdditionalWaitingStacks)
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}
167
168
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
188 waitingStack->TransferTo(urgentStack);
189 if(numberOfAdditionalWaitingStacks>0)
190 {
191 for(G4int i=0; i<numberOfAdditionalWaitingStacks; ++i)
192 {
193 if(i==0)
194 {
195 additionalWaitingStacks[0]->TransferTo(waitingStack);
196 }
197 else
198 {
199 additionalWaitingStacks[i]->TransferTo(additionalWaitingStacks[i-1]);
200 }
201 }
202 }
203 if(userStackingAction != nullptr)
204 {
205 userStackingAction->NewStage();
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 nullptr;
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}
235
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;
267 if(i<1||i>numberOfAdditionalWaitingStacks)
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}
282
284{
285 if(userStackingAction != nullptr)
286 {
287 userStackingAction->PrepareNewEvent();
288 }
289
290 // Set the urgentStack in a defined state. Not doing it would
291 // affect reproducibility
292 //
293 urgentStack->clearAndDestroy();
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;
318 if(userStackingAction != nullptr)
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;
348 if(i<1||i>numberOfAdditionalWaitingStacks)
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}
366
368{
369 if(iAdd > numberOfAdditionalWaitingStacks)
370 {
371 for(G4int i=numberOfAdditionalWaitingStacks; i<iAdd; ++i)
372 {
373 auto* newStack = new G4TrackStack;
374 additionalWaitingStacks.push_back(newStack);
375 }
376 numberOfAdditionalWaitingStacks = iAdd;
377 }
378 else if (iAdd < numberOfAdditionalWaitingStacks)
379 {
380 for(G4int i=numberOfAdditionalWaitingStacks; i>iAdd; --i)
381 {
382 delete additionalWaitingStacks[i];
383 }
384 }
385}
386
389 G4ClassificationOfNewTrack destination)
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;
407 if(i<=numberOfAdditionalWaitingStacks)
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 {
422 urgentStack->clearAndDestroy();
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;
441 if(i<=numberOfAdditionalWaitingStacks)
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}
465
468 G4ClassificationOfNewTrack destination)
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;
486 if(i<=numberOfAdditionalWaitingStacks)
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() != 0u) )
497 {
498 aStackedTrack = originStack->PopFromStack();
499 delete aStackedTrack.GetTrack();
500 delete aStackedTrack.GetTrajectory();
501 }
502 else if (urgentStack->GetNTrack() != 0u )
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;
525 if(i<=numberOfAdditionalWaitingStacks)
526 {
527 targetStack = additionalWaitingStacks[i-1];
528 }
529 break;
530 }
531 if((originStack != nullptr) && (originStack->GetNTrack() != 0u))
532 {
533 aStackedTrack = originStack->PopFromStack();
534 if(targetStack != nullptr) { targetStack->PushToStack(aStackedTrack); }
535 else { urgentStack->PushToStack(aStackedTrack); }
536 }
537 else if(urgentStack->GetNTrack() != 0u)
538 {
539 aStackedTrack = urgentStack->PopFromStack();
540 if(targetStack != nullptr) { targetStack->PushToStack(aStackedTrack); }
541 else { urgentStack->PushToStack(aStackedTrack); }
542 }
543 }
544 return;
545}
546
548{
551 for(G4int i=1; i<=numberOfAdditionalWaitingStacks; ++i)
552 {
554 }
555}
556
558{
559 urgentStack->clearAndDestroy();
560}
561
563{
564 if(i==0)
565 {
566 waitingStack->clearAndDestroy();
567 }
568 else
569 {
570 if(i<=numberOfAdditionalWaitingStacks)
571 {
572 additionalWaitingStacks[i-1]->clearAndDestroy();
573 }
574 }
575}
576
578{
579 postponeStack->clearAndDestroy();
580}
581
583{
584 std::size_t n = urgentStack->GetNTrack()
585 + waitingStack->GetNTrack()
586 + postponeStack->GetNTrack();
587 for(G4int i=1; i<=numberOfAdditionalWaitingStacks; ++i)
588 {
589 n += additionalWaitingStacks[i-1]->GetNTrack();
590 }
591 return G4int(n);
592}
593
595{
596 return (G4int)urgentStack->GetNTrack();
597}
598
600{
601 if(i==0)
602 {
603 return (G4int)waitingStack->GetNTrack();
604 }
605
606 if(i<=numberOfAdditionalWaitingStacks)
607 {
608 return (G4int)additionalWaitingStacks[i-1]->GetNTrack();
609 }
610
611 return 0;
612}
613
615{
616 return (G4int)postponeStack->GetNTrack();
617}
618
620{
621 verboseLevel = value;
622}
623
625{
626 userStackingAction = value;
627 if(userStackingAction != nullptr)
628 {
629 userStackingAction->SetStackManager(this);
630 }
631}
632
633G4ClassificationOfNewTrack G4StackManager::
634DefaultClassification(G4Track* aTrack)
635{
636 G4ClassificationOfNewTrack classification = fUrgent;
637 if( aTrack->GetTrackStatus() == fPostponeToNextEvent )
638 {
639 classification = fPostpone;
640 }
641 return classification;
642}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
@ fPostponeToNextEvent
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4int GetParticleDefinitionID() const
const G4String & GetParticleName() const
G4int GetNTotalTrack() const
void TransferOneStackedTrack(G4ClassificationOfNewTrack origin, G4ClassificationOfNewTrack destination)
G4int GetNUrgentTrack() const
G4int GetNPostponedTrack() const
void SetNumberOfAdditionalWaitingStacks(G4int iAdd)
void SetVerboseLevel(G4int const value)
G4Track * PopNextTrack(G4VTrajectory **newTrajectory)
G4int PushOneTrack(G4Track *newTrack, G4VTrajectory *newTrajectory=nullptr)
void ClearWaitingStack(G4int i=0)
G4int PrepareNewEvent()
void TransferStackedTracks(G4ClassificationOfNewTrack origin, G4ClassificationOfNewTrack destination)
void SetUserStackingAction(G4UserStackingAction *value)
void ClearPostponeStack()
G4int GetNWaitingTrack(G4int i=0) const
G4Track * GetTrack() const
G4VTrajectory * GetTrajectory() const
void PushToStack(const G4StackedTrack &aStackedTrack)
Definition: G4TrackStack.hh:58
void clearAndDestroy()
Definition: G4TrackStack.cc:41
void TransferTo(G4TrackStack *aStack)
Definition: G4TrackStack.cc:51
std::size_t GetNTrack() const
Definition: G4TrackStack.hh:67
G4StackedTrack PopFromStack()
Definition: G4TrackStack.hh:60
std::size_t GetMaxNTrack() const
Definition: G4TrackStack.hh:68
G4TrackStatus GetTrackStatus() const
G4int GetTrackID() const
const G4ParticleDefinition * GetParticleDefinition() const
const G4VProcess * GetCreatorProcess() const
void SetTrackID(const G4int aValue)
G4int GetParentID() const
void SetParentID(const G4int aValue)
virtual G4ClassificationOfNewTrack ClassifyNewTrack(const G4Track *aTrack)
virtual void PrepareNewEvent()
void SetStackManager(G4StackManager *value)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386