ATLAS Offline Software
Loading...
Searching...
No Matches
TRT_StrawStatus.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// 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
32
33const size_t nBarrelStraws { 1642 };
34const size_t nEndcapStraws { 3840 };
36
37//================ Constructor =================================================
38
39InDet::TRT_StrawStatus::TRT_StrawStatus(const std::string& name, ISvcLocator* pSvcLocator)
40:
41AthAlgorithm(name,pSvcLocator),
43m_TRTHelper(nullptr)
44{}
45
46//================ Destructor =================================================
47
50
51
52//================ Initialisation =================================================
53
55{
56
57 m_accumulateHits = std::make_unique<ACCHITS_t>();
58 assert( (*m_accumulateHits)[0][0].size() == nAllStraws );
59 clear();
60
61 // Code entered here will be executed once at program start.
62 // Initialize ReadHandleKey
63 ATH_CHECK(m_eventInfoKey.initialize());
64 ATH_CHECK(m_rdoContainerKey.initialize());
65 ATH_CHECK(m_tracksName.initialize());
66 ATH_CHECK(m_vxContainerKey.initialize());
67 ATH_CHECK( detStore()->retrieve(m_TRTHelper, "TRT_ID"));
69 ATH_CHECK( m_trt_hole_finder.retrieve() );
70 ATH_CHECK(m_updator.retrieve());
71 ATH_CHECK( m_mapSvc.retrieve());
72 ATH_CHECK( m_DCSSvc.retrieve());
73 ATH_MSG_DEBUG( "initialize() successful in " << name() << ", retrieved: ..., locR_cut = " << m_locR_cut << " mm." );
74 return StatusCode::SUCCESS;
75}
76
77//================ Finalisation =================================================
78
82 // Code entered here will be executed once at the end of the program run.
83 m_accumulateHits.reset();
84 return StatusCode::SUCCESS;
85}
86
87//================ Execution ====================================================
88
90
92 StatusCode sc = StatusCode::SUCCESS;
93 if (not eventInfo.isValid()) {
94 ATH_MSG_ERROR( "Unable to retrieve Event Info " );
95 return StatusCode::FAILURE;
96 }
97 int runNumber = (int) eventInfo->runNumber();
98 if (runNumber != m_runNumber) {
99 if (m_nEvents) { ATH_CHECK( reportResults() ); clear(); }
100 m_runNumber = runNumber;
101 }
102 int lumiBlock0 =eventInfo->lumiBlock();
104
105 if (not rdoContainer.isValid()) {
106 ATH_MSG_ERROR( "no TRT_RDO container available " );
107 return StatusCode::FAILURE;
108 }
110 if (not trkCollection.isValid()) {
111 ATH_MSG_ERROR( "Could not find Tracks Collection: " << m_tracksName );
112 return StatusCode::FAILURE;
113 }
114
115 //================ Event selection
116
118 if (not vertices.isValid()) {
119 ATH_MSG_DEBUG ("Couldn't retrieve VertexContainer with key: PrimaryVertices");
120 return StatusCode::SUCCESS; // just skip to next event in case of no vertexcontainer
121 }
122
123 int countVertices(0);
124 for (const xAOD::Vertex* vx : *(vertices.cptr()) ) {
125 if (vx->vertexType() == xAOD::VxType::PriVtx) {
126 if ( vx-> nTrackParticles() >= 3) countVertices++;
127 }
128 }
129 if (countVertices < 1) {
130 ATH_MSG_INFO( "no vertices with at least 3 tracks found" );
131 return StatusCode::SUCCESS;
132 }
133
134 if (m_skipBusyEvents) { // cosmic running
135 int countRDOhitsInEvent(0);
136 for (TRT_RDO_Container::const_iterator rdoIt = rdoContainer->begin(); rdoIt != rdoContainer->end(); ++rdoIt) {
137 const InDetRawDataCollection<TRT_RDORawData>* TRTCollection(*rdoIt);
138 if (not TRTCollection) continue;
139 for (DataVector<TRT_RDORawData>::const_iterator trtIt = TRTCollection->begin(); trtIt != TRTCollection->end(); ++trtIt) {
140 countRDOhitsInEvent++;
141 }
142 }
143 if (countRDOhitsInEvent>100000) {
144 ATH_MSG_INFO( "N RDO hits in event greater than 100000: " << countRDOhitsInEvent << ", exiting" );
145 return sc;
146 }
147
148 if (trkCollection->size() > 10) {
149 ATH_MSG_INFO( "N tracks greater than 10: " << trkCollection->size() << ", exiting" );
150 return sc;
151 }
152 }
153
154 //================ End event selection
155
156 //================ Loop over all tracks, accumulate hits on track, also find holes on track
157
158 std::vector<Identifier> holeIdentifiers;
159 std::vector<Identifier> holeIdentifiersWithHits; // holes on straws that have hits, it is just that the hit was not associalted to a track
160 for ( DataVector<Trk::Track>::const_iterator trackIt = trkCollection->begin(); trackIt != trkCollection->end(); ++trackIt) {
161 const Trk::Track *track = *trackIt;
162 //=== select track
163 const Trk::Perigee* perigee = (*trackIt)->perigeeParameters();
164 if ( not perigee ) { ATH_MSG_ERROR( "Trk::Perigee missing" ); continue; }
165 if ( std::abs(perigee->pT())/CLHEP::GeV < 1. ) continue; // 1 GeV pT cut
166
167 const DataVector<const Trk::TrackStateOnSurface>* trackStates = (**trackIt).trackStateOnSurfaces();
168 if ( not trackStates ) { ATH_MSG_ERROR( "Trk::TrackStateOnSurface empty" ); continue; }
169
170 int n_pixel_hits(0), n_sct_hits(0), n_trt_hits(0); // count hits, require minimal number of all hits
171 for ( DataVector<const Trk::TrackStateOnSurface>::const_iterator trackStatesIt = trackStates->begin(); trackStatesIt != trackStates->end(); ++trackStatesIt ) {
172 if ( *trackStatesIt == nullptr ) { ATH_MSG_ERROR( "*trackStatesIt == 0" ); continue; }
173
174 if ( !((*trackStatesIt)->type(Trk::TrackStateOnSurface::Measurement)) ) continue; // this skips outliers
175
176 if ( dynamic_cast<const InDet::TRT_DriftCircleOnTrack*> ( (*trackStatesIt)->measurementOnTrack() ) ) n_trt_hits++;
177 else if ( dynamic_cast<const InDet::SCT_ClusterOnTrack*> ( (*trackStatesIt)->measurementOnTrack() ) ) n_sct_hits++;
178 else if( dynamic_cast<const InDet::PixelClusterOnTrack*> ( (*trackStatesIt)->measurementOnTrack() ) ) n_pixel_hits++;
179 }
180 if (n_pixel_hits<2 || n_sct_hits < 6 || n_trt_hits<15) continue; // end count hits
181
182 //=== end select track
183
184 //=== loop over all hits on track, accumulate them
185
186 for ( DataVector<const Trk::TrackStateOnSurface>::const_iterator trackStatesIt = trackStates->begin(); trackStatesIt != trackStates->end(); ++trackStatesIt ) {
187
188 if ( *trackStatesIt == nullptr ) { ATH_MSG_ERROR( "*trackStatesIt == 0" ); continue; }
189
190 if ( !((*trackStatesIt)->type(Trk::TrackStateOnSurface::Measurement)) ) continue; // this skips outliers
191
192 const InDet::TRT_DriftCircleOnTrack *driftCircleOnTrack = dynamic_cast<const InDet::TRT_DriftCircleOnTrack *>( (*trackStatesIt)->measurementOnTrack() );
193 if ( not driftCircleOnTrack ) continue; // not TRT measurement - this way, keep both hits and outliers
194
195 const Trk::TrackStateOnSurface& hit = **trackStatesIt;
196
197 const Trk::TrackParameters* unbiased_track_parameters = m_updator->removeFromState( *(hit.trackParameters()),
199 hit.measurementOnTrack()->localCovariance()).release();
200
201 double unbiased_locR = unbiased_track_parameters->parameters()[Trk::locR];
202 if ( std::abs(unbiased_locR) > m_locR_cut ) continue; // same cut as the default hole search cut
203
204 const InDet::TRT_DriftCircle *driftCircle = driftCircleOnTrack->prepRawData();
205 if ( driftCircle == nullptr ) { ATH_MSG_ERROR( "driftCircle == 0" ); continue; }
206
207 Identifier id = driftCircle->identify();
208 int index[6]; myStrawIndex(id, index); // side, layer, phi, straw_layer, straw_within_layer, straw_index
209 (*m_accumulateHits)[(index[0]>0)?0:1][index[2]][index[5]][1]++; // accumulate hits on track
210 if (driftCircle->highLevel()) (*m_accumulateHits)[(index[0]>0)?0:1][index[2]][index[5]][3]++; // accumulate hits on track
211
212 } // end trackStatesIt loop
213
214 // add holeIdentifiers - fill vector
215
216 std::unique_ptr<const Trk::TrackStates> holes (m_trt_hole_finder->getHolesOnTrack( *track ));
217 if ( holes==nullptr ) continue; // no holes found
218 for (const Trk::TrackStateOnSurface* trackStates : *holes) {
219
220 if ( !trackStates->type( Trk::TrackStateOnSurface::Hole ) ) { ATH_MSG_ERROR( "m_trt_hole_finder returned something that is not a hole" ); continue; }
221
222 const Trk::TrackParameters* track_parameters = trackStates->trackParameters();
223 if (!track_parameters) { ATH_MSG_WARNING( "m_trt_hole_finder track_parameters missing" ); continue; }
224
226 if ( !(m_TRTHelper->is_trt(id)) ) { ATH_MSG_ERROR( "m_trt_hole_finder returned something that is not a TRT hole" ); continue; }
227
228 // add se same 1.4 mm locR selection, in case it is not on by default
229 if ( std::abs( track_parameters->parameters()[Trk::locR] ) > m_locR_cut ) continue;
230
231 holeIdentifiers.push_back( id );
232 } // end add holeIdentifiers
233
234 } // end trackIt loop
235
236 //================ End loop over all tracks
237
238 //================ Loop over all hits - it includes hits from dead straws that are masked off in drift circle creation
239
240 for (TRT_RDO_Container::const_iterator rdoIt = rdoContainer->begin(); rdoIt != rdoContainer->end(); ++rdoIt) {
241 const InDetRawDataCollection<TRT_RDORawData>* TRTCollection(*rdoIt);
242 if (TRTCollection==nullptr) continue;
243 for (DataVector<TRT_RDORawData>::const_iterator trtIt = TRTCollection->begin(); trtIt != TRTCollection->end(); ++trtIt) {
244 Identifier id = (*trtIt)->identify();
245 int index[6]; myStrawIndex(id, index); // side, layer, phi, straw_layer, straw_within_layer, straw_index
246 (*m_accumulateHits)[(index[0]>0)?0:1][index[2]][index[5]][0]++; // accumulate all hits
247 if ((*trtIt)->highLevel()) (*m_accumulateHits)[(index[0]>0)?0:1][index[2]][index[5]][2]++; // accumulate TR hits
248
249 if (std::find(holeIdentifiers.begin(), holeIdentifiers.end(), id) != holeIdentifiers.end()) // a hole was found on the same straw, but hits is there
250 holeIdentifiersWithHits.push_back( id );
251 }
252 }
253
254 //================ End loop over all hits
255
256 //================ End loop over all holes, each time also save whether the straw with a hole had a hit
257
258 for (unsigned int i=0; i<holeIdentifiers.size(); i++) {
259
260 Identifier id = holeIdentifiers[i];
261
262 int index[6]; myStrawIndex(id, index); // side, layer, phi, straw_layer, straw_within_layer, straw_index
263
264 (*m_accumulateHits)[(index[0]>0)?0:1][index[2]][index[5]][4]++;
265
266 if (std::find(holeIdentifiersWithHits.begin(), holeIdentifiersWithHits.end(), id) != holeIdentifiersWithHits.end())
267 (*m_accumulateHits)[(index[0]>0)?0:1][index[2]][index[5]][5]++;
268 }
269
270 //================ End loop over all hits
271
272 m_nEvents++;
273 last_lumiBlock0 = lumiBlock0;
274 if (m_nEvents%1000==0 && msgLvl(MSG::DEBUG)) ATH_CHECK( reportResults() );
275 return sc;
276}
277
278//============================================================================================
279
281 m_nEvents = 0;
282 *m_accumulateHits = {};
283 return;
284}
285
287 ATH_MSG_INFO( "InDet::TRT_StrawStatus::reportResults() for " << m_nEvents << " events." );
288 char fileName[300];
289 snprintf(fileName, 299,"%s.%07d_newFormat.txt", m_fileName.value().c_str(), m_runNumber);
290 FILE *f = fopen(fileName, "w");
291 if (!f) {
292 ATH_MSG_ERROR( "InDet::TRT_StrawStatus::reportResults: Cannot open " << fileName << " for write" );
293 return StatusCode::FAILURE;
294 }
295 fprintf(f, "%d %d %d %d %d %d %d %d %d \n", 0, 0, 0, 0, 0, 0, 0, 0, m_nEvents);
296 for (size_t i=0; i<2; i++) for (size_t j=0; j<32; j++) for (size_t k=0; k<nAllStraws; k++) {
297 int side = (i>0)?-1:1;
298 if (k>=1642) side *= 2;
299 fprintf(f, "%d %zu %zu", side, j, k);
300 for (int m=0; m<6; m++) fprintf(f, " %d", (*m_accumulateHits)[i][j][k][m]);
301 fprintf(f, "\n");
302 }
303 fclose(f);
304 return StatusCode::SUCCESS;
305}
306
308 ATH_MSG_INFO( "InDet::TRT_StrawStatus::printDetailedInformation() " );
309 char fileName[300];
310 snprintf(fileName, 299,"%s.%07d_printDetailedInformation.txt", m_fileName.value().c_str(), m_runNumber);
311 FILE *f = fopen(fileName, "w");
312 for (std::vector<Identifier>::const_iterator it = m_TRTHelper->straw_layer_begin(); it != m_TRTHelper->straw_layer_end(); ++it ) {
313 for (int i=0; i<=m_TRTHelper->straw_max( *it); i++) {
314 Identifier id = m_TRTHelper->straw_id( *it, i);
315 int index[6];
316 myStrawIndex(id, index);
317 int chip, HVpad;
318 m_TRTStrawNeighbourSvc->getChip(id, chip);
319 m_TRTStrawNeighbourSvc->getPad(id, HVpad);
320 if (!m_printStatusCount) {
321 ATH_MSG_INFO( "if the code crashes on the next line, there is a problem with m_TRTStrawStatusSummarySvc not being loaded " );
322 ATH_MSG_INFO( "in that case, running with reco turned on normally solves the problem, know of no better solution at the moment" );
323 ATH_MSG_INFO( "if you do not need the detailed print information, you can also just set printDetailedInformation to 0 to avoid this crash" );
325 }
326 int status = m_TRTStrawStatusSummaryTool->get_status( id, Gaudi::Hive::currentContext() );
327 int statusTemporary = m_TRTStrawStatusSummaryTool->getStatus( id, Gaudi::Hive::currentContext() );
328 int statusPermanent = m_TRTStrawStatusSummaryTool->getStatusPermanent( id, Gaudi::Hive::currentContext() );
329 for (int j=0; j<6; j++) fprintf(f, "%d ", index[j]);
330 fprintf(f, "%d %d %d %d %d\n", chip, HVpad, status, statusTemporary, statusPermanent);
331 }
332 }
333 fclose(f);
334 return;
335}
336
338 int side = m_TRTHelper->barrel_ec(id);
339 int layerNumber = m_TRTHelper->layer_or_wheel(id);
340 int strawLayerNumber = m_TRTHelper->straw_layer(id);
341 int strawNumber = m_TRTHelper->straw(id);
342 int straw(0);
343
344 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
345 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
346 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
347
348 if (abs(side)==1) { // barrel unique straw number
349 if (layerNumber==1) strawLayerNumber+= 19;
350 else if (layerNumber==2) strawLayerNumber+= 43;
351 straw = ( numberOfStraws[strawLayerNumber+1] - strawNumber -1 );
352 } else { // end-cap unique straw number
353 int board = layerNumber;
354 if (board<6) { board *= 2; if (strawLayerNumber>7) board++; }
355 else { board += 6; }
356 straw = board * 192 + strawNumber * 8 + strawLayerNumber % 8 ;
357 straw += 1642;
358 }
359 index[0] = side;
360 index[1] = layerNumber;
361 index[2] = m_TRTHelper->phi_module(id);
362 index[3] = strawLayerNumber;
363 index[4] = strawNumber;
364 index[5] = straw;
365 return;
366}
#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.