Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4SmartTrackStack.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// $Id$
28//
29
30#include "G4SmartTrackStack.hh"
31#include "G4VTrajectory.hh"
32#include "G4Track.hh"
33
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}
44
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}
54
56{
57 for (int i = 0; i < nTurn; i++) {
58 delete stacks[i];
59 }
60}
61
63G4SmartTrackStack::operator=(const G4SmartTrackStack &) {
64 return *this;
65}
66
67int G4SmartTrackStack::operator==(const G4SmartTrackStack &right) const {
68 return (this==&right);
69}
70
71int G4SmartTrackStack::operator!=(const G4SmartTrackStack &right) const {
72 return (this!=&right);
73}
74
76{
77 for (int i = 0; i < nTurn; i++) {
78 stacks[i]->TransferTo(aStack);
79 }
80 nTracks = 0;
81}
82
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}
103
104enum {
107
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}
141
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}
151
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}
@ electronCode
@ neutronCode
@ gammaCode
@ positronCode
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
G4double GetTotalEnergy() const
G4int GetPDGcode() const
void TransferTo(G4TrackStack *aStack)
void PushToStack(const G4StackedTrack &aStackedTrack)
G4int GetNTrack() const
G4StackedTrack PopFromStack()
G4Track * GetTrack() const
void PushToStack(const G4StackedTrack &aStackedTrack)
Definition: G4TrackStack.hh:62
void clearAndDestroy()
Definition: G4TrackStack.cc:40
G4int GetNTrack() const
Definition: G4TrackStack.hh:74
void TransferTo(G4TrackStack *aStack)
Definition: G4TrackStack.cc:49
G4double getTotalEnergy(void) const
Definition: G4TrackStack.cc:61
G4int GetSafetyValve2() const
Definition: G4TrackStack.hh:77
G4int GetSafetyValve1() const
Definition: G4TrackStack.hh:76
G4StackedTrack PopFromStack()
Definition: G4TrackStack.hh:63
const G4DynamicParticle * GetDynamicParticle() const
G4int GetParentID() const