Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CascadeHistory Class Reference

#include <G4CascadeHistory.hh>

Classes

struct  HistoryEntry
 

Public Member Functions

 G4CascadeHistory (G4int verbose=0)
 
 ~G4CascadeHistory ()
 
void setVerboseLevel (G4int verbose=0)
 
void Clear ()
 
G4int AddEntry (G4CascadParticle &cpart)
 
G4int AddVertex (G4CascadParticle &cpart, std::vector< G4CascadParticle > &daug)
 
void DropEntry (const G4CascadParticle &cpart)
 
void Print (std::ostream &os) const
 

Protected Member Functions

void AssignHistoryID (G4CascadParticle &cpart)
 
void FillDaughters (G4int iEntry, std::vector< G4CascadParticle > &daug)
 
void PrintEntry (std::ostream &os, G4int iEntry) const
 
void PrintParticle (std::ostream &os, const G4CascadParticle &cpart) const
 
G4bool PrintingDone (G4int iEntry) const
 
const char * GuessTarget (const HistoryEntry &entry) const
 
G4int size () const
 

Detailed Description

Definition at line 40 of file G4CascadeHistory.hh.

Constructor & Destructor Documentation

◆ G4CascadeHistory()

G4CascadeHistory::G4CascadeHistory ( G4int verbose = 0)
inline

Definition at line 42 of file G4CascadeHistory.hh.

42: verboseLevel(verbose) {;}

◆ ~G4CascadeHistory()

G4CascadeHistory::~G4CascadeHistory ( )
inline

Definition at line 43 of file G4CascadeHistory.hh.

43{;} // *** Do not want subclasses ***

Member Function Documentation

◆ AddEntry()

G4int G4CascadeHistory::AddEntry ( G4CascadParticle & cpart)

Definition at line 105 of file G4CascadeHistory.cc.

105 {
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}
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4int getHistoryId() const
void AssignHistoryID(G4CascadParticle &cpart)

Referenced by AddVertex(), FillDaughters(), and G4IntraNucleiCascader::generateCascade().

◆ AddVertex()

G4int G4CascadeHistory::AddVertex ( G4CascadParticle & cpart,
std::vector< G4CascadParticle > & daug )

Definition at line 58 of file G4CascadeHistory.cc.

59 {
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}
void FillDaughters(G4int iEntry, std::vector< G4CascadParticle > &daug)
G4int AddEntry(G4CascadParticle &cpart)

Referenced by G4IntraNucleiCascader::generateCascade().

◆ AssignHistoryID()

void G4CascadeHistory::AssignHistoryID ( G4CascadParticle & cpart)
protected

Definition at line 136 of file G4CascadeHistory.cc.

136 {
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}
void setHistoryId(G4int id)

Referenced by AddEntry().

◆ Clear()

void G4CascadeHistory::Clear ( )

Definition at line 49 of file G4CascadeHistory.cc.

49 {
50 if (verboseLevel>1) G4cout << " >>> G4CascadeHistory::Clear" << G4endl;
51 theHistory.clear();
52 entryPrinted.clear();
53}

Referenced by G4IntraNucleiCascader::newCascade().

◆ DropEntry()

void G4CascadeHistory::DropEntry ( const G4CascadParticle & cpart)

Definition at line 126 of file G4CascadeHistory.cc.

126 {
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}

Referenced by G4IntraNucleiCascader::processTrappedParticle().

◆ FillDaughters()

void G4CascadeHistory::FillDaughters ( G4int iEntry,
std::vector< G4CascadParticle > & daug )
protected

Definition at line 79 of file G4CascadeHistory.cc.

79 {
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}

Referenced by AddVertex().

◆ GuessTarget()

const char * G4CascadeHistory::GuessTarget ( const HistoryEntry & entry) const
protected

Definition at line 203 of file G4CascadeHistory.cc.

204 {
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}
const G4InuclElementaryParticle & getParticle() const
G4double getCharge() const

Referenced by PrintEntry().

◆ Print()

void G4CascadeHistory::Print ( std::ostream & os) const

Definition at line 155 of file G4CascadeHistory.cc.

155 {
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}
void PrintEntry(std::ostream &os, G4int iEntry) const
G4bool PrintingDone(G4int iEntry) const

Referenced by G4IntraNucleiCascader::collide(), operator<<(), and G4IntraNucleiCascader::rescatter().

◆ PrintEntry()

void G4CascadeHistory::PrintEntry ( std::ostream & os,
G4int iEntry ) const
protected

Definition at line 168 of file G4CascadeHistory.cc.

168 {
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}
double z() const
Hep3Vector unit() const
Hep3Vector vect() const
G4int getGeneration() const
G4LorentzVector getMomentum() const
G4int getCurrentZone() const
const G4ThreeVector & getPosition() const
const char * GuessTarget(const HistoryEntry &entry) const
const G4ParticleDefinition * getDefinition() const
const G4String & GetParticleName() const

Referenced by Print(), and PrintEntry().

◆ PrintingDone()

G4bool G4CascadeHistory::PrintingDone ( G4int iEntry) const
inlineprotected

Definition at line 82 of file G4CascadeHistory.hh.

82 {
83 return (entryPrinted.find(iEntry) != entryPrinted.end());
84 }

Referenced by Print(), and PrintEntry().

◆ PrintParticle()

void G4CascadeHistory::PrintParticle ( std::ostream & os,
const G4CascadParticle & cpart ) const
protected

◆ setVerboseLevel()

void G4CascadeHistory::setVerboseLevel ( G4int verbose = 0)
inline

Definition at line 45 of file G4CascadeHistory.hh.

45{ verboseLevel = verbose; }

Referenced by G4IntraNucleiCascader::setVerboseLevel().

◆ size()

G4int G4CascadeHistory::size ( ) const
inlineprotected

Definition at line 89 of file G4CascadeHistory.hh.

89{ return (G4int)theHistory.size(); }

Referenced by AddEntry(), AssignHistoryID(), Print(), and PrintEntry().


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