Geant4 11.2.2
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 ()
 
G4SmartTrackStackoperator= (const G4SmartTrackStack &)=delete
 
G4bool operator== (const G4SmartTrackStack &) const =delete
 
G4bool operator!= (const G4SmartTrackStack &) const =delete
 
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 44 of file G4SmartTrackStack.hh.

Constructor & Destructor Documentation

◆ G4SmartTrackStack()

G4SmartTrackStack::G4SmartTrackStack ( )

Definition at line 47 of file G4SmartTrackStack.cc.

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

◆ ~G4SmartTrackStack()

G4SmartTrackStack::~G4SmartTrackStack ( )

Definition at line 56 of file G4SmartTrackStack.cc.

57{
58 for (auto* sp : stacks)
59 {
60 delete sp;
61 }
62}

Member Function Documentation

◆ clear()

void G4SmartTrackStack::clear ( )

Definition at line 139 of file G4SmartTrackStack.cc.

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

◆ clearAndDestroy()

void G4SmartTrackStack::clearAndDestroy ( )

Definition at line 150 of file G4SmartTrackStack.cc.

151{
152 for (G4int i = 0; i < nTurn; ++i)
153 {
154 stacks[i]->clearAndDestroy();
155 energies[i] = 0.0;
156 fTurn = 0;
157 }
158 nTracks = 0;
159}

◆ 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:67

◆ getEnergyOfStack()

G4double G4SmartTrackStack::getEnergyOfStack ( G4TrackStack * aTrackStack)

◆ GetMaxNTrack()

G4int G4SmartTrackStack::GetMaxNTrack ( ) const
inline

Definition at line 64 of file G4SmartTrackStack.hh.

64{ return maxNTracks; }

◆ GetNTrack()

G4int G4SmartTrackStack::GetNTrack ( ) const
inline

Definition at line 63 of file G4SmartTrackStack.hh.

63{ return nTracks; }

Referenced by PopFromStack(), and PushToStack().

◆ 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 73 of file G4SmartTrackStack.cc.

74{
75 G4StackedTrack aStackedTrack;
76
77 if (nTracks != 0)
78 {
79 while (true)
80 {
81 if (stacks[fTurn]->GetNTrack() != 0u)
82 {
83 aStackedTrack = stacks[fTurn]->PopFromStack();
84 energies[fTurn] -= aStackedTrack.GetTrack()->GetDynamicParticle()->GetTotalEnergy();
85 --nTracks;
86 break;
87 }
88
89 fTurn = (fTurn+1) % nTurn;
90 }
91 }
92
93 return aStackedTrack;
94}
G4double GetTotalEnergy() const
G4int GetNTrack() const
G4Track * GetTrack() const
const G4DynamicParticle * GetDynamicParticle() const

◆ PushToStack()

void G4SmartTrackStack::PushToStack ( const G4StackedTrack & aStackedTrack)

Definition at line 101 of file G4SmartTrackStack.cc.

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

Referenced by G4TrackStack::TransferTo().

◆ TransferTo()

void G4SmartTrackStack::TransferTo ( G4TrackStack * aStack)

Definition at line 64 of file G4SmartTrackStack.cc.

65{
66 for (auto* sp : stacks)
67 {
68 sp->TransferTo(aStack);
69 }
70 nTracks = 0;
71}

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