ATLAS Offline Software
Loading...
Searching...
No Matches
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>
15
16#include <G4EventManager.hh>
17#include "G4VUserEventInformation.hh"
18
19#include <memory>
29class AtlasG4EventUserInfo: public G4VUserEventInformation {
30public:
31 AtlasG4EventUserInfo(const EventContext& ctx)
32 : G4VUserEventInformation()
33 , m_eventContext(ctx)
34 , m_eventStore(Atlas::getExtendedEventContext(ctx).proxy())
35 {}
36
41 HepMC::GenEvent* GetHepMCEvent() ;
46 void SetHepMCEvent(HepMC::GenEvent*);
47
64
78
93 void SetLastProcessedTrackID(int trackID) { m_lastProcessedTrackID = trackID; }
94
110 void SetLastProcessedStep(int stepNumber) { m_lastProcessedStep = stepNumber; }
111
116 std::shared_ptr<HitCollectionMap> GetHitCollectionMap() const { return m_hitCollectionMap; }
117
121 void SetHitCollectionMap(std::shared_ptr<HitCollectionMap> hitCollections) { m_hitCollectionMap = hitCollections; }
122
123 const EventContext& GetEventContext() const { return m_eventContext; }
124
126
127 void Print() const {}
128
129 // Static helper method to get the AtlasG4EventUserInfo from G4EventManager
131 {
132 // Event manager may be null in unit tests.
133 G4EventManager* eventManager = G4EventManager::GetEventManager();
134 if (!eventManager) {
135 return nullptr;
136 }
137 // User info may be null
138 return static_cast<AtlasG4EventUserInfo*>(eventManager->GetUserInformation());
139 }
140
141
142private:
143 const EventContext& m_eventContext;
145 HepMC::GenEvent *m_theEvent{};
148
149 std::shared_ptr<HitCollectionMap> m_hitCollectionMap{std::make_shared<HitCollectionMap>()};
150
151 // These next two variables are used by the CaloCalibrationHit
152 // recording code as event-level flags They correspond to the Track
153 // ID and step number of the last G4Step processed by a
154 // CaloCalibrationHit SD Both are needed, because a particle might
155 // have only one step
158};
159
160#endif // MCTRUTH_ATLASG4EVENTUSERINFO_H
void SetHepMCEvent(HepMC::GenEvent *)
set m_theEvent, the pointer to the HepMC::GenEvent used to create the G4Event.
void SetLastProcessedStep(int stepNumber)
record value of the G4Track::GetCurrentStepNumber() for the current G4Step.
void SetCurrentPrimaryGenParticle(HepMC::ConstGenParticlePtr p)
set m_currentPrimaryGenParticle, the pointer to the HepMC::GenParticle used to create the current G4P...
HepMC::GenParticlePtr GetCurrentGenParticle()
return a pointer to the GenParticle corresponding to the current G4Track (if there is one).
std::shared_ptr< HitCollectionMap > m_hitCollectionMap
HepMC::ConstGenParticlePtr GetCurrentPrimaryGenParticle() const
return a pointer to the HepMC::GenParticle used to create the current G4PrimaryParticle.
std::shared_ptr< HitCollectionMap > GetHitCollectionMap() const
Get the HitCollectionMap object with shared ownership.
HepMC::GenParticlePtr m_currentGenParticle
static AtlasG4EventUserInfo * GetEventUserInfo()
int GetLastProcessedTrackID() const
return the value of G4Track::GetTrackID() for the last G4Step processed by a CaloCalibrationHit Sensi...
int GetLastProcessedStep() const
return the value of the G4Track::GetCurrentStepNumber() for the last G4Step processed by a CaloCalibr...
HepMC::GenEvent * m_theEvent
HepMC::ConstGenParticlePtr GetCurrentGenParticle() const
AtlasG4EventUserInfo(const EventContext &ctx)
void SetCurrentGenParticle(HepMC::GenParticlePtr p)
set m_currentGenParticle, the pointer to the GenParticle corresponding to the current G4Track.
void SetHitCollectionMap(std::shared_ptr< HitCollectionMap > hitCollections)
Set the HitCollectionMap object.
const EventContext & m_eventContext
HepMC::GenEvent * GetHepMCEvent()
return a pointer to the HepMC::GenEvent used to create the G4Event.
void SetLastProcessedTrackID(int trackID)
record the value of G4Track::GetTrackID() for the current G4Step.
HepMC::ConstGenParticlePtr m_currentPrimaryGenParticle
const EventContext & GetEventContext() const
GenParticle * GenParticlePtr
Definition GenParticle.h:37
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38