G4SmartTrackStack Class Reference

#include <G4SmartTrackStack.hh>


Public Member Functions

 G4SmartTrackStack ()
 ~G4SmartTrackStack ()
void PushToStack (const G4StackedTrack &aStackedTrack)
G4StackedTrack PopFromStack ()
void clear ()
void clearAndDestroy ()
void TransferTo (G4TrackStack *aStack)
G4double getEnergyOfStack (G4TrackStack *aTrackStack)
void dumpStatistics ()
G4int GetNTrack () const
G4int GetMaxNTrack () const


Detailed Description

Definition at line 46 of file G4SmartTrackStack.hh.


Constructor & Destructor Documentation

G4SmartTrackStack::G4SmartTrackStack (  ) 

Definition at line 45 of file G4SmartTrackStack.cc.

00046         :fTurn(0), nTurn(5), maxNTracks(0), nTracks(0)
00047 {
00048         for(int i=0;i<nTurn;i++)
00049         {
00050           stacks[i] = new G4TrackStack(5000);
00051           energies[i] = 0.;
00052         }
00053 }

G4SmartTrackStack::~G4SmartTrackStack (  ) 

Definition at line 55 of file G4SmartTrackStack.cc.

00056 {
00057         for (int i  = 0; i < nTurn; i++) {
00058                 delete stacks[i];
00059         }
00060 }


Member Function Documentation

void G4SmartTrackStack::clear (  ) 

Definition at line 142 of file G4SmartTrackStack.cc.

00143 {
00144   for (int i = 0; i < nTurn; i++) {
00145     stacks[i]->clear();
00146     energies[i] = 0.0;
00147     fTurn = 0;
00148   }
00149   nTracks = 0;
00150 }

void G4SmartTrackStack::clearAndDestroy (  ) 

Definition at line 152 of file G4SmartTrackStack.cc.

00153 {
00154   for (int i = 0; i < nTurn; i++) {
00155     stacks[i]->clearAndDestroy();
00156     energies[i] = 0.0;
00157     fTurn = 0;
00158   }
00159   nTracks = 0;
00160 }

void G4SmartTrackStack::dumpStatistics (  ) 

Definition at line 34 of file G4SmartTrackStack.cc.

References G4cerr, G4endl, GetNTrack(), and G4TrackStack::getTotalEnergy().

00035 {
00036         // Print to stderr so that we can split stats output from normal
00037         // output of Geant4 which is typically being printed to stdout
00038         for (int i = 0; i < nTurn; i++) {
00039                 G4cerr << stacks[i]->GetNTrack() << " ";
00040                 G4cerr << stacks[i]->getTotalEnergy() << " ";
00041         }
00042         G4cerr << G4endl;
00043 }

G4double G4SmartTrackStack::getEnergyOfStack ( G4TrackStack aTrackStack  ) 

G4int G4SmartTrackStack::GetMaxNTrack (  )  const [inline]

Definition at line 82 of file G4SmartTrackStack.hh.

00082 { return maxNTracks; }

G4int G4SmartTrackStack::GetNTrack (  )  const [inline]

Definition at line 81 of file G4SmartTrackStack.hh.

Referenced by dumpStatistics(), PopFromStack(), and PushToStack().

00081 { return nTracks; }

G4StackedTrack G4SmartTrackStack::PopFromStack (  ) 

Definition at line 83 of file G4SmartTrackStack.cc.

References G4Track::GetDynamicParticle(), GetNTrack(), G4DynamicParticle::GetTotalEnergy(), G4StackedTrack::GetTrack(), and G4TrackStack::PopFromStack().

00084 {
00085         G4StackedTrack aStackedTrack;
00086 
00087         if (nTracks) {
00088                 while (true) {
00089                         if (stacks[fTurn]->GetNTrack()) {
00090                                 aStackedTrack = stacks[fTurn]->PopFromStack();
00091                                 energies[fTurn] -= aStackedTrack.GetTrack()->GetDynamicParticle()->GetTotalEnergy();
00092         nTracks--;
00093         break;
00094                         } else {
00095                                 fTurn = (fTurn+1) % nTurn;
00096                         }
00097                 }
00098         }
00099 
00100         // dumpStatistics();
00101         return aStackedTrack;
00102 }

void G4SmartTrackStack::PushToStack ( const G4StackedTrack aStackedTrack  ) 

Definition at line 108 of file G4SmartTrackStack.cc.

References electronCode, gammaCode, G4Track::GetDynamicParticle(), GetNTrack(), G4TrackStack::GetNTrack(), G4Track::GetParentID(), G4DynamicParticle::GetPDGcode(), G4TrackStack::GetSafetyValve1(), G4TrackStack::GetSafetyValve2(), G4DynamicParticle::GetTotalEnergy(), G4StackedTrack::GetTrack(), neutronCode, positronCode, and G4TrackStack::PushToStack().

Referenced by G4TrackStack::TransferTo().

00109 {
00110 
00111   G4int iDest = 0;
00112   if (aStackedTrack.GetTrack()->GetParentID()) {
00113     G4int code = aStackedTrack.GetTrack()->GetDynamicParticle()->GetPDGcode();
00114     if (code == electronCode)
00115       iDest = 2;
00116     else if (code == gammaCode)
00117       iDest = 3;
00118     else if (code == positronCode)
00119       iDest = 4;
00120     else if (code == neutronCode)
00121       iDest = 1;
00122   } else {
00123     // We have a primary track, which should go first.
00124     fTurn = 0; // reseting the turn
00125   }
00126   stacks[iDest]->PushToStack(aStackedTrack);
00127   energies[iDest] += aStackedTrack.GetTrack()->GetDynamicParticle()->GetTotalEnergy();
00128   nTracks++;
00129   
00130   G4int dy1 = stacks[iDest]->GetNTrack() - stacks[iDest]->GetSafetyValve1();
00131   G4int dy2 = stacks[fTurn]->GetNTrack() - stacks[fTurn]->GetSafetyValve2();
00132   
00133   if (dy1 > 0 || dy1 > dy2 ||
00134       (iDest == 2 &&
00135        stacks[iDest]->GetNTrack() < 50 && energies[iDest] < energies[fTurn])) {
00136         fTurn = iDest;
00137       }
00138   
00139   if (nTracks > maxNTracks) maxNTracks = nTracks;
00140 }

void G4SmartTrackStack::TransferTo ( G4TrackStack aStack  ) 

Definition at line 75 of file G4SmartTrackStack.cc.

00076 {
00077         for (int i = 0; i < nTurn; i++) {
00078                 stacks[i]->TransferTo(aStack);
00079         }
00080   nTracks = 0;
00081 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:23 2013 for Geant4 by  doxygen 1.4.7