ATLAS Offline Software
AtlasG4EventUserInfo.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 #ifndef MCTRUTH_ATLASG4EVENTUSERINFO_H
6 #define MCTRUTH_ATLASG4EVENTUSERINFO_H
7 
8 
9 #include <GaudiKernel/EventContext.h>
11 #include "AtlasHepMC/GenParticle.h"
12 #include "G4VUserEventInformation.hh"
14 #include <memory>
24 class AtlasG4EventUserInfo: public G4VUserEventInformation {
25 public:
26  AtlasG4EventUserInfo(const EventContext& ctx)
27  : G4VUserEventInformation()
28  , m_eventContext(ctx)
29  {}
30 
35  HepMC::GenEvent* GetHepMCEvent() ;
40  void SetHepMCEvent(HepMC::GenEvent*);
41 
58 
72 
87  void SetLastProcessedTrackID(int trackID) { m_lastProcessedTrackID = trackID; }
88 
97  int GetLastProcessedStep() const { return m_lastProcessedStep; }
104  void SetLastProcessedStep(int stepNumber) { m_lastProcessedStep = stepNumber; }
105 
110  std::shared_ptr<HitCollectionMap> GetHitCollectionMap() const { return m_hitCollectionMap; }
111 
115  void SetHitCollectionMap(std::shared_ptr<HitCollectionMap> hitCollections) { m_hitCollectionMap = hitCollections; }
116 
117  const EventContext& GetEventContext() const { return m_eventContext; }
118 
119  void Print() const {}
120 
121 private:
122  const EventContext& m_eventContext;
123  HepMC::GenEvent *m_theEvent{};
126 
127  std::shared_ptr<HitCollectionMap> m_hitCollectionMap{std::make_shared<HitCollectionMap>()};
128 
129  // These next two variables are used by the CaloCalibrationHit
130  // recording code as event-level flags They correspond to the Track
131  // ID and step number of the last G4Step processed by a
132  // CaloCalibrationHit SD Both are needed, because a particle might
133  // have only one step
136 };
137 
138 #endif // MCTRUTH_ATLASG4EVENTUSERINFO_H
AtlasG4EventUserInfo::m_lastProcessedTrackID
int m_lastProcessedTrackID
Definition: AtlasG4EventUserInfo.h:134
AtlasG4EventUserInfo::SetCurrentPrimaryGenParticle
void SetCurrentPrimaryGenParticle(HepMC::ConstGenParticlePtr p)
set m_currentPrimaryGenParticle, the pointer to the HepMC::GenParticle used to create the current G4P...
Definition: AtlasG4EventUserInfo.h:57
AtlasG4EventUserInfo::GetCurrentGenParticle
HepMC::ConstGenParticlePtr GetCurrentGenParticle() const
Definition: AtlasG4EventUserInfo.h:64
AtlasG4EventUserInfo::SetLastProcessedTrackID
void SetLastProcessedTrackID(int trackID)
record the value of G4Track::GetTrackID() for the current G4Step.
Definition: AtlasG4EventUserInfo.h:87
AtlasG4EventUserInfo
This class is attached to G4Event objects as UserInformation. It holds a pointer to the HepMC::GenEve...
Definition: AtlasG4EventUserInfo.h:24
AtlasG4EventUserInfo::m_currentGenParticle
HepMC::GenParticlePtr m_currentGenParticle
Definition: AtlasG4EventUserInfo.h:125
AtlasG4EventUserInfo::GetHepMCEvent
HepMC::GenEvent * GetHepMCEvent()
return a pointer to the HepMC::GenEvent used to create the G4Event.
Definition: AtlasG4EventUserInfo.cxx:8
HitCollectionMap.h
AtlasG4EventUserInfo::SetLastProcessedStep
void SetLastProcessedStep(int stepNumber)
record value of the G4Track::GetCurrentStepNumber() for the current G4Step.
Definition: AtlasG4EventUserInfo.h:104
AtlasG4EventUserInfo::SetHepMCEvent
void SetHepMCEvent(HepMC::GenEvent *)
set m_theEvent, the pointer to the HepMC::GenEvent used to create the G4Event.
Definition: AtlasG4EventUserInfo.cxx:13
HepMC::GenParticlePtr
GenParticle * GenParticlePtr
Definition: GenParticle.h:37
AtlasG4EventUserInfo::SetHitCollectionMap
void SetHitCollectionMap(std::shared_ptr< HitCollectionMap > hitCollections)
Set the HitCollectionMap object.
Definition: AtlasG4EventUserInfo.h:115
GenParticle.h
AtlasG4EventUserInfo::GetLastProcessedStep
int GetLastProcessedStep() const
return the value of the G4Track::GetCurrentStepNumber() for the last G4Step processed by a CaloCalibr...
Definition: AtlasG4EventUserInfo.h:97
AtlasG4EventUserInfo::GetLastProcessedTrackID
int GetLastProcessedTrackID() const
return the value of G4Track::GetTrackID() for the last G4Step processed by a CaloCalibrationHit Sensi...
Definition: AtlasG4EventUserInfo.h:80
AtlasG4EventUserInfo::GetCurrentGenParticle
HepMC::GenParticlePtr GetCurrentGenParticle()
return a pointer to the GenParticle corresponding to the current G4Track (if there is one).
Definition: AtlasG4EventUserInfo.h:63
python.utils.AtlRunQueryDQUtils.p
p
Definition: AtlRunQueryDQUtils.py:209
AtlasG4EventUserInfo::Print
void Print() const
Definition: AtlasG4EventUserInfo.h:119
AtlasG4EventUserInfo::m_currentPrimaryGenParticle
HepMC::ConstGenParticlePtr m_currentPrimaryGenParticle
Definition: AtlasG4EventUserInfo.h:124
GenEvent_fwd.h
AtlasG4EventUserInfo::m_hitCollectionMap
std::shared_ptr< HitCollectionMap > m_hitCollectionMap
Definition: AtlasG4EventUserInfo.h:127
HepMC::ConstGenParticlePtr
const GenParticle * ConstGenParticlePtr
Definition: GenParticle.h:38
AtlasG4EventUserInfo::GetCurrentPrimaryGenParticle
HepMC::ConstGenParticlePtr GetCurrentPrimaryGenParticle() const
return a pointer to the HepMC::GenParticle used to create the current G4PrimaryParticle.
Definition: AtlasG4EventUserInfo.h:48
AtlasG4EventUserInfo::m_theEvent
HepMC::GenEvent * m_theEvent
Definition: AtlasG4EventUserInfo.h:123
AtlasG4EventUserInfo::GetEventContext
const EventContext & GetEventContext() const
Definition: AtlasG4EventUserInfo.h:117
AtlasG4EventUserInfo::m_eventContext
const EventContext & m_eventContext
Definition: AtlasG4EventUserInfo.h:122
AtlasG4EventUserInfo::SetCurrentGenParticle
void SetCurrentGenParticle(HepMC::GenParticlePtr p)
set m_currentGenParticle, the pointer to the GenParticle corresponding to the current G4Track.
Definition: AtlasG4EventUserInfo.h:71
AtlasG4EventUserInfo::m_lastProcessedStep
int m_lastProcessedStep
Definition: AtlasG4EventUserInfo.h:135
AtlasG4EventUserInfo::AtlasG4EventUserInfo
AtlasG4EventUserInfo(const EventContext &ctx)
Definition: AtlasG4EventUserInfo.h:26
AtlasG4EventUserInfo::GetHitCollectionMap
std::shared_ptr< HitCollectionMap > GetHitCollectionMap() const
Get the HitCollectionMap object with shared ownership.
Definition: AtlasG4EventUserInfo.h:110