Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4CascadeHistory.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// G4CascadeHistory: Container to record all particles produced during
28// cascade, with daughters; printout is formatted hierarchically.
29
30#include "G4CascadeHistory.hh"
31#include "G4CascadParticle.hh"
33#include "G4InuclParticle.hh"
34#include <iostream>
35#include <iomanip>
36#include <ios>
37
38
39// Constructors for individual entries (vertices)
40
42 for (G4int i=0; i<10; i++) dId[i] = -1;
43 n = 0;
44}
45
46
47// Initialize buffers for new event
48
50 if (verboseLevel>1) G4cout << " >>> G4CascadeHistory::Clear" << G4endl;
51 theHistory.clear();
52 entryPrinted.clear();
53}
54
55
56// Record a new cascade vertex (particle and daughters)
57
59 std::vector<G4CascadParticle>& daug) {
60 if (verboseLevel>1) G4cout << " >>> G4CascadeHistory::AddVertex" << G4endl;
61
62 // Create new entry for vertex or update particle kinematics
63 G4int id = AddEntry(cpart);
64 FillDaughters(id, daug);
65
66 if (verboseLevel>3) {
67 G4cout << " entry " << id << " " << &theHistory[id] << " got "
68 << theHistory[id].n << " daughters:";
69 for (G4int i=0; i<theHistory[id].n; i++) {
70 G4cout << " " << theHistory[id].dId[i];
71 }
72 G4cout << G4endl;
73 }
74
75 return id;
76}
77
78void
79G4CascadeHistory::FillDaughters(G4int iEntry, std::vector<G4CascadParticle>& daug) {
80 G4int nDaug = (G4int)daug.size();
81
82 if (verboseLevel>1)
83 G4cout << " >>> G4CascadeHistory::FillDaughters " << iEntry << G4endl;
84
85 // NOTE: Cannot use reference to element, as push_back can invalidate refs!
86 theHistory[iEntry].clear();
87
88 theHistory[iEntry].n = nDaug;
89 for (G4int i=0; i<nDaug; i++) {
90 G4int id = AddEntry(daug[i]);
91 theHistory[iEntry].dId[i] = id;
92 }
93
94 if (verboseLevel>3) {
95 G4cout << " got " << theHistory[iEntry].n << " daughters:";
96 for (G4int i=0; i<theHistory[iEntry].n; i++) {
97 G4cout << " " << theHistory[iEntry].dId[i];
98 }
99 G4cout << G4endl;
100 }
101}
102
103// Add particle to history list, or update if already recorded
104
106 AssignHistoryID(cpart); // Make sure particle has index
107
108 G4int id = cpart.getHistoryId();
109 if (id < size()) {
110 if (verboseLevel>2)
111 G4cout << " AddEntry updating " << id << " " << &theHistory[id] << G4endl;
112 theHistory[id].cpart = cpart; // Copies kinematics
113 } else {
114 theHistory.push_back(HistoryEntry(cpart));
115 if (verboseLevel>2)
116 G4cout << " AddEntry creating " << id << " " << &theHistory.back() << G4endl;
117 }
118
119 if (verboseLevel>3) G4cout << theHistory[id].cpart << G4endl; // Sanity check
120
121 return id;
122}
123
124// Discard particle reabsorbed during cascade
125
127 if (verboseLevel>1) G4cout << " >>> G4CascadeHistory::DropEntry" << G4endl;
128
129
130 G4int id = cpart.getHistoryId(); // Particle must appear in history
131 if (id>=0) theHistory[id].n = -1; // Special flag for absorbed particle
132}
133
134// Check if particle already in history, assign ID if not there
135
137 if (cpart.getHistoryId() >= 0) return; // ID already assigned
138
139 if (verboseLevel>2) {
140 G4cout << " >>> G4CascadeHistory::NewHistoryID assigning ID "
141 << size() << G4endl;
142 }
143
144 cpart.setHistoryId(size());
145}
146
147
148// Generate hierarchical (indented) report of history
149
150std::ostream& operator<<(std::ostream& os, const G4CascadeHistory& history) {
151 history.Print(os);
152 return os;
153}
154
155void G4CascadeHistory::Print(std::ostream& os) const {
156 if (verboseLevel) os << " >>> G4CascadeHistory::Print" << std::endl;
157
158 os << " Cascade structure: vertices, (-O-) exciton, (***) outgoing"
159 << std::endl;
160
161 for (G4int i=0; i<size(); i++) {
162 if (!PrintingDone(i)) PrintEntry(os, i);
163 }
164}
165
166// Add single-line report for particle, along with daughters
167
168void G4CascadeHistory::PrintEntry(std::ostream& os, G4int iEntry) const {
169 if (iEntry >= size()) return; // Skip nonexistent entry
170 if (PrintingDone(iEntry)) return; // Skip entry already reported
171
172 entryPrinted.insert(iEntry);
173
174 const HistoryEntry& entry = theHistory[iEntry]; // For convenience
175 const G4CascadParticle& cpart = entry.cpart;
176
177 G4int indent = cpart.getGeneration()*2;
178
179 // Index and indentation of cascade vertex
180 std::ios::fmtflags osFlags = os.flags();
181 os.setf(std::ios::left); // Pushes all blanks to right end of output
182 os << "#" << std::setw(3+indent) << iEntry;
183 os.flags(osFlags);
184
185 os << cpart.getParticle().getDefinition()->GetParticleName()
186 << " p " << cpart.getMomentum() << " (cosTh "
187 << cpart.getMomentum().vect().unit().z() << ")"
188 << " @ " << cpart.getPosition()
189 << " zone " << cpart.getCurrentZone();
190
191 // Flag as final-state particle or report daughters iteratively
192 os << " (" << GuessTarget(entry) << ")";
193 if (entry.n > 0) {
194 os << " -> N=" << entry.n << std::endl;
195 for (G4int i=0; i<entry.n; i++) {
196 PrintEntry(os, entry.dId[i]);
197 }
198 } else os << std::endl;
199}
200
201// Derive target of cascade step from particle and daughters
202
205 if (verboseLevel>2) G4cout << " >>> G4CascadeHistory::GuessTarget" << G4endl;
206
207 if (entry.n < 0) return "-O-"; // Exciton or trapped-decay
208 if (entry.n == 0) return "***"; // Outgoing (final state) particle
209
210 const G4CascadParticle& cpart = entry.cpart; // For convenience
211 if (verboseLevel>3) G4cout << "cpart: " << cpart;
212
213 // Compute baryon number and charge from daughters minus projectile
214 G4int targetB = -cpart.getParticle().baryon();
215 G4int targetQ = (G4int)-cpart.getParticle().getCharge();
216
217 for (G4int i=0; i<entry.n; i++) {
218 const G4CascadParticle& cdaug = theHistory[entry.dId[i]].cpart;
219 if (verboseLevel>3)
220 G4cout << "cdaug " << i << " ID " << entry.dId[i] << ": " << cdaug;
221
222 targetB += cdaug.getParticle().baryon();
223 targetQ += (G4int)cdaug.getParticle().getCharge();
224 }
225
226 // Target possibilities are proton, neutron or dibaryon (pp, nn, pn)
227 if (targetB==1 && targetQ==0) return "n";
228 if (targetB==1 && targetQ==1) return "p";
229 if (targetB==2 && targetQ==0) return "nn";
230 if (targetB==2 && targetQ==1) return "pn";
231 if (targetB==2 && targetQ==2) return "pp";
232
233 if (verboseLevel>2) {
234 G4cout << " ERROR identifying target: deltaB " << targetB
235 << " deltaQ " << targetQ << " from\n" << cpart << " to" << G4endl;
236 for (G4int j=0; j<entry.n; j++) {
237 G4cout << theHistory[entry.dId[j]].cpart;
238 }
239 }
240
241 return "BAD TARGET"; // Should not get here if EPCollider worked right
242}
std::ostream & operator<<(std::ostream &os, const G4CascadeHistory &history)
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
Hep3Vector unit() const
Hep3Vector vect() const
G4int getHistoryId() const
G4int getGeneration() const
const G4InuclElementaryParticle & getParticle() const
void setHistoryId(G4int id)
G4LorentzVector getMomentum() const
G4int getCurrentZone() const
const G4ThreeVector & getPosition() const
void PrintEntry(std::ostream &os, G4int iEntry) const
void FillDaughters(G4int iEntry, std::vector< G4CascadParticle > &daug)
void Print(std::ostream &os) const
G4int AddVertex(G4CascadParticle &cpart, std::vector< G4CascadParticle > &daug)
G4int size() const
void DropEntry(const G4CascadParticle &cpart)
G4int AddEntry(G4CascadParticle &cpart)
void AssignHistoryID(G4CascadParticle &cpart)
G4bool PrintingDone(G4int iEntry) const
const char * GuessTarget(const HistoryEntry &entry) const
const G4ParticleDefinition * getDefinition() const
G4double getCharge() const
const G4String & GetParticleName() const