ATLAS Offline Software
Loading...
Searching...
No Matches
MCTruthSteppingAction.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
8
9#include "G4Event.hh"
10#include "G4Step.hh"
11#include "G4StepPoint.hh"
12#include "G4TouchableHistory.hh"
13
14
15namespace G4UA
16{
17
18 //---------------------------------------------------------------------------
19 // Constructor
20 //---------------------------------------------------------------------------
23 IMessageSvc* msgSvc, MSG::Level level)
24 : AthMessaging(msgSvc, "MCTruthSteppingAction"),
25 m_isInitialized(false),
26 m_volumeCollectionMap(volCollMap)
27 {
28 msg().setLevel(level);
29 }
30
31 //---------------------------------------------------------------------------
32 // Setup the recording envelopes
33 //---------------------------------------------------------------------------
35 {
36 ATH_MSG_DEBUG("Setting up " << m_volumeCollectionMap.size() <<
37 " recording envelopes:");
39 for(const auto& volCollPair : m_volumeCollectionMap) {
40 ATH_MSG_DEBUG(" " << volCollPair.first << ", " << volCollPair.second);
41
42 // Construct the helper in place on the vector
43 m_recordingEnvelopes.emplace_back(volCollPair.first, volCollPair.second);
44 RecordingEnvelope& recEnvelope = m_recordingEnvelopes.back();
45
46 // Make sure the RecEnvelope can initialize properly
47 if(!recEnvelope.Initialize()) {
48 //FIXME - should this be an error?
49 ATH_MSG_WARNING("Envelope volume " << recEnvelope.GetVolumeName() <<
50 " not found in geometry!");
51 ATH_MSG_WARNING("TrackRecordCollection " <<
52 recEnvelope.GetTrackRecordCollectionName() <<
53 " will NOT be recorded");
54 // Throw away uninitialized RecordingEnvelope
55 m_recordingEnvelopes.pop_back();
56 }
57 }
58 if (m_recordingEnvelopes.size() == 0) {
59 ATH_MSG_WARNING("No recording envelopes found!");
60 }
61 }
62
63 //---------------------------------------------------------------------------
64 // Beginning of event
65 //---------------------------------------------------------------------------
67 {
68 // First time initialization
69 if(!m_isInitialized) {
71 if (m_recordingEnvelopes.size() == 0) {
72 ATH_MSG_WARNING("No recording envelopes found!");
73 }
74 m_isInitialized = true;
75 }
76 // Every event initialization
77 for (auto& recEnvelope : m_recordingEnvelopes) {
78 if(auto* eventInfo = static_cast<AtlasG4EventUserInfo*>( event->GetUserInformation())){
79 recEnvelope.BeginOfEvent(eventInfo->GetHitCollectionMap()->Find<TrackRecordCollection>(recEnvelope.GetTrackRecordCollectionName()));
80 }
81 }
82 }
83
84 //---------------------------------------------------------------------------
85 // Process one tracking step
86 //---------------------------------------------------------------------------
88 {
89 if (m_recordingEnvelopes.size() == 0) return;
90 TrackHelper trackHelper(aStep->GetTrack());
91 if (trackHelper.IsSecondary()) return;
92
93 G4StepPoint* preStep = aStep->GetPreStepPoint();
94 G4StepPoint* postStep = aStep->GetPostStepPoint();
95
96 G4VPhysicalVolume* preVol = preStep->GetPhysicalVolume();
97 G4VPhysicalVolume* postVol = postStep->GetPhysicalVolume();
98
99 if (preVol == postVol) return;
100
101 const G4TouchableHistory* preTHist = static_cast<const G4TouchableHistory*>(preStep->GetTouchable());
102 const G4TouchableHistory* postTHist = static_cast<const G4TouchableHistory*>(postStep->GetTouchable());
103 const int preStepVolDepth = preTHist->GetHistoryDepth();
104 const int postStepVolDepth = postTHist->GetHistoryDepth();
105
106 for (auto& recEnvelope : m_recordingEnvelopes)
107 {
108 const int envelopeLevel = recEnvelope.GetLevel();
109 if (envelopeLevel <= preStepVolDepth)
110 {
111 //NB preTHist->GetVolume(preStepVolDepth) would just give us the World volume.
112 const G4LogicalVolume* logicalVolume1 =
113 preTHist->GetVolume(preStepVolDepth-envelopeLevel)->GetLogicalVolume();
114 if (logicalVolume1 != recEnvelope.GetLogicalVolume()) continue;
115
116 if (envelopeLevel <= postStepVolDepth &&
117 logicalVolume1 == postTHist->GetVolume(postStepVolDepth-envelopeLevel)
118 ->GetLogicalVolume())
119 {
120 continue;
121 }
122
123 // We have a track crossing a recording envelope
124 // volume boundary, so make a TrackRecord
125 recEnvelope.AddTrackRecord(aStep);
126
127 // Done with this volume.
128 break;
129 }
130 }
131 }
132
133} // namespace G4UA
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
AtlasHitsVector< TrackRecord > TrackRecordCollection
MsgStream & msg() const
The standard message stream.
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
This class is attached to G4Event objects as UserInformation.
std::map< std::string, std::string > VolumeCollectionMap_t
Map of volume name to output TrackRecordCollection name.
bool m_isInitialized
Used to delay initialization until the event loop, after geo is ready.
void setupRecEnvelopes()
Setup the list of RecordingEnvelope helpers.
std::vector< RecordingEnvelope > m_recordingEnvelopes
List of RecordingEnvelope helpers to invoke.
MCTruthSteppingAction(const VolumeCollectionMap_t &collMap, IMessageSvc *msgSvc, MSG::Level level)
Construct the action with specified volumes and output collections.
VolumeCollectionMap_t m_volumeCollectionMap
Map of volume name to output collection name.
virtual void BeginOfEventAction(const G4Event *) override final
Called at the start of each G4 event.
virtual void UserSteppingAction(const G4Step *) override final
Process one particle step.
Responsible for finding the G4LogicalVolume pointer for each recording envelope and for creating and ...
const std::string & GetTrackRecordCollectionName() const
Returns the name of the TrackRecordCollection to which tracks crossing this recording envelope should...
const std::string & GetVolumeName() const
Returns the name of the recording envelope volume.
bool Initialize()
Finds the pointer to the G4LogicalVolume called m_envelopeVolumeName and the number of levels beneath...
bool IsSecondary() const