ATLAS Offline Software
Loading...
Searching...
No Matches
MaterialTrackReader.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#include "GaudiKernel/IInterface.h"
7
8#include "TTree.h"
9#include "TFile.h"
10
11
12ActsTrk::MaterialTrackReader::MaterialTrackReader(const std::string& name, ISvcLocator* pSvcLocator) :
13 AthReentrantAlgorithm(name, pSvcLocator),
15{}
16
22
24{
26
27 if (m_fileNames.empty()) {
28 ATH_MSG_ERROR("No input files given... Please check!!");
29 return StatusCode::FAILURE;
30 }
31
32 // set up the input chain
33 m_inputChain = std::make_unique<TChain>(m_treeName.value().c_str());
34
35 // loop over the input files
36 for (const auto& inputFile : m_fileNames) {
37 // add file to the input chain
38 m_inputChain->Add(inputFile.c_str());
39 ATH_MSG_DEBUG("Adding File " << inputFile << " to tree '" << m_treeName << "'.");
40 }
41
42 // Connect the branches
43 m_accessor.connectForRead(*m_inputChain);
44
45 // get the number of events, which also loads the tree
46 std::size_t nentries = m_inputChain->GetEntries();
47 m_events = static_cast<std::size_t>(m_inputChain->GetMaximum("event_id") + 1);
48 m_batchSize = nentries / m_events;
49
50 ATH_MSG_DEBUG("The full chain has " << nentries << " entries for " << m_events
51 << " events this corresponds to a batch size of: " << m_batchSize);
52
53 return StatusCode::SUCCESS;
54}
55
56StatusCode
57ActsTrk::MaterialTrackReader::execute (const EventContext& ctx) const
58{
59 if (m_inputChain == nullptr) {
60 ATH_MSG_DEBUG("Invalid pointer to input chain");
61 return StatusCode::SUCCESS;
62 }
63
64 if (ctx.evt() >= m_events) {
65 ATH_MSG_INFO("Maximum number of events in file reached, nothing to be done...");
66 return StatusCode::SUCCESS;
67 }
68
69 // Write to the collection to the EventStore
71
72 // Record the collection once per event if not already there
73 if (!materialTracks.isPresent()) {
74 auto coll = std::make_unique<ActsTrk::RecordedMaterialTrackCollection>();
75 ATH_CHECK(materialTracks.record(std::move(coll)));
76 }
77
78 // Add the track to the recorded collection
79 auto* coll = materialTracks.ptr();
80 if (!coll) {
81 ATH_MSG_ERROR("RecordedMaterialTrackCollection ptr() is null for key "
83 return StatusCode::FAILURE;
84 }
85
86 // lock the mutex
87 std::lock_guard<std::mutex> lock(m_readMutex);
88 for (std::size_t index = 0; index < m_batchSize; index++) {
89 // evaluate the entry number based on the batch size
90 auto entry = m_batchSize * ctx.evt() + index;
91
92 ATH_MSG_DEBUG("Reading event / entry (in batch) : " << ctx.evt() << " / " << entry << "(" << index << ")");
93
94 // get the correspoing entry and read it
95 m_inputChain->GetEntry(entry);
96 Acts::RecordedMaterialTrack rmTrack = m_accessor.read();
97
98 ATH_MSG_DEBUG("Track vertex: " << rmTrack.first.first);
99 ATH_MSG_DEBUG("Track momentum:" << rmTrack.first.second);
100
101 // filling the collection
102 coll->push_back(std::move(rmTrack));
103 }
104
105 return StatusCode::SUCCESS;
106
107}
108
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
virtual StatusCode execute(const EventContext &ctx) const override
ActsPlugins::RootMaterialTrackIo m_accessor
The read - write payload.
virtual StatusCode initialize() override
MaterialTrackReader(const std::string &name, ISvcLocator *pSvcLocator)
std::mutex m_readMutex
mutex used to protect multi-threaded reads
Gaudi::Property< std::vector< std::string > > m_fileNames
The list of input filenames.
Gaudi::Property< std::string > m_treeName
The name of the input tree.
std::size_t m_events
The number of events.
std::size_t m_batchSize
The batch size (number of track per events)
SG::WriteHandleKey< RecordedMaterialTrackCollection > m_materialTrackCollectionKey
The RecordedMaterialTrackCollection to write.
Gaudi::Property< bool > m_readCachedSurfaceInformation
Read surface information for the root file.
An algorithm that can be simultaneously executed in multiple threads.
bool isPresent() const
Is the referenced object present in SG?
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
pointer_type ptr()
Dereference the pointer.
Definition index.py:1