|
DD4hep
1.30.0
Detector Description Toolkit for High Energy Physics
|
Go to the documentation of this file.
53 print(
"+++ ******** MC Particle Printout for event ID:%d ********",e->GetEventID());
58 except(
"No Geant4MonteCarloTruth object present in the event structure to access the particle information!",
c_name());
86 ::snprintf(equiv,
sizeof(equiv),
"/%d",p->
g4Parent);
88 print(
"+++ %s ID:%7d %12s %6d%-7s %7s %3s %5d %3s %+.3e %-4s %-7s %-3s %-3s %2d [%s%s%s] %c%c%c%c",
113 G4HCofThisEvent* hc = e->GetHCofThisEvent();
114 for (
int ihc=0, nhc=hc->GetNumberOfCollections(); ihc<nhc; ++ihc) {
118 size_t nhits = coll->
GetSize();
119 for(
size_t i=0; i<nhits; ++i) {
122 if ( 0 != trk_hit ) {
126 if ( trueID == p->
id ) {
127 print(
"+++ %20s %s[%d] (%+.2e,%+.2e,%+.2e)[mm]",
"",c->GetName().c_str(),i,
132 if ( 0 != cal_hit ) {
134 for(Geant4HitData::Contributions::iterator j=contrib.begin(); j!=contrib.end(); ++j) {
138 if ( trueID == p->
id ) {
139 print(
"+++ %20s %s[%d] (%+.2e,%+.2e,%+.2e)[mm]",
"",c->GetName().c_str(),i,
147 print(
"+++ Hit unknown hit collection type: %s --> %s",
148 c->GetName().c_str(),typeName(
typeid(*c)).c_str());
160 int num_secondaries = 0;
161 int num_calo_hits = 0;
162 int num_tracker_hits = 0;
164 print(
"+++ MC Particles #Tracks:%7d ParticleType Parent/Geant4 "
165 "Primary Secondary Energy in [MeV] Calo Tracker Process/Par Details",
166 int(particles.size()));
167 for(ParticleMap::const_iterator i=particles.begin(); i!=particles.end(); ++i) {
171 num_secondaries += int(p->
daughters.size());
179 print(
"+++ MC Particles #Tracks:%7d ParticleType Parent/Geant4 "
180 "Primary Secondary Energy Calo Tracker Process/Par",
181 int(particles.size()));
182 print(
"+++ MC Particle Summary: %7d %10d %7d %7d %9d %5d %6d",
183 num_primary, num_secondaries, num_energy,
184 num_calo_hits,num_tracker_hits,num_process,num_parent);
193 size_t len =
sizeof(txt)-33;
195 if ( level>
int(len)-3 ) level = len-3;
197 ::snprintf(txt,
sizeof(txt),
"%5d ",level);
198 ::memset(txt+6,
' ',len-6);
201 if (
size_t(level + 6) < len ) {
203 ::memset(txt+level+6+1,
'-',len-level-3-6);
209 auto iter = particles.find(id_dau);
210 if ( iter != particles.end() ) {
219 print(
"+++ MC Particle Parent daughter relationships. [%d particles]",
int(particles.size()));
220 print(
"+++ MC Particles %12s #Tracks:%7d %-12s Parent%-7s "
221 "Primary Secondary Energy %-8s Calo Tracker Process/Par Details",
222 "",
int(particles.size()),
"ParticleType",
"",
"in [MeV]");
223 for(ParticleMap::const_iterator i=particles.begin(); i!=particles.end(); ++i) {
@ G4PARTICLE_CREATED_TRACKER_HIT
Position position
Hit position.
virtual ~Geant4ParticlePrint()
Default destructor.
virtual size_t GetSize() const override
Access the collection size.
@ G4PARTICLE_ABOVE_ENERGY_THRESHOLD
bool m_printBegin
Property: Flag to indicate output type at begin of event.
std::string particleName() const
Access to the Geant4 particle name.
bool m_printHits
Property: Flag to indicate output of hit data in tree.
const ParticleMap & particles() const
Access the particle map.
@ G4PARTICLE_KEEP_PROCESS
Generic hit container class using Geant4HitWrapper objects.
Geant4ParticleMap::ParticleMap ParticleMap
std::string processName() const
Access to the creator process name.
Geant4HitWrapper & hit(size_t which)
Access the hit wrapper.
Class of the Geant4 toolkit. See http://www-geant4.kek.jp/Reference.
int g4Parent
not persistent
static void increment(T *)
Increment count according to type information.
Geant4Event & event() const
Access the geant4 event – valid only between BeginEvent() and EndEvent()!
Geant4ParticlePrint(Geant4Context *context, const std::string &nam)
Standard constructor.
virtual void end(const G4Event *event) override
Post-event action callback.
virtual void begin(const G4Event *event) override
Pre-event action callback.
@ G4PARTICLE_GEN_DOCUMENTATION
Data structure to map particles produced during the generation and the simulation.
T * extension(bool alert=true)
Access to type safe extension object. Exception is thrown if the object is invalid.
void except(const char *fmt,...) const
Support of exceptions: Print fatal message and throw runtime_error.
std::vector< MonteCarloContrib > Contributions
Concrete basic implementation of the Geant4 event action.
void printParticles(const G4Event *e, const ParticleMap &particles) const
Print record of kept particles.
DDG4 tracker hit class used by the generic DDG4 tracker sensitive detector.
size_t numParent() const
Accessor to the number of particle parents.
Base class for geant4 hit structures used by the default DDG4 sensitive detector implementations.
Geant4Action & declareProperty(const std::string &nam, T &val)
Declare property.
static void decrement(T *)
Decrement count according to type information.
bool m_printGeneration
Property: Flag to indicate output type as part of the generator action.
Particles daughters
The list of daughters of this MC particle.
Contribution truth
Monte Carlo / Geant4 information.
std::string processTypeName() const
Access to the creator process type name.
void printParticleTree(const G4Event *e, const ParticleMap &particles, int level, Geant4ParticleHandle p) const
Print tree of kept particles.
const G4VProcess * process
Reference to the G4VProcess, which created this track.
DDG4 calorimeter hit class used by the generic DDG4 calorimeter sensitive detector.
void makePrintout(const G4Event *e) const
Print particle table.
int particleID(int track, bool throw_if_not_found=true) const
Access the equivalent track id (shortcut to the usage of TrackEquivalents)
void print(const char *fmt,...) const
Support for messages with variable output level using output level.
int trackID
Geant 4 Track identifier.
@ G4PARTICLE_HAS_SECONDARIES
Namespace for the Geant4 based simulation part of the AIDA detector description toolkit.
dd4hep::detail::ReferenceBitMask< int > PropertyMask
const char * c_name() const
Access name of the action.
Contributions truth
Hit contributions by individual particles.
void printParticle(const std::string &prefix, const G4Event *e, Geant4ParticleHandle p) const
Data structure to access derived MC particle information.
Particles parents
The list of parents of this MC particle.
virtual void operator()(G4Event *event)
Generation action callback.
@ G4PARTICLE_CREATED_CALORIMETER_HIT
Utility class describing the monte carlo contribution of a given particle to a hit.
double energy() const
Scalar particle energy.
int m_outputType
Property: Flag to indicate output type: 1: TABLE, 2:TREE, 3:BOTH (default)
Position position
Hit position.
bool m_printEnd
Property: Flag to indicate output type at end of event.
Generic context to extend user, run and event information.
Geant4Context * context() const
Access the context.