Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 //
27 // $Id: G4StackManager.cc 73760 2013-09-10 12:49:34Z gcosmo $
28 //
29 //
30 // Last Modification : 09/Dec/96 M.Asai
31 //
32 
33 #include "G4StackManager.hh"
34 #include "G4StackingMessenger.hh"
35 #include "G4VTrajectory.hh"
36 #include "evmandefs.hh"
37 #include "G4ios.hh"
38 
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 }
53 
55 {
56  if(userStackingAction) delete userStackingAction;
57 
58 #ifdef G4VERBOSE
59  if(verboseLevel>0)
60  {
61  G4cout << "++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++" << G4endl;
62  G4cout << " Maximum number of tracks in the urgent stack : " << urgentStack->GetMaxNTrack() << G4endl;
63  G4cout << "++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++" << G4endl;
64  }
65 #endif
66  delete urgentStack;
67  delete waitingStack;
68  delete postponeStack;
69  delete theMessenger;
70  if(numberOfAdditionalWaitingStacks>0) {
71  for(int i=0;i<numberOfAdditionalWaitingStacks;i++) {
72  delete additionalWaitingStacks[i];
73  }
74  }
75 }
76 
77 const G4StackManager & G4StackManager::operator=
78 (const G4StackManager &) { return *this; }
79 G4int G4StackManager::operator==(const G4StackManager &)
80 const{ return false; }
81 G4int G4StackManager::operator!=(const G4StackManager &)
82 const{ return true; }
83 
84 #include "G4ParticleDefinition.hh"
85 #include "G4VProcess.hh"
86 
87 //Needed for temporal service
88 #include "G4ParticleTable.hh"
89 #include "G4ProcessManager.hh"
90 #include "G4ProcessVector.hh"
91 
93 {
94  const G4ParticleDefinition* pd = newTrack->GetParticleDefinition();
95  if(pd->GetParticleDefinitionID() < 0)
96  {
97 #ifdef G4VERBOSE
99  if(verboseLevel>0) {
100  ED << "A track without proper process manager is pushed into the track stack.\n"
101  << " Particle name : " << pd->GetParticleName() << " -- ";
102  if(newTrack->GetParentID()<0)
103  { ED << "created by a primary particle generator."; }
104  else
105  {
106  const G4VProcess* vp = newTrack->GetCreatorProcess();
107  if(vp)
108  { ED << "created by " << vp->GetProcessName() << "."; }
109  else
110  { ED << "creaded by unknown process."; }
111  }
112  }
113 #endif
114 //// Temporal care of setting process manager for general ion.
115  if(pd->IsGeneralIon())
116  {
117 #ifdef G4VERBOSE
118  if( verboseLevel > 0 ) {
119  ED << "\n Process manager is temporally set, but this operation is thread-unsafe\n"
120  << "and will be replaced with other methods at version 10.0.";
121  G4Exception("G4StackManager::PushOneTrack","Event10051",JustWarning,ED);
122  }
123 #endif
125  G4ProcessManager* pman=0;
126  if (genericIon!=0) pman = genericIon->GetProcessManager();
127  if ((genericIon ==0) || (pman==0)){
128  G4Exception( "G4IonTable::AddProcessManager()","PART10052", FatalException,
129  "Can not define process manager. GenericIon is not available.");
130  }
131  G4ParticleDefinition* ion = const_cast<G4ParticleDefinition*>(pd);
133 #ifdef G4VERBOSE
134  if( verboseLevel > 1 )
135  {
136  G4ProcessManager* ionPman = ion->GetProcessManager();
137  G4cout << "Now " << ion->GetParticleName() << " has a process manaegr at " << ionPman
138  << " that is equivalent to " << pman << G4endl;
139  G4ProcessVector* ionPvec = ionPman->GetProcessList();
140  for(G4int ip1=0;ip1<ionPvec->size();ip1++)
141  {
142  G4cout << " " << ip1 << " - " << (*ionPvec)[ip1]->GetProcessName()
143  << " AtRest " << ionPman->GetAtRestIndex((*ionPvec)[ip1])
144  << ", AlongStep " << ionPman->GetAlongStepIndex((*ionPvec)[ip1])
145  << ", PostStep " << ionPman->GetPostStepIndex((*ionPvec)[ip1])
146  << G4endl;
147  }
148  }
149 #endif
150  }
151 //// End of temporal care of setting process manager
152  else
153  {
154 #ifdef G4VERBOSE
155  if( verboseLevel > 0 ) {
156  ED << "\nThis track is deleted.";
157  G4Exception("G4StackManager::PushOneTrack","Event10051",
158  JustWarning,ED);
159  }
160 #endif
161  delete newTrack;
162  return GetNUrgentTrack();
163  }
164  }
165 
166  G4ClassificationOfNewTrack classification = DefaultClassification( newTrack );
167  if(userStackingAction)
168  { classification = userStackingAction->ClassifyNewTrack( newTrack ); }
169 
170  if(classification==fKill) // delete newTrack without stacking
171  {
172 #ifdef G4VERBOSE
173  if( verboseLevel > 1 )
174  {
175  G4cout << " ---> G4Track " << newTrack << " (trackID "
176  << newTrack->GetTrackID() << ", parentID "
177  << newTrack->GetParentID() << ") is not to be stored." << G4endl;
178  }
179 #endif
180  delete newTrack;
181  delete newTrajectory;
182  }
183  else
184  {
185  G4StackedTrack newStackedTrack( newTrack, newTrajectory );
186  switch (classification)
187  {
188  case fUrgent:
189  urgentStack->PushToStack( newStackedTrack );
190  break;
191  case fWaiting:
192  waitingStack->PushToStack( newStackedTrack );
193  break;
194  case fPostpone:
195  postponeStack->PushToStack( newStackedTrack );
196  break;
197  default:
198  G4int i = classification - 10;
199  if(i<1||i>numberOfAdditionalWaitingStacks) {
201  ED << "invalid classification " << classification << G4endl;
202  G4Exception("G4StackManager::PushOneTrack","Event0051",
203  FatalException,ED);
204  } else {
205  additionalWaitingStacks[i-1]->PushToStack( newStackedTrack );
206  }
207  break;
208  }
209  }
210 
211  return GetNUrgentTrack();
212 }
213 
214 
216 {
217 #ifdef G4VERBOSE
218  if( verboseLevel > 1 )
219  {
220  G4cout << "### pop requested out of "
221  << GetNUrgentTrack() << " stacked tracks." << G4endl;
222  }
223 #endif
224 
225  while( GetNUrgentTrack() == 0 )
226  {
227 #ifdef G4VERBOSE
228  if( verboseLevel > 1 ) G4cout << "### " << GetNWaitingTrack()
229  << " waiting tracks are re-classified to" << G4endl;
230 #endif
231  waitingStack->TransferTo(urgentStack);
232  if(numberOfAdditionalWaitingStacks>0) {
233  for(int i=0;i<numberOfAdditionalWaitingStacks;i++) {
234  if(i==0) {
235  additionalWaitingStacks[0]->TransferTo(waitingStack);
236  } else {
237  additionalWaitingStacks[i]->TransferTo(additionalWaitingStacks[i-1]);
238  }
239  }
240  }
241  if(userStackingAction) userStackingAction->NewStage();
242 #ifdef G4VERBOSE
243  if( verboseLevel > 1 ) G4cout << " " << GetNUrgentTrack()
244  << " urgent tracks and " << GetNWaitingTrack()
245  << " waiting tracks." << G4endl;
246 #endif
247  if( ( GetNUrgentTrack()==0 ) && ( GetNWaitingTrack()==0 ) ) return 0;
248  }
249 
250  G4StackedTrack selectedStackedTrack = urgentStack->PopFromStack();
251  G4Track * selectedTrack = selectedStackedTrack.GetTrack();
252  *newTrajectory = selectedStackedTrack.GetTrajectory();
253 
254 #ifdef G4VERBOSE
255  if( verboseLevel > 2 )
256  {
257  G4cout << "Selected G4StackedTrack : " << &selectedStackedTrack
258  << " with G4Track " << selectedStackedTrack.GetTrack()
259  << " (trackID " << selectedStackedTrack.GetTrack()->GetTrackID()
260  << ", parentID " << selectedStackedTrack.GetTrack()->GetParentID()
261  << ")" << G4endl;
262  }
263 #endif
264 
265  return selectedTrack;
266 }
267 
269 {
270  G4StackedTrack aStackedTrack;
271  G4TrackStack tmpStack;
272 
273  if( !userStackingAction ) return;
274  if( GetNUrgentTrack() == 0 ) return;
275 
276  urgentStack->TransferTo(&tmpStack);
277  while( tmpStack.GetNTrack() > 0 )
278  {
279  aStackedTrack=tmpStack.PopFromStack();
280  G4ClassificationOfNewTrack classification =
281  userStackingAction->ClassifyNewTrack( aStackedTrack.GetTrack() );
282  switch (classification)
283  {
284  case fKill:
285  delete aStackedTrack.GetTrack();
286  delete aStackedTrack.GetTrajectory();
287  break;
288  case fUrgent:
289  urgentStack->PushToStack( aStackedTrack );
290  break;
291  case fWaiting:
292  waitingStack->PushToStack( aStackedTrack );
293  break;
294  case fPostpone:
295  postponeStack->PushToStack( aStackedTrack );
296  break;
297  default:
298  G4int i = classification - 10;
299  if(i<1||i>numberOfAdditionalWaitingStacks) {
301  ED << "invalid classification " << classification << G4endl;
302  G4Exception("G4StackManager::ReClassify","Event0052",
303  FatalException,ED);
304  } else {
305  additionalWaitingStacks[i-1]->PushToStack( aStackedTrack );
306  }
307  break;
308  }
309  }
310 }
311 
313 {
314  if(userStackingAction) userStackingAction->PrepareNewEvent();
315 
316  urgentStack->clearAndDestroy(); // Set the urgentStack in a defined state. Not doing it would affect reproducibility.
317 
318  G4int n_passedFromPrevious = 0;
319 
320  if( GetNPostponedTrack() > 0 )
321  {
322 #ifdef G4VERBOSE
323  if( verboseLevel > 1 )
324  {
326  << " postponed tracked are now shifted to the stack." << G4endl;
327  }
328 #endif
329 
330  G4StackedTrack aStackedTrack;
331  G4TrackStack tmpStack;
332 
333  postponeStack->TransferTo(&tmpStack);
334 
335  while( tmpStack.GetNTrack() > 0 )
336  {
337  aStackedTrack=tmpStack.PopFromStack();
338  G4Track* aTrack = aStackedTrack.GetTrack();
339  aTrack->SetParentID(-1);
340  G4ClassificationOfNewTrack classification;
341  if(userStackingAction)
342  { classification = userStackingAction->ClassifyNewTrack( aTrack ); }
343  else
344  { classification = DefaultClassification( aTrack ); }
345 
346  if(classification==fKill)
347  {
348  delete aTrack;
349  delete aStackedTrack.GetTrajectory();
350  }
351  else
352  {
353  aTrack->SetTrackID(-(++n_passedFromPrevious));
354  switch (classification)
355  {
356  case fUrgent:
357  urgentStack->PushToStack( aStackedTrack );
358  break;
359  case fWaiting:
360  waitingStack->PushToStack( aStackedTrack );
361  break;
362  case fPostpone:
363  postponeStack->PushToStack( aStackedTrack );
364  break;
365  default:
366  G4int i = classification - 10;
367  if(i<1||i>numberOfAdditionalWaitingStacks) {
369  ED << "invalid classification " << classification << G4endl;
370  G4Exception("G4StackManager::PrepareNewEvent","Event0053",
371  FatalException,ED);
372  } else {
373  additionalWaitingStacks[i-1]->PushToStack( aStackedTrack );
374  }
375  break;
376  }
377  }
378  }
379  }
380 
381  return n_passedFromPrevious;
382 }
383 
385 {
386  if(iAdd > numberOfAdditionalWaitingStacks)
387  {
388  for(int i=numberOfAdditionalWaitingStacks;i<iAdd;i++)
389  {
390  G4TrackStack* newStack = new G4TrackStack;
391  additionalWaitingStacks.push_back(newStack);
392  }
393  numberOfAdditionalWaitingStacks = iAdd;
394  }
395  else if (iAdd < numberOfAdditionalWaitingStacks)
396  {
397  for(int i=numberOfAdditionalWaitingStacks;i>iAdd;i--)
398  {
399  delete additionalWaitingStacks[i];
400  }
401  }
402 }
403 
405 {
406  if(origin==destination) return;
407  if(origin==fKill) return;
408  G4TrackStack* originStack = 0;
409  switch(origin)
410  {
411  case fUrgent:
412  originStack = 0;
413  break;
414  case fWaiting:
415  originStack = waitingStack;
416  break;
417  case fPostpone:
418  originStack = postponeStack;
419  break;
420  default:
421  int i = origin - 10;
422  if(i<=numberOfAdditionalWaitingStacks) originStack = additionalWaitingStacks[i-1];
423  break;
424  }
425 
426  if(destination==fKill)
427  {
428  if(originStack)
429  { originStack->clearAndDestroy(); }
430  else
431  { urgentStack->clearAndDestroy(); }
432  }
433  else
434  {
435  G4TrackStack* targetStack = 0;
436  switch(destination)
437  {
438  case fUrgent:
439  targetStack = 0;
440  break;
441  case fWaiting:
442  targetStack = waitingStack;
443  break;
444  case fPostpone:
445  targetStack = postponeStack;
446  break;
447  default:
448  int i = destination - 10;
449  if(i<=numberOfAdditionalWaitingStacks) targetStack = additionalWaitingStacks[i-1];
450  break;
451  }
452  if(originStack)
453  {
454  if(targetStack)
455  { originStack->TransferTo(targetStack); }
456  else
457  { originStack->TransferTo(urgentStack); }
458  }
459  else
460  { urgentStack->TransferTo(targetStack); }
461  }
462  return;
463 }
464 
466 {
467  if(origin==destination) return;
468  if(origin==fKill) return;
469  G4TrackStack* originStack = 0;
470  switch(origin)
471  {
472  case fUrgent:
473  originStack = 0;
474  break;
475  case fWaiting:
476  originStack = waitingStack;
477  break;
478  case fPostpone:
479  originStack = postponeStack;
480  break;
481  default:
482  int i = origin - 10;
483  if(i<=numberOfAdditionalWaitingStacks) originStack = additionalWaitingStacks[i-1];
484  break;
485  }
486 
487  G4StackedTrack aStackedTrack;
488  if(destination==fKill)
489  {
490  if( originStack && originStack->GetNTrack() ) {
491  aStackedTrack = originStack->PopFromStack();
492  delete aStackedTrack.GetTrack();
493  delete aStackedTrack.GetTrajectory();
494  }
495  else if (urgentStack->GetNTrack() ) {
496  aStackedTrack = urgentStack->PopFromStack();
497  delete aStackedTrack.GetTrack();
498  delete aStackedTrack.GetTrajectory();
499  }
500  }
501  else
502  {
503  G4TrackStack* targetStack = 0;
504  switch(destination)
505  {
506  case fUrgent:
507  targetStack = 0;
508  break;
509  case fWaiting:
510  targetStack = waitingStack;
511  break;
512  case fPostpone:
513  targetStack = postponeStack;
514  break;
515  default:
516  int i = destination - 10;
517  if(i<=numberOfAdditionalWaitingStacks) targetStack = additionalWaitingStacks[i-1];
518  break;
519  }
520  if(originStack && originStack->GetNTrack()) {
521  aStackedTrack = originStack->PopFromStack();
522  if(targetStack) { targetStack->PushToStack(aStackedTrack); }
523  else { urgentStack->PushToStack(aStackedTrack); }
524  }
525  else if(urgentStack->GetNTrack()) {
526  aStackedTrack = urgentStack->PopFromStack();
527  if(targetStack) { targetStack->PushToStack(aStackedTrack); }
528  else { urgentStack->PushToStack(aStackedTrack); }
529  }
530  }
531  return;
532 }
533 
535 {
538  for(int i=1;i<=numberOfAdditionalWaitingStacks;i++) {ClearWaitingStack(i);}
539 }
540 
542 {
543  urgentStack->clearAndDestroy();
544 }
545 
547 {
548  if(i==0) {
549  waitingStack->clearAndDestroy();
550  } else {
551  if(i<=numberOfAdditionalWaitingStacks) additionalWaitingStacks[i-1]->clearAndDestroy();
552  }
553 }
554 
556 {
557  postponeStack->clearAndDestroy();
558 }
559 
561 {
562  int n = urgentStack->GetNTrack() + waitingStack->GetNTrack() + postponeStack->GetNTrack();
563  for(int i=1;i<=numberOfAdditionalWaitingStacks;i++) {n += additionalWaitingStacks[i-1]->GetNTrack();}
564  return n;
565 }
566 
568 {
569  return urgentStack->GetNTrack();
570 }
571 
573 {
574  if(i==0) { return waitingStack->GetNTrack(); }
575  else {
576  if(i<=numberOfAdditionalWaitingStacks) { return additionalWaitingStacks[i-1]->GetNTrack();}
577  }
578  return 0;
579 }
580 
582 {
583  return postponeStack->GetNTrack();
584 }
585 
587 {
588  verboseLevel = value;
589 }
590 
592 {
593  userStackingAction = value;
594  if(userStackingAction) userStackingAction->SetStackManager(this);
595 }
596 
597 G4ClassificationOfNewTrack G4StackManager::DefaultClassification(G4Track *aTrack)
598 {
599  G4ClassificationOfNewTrack classification = fUrgent;
600  if( aTrack->GetTrackStatus() == fPostponeToNextEvent )
601  { classification = fPostpone; }
602  return classification;
603 }
604 
605 
606 
607 
G4Track * GetTrack() const
G4int GetParticleDefinitionID() const
G4int GetParentID() const
void SetParticleDefinitionID(G4int id=-1)
G4int GetAlongStepIndex(G4VProcess *aProcess, G4ProcessVectorTypeIndex typ=typeGPIL) const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void SetStackManager(G4StackManager *value)
void PushToStack(const G4StackedTrack &aStackedTrack)
Definition: G4TrackStack.hh:62
G4int GetNTotalTrack() const
G4TrackStatus GetTrackStatus() const
G4Track * PopNextTrack(G4VTrajectory **newTrajectory)
G4ParticleDefinition * GetGenericIon() const
G4ProcessManager * GetProcessManager() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4int GetAtRestIndex(G4VProcess *aProcess, G4ProcessVectorTypeIndex typ=typeGPIL) const
G4int GetNUrgentTrack() const
G4bool IsGeneralIon() const
G4int GetNTrack() const
Definition: G4TrackStack.hh:74
void TransferTo(G4TrackStack *aStack)
Definition: G4TrackStack.cc:49
const G4VProcess * GetCreatorProcess() const
void SetUserStackingAction(G4UserStackingAction *value)
G4GLOB_DLL std::ostream G4cout
G4VTrajectory * GetTrajectory() const
G4StackedTrack PopFromStack()
Definition: G4TrackStack.hh:63
G4int PushOneTrack(G4Track *newTrack, G4VTrajectory *newTrajectory=0)
const G4ParticleDefinition * GetParticleDefinition() const
G4int PrepareNewEvent()
virtual void PrepareNewEvent()
G4int GetTrackID() const
const G4int n
void SetVerboseLevel(G4int const value)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void TransferOneStackedTrack(G4ClassificationOfNewTrack origin, G4ClassificationOfNewTrack destination)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4int size() const
void TransferStackedTracks(G4ClassificationOfNewTrack origin, G4ClassificationOfNewTrack destination)
G4int GetNWaitingTrack(int i=0) const
void ClearPostponeStack()
virtual G4ClassificationOfNewTrack ClassifyNewTrack(const G4Track *aTrack)
static G4ParticleTable * GetParticleTable()
void SetParentID(const G4int aValue)
G4int GetPostStepIndex(G4VProcess *aProcess, G4ProcessVectorTypeIndex typ=typeGPIL) const
void clearAndDestroy()
Definition: G4TrackStack.cc:40
void SetNumberOfAdditionalWaitingStacks(G4int iAdd)
const XML_Char int const XML_Char * value
G4int GetNPostponedTrack() const
#define G4endl
Definition: G4ios.hh:61
G4int GetMaxNTrack() const
Definition: G4TrackStack.hh:75
void SetTrackID(const G4int aValue)
G4ProcessVector * GetProcessList() const
void ClearWaitingStack(int i=0)