Geant4  10.02.p03
G4SmartTrackStack Class Reference

#include <G4SmartTrackStack.hh>

Collaboration diagram for G4SmartTrackStack:

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
 

Private Member Functions

const G4SmartTrackStackoperator= (const G4SmartTrackStack &right)
 
G4int operator== (const G4SmartTrackStack &right) const
 
G4int operator!= (const G4SmartTrackStack &right) const
 
G4int n_stackedTrack () const
 

Private Attributes

G4int fTurn
 
G4int nTurn
 
G4double energies [5]
 
G4TrackStackstacks [5]
 
G4int maxNTracks
 
G4int nTracks
 

Detailed Description

Definition at line 46 of file G4SmartTrackStack.hh.

Constructor & Destructor Documentation

◆ G4SmartTrackStack()

G4SmartTrackStack::G4SmartTrackStack ( )

Definition at line 45 of file G4SmartTrackStack.cc.

46  :fTurn(0), nTurn(5), maxNTracks(0), nTracks(0)
47 {
48  for(int i=0;i<nTurn;i++)
49  {
50  stacks[i] = new G4TrackStack(5000);
51  energies[i] = 0.;
52  }
53 }
G4TrackStack * stacks[5]

◆ ~G4SmartTrackStack()

G4SmartTrackStack::~G4SmartTrackStack ( )

Definition at line 55 of file G4SmartTrackStack.cc.

56 {
57  for (int i = 0; i < nTurn; i++) {
58  delete stacks[i];
59  }
60 }
G4TrackStack * stacks[5]

Member Function Documentation

◆ clear()

void G4SmartTrackStack::clear ( void  )

Definition at line 142 of file G4SmartTrackStack.cc.

143 {
144  for (int i = 0; i < nTurn; i++) {
145  stacks[i]->clear();
146  energies[i] = 0.0;
147  fTurn = 0;
148  }
149  nTracks = 0;
150 }
G4TrackStack * stacks[5]

◆ clearAndDestroy()

void G4SmartTrackStack::clearAndDestroy ( )

Definition at line 152 of file G4SmartTrackStack.cc.

153 {
154  for (int i = 0; i < nTurn; i++) {
155  stacks[i]->clearAndDestroy();
156  energies[i] = 0.0;
157  fTurn = 0;
158  }
159  nTracks = 0;
160 }
G4TrackStack * stacks[5]
void clearAndDestroy()
Definition: G4TrackStack.cc:40
Here is the call graph for this function:

◆ dumpStatistics()

void G4SmartTrackStack::dumpStatistics ( )

Definition at line 34 of file G4SmartTrackStack.cc.

35 {
36  // Print to stderr so that we can split stats output from normal
37  // output of Geant4 which is typically being printed to stdout
38  for (int i = 0; i < nTurn; i++) {
39  G4cerr << stacks[i]->GetNTrack() << " ";
40  G4cerr << stacks[i]->getTotalEnergy() << " ";
41  }
42  G4cerr << G4endl;
43 }
G4TrackStack * stacks[5]
G4double getTotalEnergy(void) const
Definition: G4TrackStack.cc:61
G4int GetNTrack() const
Definition: G4TrackStack.hh:74
#define G4endl
Definition: G4ios.hh:61
G4GLOB_DLL std::ostream G4cerr
Here is the call graph for this function:

◆ getEnergyOfStack()

G4double G4SmartTrackStack::getEnergyOfStack ( G4TrackStack aTrackStack)

◆ GetMaxNTrack()

G4int G4SmartTrackStack::GetMaxNTrack ( ) const
inline

Definition at line 82 of file G4SmartTrackStack.hh.

82 { return maxNTracks; }

◆ GetNTrack()

G4int G4SmartTrackStack::GetNTrack ( ) const
inline

Definition at line 81 of file G4SmartTrackStack.hh.

81 { return nTracks; }
Here is the caller graph for this function:

◆ n_stackedTrack()

G4int G4SmartTrackStack::n_stackedTrack ( ) const
inlineprivate

Definition at line 85 of file G4SmartTrackStack.hh.

86  {
87  return stacks[0]->GetNTrack() +
88  stacks[1]->GetNTrack() +
89  stacks[2]->GetNTrack() +
90  stacks[3]->GetNTrack() +
91  stacks[4]->GetNTrack();
92  }
G4TrackStack * stacks[5]
G4int GetNTrack() const
Definition: G4TrackStack.hh:74
Here is the call graph for this function:

◆ operator!=()

int G4SmartTrackStack::operator!= ( const G4SmartTrackStack right) const
private

Definition at line 71 of file G4SmartTrackStack.cc.

