|
DD4hep
1.30.0
Detector Description Toolkit for High Energy Physics
|
Go to the documentation of this file.
13 #ifndef DD4HEP_DDG4_GEANT4HITTRUTHHANDLER_H
14 #define DD4HEP_DDG4_GEANT4HITTRUTHHANDLER_H
51 virtual void begin(
const G4Event* event)
override;
53 virtual void end(
const G4Event* event)
override;
80 #include <G4HCofThisEvent.hh>
106 size_t nhits = coll->
GetSize();
107 for(
size_t i=0; i<nhits; ++i) {
110 if ( 0 != trk_hit ) {
118 if ( 0 != cal_hit ) {
121 for(Geant4HitData::Contributions::iterator j=c.begin(); j!=c.end(); ++j) {
134 G4HCofThisEvent* hce =
event->GetHCofThisEvent();
136 int nCol = hce->GetNumberOfCollections();
138 if ( truth && !truth->
isValid() ) {
140 printout(WARNING,
name(),
"+++ [Event:%d] No valid MC truth info present. "
141 "Is a Particle handler installed ?",event->GetEventID());
143 for (
int i = 0; i < nCol; ++i) {
149 warning(
"+++ [Event:%d] The value of G4HCofThisEvent is NULL. No collections saved!",
150 event->GetEventID());
virtual ~Geant4HitTruthHandler()
Default destructor.
bool m_needsControl
Default property: Flag to create control instance.
virtual size_t GetSize() const override
Access the collection size.
Generic hit container class using Geant4HitWrapper objects.
Geant4HitWrapper & hit(size_t which)
Access the hit wrapper.
Class of the Geant4 toolkit. See http://www-geant4.kek.jp/Reference.
static void increment(T *)
Increment count according to type information.
Geant4Event & event() const
Access the geant4 event – valid only between BeginEvent() and EndEvent()!
#define DECLARE_GEANT4ACTION(name)
Plugin defintion to create Geant4Action objects.
void warning(const char *fmt,...) const
Support of warning messages.
virtual void begin(const G4Event *event) override
Geant4EventAction interface: Begin-of-event callback.
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.
std::vector< MonteCarloContrib > Contributions
Concrete basic implementation of the Geant4 event action.
std::vector< std::string > CollectionNames
virtual void end(const G4Event *event) override
Geant4EventAction interface: End-of-event callback.
DDG4 tracker hit class used by the generic DDG4 tracker sensitive detector.
Class to measure the energy of escaping tracks.
Base class for geant4 hit structures used by the default DDG4 sensitive detector implementations.
static void decrement(T *)
Decrement count according to type information.
Contribution truth
Monte Carlo / Geant4 information.
const std::string & name() const
Access name of the action.
DDG4 calorimeter hit class used by the generic DDG4 calorimeter sensitive detector.
int particleID(int track, bool throw_if_not_found=true) const
Access the equivalent track id (shortcut to the usage of TrackEquivalents)
int trackID
Geant 4 Track identifier.
Namespace for the Geant4 based simulation part of the AIDA detector description toolkit.
Contributions truth
Hit contributions by individual particles.
Namespace for the AIDA detector description toolkit.
Utility class describing the monte carlo contribution of a given particle to a hit.
void handleCollection(Geant4ParticleMap *truth, G4VHitsCollection *hc)
Dump single container of hits.
bool isValid() const
Check if the particle map was ever filled (ie. some particle handler was present)
Geant4HitTruthHandler(Geant4Context *context, const std::string &nam)
Standard constructor.
Generic context to extend user, run and event information.
Geant4Context * context() const
Access the context.