Geant4-11
Public Member Functions | Private Member Functions | Private Attributes
G4SmartTrackStack Class Reference

#include <G4SmartTrackStack.hh>

Public Member Functions

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

Private Member Functions

G4int n_stackedTrack () const
 

Private Attributes

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

Detailed Description

Definition at line 42 of file G4SmartTrackStack.hh.

Constructor & Destructor Documentation

◆ G4SmartTrackStack()

G4SmartTrackStack::G4SmartTrackStack ( )

Definition at line 47 of file G4SmartTrackStack.cc.

48 : fTurn(0), nTurn(5), maxNTracks(0), nTracks(0)
49{
50 for(G4int i=0; i<nTurn; ++i)
51 {
52 stacks[i] = new G4TrackStack(5000);
53 energies[i] = 0.;
54 }
55}
int G4int
Definition: G4Types.hh:85
G4TrackStack * stacks[5]

References energies, nTurn, and stacks.

◆ ~G4SmartTrackStack()

G4SmartTrackStack::~G4SmartTrackStack ( )

Definition at line 57 of file G4SmartTrackStack.cc.

58{
59 for (G4int i=0; i<nTurn; ++i)
60 {
61 delete stacks[i];
62 }
63}

References nTurn, and stacks.

Member Function Documentation

◆ clear()

void G4SmartTrackStack::clear ( )

Definition at line 142 of file G4SmartTrackStack.cc.

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

References energies, fTurn, nTracks, nTurn, and stacks.

◆ clearAndDestroy()

void G4SmartTrackStack::clearAndDestroy ( )

Definition at line 153 of file G4SmartTrackStack.cc.

154{
155 for (G4int i = 0; i < nTurn; ++i)
156 {
158 energies[i] = 0.0;
159 fTurn = 0;
160 }
161 nTracks = 0;
162}
void clearAndDestroy()
Definition: G4TrackStack.cc:41

References G4TrackStack::clearAndDestroy(), energies, fTurn, nTracks, nTurn, and stacks.

◆ dumpStatistics()

void G4SmartTrackStack::dumpStatistics ( )

Definition at line 35 of file G4SmartTrackStack.cc.

36{
37 // Print to stderr so that we can split stats output from normal
38 // output of Geant4 which is typically being printed to stdout
39 for (G4int i=0; i<nTurn; ++i)
40 {
41 G4cerr << stacks[i]->GetNTrack() << " ";
42 G4cerr << stacks[i]->getTotalEnergy() << " ";
43 }
44 G4cerr << G4endl;
45}
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4double getTotalEnergy(void) const
Definition: G4TrackStack.cc:68
std::size_t GetNTrack() const
Definition: G4TrackStack.hh:68

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

◆ getEnergyOfStack()

G4double G4SmartTrackStack::getEnergyOfStack ( G4TrackStack aTrackStack)

◆ GetMaxNTrack()

G4int G4SmartTrackStack::GetMaxNTrack ( ) const
inline

Definition at line 62 of file G4SmartTrackStack.hh.

62{ return maxNTracks; }

References maxNTracks.

◆ GetNTrack()

G4int G4SmartTrackStack::GetNTrack ( ) const
inline

Definition at line 61 of file G4SmartTrackStack.hh.

61{ return nTracks; }

References nTracks.

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

◆ n_stackedTrack()

G4int G4SmartTrackStack::n_stackedTrack ( ) const
inlineprivate

Definition at line 66 of file G4SmartTrackStack.hh.

67 {
68 return G4int(stacks[0]->GetNTrack() +
69 stacks[1]->GetNTrack() +
70 stacks[2]->GetNTrack() +
71 stacks[3]->GetNTrack() +
72 stacks[4]->GetNTrack());
73 }
G4int GetNTrack() const

References GetNTrack(), and stacks.

◆ operator!=()

G4bool G4SmartTrackStack::operator!= ( const G4SmartTrackStack ) const
delete

◆ operator=()

G4SmartTrackStack & G4SmartTrackStack::operator= ( const G4SmartTrackStack )
delete

◆ operator==()

G4bool G4SmartTrackStack::operator== ( const G4SmartTrackStack ) const
delete

◆ PopFromStack()

G4StackedTrack G4SmartTrackStack::PopFromStack ( )

Definition at line 74 of file G4SmartTrackStack.cc.

