ATLAS Offline Software
Loading...
Searching...
No Matches
TRT_StrawStatus.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// TRT_StrawStatus.cxx, (c) ATLAS Detector software
10#include "Identifier/Identifier.h"
12
16
17#include "TrkTrack/Track.h"
19
21
23
25
28
30
31#include <fstream>
32#include <format>
33#include <system_error>
34
36
37const size_t nBarrelStraws { 1642 };
38const size_t nEndcapStraws { 3840 };
40
41//================ Constructor =================================================
42
43InDet::TRT_StrawStatus::TRT_StrawStatus(const std::string& name, ISvcLocator* pSvcLocator)
44:
45AthAlgorithm(name,pSvcLocator),
47m_TRTHelper(nullptr)
48{}
49
50//================ Destructor =================================================
51
54
55
56//================ Initialisation =================================================
57
59{
60
61 m_accumulateHits = std::make_unique<ACCHITS_t>();
62 assert( (*m_accumulateHits)[0][0].size() == nAllStraws );
63 clear();
64
65 // Code entered here will be executed once at program start.
66 // Initialize ReadHandleKey
67 ATH_CHECK(m_eventInfoKey.initialize());
68 ATH_CHECK(m_rdoContainerKey.initialize());
69 ATH_CHECK(m_tracksName.initialize());
70 ATH_CHECK(m_vxContainerKey.initialize());
71 ATH_CHECK( detStore()->retrieve(m_TRTHelper, "TRT_ID"));
73 ATH_CHECK( m_trt_hole_finder.retrieve() );
74 ATH_CHECK(m_updator.retrieve());
75 ATH_CHECK( m_mapSvc.retrieve());
76 ATH_CHECK( m_DCSSvc.retrieve());
77 ATH_MSG_DEBUG( "initialize() successful in " << name() << ", retrieved: ..., locR_cut = " << m_locR_cut << " mm." );
78 return StatusCode::SUCCESS;
79}
80
81//================ Finalisation =================================================
82
86 // Code entered here will be executed once at the end of the program run.
87 m_accumulateHits.reset();
88 return StatusCode::SUCCESS;
89}
90
91//================ Execution ====================================================
92
94
96 StatusCode sc = StatusCode::SUCCESS;
97 if (not eventInfo.isValid()) {
98 ATH_MSG_ERROR( "Unable to retrieve Event Info " );
99 return StatusCode::FAILURE;
100 }
101 int runNumber = (int) eventInfo->runNumber();
102 if (runNumber != m_runNumber) {
103 if (m_nEvents) { ATH_CHECK( reportResults() ); clear(); }
104 m_runNumber = runNumber;
105 }
106 int lumiBlock0 =eventInfo->lumiBlock();
108
109 if (not rdoContainer.isValid()) {
110 ATH_MSG_ERROR( "no TRT_RDO container available " );
111 return StatusCode::FAILURE;
112 }
114 if (not trkCollection.isValid()) {
115 ATH_MSG_ERROR( "Could not find Tracks Collection: " << m_tracksName );
116 return StatusCode::FAILURE;
117 }
118
119 //================ Event selection
120
122 if (not vertices.isValid()) {
123 ATH_MSG_DEBUG ("Couldn't retrieve VertexContainer with key: PrimaryVertices");
124 return StatusCode::SUCCESS; // just skip to next event in case of no vertexcontainer
125 }
126
127 int countVertices(0);
128 for (const xAOD::Vertex* vx : *(vertices.cptr()) ) {
129 if (vx->vertexType() == xAOD::VxType::PriVtx) {
130 if ( vx-> nTrackParticles() >= 3) countVertices++;
131 }
132 }
133 if (countVertices < 1) {
134 ATH_MSG_INFO( "no vertices with at least 3 tracks found" );
135 return StatusCode::SUCCESS;
136 }
137
138 if (m_skipBusyEvents) { // cosmic running
139 int countRDOhitsInEvent(0);
140 for (TRT_RDO_Container::const_iterator rdoIt = rdoContainer->begin(); rdoIt != rdoContainer->end(); ++rdoIt) {
141 const InDetRawDataCollection<TRT_RDORawData>* TRTCollection(*rdoIt);
142 if (not TRTCollection) continue;
143 for (DataVector<TRT_RDORawData>::const_iterator trtIt = TRTCollection->begin(); trtIt != TRTCollection->end(); ++trtIt) {
144 countRDOhitsInEvent++;
145 }
146 }
147 if (countRDOhitsInEvent>100000) {
148 ATH_MSG_INFO( "N RDO hits in event greater than 100000: " << countRDOhitsInEvent << ", exiting" );
149 return sc;
150 }
151
152 if (trkCollection->size() > 10) {
153 ATH_MSG_INFO( "N tracks greater than 10: " << trkCollection->size() << ", exiting" );
154 return sc;
155 }
156 }
157
158 //================ End event selection
159
160 //================ Loop over all tracks, accumulate hits on track, also find holes on track
161
162 std::vector<Identifier> holeIdentifiers;
163 std::vector<Identifier> holeIdentifiersWithHits; // holes on straws that have hits, it is just that the hit was not associalted to a track
164 for ( DataVector<Trk::Track>::const_iterator trackIt = trkCollection->begin(); trackIt != trkCollection->end(); ++trackIt) {
165 const Trk::Track *track = *trackIt;
166 //=== select track
167 const Trk::Perigee* perigee = (*trackIt)->perigeeParameters();
168 if ( not perigee ) { ATH_MSG_ERROR( "Trk::Perigee missing" ); continue; }
169 if ( std::abs(perigee->pT())/CLHEP::GeV < 1. ) continue; // 1 GeV pT cut
170
171 const DataVector<const Trk::TrackStateOnSurface>* trackStates = (**trackIt).trackStateOnSurfaces();
172 if ( not trackStates ) { ATH_MSG_ERROR( "Trk::TrackStateOnSurface empty" ); continue; }
173
174 int n_pixel_hits(0), n_sct_hits(0), n_trt_hits(0); // count hits, require minimal number of all hits
175 for ( DataVector<const Trk::TrackStateOnSurface>::const_iterator trackStatesIt = trackStates->begin(); trackStatesIt != trackStates->end(); ++trackStatesIt ) {
176 if ( *trackStatesIt == nullptr ) { ATH_MSG_ERROR( "*trackStatesIt == 0" ); continue; }
177
178 if ( !((*trackStatesIt)->type(Trk::TrackStateOnSurface::Measurement)) ) continue; // this skips outliers
179
180 if ( dynamic_cast<const InDet::TRT_DriftCircleOnTrack*> ( (*trackStatesIt)->measurementOnTrack() ) ) n_trt_hits++;
181 else if ( dynamic_cast<const InDet::SCT_ClusterOnTrack*> ( (*trackStatesIt)->measurementOnTrack() ) ) n_sct_hits++;
182 else if( dynamic_cast<const InDet::PixelClusterOnTrack*> ( (*trackStatesIt)->measurementOnTrack() ) ) n_pixel_hits++;
183 }
184 if (n_pixel_hits<2 || n_sct_hits < 6 || n_trt_hits<15) continue; // end count hits
185
186 //=== end select track
187
188 //=== loop over all hits on track, accumulate them
189
190 for ( DataVector<const Trk::TrackStateOnSurface>::const_iterator trackStatesIt = trackStates->begin(); trackStatesIt != trackStates->end(); ++trackStatesIt ) {
191
192 if ( *trackStatesIt == nullptr ) { ATH_MSG_ERROR( "*trackStatesIt == 0" ); continue; }
193
194 if ( !((*trackStatesIt)->type(Trk::TrackStateOnSurface::Measurement)) ) continue; // this skips outliers
195
196 const InDet::TRT_DriftCircleOnTrack *driftCircleOnTrack = dynamic_cast<const InDet::TRT_DriftCircleOnTrack *>( (*trackStatesIt)->measurementOnTrack() );
197 if ( not driftCircleOnTrack ) continue; // not TRT measurement - this way, keep both hits and outliers
198
199 const Trk::TrackStateOnSurface& hit = **trackStatesIt;
200
201 const Trk::TrackParameters* unbiased_track_parameters = m_updator->removeFromState( *(hit.trackParameters()),
203 hit.measurementOnTrack()->localCovariance()).release();
204
205 double unbiased_locR = unbiased_track_parameters->parameters()[Trk::locR];
206 if ( std::abs(unbiased_locR) > m_locR_cut ) continue; // same cut as the default hole search cut
207
208 const InDet::TRT_DriftCircle *driftCircle = driftCircleOnTrack->prepRawData();
209 if ( driftCircle == nullptr ) { ATH_MSG_ERROR( "driftCircle == 0" ); continue; }
210
211 Identifier id = driftCircle->identify();
212 int index[6]; myStrawIndex(id, index); // side, layer, phi, straw_layer, straw_within_layer, straw_index
213 (*m_accumulateHits)[(index[0]>0)?0:1][index[2]][index[5]][1]++; // accumulate hits on track
214 if (driftCircle->highLevel()) (*m_accumulateHits)[(index[0]>0)?0:1][index[2]][index[5]][3]++; // accumulate hits on track
215
216 } // end trackStatesIt loop
217
218 // add holeIdentifiers - fill vector
219
220 std::unique_ptr<const Trk::TrackStates> holes (m_trt_hole_finder->getHolesOnTrack( *track ));
221 if ( holes==nullptr ) continue; // no holes found
222 for (const Trk::TrackStateOnSurface* trackStates : *holes) {
223
224 if ( !trackStates->type( Trk::TrackStateOnSurface::Hole ) ) { ATH_MSG_ERROR( "m_trt_hole_finder returned something that is not a hole" ); continue; }
225
226 const Trk::TrackParameters* track_parameters = trackStates->trackParameters();
227 if (!track_parameters) { ATH_MSG_WARNING( "m_trt_hole_finder track_parameters missing" ); continue; }
228
230 if ( !(m_TRTHelper->is_trt(id)) ) { ATH_MSG_ERROR( "m_trt_hole_finder returned something that is not a TRT hole" ); continue; }
231
232 // add se same 1.4 mm locR selection, in case it is not on by default
233 if ( std::abs( track_parameters->parameters()[Trk::locR] ) > m_locR_cut ) continue;
234
235 holeIdentifiers.push_back( id );
236 } // end add holeIdentifiers
237
238 } // end trackIt loop
239
240 //================ End loop over all tracks
241
242 //================ Loop over all hits - it includes hits from dead straws that are masked off in drift circle creation
243
244 for (TRT_RDO_Container::const_iterator rdoIt = rdoContainer->begin(); rdoIt != rdoContainer->end(); ++rdoIt) {
245 const InDetRawDataCollection<TRT_RDORawData>* TRTCollection(*rdoIt);
246 if (TRTCollection==nullptr) continue;
247 for (DataVector<TRT_RDORawData>::const_iterator trtIt = TRTCollection->begin(); trtIt != TRTCollection->end(); ++trtIt) {
248 Identifier id = (*trtIt)->identify();
249 int index[6]; myStrawIndex(id, index); // side, layer, phi, straw_layer, straw_within_layer, straw_index
250 (*m_accumulateHits)[(index[0]>0)?0:1][index[2]][index[5]][0]++; // accumulate all hits
251 if ((*trtIt)->highLevel()) (*m_accumulateHits)[(index[0]>0)?0:1][index[2]][index[5]][2]++; // accumulate TR hits
252
253 if (std::find(holeIdentifiers.begin(), holeIdentifiers.end(), id) != holeIdentifiers.end()) // a hole was found on the same straw, but hits is there
254 holeIdentifiersWithHits.push_back( id );
255 }
256 }
257
258 //================ End loop over all hits
259
260 //================ End loop over all holes, each time also save whether the straw with a hole had a hit
261
262 for (unsigned int i=0; i<holeIdentifiers.size(); i++) {
263
264 Identifier id = holeIdentifiers[i];
265
266 int index[6]; myStrawIndex(id, index); // side, layer, phi, straw_layer, straw_within_layer, straw_index
267
268 (*m_accumulateHits)[(index[0]>0)?0:1][index[2]][index[5]][4]++;
269
270 if (std::find(holeIdentifiersWithHits.begin(), holeIdentifiersWithHits.end(), id) != holeIdentifiersWithHits.end())
271 (*m_accumulateHits)[(index[0]>0)?0:1][index[2]][index[5]][5]++;
272 }
273
274 //================ End loop over all hits
275
276 m_nEvents++;
277 last_lumiBlock0 = lumiBlock0;
278 if (m_nEvents%1000==0 && msgLvl(MSG::DEBUG)) ATH_CHECK( reportResults() );
279 return sc;
280}
281
282//============================================================================================
283
285 m_nEvents = 0;
286 *m_accumulateHits = {};
287 return;
288}
289
291 ATH_MSG_INFO( "InDet::TRT_StrawStatus::reportResults() for " << m_nEvents << " events." );
292 char fileName[300];
293 snprintf(fileName, 299,"%s.%07d_newFormat.txt", m_fileName.value().c_str(), m_runNumber);
294 FILE *f = fopen(fileName, "w");
295 if (!f) {
296 ATH_MSG_ERROR( "InDet::TRT_StrawStatus::reportResults: Cannot open " << fileName << " for write" );
297 return StatusCode::FAILURE;
298 }
299 fprintf(f, "%d %d %d %d %d %d %d %d %d \n", 0, 0, 0, 0, 0, 0, 0, 0, m_nEvents);
300 for (size_t i=0; i<2; i++) for (size_t j=0; j<32; j++) for (size_t k=0; k<nAllStraws; k++) {
301 int side = (i>0)?-1:1;
302 if (k>=1642) side *= 2;
303 fprintf(f, "%d %zu %zu", side, j, k);
304 for (int m=0; m<6; m++) fprintf(f, " %d", (*m_accumulateHits)[i][j][k][m]);
305 fprintf(f, "\n");
306 }
307 fclose(f);
308 return StatusCode::SUCCESS;
309}
310
311void
313 ATH_MSG_INFO("InDet::TRT_StrawStatus::printDetailedInformation()");
314 const std::string fileName =
315 std::format("{}.{}{}_printDetailedInformation.txt",
316 m_fileName.value(),
317 std::string(7 - std::min<std::size_t>(7, std::to_string(m_runNumber).size()), '0'),
319 std::ofstream out(fileName, std::ios::out | std::ios::trunc);
320 if (!out) {
321 ATH_MSG_ERROR(std::format("Failed to open '{}' for writing", fileName));
322 return;
323 }
324
325 // Iterate over straw layers
326 for (auto it = m_TRTHelper->straw_layer_begin(); it != m_TRTHelper->straw_layer_end(); ++it) {
327 const Identifier layerId = *it;
328 const int maxStraw = m_TRTHelper->straw_max(layerId);
329
330 for (int i = 0; i <= maxStraw; ++i) {
331 const Identifier id = m_TRTHelper->straw_id(layerId, i);
332 std::array<int, 6> index{};
333 myStrawIndex(id, index.data());
334 int chip = 0;
335 int HVpad = 0;
336 m_TRTStrawNeighbourSvc->getChip(id, chip);
337 m_TRTStrawNeighbourSvc->getPad(id, HVpad);
338 if (!m_printStatusCount) {
339 ATH_MSG_INFO("if the code crashes on the next line, there is a problem with "
340 "m_TRTStrawStatusSummarySvc not being loaded ");
341 ATH_MSG_INFO("in that case, running with reco turned on normally solves the problem, "
342 "know of no better solution at the moment");
343 ATH_MSG_INFO("if you do not need the detailed print information, you can also just set "
344 "printDetailedInformation to 0 to avoid this crash");
346 }
347 const auto& ctx = Gaudi::Hive::currentContext();
348 const int status = m_TRTStrawStatusSummaryTool->get_status(id, ctx);
349 const int statusTemporary = m_TRTStrawStatusSummaryTool->getStatus(id, ctx);
350 const int statusPermanent = m_TRTStrawStatusSummaryTool->getStatusPermanent(id, ctx);
351 // Write: six indices, then chip/HVpad/status trio
352 for (int j = 0; j < 6; ++j) out << index[j] << ' ';
353 out << chip << ' ' << HVpad << ' ' << status << ' '
354 << statusTemporary << ' ' << statusPermanent << '\n';
355 }
356 }
357}
358
359void
361 int side = m_TRTHelper->barrel_ec(id);
362 int layerNumber = m_TRTHelper->layer_or_wheel(id);
363 int strawLayerNumber = m_TRTHelper->straw_layer(id);
364 int strawNumber = m_TRTHelper->straw(id);
365 int straw(0);
366
367 const int numberOfStraws[74] = { 0, 15, 31, 47, 63, 79, 96, 113, 130, 147, 164, 182, 200, 218, 236, 254, 273, 292, 311, 329, // layer 0, 329 straws, strawlayers 0-18
368 348, 368, 388, 408, 428, 448, 469, 490, 511, 532, 553, 575, 597, 619, 641, 663, 686, 709, 732, 755, 778, 802, 826, 849, // layer 1, 520 straws, strawLayers 0-23
369 872, 896, 920, 944, 968, 993, 1018, 1043, 1068, 1093, 1119, 1145, 1171, 1197, 1223, 1250, 1277, 1304, 1331, 1358, 1386, 1414, 1442, 1470, 1498, 1527, 1556, 1585, 1614, 1642 }; // layer 2
370
371 if (abs(side)==1) { // barrel unique straw number
372 if (layerNumber==1) strawLayerNumber+= 19;
373 else if (layerNumber==2) strawLayerNumber+= 43;
374 straw = ( numberOfStraws[strawLayerNumber+1] - strawNumber -1 );
375 } else { // end-cap unique straw number
376 int board = layerNumber;
377 if (board<6) { board *= 2; if (strawLayerNumber>7) board++; }
378 else { board += 6; }
379 straw = board * 192 + strawNumber * 8 + strawLayerNumber % 8 ;
380 straw += 1642;
381 }
382 index[0] = side;
383 index[1] = layerNumber;
384 index[2] = m_TRTHelper->phi_module(id);
385 index[3] = strawLayerNumber;
386 index[4] = strawNumber;
387 index[5] = straw;
388 return;
389}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
An STL vector of pointers that by default owns its pointed-to elements.
static Double_t sc
Handle class for reading from StoreGate.
This is an Identifier helper class for the TRT subdetector.
const size_t nEndcapStraws
const size_t nBarrelStraws
int last_lumiBlock0
const size_t nAllStraws
AthAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
const ServiceHandle< StoreGateSvc > & detStore() const
bool msgLvl(const MSG::Level lvl) const
Derived DataVector<T>.
Definition DataVector.h:795
DataModel_detail::const_iterator< DataVector > const_iterator
Standard const_iterator.
Definition DataVector.h:838
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
Specific class to represent the pixel measurements.
Specific class to represent the SCT measurements.
Represents 'corrected' measurements from the TRT (for example, corrected for wire sag).
virtual const TRT_DriftCircle * prepRawData() const override final
returns the PrepRawData - is a TRT_DriftCircle in this scope
std::unique_ptr< ACCHITS_t > m_accumulateHits
SG::ReadHandleKey< xAOD::VertexContainer > m_vxContainerKey
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfoKey
StatusCode finalize()
standard Athena-Algorithm method
ToolHandle< Trk::ITrackHoleSearchTool > m_trt_hole_finder
Gaudi::Property< int > m_printDetailedInformation
ServiceHandle< ITRT_DCS_ConditionsSvc > m_DCSSvc
void myStrawIndex(Identifier id, int *index)
function that returns straw index (in range 0-5481; 0-1641 for barrel, the rest for endcap) same conv...
TRT_StrawStatus(const std::string &name, ISvcLocator *pSvcLocator)
Standard Athena-Algorithm Constructor.
StatusCode initialize()
standard Athena-Algorithm method
ServiceHandle< ITRT_HWMappingSvc > m_mapSvc
ToolHandle< ITRT_StrawStatusSummaryTool > m_TRTStrawStatusSummaryTool
Gaudi::Property< double > m_locR_cut
std::atomic< int > m_printStatusCount
int m_nEvents
returns index of hardware units: board, chip, pad private fix for now, will call TRTStrawNeighbourSvc...
SG::ReadHandleKey< DataVector< Trk::Track > > m_tracksName
~TRT_StrawStatus()
Default Destructor.
PublicToolHandle< Trk::IUpdator > m_updator
Gaudi::Property< int > m_skipBusyEvents
Gaudi::Property< std::string > m_fileName
ServiceHandle< ITRT_StrawNeighbourSvc > m_TRTStrawNeighbourSvc
StatusCode execute()
standard Athena-Algorithm method
SG::ReadHandleKey< TRT_RDO_Container > m_rdoContainerKey
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type cptr()
Dereference the pointer.
const LocalParameters & localParameters() const
Interface method to get the LocalParameters.
const Amg::MatrixX & localCovariance() const
Interface method to get the localError.
virtual const Surface & associatedSurface() const override=0
Access to the Surface associated to the Parameters.
double pT() const
Access method for transverse momentum.
Identifier associatedDetectorElementIdentifier() const
return Identifier of the associated Detector Element
represents the track state (measurement, material, fit parameters and quality) at a surface.
const MeasurementBase * measurementOnTrack() const
returns MeasurementBase const overload
const TrackParameters * trackParameters() const
return ptr to trackparameters const overload
@ Measurement
This is a measurement, and will at least contain a Trk::MeasurementBase.
@ Hole
A hole on the track - this is defined in the following way.
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ locR
Definition ParamDefs.h:44
ParametersBase< TrackParametersDim, Charged > TrackParameters
Definition index.py:1
@ PriVtx
Primary vertex.
Vertex_v1 Vertex
Define the latest version of the vertex class.