ATLAS Offline Software
Loading...
Searching...
No Matches
MaterialTrackRecorder.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7#include "GaudiKernel/Bootstrap.h"
8#include "GaudiKernel/ISvcLocator.h"
9#include "GaudiKernel/IMessageSvc.h"
10#include "GaudiKernel/ThreadLocalContext.h"
11
13
14#include "G4Event.hh"
15#include "G4Run.hh"
16#include "G4Step.hh"
17
18namespace ActsTrk
19{
21 AthMessaging(Gaudi::svcLocator()->service< IMessageSvc >( "MessageSvc" ),"MaterialTrackRecorder"),
23 {}
24
26 {
27 ATH_MSG_DEBUG(" BeginOfEventAction");
28
29 // Get the event context for the write handle below
30 const EventContext& ctx = Gaudi::Hive::currentContext();
31
32 SG::WriteHandle<RecordedMaterialTrackCollection> materialTracks(m_cfg.materialTrackCollectionName, ctx);
33
34 // Only record if it does NOT already exist
35 if (!materialTracks.isPresent()) {
36 auto coll = std::make_unique<RecordedMaterialTrackCollection>();
37 if (materialTracks.record(std::move(coll)).isFailure()) {
38 ATH_MSG_ERROR("Failed to record RecordedMaterialTrackCollection with key "
39 << m_cfg.materialTrackCollectionName);
40 m_rmtCollection = nullptr;
41 return;
42 }
43 } else {
44 // Sanity check (should normally not happen unless someone else records it)
45 ATH_MSG_DEBUG("Collection already present: " << m_cfg.materialTrackCollectionName);
46 }
47
48 m_rmtCollection = materialTracks.ptr();
49
50 if (!m_rmtCollection) {
51 ATH_MSG_ERROR("Recorded collection pointer is null right after record/present for key "
52 << m_cfg.materialTrackCollectionName);
53 }
54
55 ATH_MSG_DEBUG("ctx event=" << ctx.eventID().event_number()
56 << " slot=" << ctx.slot() << " key=" << m_cfg.materialTrackCollectionName);
57 ATH_MSG_DEBUG("recorded? present=" << materialTracks.isPresent() << " valid=" << materialTracks.isValid());
58 }
59
61 {
62 ATH_MSG_DEBUG(" EndOfEventAction");
63 }
64
66 {
67 ATH_MSG_DEBUG(" BeginOfRunAction");
68 m_rmtCollection = nullptr;
69 }
70
72 {
73 ATH_MSG_DEBUG(" UserSteppingAction");
74
75 if (!m_rmtCollection) {
76 ATH_MSG_ERROR("No per-event collection pointer cached (m_coll is null).");
77 return;
78 }
79
80 // Get the material & check if it is present
81 G4Material* material = step->GetPreStepPoint()->GetMaterial();
82 if (material == nullptr) {
83 return;
84 }
85
86 // First check for exclusion
87 std::string materialName = material->GetName();
88 for (const auto& emat : m_cfg.excludeMaterials) {
89 if (emat == materialName) {
90 ATH_MSG_DEBUG("Exclude step in material '" << materialName << "'.");
91 return;
92 }
93 }
94
95 ATH_MSG_DEBUG("Performing a step with step size = "
96 << s_convertLength * step->GetStepLength());
97
98 // Quantities valid for elemental materials and mixtures
99 double X0 = s_convertLength * material->GetRadlen();
100 double L0 = s_convertLength * material->GetNuclearInterLength();
101 double rho = s_convertDensity * material->GetDensity();
102
103 // Get{A,Z} is only meaningful for single-element materials (according to
104 // the Geant4 docs). Need to compute average manually.
105 const G4ElementVector* elements = material->GetElementVector();
106 const G4double* fraction = material->GetFractionVector();
107 std::size_t nElements = material->GetNumberOfElements();
108
109 double Ar = 0.;
110 double Z = 0.;
111 if (nElements == 1) {
112 Ar = material->GetA() / (CLHEP::gram / CLHEP::mole);
113 Z = material->GetZ();
114 } else {
115 for (std::size_t i = 0; i < nElements; i++) {
116 Ar += elements->at(i)->GetA() * fraction[i] / (CLHEP::gram / CLHEP::mole);
117 Z += elements->at(i)->GetZ() * fraction[i];
118 }
119 }
120
121 // Construct passed material slab for the step
122 const auto slab = Acts::MaterialSlab(Acts::Material::fromMassDensity(X0, L0, Ar, Z, rho), s_convertLength * step->GetStepLength());
123
124 // Create the RecordedMaterialSlab
125 Acts::MaterialInteraction mInteraction;
126 mInteraction.position = convertPosition(step->GetPreStepPoint()->GetPosition());
127 mInteraction.direction = convertDirection(step->GetPreStepPoint()->GetMomentum()).normalized();
128 mInteraction.materialSlab = slab;
129 mInteraction.pathCorrection = (step->GetStepLength() / CLHEP::mm);
130
131 G4Track* g4Track = step->GetTrack();
132
134 Acts::Vector3 vertex = convertPosition(g4Track->GetVertexPosition());
135 Acts::Vector3 direction = convertDirection(g4Track->GetMomentumDirection());
136 rmTrack.first = {vertex, direction};
137 rmTrack.second.materialInteractions.push_back(mInteraction);
138 m_rmtCollection->push_back(rmTrack);
139 }
140}
141
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
Handle class for recording to StoreGate.
virtual void UserSteppingAction(const G4Step *) override
static constexpr double s_convertDensity
Acts::Vector3 convertPosition(const G4ThreeVector &g4vec)
Useful fuction for position coversion.
static constexpr double s_convertLength
Conversion constants.
virtual void EndOfEventAction(const G4Event *) override
Acts::Vector3 convertDirection(const G4ThreeVector &g4vec)
Useful fuction for direction coversion.
virtual void BeginOfRunAction(const G4Run *) override
virtual void BeginOfEventAction(const G4Event *) override
MaterialTrackRecorder(const Config &config)
RecordedMaterialTrackCollection * m_rmtCollection
Pointer to the collection, non-owning.
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
bool isPresent() const
Is the referenced object present in SG?
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
pointer_type ptr()
Dereference the pointer.
The AlignStoreProviderAlg loads the rigid alignment corrections and pipes them through the readout ge...
std::pair< std::pair< Acts::Vector3, Acts::Vector3 >, RecordedMaterial > RecordedMaterialTrack
Recorded material track.
=============================================================================