71  {
72  return (this!=&right);
73 }

◆ operator=()

const G4SmartTrackStack & G4SmartTrackStack::operator= ( const G4SmartTrackStack right)
private

Definition at line 63 of file G4SmartTrackStack.cc.

63  {
64  return *this;
65 }

◆ operator==()

int G4SmartTrackStack::operator== ( const G4SmartTrackStack right) const
private

Definition at line 67 of file G4SmartTrackStack.cc.

67  {
68  return (this==&right);
69 }

◆ PopFromStack()

G4StackedTrack G4SmartTrackStack::PopFromStack ( )

Definition at line 83 of file G4SmartTrackStack.cc.

84 {
85  G4StackedTrack aStackedTrack;
86 
87  if (nTracks) {
88  while (true) {
89  if (stacks[fTurn]->GetNTrack()) {
90  aStackedTrack = stacks[fTurn]->PopFromStack();
91  energies[fTurn] -= aStackedTrack.GetTrack()->GetDynamicParticle()->GetTotalEnergy();
92  nTracks--;
93  break;
94  } else {
95  fTurn = (fTurn+1) % nTurn;
96  }
97  }
98  }
99 
100  // dumpStatistics();
101  return aStackedTrack;
102 }
G4TrackStack * stacks[5]
G4Track * GetTrack() const
G4StackedTrack PopFromStack()
Definition: G4TrackStack.hh:63
G4int GetNTrack() const
Here is the call graph for this function:

◆ PushToStack()

void G4SmartTrackStack::PushToStack ( const G4StackedTrack aStackedTrack)

Definition at line 108 of file G4SmartTrackStack.cc.

109 {
110 
111  G4int iDest = 0;
112  if (aStackedTrack.GetTrack()->GetParentID()) {
113  G4int code = aStackedTrack.GetTrack()->GetDynamicParticle()->GetPDGcode();
114  if (code == electronCode)
115  iDest = 2;
116  else if (code == gammaCode)
117  iDest = 3;
118  else if (code == positronCode)
119  iDest = 4;
120  else if (code == neutronCode)
121  iDest = 1;
122  } else {
123  // We have a primary track, which should go first.
124  fTurn = 0; // reseting the turn
125  }
126  stacks[iDest]->PushToStack(aStackedTrack);
127  energies[iDest] += aStackedTrack.GetTrack()->GetDynamicParticle()->GetTotalEnergy();
128  nTracks++;
129 
130  G4int dy1 = stacks[iDest]->GetNTrack() - stacks[iDest]->GetSafetyValve1();
132 
133  if (dy1 > 0 || dy1 > dy2 ||
134  (iDest == 2 &&
135  stacks[iDest]->GetNTrack() < 50 && energies[iDest] < energies[fTurn])) {
136  fTurn = iDest;
137  }
138 
140 }
G4TrackStack * stacks[5]
void PushToStack(const G4StackedTrack &aStackedTrack)
Definition: G4TrackStack.hh:62
G4int GetSafetyValve1() const
Definition: G4TrackStack.hh:76
G4Track * GetTrack() const
int G4int
Definition: G4Types.hh:78
G4int GetNTrack() const
Definition: G4TrackStack.hh:74
Definition: inftrees.h:24
G4int GetSafetyValve2() const
Definition: G4TrackStack.hh:77
G4int GetNTrack() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TransferTo()

void G4SmartTrackStack::TransferTo ( G4TrackStack aStack)

Definition at line 75 of file G4SmartTrackStack.cc.

76 {
77  for (int i = 0; i < nTurn; i++) {
78  stacks[i]->TransferTo(aStack);
79  }
80  nTracks = 0;
81 }
G4TrackStack * stacks[5]
void TransferTo(G4TrackStack *aStack)
Definition: G4TrackStack.cc:49
Here is the call graph for this function:

Member Data Documentation

◆ energies

G4double G4SmartTrackStack::energies[5]
private

Definition at line 70 of file G4SmartTrackStack.hh.

◆ fTurn

G4int G4SmartTrackStack::fTurn
private

Definition at line 68 of file G4SmartTrackStack.hh.

◆ maxNTracks

G4int G4SmartTrackStack::maxNTracks
private

Definition at line 77 of file G4SmartTrackStack.hh.

◆ nTracks

G4int G4SmartTrackStack::nTracks
private

Definition at line 78 of file G4SmartTrackStack.hh.

◆ nTurn

G4int G4SmartTrackStack::nTurn
private

Definition at line 69 of file G4SmartTrackStack.hh.

◆ stacks

G4TrackStack* G4SmartTrackStack::stacks[5]
private

Definition at line 71 of file G4SmartTrackStack.hh.


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