|
||||
File indexing completed on 2025-01-18 09:58:00
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // 0027 // G4CascadeHistory: Container to record all particles produced during 0028 // cascade, with daughters; printout is formatted hierarchically. 0029 0030 #ifndef G4CASCADE_HISTORY_HH 0031 #define G4CASCADE_HISTORY_HH 0032 0033 #include "globals.hh" 0034 #include "G4CascadParticle.hh" 0035 #include <iosfwd> 0036 #include <set> 0037 #include <vector> 0038 0039 0040 class G4CascadeHistory { 0041 public: 0042 G4CascadeHistory(G4int verbose=0) : verboseLevel(verbose) {;} 0043 ~G4CascadeHistory() {;} // *** Do not want subclasses *** 0044 0045 void setVerboseLevel(G4int verbose=0) { verboseLevel = verbose; } 0046 0047 // Reset buffers for new event 0048 void Clear(); 0049 0050 // Add particle to history list, assigning ID number (non-const input) 0051 G4int AddEntry(G4CascadParticle& cpart); 0052 0053 // Record full interaction vertex (non-const input) 0054 G4int AddVertex(G4CascadParticle& cpart, std::vector<G4CascadParticle>& daug); 0055 0056 // Discard particle reabsorbed during cascade 0057 void DropEntry(const G4CascadParticle& cpart); 0058 0059 // Report cascade structure hierarchically 0060 void Print(std::ostream& os) const; 0061 0062 protected: 0063 struct HistoryEntry { 0064 G4CascadParticle cpart; 0065 G4int n; 0066 G4int dId[10]; // Must be fixed size for allocation in vector 0067 0068 HistoryEntry() { clear(); }; 0069 HistoryEntry(const G4CascadParticle& cp) : cpart(cp) { clear(); } 0070 void clear(); 0071 }; 0072 0073 // Assign ID number to particle (non-const input) 0074 void AssignHistoryID(G4CascadParticle& cpart); 0075 0076 // Populate list of daughters in interaction, adding them to history list 0077 void FillDaughters(G4int iEntry, std::vector<G4CascadParticle>& daug); 0078 0079 // Add single-line report for particle, along with daughters 0080 void PrintEntry(std::ostream& os, G4int iEntry) const; 0081 void PrintParticle(std::ostream& os, const G4CascadParticle& cpart) const; 0082 G4bool PrintingDone(G4int iEntry) const { 0083 return (entryPrinted.find(iEntry) != entryPrinted.end()); 0084 } 0085 0086 // Derive target of cascade step from particle and daughters 0087 const char* GuessTarget(const HistoryEntry& entry) const; 0088 0089 G4int size() const { return (G4int)theHistory.size(); } 0090 0091 private: 0092 G4int verboseLevel; 0093 0094 std::vector<HistoryEntry> theHistory; // List of particles and daughters 0095 mutable std::set<G4int> entryPrinted; // Particle indices already reported 0096 0097 // No copying allowed 0098 G4CascadeHistory(const G4CascadeHistory& rhs); 0099 G4CascadeHistory& operator=(const G4CascadeHistory& rhs); 0100 }; 0101 0102 // Inline reporting 0103 std::ostream& operator<<(std::ostream& os, const G4CascadeHistory& history); 0104 0105 #endif /* G4CASCADE_HISTORY_HH */
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |