ATLAS Offline Software
Loading...
Searching...
No Matches
AthenaTrackingAction.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7#include <iostream>
8
9#include "G4Event.hh"
10#include "G4EventManager.hh"
11
14#include "MCTruth/TrackHelper.h"
17
18namespace G4UA
19{
20
21 //---------------------------------------------------------------------------
22 // Constructor
23 //---------------------------------------------------------------------------
25 int secondarySavingLevel, int subDetVolLevel)
26 : AthMessaging("AthenaTrackingAction")
27 , m_secondarySavingLevel(secondarySavingLevel)
28 , m_subDetVolLevel(subDetVolLevel)
29 {
30 setLevel(lvl);
31 }
32
33 //---------------------------------------------------------------------------
34 // Pre-tracking action.
35 //---------------------------------------------------------------------------
37 {
38 ATH_MSG_DEBUG("Starting to track a new particle");
39
40 // Use the TrackHelper code to identify the kind of particle.
41 TrackHelper trackHelper(track);
42
43 // Condition for storing the GenParticle in the AtlasG4EventUserInfo for later.
44 if (trackHelper.IsPrimary() || trackHelper.IsRegisteredSecondary())
45 {
46 HepMC::GenParticlePtr currentGenParticle = trackHelper.GetTrackInformation()->GetCurrentGenParticle();
47
48 // Assign the GenParticle to the AtlasG4EventUserInfo.
49 AtlasG4EventUserInfo* atlasG4EvtUserInfo = static_cast<AtlasG4EventUserInfo*>
50 (G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetUserInformation());
51 if (trackHelper.IsPrimary()) atlasG4EvtUserInfo->SetCurrentPrimaryGenParticle(currentGenParticle);
52 atlasG4EvtUserInfo->SetCurrentGenParticle(std::move(currentGenParticle));
53 }
54
55 // Condition for creating a trajectory object to store truth.
56 if (trackHelper.IsPrimary() ||
57 (trackHelper.IsRegisteredSecondary() && m_secondarySavingLevel>1) ||
58 (trackHelper.IsSecondary() && m_secondarySavingLevel>2))
59 {
60 ATH_MSG_DEBUG("Preparing an AtlasTrajectory for saving truth");
61
62 // Create a new AtlasTrajectory for this particle
63 AtlasTrajectory* trajectory = new AtlasTrajectory(track, m_subDetVolLevel);
64
65 // Assign the trajectory to the tracking manager.
66 // TODO: consider caching the tracking manager once to reduce overhead.
67 auto trkMgr = G4EventManager::GetEventManager()->GetTrackingManager();
68 //trajectory->setTrackingManager(trkMgr);
69 trkMgr->SetStoreTrajectory(true);
70 trkMgr->SetTrajectory(trajectory);
71 }
72 }
73
74 //---------------------------------------------------------------------------
75 // Post-tracking action.
76 //---------------------------------------------------------------------------
77 void AthenaTrackingAction::PostUserTrackingAction(const G4Track* /*track*/)
78 {
79 ATH_MSG_DEBUG("Finished tracking a particle");
80
81 // We are done tracking this particle, so reset the trajectory.
82 // TODO: consider caching the tracking manager once to reduce overhead.
83 G4EventManager::GetEventManager()->GetTrackingManager()->
84 SetStoreTrajectory(false);
85 }
86
87} // namespace G4UA
#define ATH_MSG_DEBUG(x)
void setLevel(MSG::Level lvl)
Change the current logging level.
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
This class is attached to G4Event objects as UserInformation.
void SetCurrentPrimaryGenParticle(HepMC::ConstGenParticlePtr p)
set m_currentPrimaryGenParticle, the pointer to the HepMC::GenParticle used to create the current G4P...
void SetCurrentGenParticle(HepMC::GenParticlePtr p)
set m_currentGenParticle, the pointer to the GenParticle corresponding to the current G4Track.
Class to store G4 trajectory information.
AthenaTrackingAction(MSG::Level lvl, int secondarySavingLevel, int subDetVolLevel)
Constructor.
int m_subDetVolLevel
The level in the G4 volume hierarchy at which can we find the sub-detector name.
virtual void PreUserTrackingAction(const G4Track *) override final
Called before tracking a new particle.
virtual void PostUserTrackingAction(const G4Track *) override final
Called after tracking a particle.
int m_secondarySavingLevel
The saving level for secondaries.
TrackInformation * GetTrackInformation()
Definition TrackHelper.h:28
bool IsPrimary() const
bool IsSecondary() const
bool IsRegisteredSecondary() const
virtual HepMC::ConstGenParticlePtr GetCurrentGenParticle() const override
return a pointer to the GenParticle corresponding to the current G4Track (if there is one).
GenParticle * GenParticlePtr
Definition GenParticle.h:37