75{
76 G4StackedTrack aStackedTrack;
77
78 if (nTracks)
79 {
80 while (true)
81 {
82 if (stacks[fTurn]->GetNTrack())
83 {
84 aStackedTrack = stacks[fTurn]->PopFromStack();
86 --nTracks;
87 break;
88 }
89 else
90 {
91 fTurn = (fTurn+1) % nTurn;
92 }
93 }
94 }
95
96 return aStackedTrack;
97}
G4double GetTotalEnergy() const
G4Track * GetTrack() const
G4StackedTrack PopFromStack()
Definition: G4TrackStack.hh:61
const G4DynamicParticle * GetDynamicParticle() const

References energies, fTurn, G4Track::GetDynamicParticle(), GetNTrack(), G4DynamicParticle::GetTotalEnergy(), G4StackedTrack::GetTrack(), nTracks, nTurn, G4TrackStack::PopFromStack(), and stacks.

◆ PushToStack()

void G4SmartTrackStack::PushToStack ( const G4StackedTrack aStackedTrack)

Definition at line 104 of file G4SmartTrackStack.cc.

105{
106
107 G4int iDest = 0;
108 if (aStackedTrack.GetTrack()->GetParentID())
109 {
110 G4int code = aStackedTrack.GetTrack()->GetDynamicParticle()->GetPDGcode();
111 if (code == electronCode)
112 iDest = 2;
113 else if (code == gammaCode)
114 iDest = 3;
115 else if (code == positronCode)
116 iDest = 4;
117 else if (code == neutronCode)
118 iDest = 1;
119 }
120 else
121 {
122 // We have a primary track, which should go first.
123 fTurn = 0; // reseting the turn
124 }
125 stacks[iDest]->PushToStack(aStackedTrack);
126 energies[iDest] += aStackedTrack.GetTrack()->GetDynamicParticle()->GetTotalEnergy();
127 ++nTracks;
128
129 G4int dy1 = stacks[iDest]->GetNTrack() - stacks[iDest]->GetSafetyValue1();
131
132 if (dy1 > 0 || dy1 > dy2 ||
133 (iDest == 2 &&
134 stacks[iDest]->GetNTrack() < 50 && energies[iDest] < energies[fTurn]))
135 {
136 fTurn = iDest;
137 }
138
140}
@ electronCode
@ neutronCode
@ gammaCode
@ positronCode
G4int GetPDGcode() const
void PushToStack(const G4StackedTrack &aStackedTrack)
Definition: G4TrackStack.hh:59
G4int GetSafetyValue2() const
Definition: G4TrackStack.hh:71
G4int GetSafetyValue1() const
Definition: G4TrackStack.hh:70
G4int GetParentID() const
Definition: inftrees.h:24

References electronCode, energies, fTurn, gammaCode, G4Track::GetDynamicParticle(), GetNTrack(), G4TrackStack::GetNTrack(), G4Track::GetParentID(), G4DynamicParticle::GetPDGcode(), G4TrackStack::GetSafetyValue1(), G4TrackStack::GetSafetyValue2(), G4DynamicParticle::GetTotalEnergy(), G4StackedTrack::GetTrack(), maxNTracks, neutronCode, nTracks, positronCode, G4TrackStack::PushToStack(), and stacks.

Referenced by G4TrackStack::TransferTo().

◆ TransferTo()

void G4SmartTrackStack::TransferTo ( G4TrackStack aStack)

Definition at line 65 of file G4SmartTrackStack.cc.

66{
67 for (G4int i=0; i<nTurn; ++i)
68 {
69 stacks[i]->TransferTo(aStack);
70 }
71 nTracks = 0;
72}
void TransferTo(G4TrackStack *aStack)
Definition: G4TrackStack.cc:51

References nTracks, nTurn, stacks, and G4TrackStack::TransferTo().

Field Documentation

◆ energies

G4double G4SmartTrackStack::energies[5]
private

◆ fTurn

G4int G4SmartTrackStack::fTurn
private

Definition at line 77 of file G4SmartTrackStack.hh.

Referenced by clear(), clearAndDestroy(), PopFromStack(), and PushToStack().

◆ maxNTracks

G4int G4SmartTrackStack::maxNTracks
private

Definition at line 86 of file G4SmartTrackStack.hh.

Referenced by GetMaxNTrack(), and PushToStack().

◆ nTracks

G4int G4SmartTrackStack::nTracks
private

◆ nTurn

G4int G4SmartTrackStack::nTurn
private

◆ stacks

G4TrackStack* G4SmartTrackStack::stacks[5]
private

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