Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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::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}

◆ ~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}

Member Function Documentation

◆ clear()

void G4SmartTrackStack::clear ( )

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}

◆ 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}
void clearAndDestroy()
Definition: G4TrackStack.cc:40

◆ 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}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
G4int GetNTrack() const
Definition: G4TrackStack.hh:74
G4double getTotalEnergy(void) const
Definition: G4TrackStack.cc:61

◆ 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; }

Referenced by PopFromStack(), and PushToStack().

◆ 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}
G4double GetTotalEnergy() const
G4int GetNTrack() const
G4Track * GetTrack() const
G4StackedTrack PopFromStack()
Definition: G4TrackStack.hh:63
const G4DynamicParticle * GetDynamicParticle() const

◆ 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();
131 G4int dy2 = stacks[fTurn]->GetNTrack() - stacks[fTurn]->GetSafetyValve2();
132
133 if (dy1 > 0 || dy1 > dy2 ||
134 (iDest == 2 &&
135 stacks[iDest]->GetNTrack() < 50 && energies[iDest] < energies[fTurn])) {
136 fTurn = iDest;
137 }
138
139 if (nTracks > maxNTracks) maxNTracks = nTracks;
140}
@ electronCode
@ neutronCode
@ gammaCode
@ positronCode
int G4int
Definition: G4Types.hh:66
G4int GetPDGcode() const
void PushToStack(const G4StackedTrack &aStackedTrack)
Definition: G4TrackStack.hh:62
G4int GetSafetyValve2() const
Definition: G4TrackStack.hh:77
G4int GetSafetyValve1() const
Definition: G4TrackStack.hh:76
G4int GetParentID() const

Referenced by G4TrackStack::TransferTo().

◆ 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}
void TransferTo(G4TrackStack *aStack)
Definition: G4TrackStack.cc:49

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