ATLAS Offline Software
Loading...
Searching...
No Matches
TRTMonitoringRun3ESD_Alg.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
14#include "TrkTrack/Track.h"
21
24
25#include <sstream>
26#include <iomanip>
27#include <memory>
28#include <cmath>
29using namespace std;
30
31//Private Static Const data member initialization
34const int TRTMonitoringRun3ESD_Alg::s_Straw_max[2] = {1642, 3840};
35const int TRTMonitoringRun3ESD_Alg::s_iStack_max[2] = {32, 64};
36const int TRTMonitoringRun3ESD_Alg::s_iChip_max[2] = {104, 240};
37const int TRTMonitoringRun3ESD_Alg::s_numberOfStacks[2] = {32, 32};
38const int TRTMonitoringRun3ESD_Alg::s_moduleNum[2] = {96, 64};
39
40TRTMonitoringRun3ESD_Alg::TRTMonitoringRun3ESD_Alg( const std::string& name, ISvcLocator* pSvcLocator )
41:AthMonitorAlgorithm(name,pSvcLocator)
42{
43}
44
46
48 using namespace Monitored;
49
50 ATH_MSG_VERBOSE("Initializing TRT Monitoring");
51
52 // Initialize superclass
54
55 // Retrieve detector manager
56 ATH_CHECK( detStore()->retrieve(m_mgr, "TRT") );
57 // Get ID helper for TRT to access various detector components like straw, straw_layer, layer_or_wheel, phi_module, etc.
58 ATH_CHECK( detStore()->retrieve(m_pTRTHelper, "TRT_ID") );
59 ATH_CHECK( detStore()->retrieve(m_idHelper, "AtlasID") );
60
61 // InDetTrackSelectionTools initialization:
62 ATH_CHECK( m_trackSelTool.retrieve() );
63
64 if (m_doExpert) {
65 // Retrieve the TRT_Straw Status Service
66 if (m_sumTool.name().empty()) {
67 ATH_MSG_WARNING("TRT_StrawStatusTool not given.");
68 } else {
69 ATH_CHECK( m_sumTool.retrieve() );
70 }
71
72 Identifier ident;
73
74 if (m_sumTool.name() != "") {
75 ATH_MSG_VERBOSE("Trying " << m_sumTool << " isGood");
76 ATH_MSG_VERBOSE("TRT_StrawStatusTool reports status = " << m_sumTool->getStatus(ident, Gaudi::Hive::currentContext()));
77 }
78 } //If do expert
79
80 // Retrieve TRT_StrawNeighbourService
81 if (m_TRTStrawNeighbourSvc.name().empty()) {
82 ATH_MSG_WARNING("TRT_StrawNeighbourSvc not given.");
83 } else {
84 if (m_TRTStrawNeighbourSvc.retrieve().isFailure()) {
85 ATH_MSG_FATAL("Could not get StrawNeighbourSvc.");
86 }
87 }
88
89 // Get Track summary tool
90 if (m_TrackSummaryTool.retrieve().isFailure())
91 ATH_MSG_ERROR("Cannot get TrackSummaryTool");
92 else
93 ATH_MSG_DEBUG("Retrieved succesfully the track summary tool" << m_TrackSummaryTool);
94
95 //Get TRTCalDbTool
96 if (m_TRTCalDbTool.name().empty()) {
97 ATH_MSG_WARNING("TRT_CalDbTool not given.");
98 } else {
99 if (m_TRTCalDbTool.retrieve().isFailure()) {
100 ATH_MSG_ERROR("Cannot get TRTCalDBTool.");
101 }
102 }
103
104 ATH_CHECK(m_drifttool.retrieve());
105
106 // Initialize arrays
107 // These arrays store information about each entry to the HitMap histograms
108
109 if (true) {
110 //loop over straw hash index to create straw number mapping for TRTViewer
111 unsigned int maxHash = m_pTRTHelper->straw_layer_hash_max();
112
113 for (int ibe = 0; ibe < 2; ibe++) { // ibe=0(barrel), ibe=1(endcap)
114 for (unsigned int index = 0; index < maxHash; index++) {
115 IdentifierHash idHash = index;
116 Identifier id = m_pTRTHelper->layer_id(idHash);
117 int idBarrelEndcap = m_pTRTHelper->barrel_ec(id);
118 int idLayerWheel = m_pTRTHelper->layer_or_wheel(id);
119 int idPhiModule = m_pTRTHelper->phi_module(id);
120 int idStrawLayer = m_pTRTHelper->straw_layer(id);
121 bool isBarrel = m_pTRTHelper->is_barrel(id);
122 int idSide;
123 int sectorflag = 0;
124 const InDetDD::TRT_BaseElement *element = nullptr;
125
126 if (ibe == 0) { // barrel
127 idSide = idBarrelEndcap ? 1 : -1;
128
129 if (isBarrel && (idBarrelEndcap == -1)) {
130 sectorflag = 1;
131 element = m_mgr->getBarrelElement(idSide, idLayerWheel, idPhiModule, idStrawLayer);
132 }
133 } else if (ibe == 1) { // endcap
134 idSide = idBarrelEndcap ? 1 : 0;
135
136 if (!isBarrel && ((idBarrelEndcap == -2) || (idBarrelEndcap == 2))) {
137 sectorflag = 1;
138 element = m_mgr->getEndcapElement(idSide, idLayerWheel, idStrawLayer, idPhiModule);//, idStrawLayer_ec);
139 }
140 }
141
142 if (sectorflag == 1) {
143 if (!element) continue;
144
145 for (unsigned int istraw = 0; istraw < element->nStraws(); istraw++) {
146 std::vector<Identifier> neighbourIDs;
147 Identifier strawID = m_pTRTHelper->straw_id(id, int(istraw));
148 int i_chip, i_pad;
149 m_TRTStrawNeighbourSvc->getChip(strawID, i_chip);
150 m_TRTStrawNeighbourSvc->getPad(id, i_pad);
151
152 if (ibe == 0) { //barrel
153 if (idLayerWheel == 1) i_chip += 21;
154
155 if (idLayerWheel == 2) i_chip += 54;
156
157 int tempStrawNumber = strawNumber(istraw, idStrawLayer, idLayerWheel);
158
159 if (tempStrawNumber < 0 || tempStrawNumber > (s_Straw_max[ibe] - 1)) {
160 ATH_MSG_WARNING("Found tempStrawNumber = " << tempStrawNumber << " out of range.");
161 } else {
162 m_mat_chip_B.at(idPhiModule).at(tempStrawNumber) = i_chip;
163 m_mat_chip_B.at(idPhiModule + 32).at(tempStrawNumber) = i_chip;
164 }
165 } else if (ibe == 1) { //endcap
166 ++i_chip -= 104;
167 int tempStrawNumber = strawNumberEndCap(istraw, idStrawLayer, idLayerWheel, idPhiModule, idSide);
168
169 if (tempStrawNumber < 0 || tempStrawNumber > (s_Straw_max[ibe] - 1)) {
170 ATH_MSG_WARNING("Found tempStrawNumber = " << tempStrawNumber << " out of range.");
171 } else {
172 m_mat_chip_E.at(idPhiModule).at(tempStrawNumber) = i_chip;
173 m_mat_chip_E.at(idPhiModule + 32).at(tempStrawNumber) = i_chip;
174 }
175 }
176 }
177 }
178 }
179 }
180 }
181
182 // Initialization of VarHandleKeys
183 ATH_CHECK( m_trackCollectionKey.initialize() );
184 ATH_CHECK( m_xAODEventInfoKey.initialize() );
185 ATH_CHECK( m_TRT_BCIDCollectionKey.initialize() );
188
189 return StatusCode::SUCCESS;
190}
191
192
193//----------------------------------------------------------------------------------//
195//----------------------------------------------------------------------------------//
196 switch (LayerNumber) {
197 case 0:
198 return strawLayerNumber;
199
200 case 1:
201 return strawLayerNumber + 19;
202
203 case 2:
204 return strawLayerNumber + 43;
205
206 default:
207 return strawLayerNumber;
208 }
209}
210
211//----------------------------------------------------------------------------------//
212int TRTMonitoringRun3ESD_Alg::strawNumber(int strawNumber, int strawlayerNumber, int LayerNumber) const {
213//----------------------------------------------------------------------------------//
214 int addToStrawNumber = 0;
215 int addToStrawNumberNext = 0;
216 int i = 0;
217 const int numberOfStraws[75] = {
218 0,
219 15,
220 16, 16, 16, 16,
221 17, 17, 17, 17, 17,
222 18, 18, 18, 18, 18,
223 19, 19, 19,
224 18,
225 19,
226 20, 20, 20, 20, 20,
227 21, 21, 21, 21, 21,
228 22, 22, 22, 22, 22,
229 23, 23, 23, 23, 23,
230 24, 24,
231 23, 23,
232 24, 24, 24, 24,
233 25, 25, 25, 25, 25,
234 26, 26, 26, 26, 26,
235 27, 27, 27, 27, 27,
236 28, 28, 28, 28, 28,
237 29, 29, 29, 29,
238 28,
239 0
240 };
241
242 do {
243 i++;
244 addToStrawNumber += numberOfStraws[i - 1];
245 addToStrawNumberNext = addToStrawNumber + numberOfStraws[i];
246 } while (strawLayerNumber(strawlayerNumber, LayerNumber) != i - 1);
247
248 strawNumber = addToStrawNumberNext - strawNumber - 1;
249
250 if (strawNumber < 0 || strawNumber > s_Straw_max[0] - 1) {
251 ATH_MSG_WARNING("strawNumber = " << strawNumber << " out of range. Will set to 0.");
252 strawNumber = 0;
253 }
254
255 return strawNumber;
256}
257
258
259//----------------------------------------------------------------------------------//
260int TRTMonitoringRun3ESD_Alg::strawNumberEndCap(int strawNumber, int strawLayerNumber, int LayerNumber, int phi_stack, int side) const {
261//----------------------------------------------------------------------------------//
262 // Before perfoming map, corrections need to be perfomed.
263 // apply special rotations for endcap mappings
264 // for eca, rotate triplets by 180 for stacks 9-16, and 25-32.
265 static const int TripletOrientation[2][32] = {
266 {
267 1, 1, 1, 1, 1, 1, 1, 1,
268 0, 0, 0, 0, 0, 0, 0, 0,
269 1, 1, 1, 1, 1, 1, 1, 1,
270 0, 0, 0, 0, 0, 0, 0, 0
271 },
272 {
273 1, 1, 1, 1, 1, 1, 1, 1,
274 0, 0, 0, 0, 0, 0, 0, 0,
275 1, 1, 1, 1, 1, 1, 1, 1,
276 0, 0, 0, 0, 0, 0, 0, 0
277 }
278 };
279 int phi1 = -1;
280
281 if (side == 2) phi1 = phi_stack, side = 1;
282 else if (side == -2) phi1 = 31 - phi_stack, side = 0;
283
284 if (phi1 > -1) {
285 if (TripletOrientation[side][phi1]) {
286 // Change straw number from 0-23 in straw layer to 0-192
288
290
291 strawNumber = (192 - 1) * TripletOrientation[side][phi1] + strawNumber * (1 - 2 * TripletOrientation[side][phi1]); //actual rotation
292
293 // Take strawNumber back to 0-23
294 if (strawLayerNumber < 8) strawLayerNumber = int(strawNumber / 24);
295
296 if (strawLayerNumber > 7) strawLayerNumber = int(strawNumber / 24) + 8;
297
299 }
300
301 // Finish rotation
302 // Flip straw in layer.
303
304 if (side == 0) strawNumber = 23 - strawNumber;
305
306 // Finish Flipping
307 }
308
309 // Done with corrections
310 // Start mapping from athena identifiers to TRTViewer maps
311 int strawNumberNew = 0;
312
313 if (LayerNumber < 6 && strawLayerNumber > 7) {
314 strawNumberNew = strawNumberNew + (384 * LayerNumber);
315 strawNumberNew = strawNumberNew + 192 + (strawLayerNumber % 8) + (strawNumber * 8);
316 } else if (LayerNumber < 6 && strawLayerNumber < 8) {
317 strawNumberNew = strawNumberNew + (384 * LayerNumber);
318 strawNumberNew = strawNumberNew + (strawLayerNumber % 8) + (strawNumber * 8);
319 } else if (LayerNumber > 5 && strawLayerNumber > 7) {
320 strawNumberNew = strawNumberNew + 2304 + 192 * (LayerNumber - 6);
321 strawNumberNew = strawNumberNew + 192 + (strawLayerNumber % 8) + (8 * strawNumber);
322 } else if (LayerNumber > 5 && strawLayerNumber < 8) {
323 strawNumberNew = strawNumberNew + 2304 + 192 * (LayerNumber - 6);
324 strawNumberNew = strawNumberNew + (strawLayerNumber % 8) + (8 * strawNumber);
325 }
326
327 strawNumber = strawNumberNew;
328
329 if (strawNumber < 0 || strawNumber > s_Straw_max[1] - 1) {
330 ATH_MSG_WARNING("strawNumber = " << strawNumber << " out of range. Will set to 0.");
331 strawNumber = 0;
332 }
333
334 return strawNumber;
335}
336
337//----------------------------------------------------------------------------------//
338float TRTMonitoringRun3ESD_Alg::radToDegrees(float radValue) const {
339//----------------------------------------------------------------------------------//
340 float degreeValue = radValue / M_PI * 180;
341
342 if (degreeValue < 0) degreeValue += 360;
343
344 return degreeValue;
345}
346
347
348// Check for EventBurst: Counts high level hits and returns true if the count is less than m_passEventBurstCut.
349// If m_EventBurstCut is less than zero, returns allways true
350//----------------------------------------------------------------------------------//
352//----------------------------------------------------------------------------------//
353 if (m_EventBurstCut <= 0) return true;
354
355 int nHLHits = 0;
356 TRT_RDO_Container::const_iterator RDO_CollectionBegin = rdoContainer.begin();
357 TRT_RDO_Container::const_iterator RDO_CollectionEnd = rdoContainer.end();
358
359 for (; RDO_CollectionBegin != RDO_CollectionEnd; ++RDO_CollectionBegin) {
360 const InDetRawDataCollection<TRT_RDORawData> *TRT_Collection(*RDO_CollectionBegin);
361
362 if (!TRT_Collection) continue;
363 else {
364 DataVector<TRT_RDORawData>::const_iterator p_rdo = TRT_Collection->begin();
365
366 for (; p_rdo != TRT_Collection->end(); ++p_rdo) {
367 const TRT_LoLumRawData *p_lolum = dynamic_cast<const TRT_LoLumRawData *>(*p_rdo);
368
369 if (!p_lolum) continue;
370
371 if (p_lolum->highLevel()) nHLHits++;
372 }
373 }
374 }
375
376 if (nHLHits > m_EventBurstCut) return false;
377 else return true;
378}
379
380
381// Fill the TRT Track level histograms
382//----------------------------------------------------------------------------------//
383StatusCode TRTMonitoringRun3ESD_Alg::fillTRTTracks(const EventContext& ctx,
384 const xAOD::TrackParticleContainer& trackCollection,
385 const xAOD::TrigDecision* trigDecision,
386 const ComTime* comTimeObject,
387 const xAOD::EventInfo& eventInfo) const {
388//----------------------------------------------------------------------------------//
389 ATH_MSG_VERBOSE("Filling TRT Tracks Histos");
390
391 // TProfile
392 auto ValidRawDriftTimeonTrkS_x = Monitored::Scalar<float>("ValidRawDriftTimeonTrkS_x", 0.0);
393 auto ValidRawDriftTimeonTrkS_y = Monitored::Scalar<float>("ValidRawDriftTimeonTrkS_y", 0.0);
394 auto ValidRawDriftTimeonTrkC_x = Monitored::Scalar<float>("ValidRawDriftTimeonTrkC_x", 0.0);
395 auto ValidRawDriftTimeonTrkC_y = Monitored::Scalar<float>("ValidRawDriftTimeonTrkC_y", 0.0);
396 auto HitTronTMapC_x = Monitored::Scalar<float>("HitTronTMapC_x", 0.0);
397 auto HitTronTMapC_y = Monitored::Scalar<float>("HitTronTMapC_y", 0.0);
398 auto HitTronTwEPCMapS_x = Monitored::Scalar<float>("HitTronTwEPCMapS_x", 0.0);
399 auto HitTronTwEPCMapS_y = Monitored::Scalar<float>("HitTronTwEPCMapS_y", 0.0);
400 auto HitTronTwEPCMapC_x = Monitored::Scalar<float>("HitTronTwEPCMapC_x", 0.0);
401 auto HitTronTwEPCMapC_y = Monitored::Scalar<float>("HitTronTwEPCMapC_y", 0.0);
402 auto AvgTroTDetPhi_B_Ar_x = Monitored::Scalar<float>("AvgTroTDetPhi_B_Ar_x", 0.0);
403 auto AvgTroTDetPhi_B_Ar_y = Monitored::Scalar<float>("AvgTroTDetPhi_B_Ar_y", 0.0);
404 auto AvgTroTDetPhi_B_x = Monitored::Scalar<float>("AvgTroTDetPhi_B_x", 0.0);
405 auto AvgTroTDetPhi_B_y = Monitored::Scalar<float>("AvgTroTDetPhi_B_y", 0.0);
406 auto AvgTroTDetPhi_E_Ar_x = Monitored::Scalar<float>("AvgTroTDetPhi_E_Ar_x", 0.0);
407 auto AvgTroTDetPhi_E_Ar_y = Monitored::Scalar<float>("AvgTroTDetPhi_E_Ar_y", 0.0);
408 auto AvgTroTDetPhi_E_x = Monitored::Scalar<float>("AvgTroTDetPhi_E_x", 0.0);
409 auto AvgTroTDetPhi_E_y = Monitored::Scalar<float>("AvgTroTDetPhi_E_y", 0.0);
410 auto NumHoTDetPhi_B_x = Monitored::Scalar<float>("NumHoTDetPhi_B_x", 0.0);
411 auto NumHoTDetPhi_B_y = Monitored::Scalar<float>("NumHoTDetPhi_B_y", 0.0);
412 auto NumHoTDetPhi_E_x = Monitored::Scalar<float>("NumHoTDetPhi_E_x", 0.0);
413 auto NumHoTDetPhi_E_y = Monitored::Scalar<float>("NumHoTDetPhi_E_y", 0.0);
414 auto EvtPhaseDetPhi_B_x = Monitored::Scalar<float>("EvtPhaseDetPhi_B_x", 0.0);
415 auto EvtPhaseDetPhi_B_y = Monitored::Scalar<float>("EvtPhaseDetPhi_B_y", 0.0);
416 auto EvtPhaseDetPhi_E_x = Monitored::Scalar<float>("EvtPhaseDetPhi_E_x", 0.0);
417 auto EvtPhaseDetPhi_E_y = Monitored::Scalar<float>("EvtPhaseDetPhi_E_y", 0.0);
418 auto NTrksperLB_x = Monitored::Scalar<float>("NTrksperLB_x", 0.0);
419 auto NTrksperLB_y = Monitored::Scalar<float>("NTrksperLB_y", 0.0);
420
421
422 // TH1F
423 auto DriftTimeonTrkDist_B = Monitored::Scalar<float>("DriftTimeonTrkDist_B", 0.0);
424 auto DriftTimeonTrkDist_B_Ar = Monitored::Scalar<float>("DriftTimeonTrkDist_B_Ar", 0.0);
425 auto DriftTimeonTrkDist_E_Ar = Monitored::Scalar<float>("DriftTimeonTrkDist_E_Ar", 0.0);
426 auto DriftTimeonTrkDist_E = Monitored::Scalar<float>("DriftTimeonTrkDist_E", 0.0);
427 auto NumTrksDetPhi_B = Monitored::Scalar<float>("NumTrksDetPhi_B", 0.0);
428 auto NumTrksDetPhi_E = Monitored::Scalar<float>("NumTrksDetPhi_E", 0.0);
429 auto Pull_Biased_Barrel = Monitored::Scalar<float>("Pull_Biased_Barrel", 0.0);
430 auto Pull_Biased_EndCap = Monitored::Scalar<float>("Pull_Biased_EndCap", 0.0);
431 auto Residual_B = Monitored::Scalar<float>("Residual_B", 0.0);
432 auto Residual_B_Ar = Monitored::Scalar<float>("Residual_B_Ar", 0.0);
433 auto Residual_B_20GeV = Monitored::Scalar<float>("Residual_B_20GeV", 0.0);
434 auto Residual_B_Ar_20GeV = Monitored::Scalar<float>("Residual_B_Ar_20GeV", 0.0);
435 auto Residual_E = Monitored::Scalar<float>("Residual_E", 0.0);
436 auto Residual_E_Ar = Monitored::Scalar<float>("Residual_E_Ar", 0.0);
437 auto Residual_E_20GeV = Monitored::Scalar<float>("Residual_E_20GeV", 0.0);
438 auto Residual_E_Ar_20GeV = Monitored::Scalar<float>("Residual_E_Ar_20GeV", 0.0);
439 auto Residual_noTubeHits_B = Monitored::Scalar<float>("Residual_noTubeHits_B", 0.0);
440 auto Residual_noTubeHits_B_Ar = Monitored::Scalar<float>("Residual_noTubeHits_B_Ar", 0.0);
441 auto Residual_noTubeHits_B_20GeV = Monitored::Scalar<float>("Residual_noTubeHits_B_20GeV", 0.0);
442 auto Residual_noTubeHits_B_Ar_20GeV = Monitored::Scalar<float>("Residual_noTubeHits_B_Ar_20GeV", 0.0);
443 auto Residual_noTubeHits_E = Monitored::Scalar<float>("Residual_noTubeHits_E", 0.0);
444 auto Residual_noTubeHits_E_Ar = Monitored::Scalar<float>("Residual_noTubeHits_E_Ar", 0.0);
445 auto Residual_noTubeHits_E_20GeV = Monitored::Scalar<float>("Residual_noTubeHits_E_20GeV", 0.0);
446 auto Residual_noTubeHits_E_Ar_20GeV = Monitored::Scalar<float>("Residual_noTubeHits_E_Ar_20GeV", 0.0);
447 auto TimeResidual_B = Monitored::Scalar<float>("TimeResidual_B", 0.0);
448 auto TimeResidual_B_Ar = Monitored::Scalar<float>("TimeResidual_B_Ar", 0.0);
449 auto TimeResidual_E = Monitored::Scalar<float>("TimeResidual_E", 0.0);
450 auto TimeResidual_E_Ar = Monitored::Scalar<float>("TimeResidual_E_Ar", 0.0);
451 auto TimeResidual_noTubeHits_B = Monitored::Scalar<float>("TimeResidual_noTubeHits_B", 0.0);
452 auto TimeResidual_noTubeHits_B_Ar = Monitored::Scalar<float>("TimeResidual_noTubeHits_B_Ar", 0.0);
453 auto TimeResidual_noTubeHits_E = Monitored::Scalar<float>("TimeResidual_noTubeHits_E", 0.0);
454 auto TimeResidual_noTubeHits_E_Ar = Monitored::Scalar<float>("TimeResidual_noTubeHits_E_Ar", 0.0);
455 auto TronTDist_E = Monitored::Scalar<float>("TronTDist_E", 0.0);
456 auto TronTDist_B = Monitored::Scalar<float>("TronTDist_B", 0.0);
457 auto TronTDist_B_Ar = Monitored::Scalar<float>("TronTDist_B_Ar", 0.0);
458 auto TronTDist_E_Ar = Monitored::Scalar<float>("TronTDist_E_Ar", 0.0);
459 auto WireToTrkPosition_B_Ar = Monitored::Scalar<float>("WireToTrkPosition_B_Ar", 0.0);
460 auto WireToTrkPosition_B = Monitored::Scalar<float>("WireToTrkPosition_B", 0.0);
461 auto WireToTrkPosition_E_Ar = Monitored::Scalar<float>("WireToTrkPosition_E_Ar", 0.0);
462 auto WireToTrkPosition_E = Monitored::Scalar<float>("WireToTrkPosition_E", 0.0);
463 auto EvtPhase = Monitored::Scalar<float>("EvtPhase", 0.0);
464 auto Summary = Monitored::Scalar<float>("Summary", 0.0);
465 auto SummaryWeight = Monitored::Scalar<float>("SummaryWeight", 0.0);
466
467 // TH2F
468 auto RtRelation_B_Ar_x = Monitored::Scalar<float>("RtRelation_B_Ar_x", 0.0);
469 auto RtRelation_B_Ar_y = Monitored::Scalar<float>("RtRelation_B_Ar_y", 0.0);
470 auto RtRelation_B_x = Monitored::Scalar<float>("RtRelation_B_x", 0.0);
471 auto RtRelation_B_y = Monitored::Scalar<float>("RtRelation_B_y", 0.0);
472 auto RtRelation_E_Ar_x = Monitored::Scalar<float>("RtRelation_E_Ar_x", 0.0);
473 auto RtRelation_E_Ar_y = Monitored::Scalar<float>("RtRelation_E_Ar_y", 0.0);
474 auto RtRelation_E_x = Monitored::Scalar<float>("RtRelation_E_x", 0.0);
475 auto RtRelation_E_y = Monitored::Scalar<float>("RtRelation_E_y", 0.0);
476 auto EvtPhaseVsTrig_x = Monitored::Scalar<float>("EvtPhaseVsTrig_x", 0.0);
477 auto EvtPhaseVsTrig_y = Monitored::Scalar<float>("EvtPhaseVsTrig_y", 0.0);
478
479 // Initialize a bunch of stuff before looping over the track collection. Fill some basic histograms.
480 const float timeCor = comTimeObject ? comTimeObject->getTime() : 0;
481
482 auto p_trk = trackCollection.begin();
483
484 const Trk::Perigee *mPer = nullptr;
485 const DataVector<const Trk::TrackParameters> *AllTrkPar(nullptr);
487
488 int ntrackstack[2][64];
489 int nTotalTracks = 0;
490 int nTracksB[2] = {0, 0};
491 int nTracksEC[2] = {0, 0};
492 int nTracksEC_B[2] = {0, 0};
493 int nTrksperLB_B = 0;
494 int nTrksperLB_E[2] = {0, 0};
495
496 for (int ibe = 0; ibe < 2; ibe++) {
497 std::fill(ntrackstack[ibe], ntrackstack[ibe] + 64, 0);
498 }
499
500for (; p_trk != trackCollection.end(); ++p_trk) {
501
502 uint8_t tempHitsVariable(0);
503 (*p_trk)->summaryValue(tempHitsVariable, xAOD::SummaryType::numberOfTRTHits);
504 int nTRTHits = unsigned(tempHitsVariable);
505
506
507 if (nTRTHits < m_minTRThits) continue;
508
509 AllTrkPar = ((*p_trk)->track())->trackParameters();
510
511 // Search of MeasuredPerigee in TrackParameters
512 // The following algorithm only finds the First perigee measurement.
513 // As there should be one and only one perigee measurement then this assumption should be valid.
514 // But no check is done to see if there is more than one perigee measurement.
515 for (p_trkpariter = AllTrkPar->begin(); p_trkpariter != AllTrkPar->end(); ++p_trkpariter) {
516 // If track parameter does have a measured perigee then the track parameter is a keeper and break out of the loop
517 if ((mPer = dynamic_cast<const Trk::Perigee *>(*p_trkpariter))) break;
518 }
519
520 if (!mPer) continue;
521
522 float theta = mPer->parameters()[Trk::theta];
523 float p = (mPer->parameters()[Trk::qOverP] != 0.) ? std::abs(1. / (mPer->parameters()[Trk::qOverP])) : 10e7;
524 float pT = (p * std::sin(theta));
525 pT = pT * 1e-3; // GeV
526
527 if (p < m_minP) continue;
528
529 const Trk::TrackStates *trackStates = ((*p_trk)->track())->trackStateOnSurfaces();
530
531 if (trackStates == nullptr) continue;
532
533 Trk::TrackStates::const_iterator TSOSItBegin0 = trackStates->begin();
534 Trk::TrackStates::const_iterator TSOSItBegin = trackStates->begin();
535 Trk::TrackStates::const_iterator TSOSItBeginTemp = trackStates->begin();
536 Trk::TrackStates::const_iterator TSOSItEnd = trackStates->end();
537
538 (*p_trk)->summaryValue(tempHitsVariable, xAOD::SummaryType::numberOfTRTHits);
539 int n_trt_hits = unsigned(tempHitsVariable);
540
541 bool is_pT_over_20GeV = false;
542
543 if (mPer->pT() > 20 * CLHEP::GeV) {
544 is_pT_over_20GeV = true;
545 } else {
546 is_pT_over_20GeV = false;
547 }
548
549 const bool cnst_is_pT_over_20GeV = is_pT_over_20GeV;
550
551 const bool passed_track_preselection = (static_cast<bool>(m_trackSelTool->accept(**p_trk)) || m_isCosmics) &&
552 n_trt_hits >= m_min_trt_hits &&
553 mPer->pT() > (m_isCosmics?m_min_pT.value() : 2.0 * CLHEP::GeV); // Hardcoded cut for pT 2.0 GeV for collision setup
554 if (!passed_track_preselection) continue;
555
556 nTotalTracks++;
557 int checkB[2] = {0, 0};
558 int checkEC[2] = {0, 0};
559 int checkEC_B[2] = {0, 0};
560 int nTRTHitsW[2][2];
561 int nTRTHits_side[2][2];
562 int nTRTHitsW_perwheel[2][18];
563 int hitontrack[2] = {0, 0};
564 int hitontrack_E_side[2] = {0, 0};
565
566 for (int ibe = 0; ibe < 2; ibe++) {
567 for (int iside = 0; iside < 2; iside++) {
568 nTRTHits_side[ibe][iside] = -1;
569 nTRTHitsW[ibe][iside] = 0;
570 }
571 std::fill(nTRTHitsW_perwheel[ibe], nTRTHitsW_perwheel[ibe] + 18, 0);
572 }
573
574 int barrel_ec = 0;
575 int layer_or_wheel = 0;
576 int phi_module = 0;
577 int straw_layer = 0;
578 int straw = 0;
579 int nearest_straw_layer[2] = {100, 100};
580 int nearest_straw[2] = {0, 0};
581 int testLayer[2] = {100, 100};
582 float phi2D[2] = {-100, -100};
583
584 for (TSOSItBeginTemp = TSOSItBegin0; TSOSItBeginTemp != TSOSItEnd; ++TSOSItBeginTemp) {
585 if ((*TSOSItBeginTemp) == nullptr) continue;
586
587 if (! ((*TSOSItBeginTemp)->type(Trk::TrackStateOnSurface::Measurement)) ) continue;
588 const InDet::TRT_DriftCircleOnTrack *trtCircle = dynamic_cast<const InDet::TRT_DriftCircleOnTrack *>((*TSOSItBeginTemp)->measurementOnTrack());
589
590 if (!trtCircle) continue;
591 const Trk::TrackParameters *aTrackParam = dynamic_cast<const Trk::TrackParameters *>((*TSOSItBeginTemp)->trackParameters());
592
593 if (!aTrackParam) continue;
594 Identifier DCoTId = trtCircle->identify();
595 barrel_ec = m_pTRTHelper->barrel_ec(DCoTId);
596 int ibe = std::abs(barrel_ec) - 1;
597 layer_or_wheel = m_pTRTHelper->layer_or_wheel (DCoTId);
598 straw_layer = m_pTRTHelper->straw_layer(DCoTId);
599 straw = m_pTRTHelper->straw(DCoTId);
600
601 // Restrict ourselves to the inner most TRT layers To get detector phi.
602 if (layer_or_wheel >= testLayer[ibe]) continue;
603 testLayer[ibe] = layer_or_wheel;
604
605 if (straw_layer < nearest_straw_layer[ibe]) {
606 nearest_straw_layer[ibe] = straw_layer;
607 nearest_straw[ibe] = straw;
608 const InDetDD::TRT_BaseElement *circleElement = nullptr;
609 circleElement = trtCircle->detectorElement();
610 phi2D[ibe] = radToDegrees(circleElement->strawCenter(nearest_straw[ibe]).phi());
611 circleElement = nullptr;
612 }
613 }
614
615 if (phi2D[0] == -999) {
616 ATH_MSG_DEBUG("Track did not go through inner layer of Barrel.");
617 } else {
618 ATH_MSG_VERBOSE("Track's closest approach is m_layer_or_wheel: " <<
619 testLayer[0] << " m_straw_layer: " <<
620 nearest_straw_layer[0] << " (in the Barrel).");
621 }
622
623 if (phi2D[1] == -999) {
624 ATH_MSG_DEBUG("Track did not go through any inner layer of EndCap A or C.");
625 } else {
626 ATH_MSG_VERBOSE("Track's closest approach is m_layer_or_wheel: " <<
627 testLayer[1] << " m_straw_layer: " <<
628 nearest_straw_layer[1] << " (in the EndCaps).");
629 }
630
631 bool trackfound[2][64];
632
633 for (int i = 0; i < 2; i++) {
634 std::fill(trackfound[i], trackfound[i] + 64, false);
635 }
636
637 for (TSOSItBegin = TSOSItBegin0; TSOSItBegin != TSOSItEnd; ++TSOSItBegin) {
638 // Select a TSOS which is non-empty, measurement type and contains both drift circle and track parameters informations
639 if ((*TSOSItBegin) == nullptr) continue;
640
641 if ( !((*TSOSItBegin)->type(Trk::TrackStateOnSurface::Measurement)) ) continue;
642
643 const InDet::TRT_DriftCircleOnTrack *trtCircle = dynamic_cast<const InDet::TRT_DriftCircleOnTrack *>((*TSOSItBegin)->measurementOnTrack());
644
645 if (!trtCircle) continue;
646
647 const Trk::TrackParameters *aTrackParam = dynamic_cast<const Trk::TrackParameters *>((*TSOSItBegin)->trackParameters());
648
649 if (!aTrackParam) continue;
650
651 Identifier DCoTId = trtCircle->identify();
652 barrel_ec = m_pTRTHelper->barrel_ec(DCoTId);
653 layer_or_wheel = m_pTRTHelper->layer_or_wheel(DCoTId);
654 phi_module = m_pTRTHelper->phi_module(DCoTId);
655 straw_layer = m_pTRTHelper->straw_layer(DCoTId);
656 straw = m_pTRTHelper->straw(DCoTId);
657 int ibe = std::abs(barrel_ec) - 1; // ibe = 0 (Barrel), ibe = 1 (Endcap)
658 int iside = barrel_ec > 0 ? 0 : 1; // iside = 0 (Side A), iside = 1 (Side C)
659 int thisStrawNumber[2] = {-1, -1};
660 int chip[2] = {0, 0};
661
662 if (ibe == 0) {
663 thisStrawNumber[ibe] = strawNumber(straw, straw_layer, layer_or_wheel);
664
665 if (thisStrawNumber[ibe] >= 0 && thisStrawNumber[ibe] < s_Straw_max[ibe]) {
666 chip[ibe] = m_mat_chip_B.at(phi_module).at(thisStrawNumber[ibe]);
667 }
668 } else if (ibe == 1) {
669 thisStrawNumber[ibe] = strawNumberEndCap(straw, straw_layer, layer_or_wheel, phi_module, barrel_ec);
670
671 if (thisStrawNumber[ibe] >= 0 && thisStrawNumber[ibe] < s_Straw_max[ibe]) {
672 chip[ibe] = m_mat_chip_E.at(phi_module).at(thisStrawNumber[ibe]);
673 }
674 } else {
675 thisStrawNumber[ibe] = -1;
676 }
677
678 if (thisStrawNumber[ibe] < 0 || thisStrawNumber[ibe] >= s_Straw_max[ibe]) continue;
679
680 if (checkB[iside] == 0 && ibe == 0) {
681 nTracksB[iside]++;
682 checkB[iside] = 1;
683 }
684
685 if (checkEC[iside] == 0 && ibe == 1) {
686 nTracksEC[iside]++;
687 checkEC[iside] = 1;
688 }
689
690 if (checkEC_B[iside] == 0 && checkB[iside] == 1 && ibe == 1 ) {
691 nTracksEC_B[iside]++;
692 checkEC_B[iside] = 1;
693 }
694
695 Identifier surfaceID;
696 const Trk::MeasurementBase *mesb = (*TSOSItBegin)->measurementOnTrack();
697 surfaceID = trtCircle->identify();
698 const bool isArgonStraw = ( Straw_Gastype( m_sumTool->getStatusHT(surfaceID, ctx) ) == GasType::Ar );
699 // Assume always Xe if m_ArgonXenonSplitter is not enabled, otherwise check the straw status (good is Xe, non-good is Ar)
700 float temp_locr = aTrackParam->parameters()[Trk::driftRadius];
701 TRTCond::RtRelation const *rtr = m_TRTCalDbTool->getRtRelation(surfaceID);
702 int iphi_module = -9999;
703
704 if (iside == 0) iphi_module = phi_module;
705 else if (iside == 1) iphi_module = phi_module + 32;
706
707 if (iphi_module >= 0 && iphi_module < 64) trackfound[ibe][iphi_module] = true;
708 else ATH_MSG_ERROR("Variable iphi_module is out of range!");
709
710 if (((ibe == 0) && (temp_locr < m_DistToStraw)) ||
711 ((ibe == 1) && ((*TSOSItBegin)->type(Trk::TrackStateOnSurface::Measurement) ||
712 (*TSOSItBegin)->type(Trk::TrackStateOnSurface::Outlier) ||
713 (*TSOSItBegin)->type(Trk::TrackStateOnSurface::Hole)) &&
714 (temp_locr < m_DistToStraw))) {
715 if (m_idHelper->is_trt(DCoTId)) {
716 if (ibe == 0) {
717 hitontrack[ibe]++;
718 } else if (ibe == 1) {
719 hitontrack[ibe]++;
720 hitontrack_E_side[iside]++;
721 }
722 }
723 }
724 const InDet::TRT_DriftCircle *RawDriftCircle = dynamic_cast<const InDet::TRT_DriftCircle *>(trtCircle->prepRawData());
725 bool isTubeHit = (mesb->localCovariance()(Trk::locX, Trk::locX) > 1.0) ? 1 : 0;
726 if (RawDriftCircle) {
727 nTRTHits_side[ibe][iside]++;
729
730 const bool driftTimeValid = RawDriftCircle->driftTimeValid();
731
732 if (driftTimeValid) {
733 const float validRawDriftTime = RawDriftCircle->rawDriftTime();
734
735 if (m_doExpert && m_doStraws) {
736 ValidRawDriftTimeonTrkS_x = thisStrawNumber[ibe];
737 ValidRawDriftTimeonTrkS_y = validRawDriftTime;
738 fill("TRTTrackHistograms"+std::to_string(ibe)+std::to_string(iphi_module), ValidRawDriftTimeonTrkS_x, ValidRawDriftTimeonTrkS_y);
739 }
740
741 if (m_doExpert && m_doChips) {
742 ValidRawDriftTimeonTrkC_x = chip[ibe] - 1;
743 ValidRawDriftTimeonTrkC_y = validRawDriftTime;
744 fill("TRTTrackHistograms"+std::to_string(ibe)+std::to_string(iphi_module), ValidRawDriftTimeonTrkC_x, ValidRawDriftTimeonTrkC_y);
745 }
746 }
747
748 if (m_doShift && m_doStraws) {
749 if (ibe == 0) {
750 if (isArgonStraw) {
751 DriftTimeonTrkDist_B_Ar = RawDriftCircle->rawDriftTime();
752 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), DriftTimeonTrkDist_B_Ar);
753 }
754 else {
755 DriftTimeonTrkDist_B = RawDriftCircle->rawDriftTime();
756 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), DriftTimeonTrkDist_B);
757 }
758 } else if (ibe == 1) {
759 if (isArgonStraw) {
760 DriftTimeonTrkDist_E_Ar = RawDriftCircle->rawDriftTime();
761 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), DriftTimeonTrkDist_E_Ar);
762 }
763 else {
764 DriftTimeonTrkDist_E = RawDriftCircle->rawDriftTime();
765 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), DriftTimeonTrkDist_E);
766 }
767 }
768 }
769
770 float locR_err = 0.0;
771 const AmgSymMatrix(5)* b_err = aTrackParam->covariance();
772
773 if (b_err) {
774 if (!Amg::hasPositiveDiagElems(*b_err)) {
775 ATH_MSG_WARNING("Some diagonal element(s) of the covariance matrix is (are) infinite or smaller than / too close to zero or above the covariance cutoff");
776 }
777 else {
778 locR_err = Amg::error(*b_err, Trk::locR);
779 }
780 } else {
781 ATH_MSG_ERROR("Track parameters have no covariance attached.");
782 }
783
784 float loc_err = Amg::error(trtCircle->localCovariance(), Trk::driftRadius);
785 float locR = aTrackParam->parameters()[Trk::driftRadius];
786 float loc = trtCircle->localParameters()[Trk::driftRadius];
787
788 if (isTubeHit) {
789 bool isOK = false;
790 loc = m_drifttool->driftRadius(RawDriftCircle->rawDriftTime(), DCoTId, t0, isOK);
791
792 if ((loc * locR) < 0) loc = -loc;
793 }
794
795 // Calculate Residuals for hit
796 if (m_doShift && m_doStraws) {
797 bool pull_b_fill;
798 double pull_b = -999.;
799 const double diff_loc_err = std::abs(loc_err-locR_err);
800 if ( diff_loc_err > 0 ) {
801 pull_b = (loc - locR) /diff_loc_err ;
802 pull_b_fill = true;
803 }
804 else pull_b_fill = false;
805 const double thist0 = m_TRTCalDbTool->getT0(surfaceID);
806 const double trkdrifttime = (!rtr) ? 0 : rtr->drifttime(std::abs(locR));
807 const double timeresidual = RawDriftCircle->rawDriftTime() - thist0 - trkdrifttime;
808
809 if (ibe == 0) {
810 if (!isTubeHit) {
811 if (pull_b_fill) {
812 Pull_Biased_Barrel = pull_b;
813 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), Pull_Biased_Barrel);
814 }
815 }
816
817 if (isArgonStraw) {
818 Residual_B_Ar = loc - locR;
819 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), Residual_B_Ar);
820 Residual_noTubeHits_B_Ar = loc - locR;
821 if (!isTubeHit) fill("ShiftTRTTrackHistograms"+std::to_string(ibe), Residual_noTubeHits_B_Ar);
822
823 if (cnst_is_pT_over_20GeV) {
824 Residual_B_Ar_20GeV = loc - locR;
825 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), Residual_B_Ar_20GeV);
826 Residual_noTubeHits_B_Ar_20GeV = loc - locR;
827 if (!isTubeHit) fill("ShiftTRTTrackHistograms"+std::to_string(ibe), Residual_noTubeHits_B_Ar_20GeV);
828 }
829 TimeResidual_B_Ar = timeresidual;
830 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), TimeResidual_B_Ar);
831 TimeResidual_noTubeHits_B_Ar = timeresidual;
832 if (!isTubeHit) fill("ShiftTRTTrackHistograms"+std::to_string(ibe), TimeResidual_noTubeHits_B_Ar);
833 } else {
834 Residual_B = loc - locR;
835 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), Residual_B);
836 Residual_noTubeHits_B = loc - locR;
837 if (!isTubeHit) fill("ShiftTRTTrackHistograms"+std::to_string(ibe), Residual_noTubeHits_B);
838 TimeResidual_B = timeresidual;
839 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), TimeResidual_B);
840 TimeResidual_noTubeHits_B = timeresidual;
841 if (!isTubeHit) fill("ShiftTRTTrackHistograms"+std::to_string(ibe), TimeResidual_noTubeHits_B);
842
843 if (cnst_is_pT_over_20GeV) {
844 Residual_B_20GeV = loc - locR;
845 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), Residual_B_20GeV);
846 Residual_noTubeHits_B_20GeV = loc - locR;
847 if (!isTubeHit) fill("ShiftTRTTrackHistograms"+std::to_string(ibe), Residual_noTubeHits_B_20GeV);
848 }
849 }
850 } else if (ibe == 1) {
851 if (!isTubeHit) {
852 if (pull_b_fill) {
853 Pull_Biased_EndCap = pull_b;
854 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), Pull_Biased_EndCap);
855 }
856 }
857
858 if (isArgonStraw) {
859 Residual_E_Ar = loc - locR;
860 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), Residual_E_Ar);
861 Residual_noTubeHits_E_Ar = loc - locR;
862 if (!isTubeHit) fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), Residual_noTubeHits_E_Ar);
863 TimeResidual_E_Ar = timeresidual;
864 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), TimeResidual_E_Ar);
865 TimeResidual_noTubeHits_E_Ar = timeresidual;
866 if (!isTubeHit) fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), TimeResidual_noTubeHits_E_Ar);
867
868 if (cnst_is_pT_over_20GeV) {
869 Residual_E_Ar_20GeV = loc - locR;
870 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), Residual_E_Ar_20GeV);
871 Residual_noTubeHits_E_Ar_20GeV = loc - locR;
872 if (!isTubeHit) fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), Residual_noTubeHits_E_Ar_20GeV);
873 }
874 } else {
875 Residual_E = loc - locR;
876 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), Residual_E);
877 Residual_noTubeHits_E = loc - locR;
878 if (!isTubeHit) fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), Residual_noTubeHits_E);
879 TimeResidual_E = timeresidual;
880 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), TimeResidual_E);
881 TimeResidual_noTubeHits_E = timeresidual;
882 if (!isTubeHit) fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), TimeResidual_noTubeHits_E);
883
884 if (cnst_is_pT_over_20GeV) {
885 Residual_E_20GeV = loc - locR;
886 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), Residual_E_20GeV);
887 Residual_noTubeHits_E_20GeV = loc - locR;
888 if (!isTubeHit) fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), Residual_noTubeHits_E_20GeV);
889 }
890 }
891 }
892 }
893
894 if (m_doShift) {
895 if (ibe == 0) {
896 if (isArgonStraw) {
897 WireToTrkPosition_B_Ar = locR;
898 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), WireToTrkPosition_B_Ar);
899 } else {
900 WireToTrkPosition_B = locR;
901 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), WireToTrkPosition_B);
902 }
903 } else if (ibe == 1) {
904 if (isArgonStraw) {
905 WireToTrkPosition_E_Ar = locR;
906 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), WireToTrkPosition_E_Ar);
907 } else {
908 WireToTrkPosition_E = locR;
909 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), WireToTrkPosition_E);
910 }
911 }
912 }
913
914 const float LE = (RawDriftCircle->driftTimeBin()) * 3.125;
915 const float EP = timeCor;
916
917 if (m_doShift && m_doStraws) {
918 if (ibe == 0) {
919 if (isArgonStraw) {
920 if (m_isCosmics) {
921 RtRelation_B_Ar_x = LE - EP - t0;
922 RtRelation_B_Ar_y = std::abs(locR);
923 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), RtRelation_B_Ar_x, RtRelation_B_Ar_y);
924 } else {
925 RtRelation_B_Ar_x = LE - t0;
926 RtRelation_B_Ar_y = std::abs(locR);
927 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), RtRelation_B_Ar_x, RtRelation_B_Ar_y);
928 }
929 } else {
930 if (m_isCosmics) {
931 RtRelation_B_x = LE - EP - t0;
932 RtRelation_B_y = std::abs(locR);
933 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), RtRelation_B_x, RtRelation_B_y);
934 } else {
935 RtRelation_B_x = LE - t0;
936 RtRelation_B_y = std::abs(locR);
937 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), RtRelation_B_x, RtRelation_B_y);
938 }
939 }
940 } else if (ibe == 1) {
941 if (isArgonStraw) {
942 if (m_isCosmics) {
943 RtRelation_E_Ar_x = LE - EP - t0;
944 RtRelation_E_Ar_y = std::abs(locR);
945 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), RtRelation_E_Ar_x, RtRelation_E_Ar_y);
946 } else {
947 RtRelation_E_Ar_x = LE - t0;
948 RtRelation_E_Ar_y = std::abs(locR);
949 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), RtRelation_E_Ar_x, RtRelation_E_Ar_y);
950 }
951 } else {
952 if (m_isCosmics) {
953 RtRelation_E_x = LE - EP - t0;
954 RtRelation_E_y = std::abs(locR);
955 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), RtRelation_E_x, RtRelation_E_y);
956 } else {
957 RtRelation_E_x = LE - t0;
958 RtRelation_E_y = std::abs(locR);
959 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), RtRelation_E_x, RtRelation_E_y);
960 }
961 }
962 }
963 }
964
965 const int driftTimeBin = RawDriftCircle->driftTimeBin();
966 const int firstBinHigh = RawDriftCircle->firstBinHigh();
967 const int lastBinHigh = RawDriftCircle->lastBinHigh();
968 const int trailingEdge = RawDriftCircle->trailingEdge();
969 float trailingEdgeScaled = (trailingEdge + 1) * 3.125;
970
971 if (firstBinHigh || lastBinHigh || driftTimeBin > 0 || trailingEdge < 23) nTRTHitsW[ibe][iside]++;
972
973 if ((trailingEdge < 23) &&
974 !(RawDriftCircle->lastBinHigh()) &&
975 !(RawDriftCircle->firstBinHigh())) {
976
977 if (m_doExpert && m_doChips) {
978 HitTronTMapC_x = chip[ibe] - 1;
979 HitTronTMapC_y = trailingEdgeScaled;
980 fill("TRTTrackHistograms"+std::to_string(ibe)+std::to_string(iphi_module), HitTronTMapC_x, HitTronTMapC_y);
981 }
982
983 if (m_doExpert && m_doStraws) {
984 HitTronTwEPCMapS_x = thisStrawNumber[ibe];
985 HitTronTwEPCMapS_y = trailingEdgeScaled - timeCor;
986 fill("TRTTrackHistograms"+std::to_string(ibe)+std::to_string(iphi_module), HitTronTwEPCMapS_x, HitTronTwEPCMapS_y);
987 }
988
989 if (m_doExpert && m_doChips) {
990 HitTronTwEPCMapC_x = chip[ibe] - 1;
991 HitTronTwEPCMapC_y = trailingEdgeScaled - timeCor;
992 fill("TRTTrackHistograms"+std::to_string(ibe)+std::to_string(iphi_module), HitTronTwEPCMapC_x, HitTronTwEPCMapC_y);
993 }
994
995 if (m_doShift && m_doStraws) {
996 if (RawDriftCircle->driftTimeValid()) {
997 if (ibe == 0) {
998 if (isArgonStraw) {
999 TronTDist_B_Ar = trailingEdgeScaled;
1000 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), TronTDist_B_Ar);
1001 AvgTroTDetPhi_B_Ar_x = phi2D[ibe];
1002 AvgTroTDetPhi_B_Ar_y = trailingEdgeScaled;
1003 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), AvgTroTDetPhi_B_Ar_x, AvgTroTDetPhi_B_Ar_y);
1004 } else {
1005 TronTDist_B = trailingEdgeScaled;
1006 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), TronTDist_B);
1007 AvgTroTDetPhi_B_x = phi2D[ibe];
1008 AvgTroTDetPhi_B_y = trailingEdgeScaled;
1009 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), AvgTroTDetPhi_B_x, AvgTroTDetPhi_B_y);
1010 }
1011 } else if (ibe == 1) {
1012 if (isArgonStraw) {
1013 TronTDist_E_Ar = trailingEdgeScaled;
1014 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), TronTDist_E_Ar);
1015 AvgTroTDetPhi_E_Ar_x = phi2D[ibe];
1016 AvgTroTDetPhi_E_Ar_y = trailingEdgeScaled;
1017 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), AvgTroTDetPhi_E_Ar_x, AvgTroTDetPhi_E_Ar_y);
1018 } else {
1019 TronTDist_E = trailingEdgeScaled;
1020 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), TronTDist_E);
1021 AvgTroTDetPhi_E_x = phi2D[ibe];
1022 AvgTroTDetPhi_E_y = trailingEdgeScaled;
1023 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), AvgTroTDetPhi_E_x, AvgTroTDetPhi_E_y);
1024 }
1025 }
1026 }
1027 }
1028 }
1029 }
1030 }
1031
1032 // ToDo: work on the part below
1033 for (int ibe = 0; ibe < 2; ibe++) {
1034 for (int i = 0; i < 64; i++)
1035 if (trackfound[ibe][i])
1036 ntrackstack[ibe][i]++;
1037
1038 if (m_doShift) {
1039 if (ibe == 0) {
1040 if (hitontrack[ibe] >= m_minTRThits) {
1041 NumHoTDetPhi_B_x = phi2D[ibe];
1042 NumHoTDetPhi_B_y = hitontrack[ibe];
1043 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), NumHoTDetPhi_B_x, NumHoTDetPhi_B_y);
1044 }
1045 }
1046
1047 if (ibe == 1) {
1048 if (hitontrack_E_side[0] >= m_minTRThits) {
1049 NumHoTDetPhi_E_x = phi2D[ibe];
1050 NumHoTDetPhi_E_y = hitontrack_E_side[0];
1051 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+"0", NumHoTDetPhi_E_x, NumHoTDetPhi_E_y);
1052 }
1053
1054 if (hitontrack_E_side[1] >= m_minTRThits) {
1055 NumHoTDetPhi_E_x = phi2D[ibe];
1056 NumHoTDetPhi_E_y = hitontrack_E_side[1];
1057 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+"1", NumHoTDetPhi_E_x, NumHoTDetPhi_E_y);
1058 }
1059 }
1060 }
1061
1062 if (phi2D[ibe] < 0) continue;
1063
1064 if (m_doShift) {
1065 if (ibe == 0) {
1066 if (nTRTHitsW[ibe][0] + nTRTHitsW[ibe][1] > 0) {
1067 NumTrksDetPhi_B = phi2D[ibe];
1068 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), NumTrksDetPhi_B);
1069 }
1070 } else if (ibe == 1) {
1071 if (nTRTHitsW[ibe][0] > 0) {
1072 NumTrksDetPhi_E = phi2D[ibe];
1073 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+"0", NumTrksDetPhi_E);
1074 }
1075
1076 if (nTRTHitsW[ibe][1] > 0) {
1077 NumTrksDetPhi_E = phi2D[ibe];
1078 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+"1", NumTrksDetPhi_E);
1079 }
1080 }
1081 }
1082
1083 if (ibe == 0) {
1084 if ((nTRTHitsW[ibe][0] + nTRTHitsW[ibe][1]) > 0) {
1085 nTrksperLB_B++;
1086 }
1087 if (comTimeObject) {
1088 if (m_doShift && (phi2D[ibe] > 0) && (std::abs(timeCor) > 1e-8)) {
1089 EvtPhaseDetPhi_B_x = phi2D[ibe];
1090 EvtPhaseDetPhi_B_y = timeCor;
1091 fill("ShiftTRTTrackHistograms"+std::to_string(ibe), EvtPhaseDetPhi_B_x, EvtPhaseDetPhi_B_y);
1092 }
1093 }
1094 } else if (ibe == 1) {
1095 for (int iside = 0; iside < 2; iside++) {
1096 if (nTRTHitsW[ibe][iside] > 0) nTrksperLB_E[iside]++;
1097 if (comTimeObject) {
1098 if (nTRTHits_side[ibe][iside] > 5 && (std::abs(timeCor)
1099 > 1e-8)) {
1100 if (m_doShift) {
1101 EvtPhaseDetPhi_E_x = phi2D[ibe];
1102 EvtPhaseDetPhi_E_y = timeCor;
1103 fill("ShiftTRTTrackHistograms"+std::to_string(ibe)+std::to_string(iside), EvtPhaseDetPhi_E_x, EvtPhaseDetPhi_E_y);
1104 }
1105 }
1106 }
1107 }
1108 }
1109 }
1110 }
1111
1112 if (comTimeObject) {
1113 if (std::abs(timeCor) > 1e-8) {
1114 if (m_doShift) {
1115 EvtPhase = timeCor;
1116 fill("ShiftTRTTrackHistograms0", EvtPhase);
1117 }
1118
1119 if (m_doShift && trigDecision) {
1120 std::vector<int> trigid;
1121 trigid.clear(); // Trigger ID
1122 // Get bits for trigger after veto
1123 std::vector<unsigned int> level1TAV = trigDecision->tav();
1124
1125 for (unsigned int j = 0; j < 8 && j < level1TAV.size(); ++j) {
1126 for (unsigned int i = 0; i < 32; ++i) {
1127 if ((level1TAV[j] >> i) & 0x1) {
1128 trigid.push_back(i + (j % 8) * 32); // Found the ID
1129 }
1130 }
1131 }
1132
1133 for (unsigned int j = 0; j < trigid.size(); ++j) {
1134 EvtPhaseVsTrig_x = timeCor;
1135 EvtPhaseVsTrig_y = trigid[j];
1136 fill("ShiftTRTTrackHistograms0", EvtPhaseVsTrig_x, EvtPhaseVsTrig_y);
1137 }
1138 }
1139 }
1140 }
1141
1142 if (m_doShift) {
1143 Summary = 0;
1144 SummaryWeight = 1.;
1145 fill("SmryHistograms", SummaryWeight, Summary);
1146
1147 if (m_doTracksMon) {
1148 Summary = 1;
1149 SummaryWeight = nTotalTracks;
1150 fill("SmryHistograms", SummaryWeight, Summary);
1151 Summary = 2;
1152 SummaryWeight = nTracksB[0];
1153 fill("SmryHistograms", SummaryWeight, Summary);
1154 Summary = 3;
1155 SummaryWeight = nTracksB[1];
1156 fill("SmryHistograms", SummaryWeight, Summary);
1157 Summary = 4;
1158 SummaryWeight = nTracksEC[0];
1159 fill("SmryHistograms", SummaryWeight, Summary);
1160 Summary = 5;
1161 SummaryWeight = nTracksEC[1];
1162 fill("SmryHistograms", SummaryWeight, Summary);
1163 Summary = 6;
1164 SummaryWeight = nTracksEC_B[0];
1165 fill("SmryHistograms", SummaryWeight, Summary);
1166 Summary = 7;
1167 SummaryWeight = nTracksEC_B[1];
1168 fill("SmryHistograms", SummaryWeight, Summary);
1169 }
1170
1171 const unsigned int lumiBlock = eventInfo.lumiBlock();
1172 ATH_MSG_VERBOSE("This is lumiblock : " << lumiBlock);
1173 int lastLumiBlock = -99; // ToDo - last lumiblock calculation is not correct
1174 if ((int)lumiBlock != lastLumiBlock) {
1175 lastLumiBlock = lumiBlock;
1176 }
1177 float evtLumiBlock = 1.;
1178 float lumiBlockScale = (evtLumiBlock > 0) ? (1. / evtLumiBlock) : 0;
1179
1180 if (m_doTracksMon && evtLumiBlock > 0) {
1181 NTrksperLB_x = lastLumiBlock;
1182 NTrksperLB_y = (float)nTrksperLB_B * lumiBlockScale;
1183 fill("ShiftTRTTrackHistograms0", NTrksperLB_x, NTrksperLB_y);
1184
1185 for (int iside = 0; iside < 2; iside++) {
1186 NTrksperLB_x = lastLumiBlock;
1187 NTrksperLB_y = (float)nTrksperLB_E[iside] * lumiBlockScale;
1188 fill("ShiftTRTTrackHistograms1"+std::to_string(iside), NTrksperLB_x, NTrksperLB_y);
1189 }
1190
1191 nTrksperLB_B = 0;
1192
1193 for (int iside = 0; iside < 2; iside++) {
1194 nTrksperLB_E[iside] = 0;
1195 }
1196 }
1197 }
1198
1199 ATH_MSG_DEBUG("end of event and lumi block");
1200 //number of events in lumiblock counter setted to zero since it is end of the run or the lumiblock
1201
1202 return StatusCode::SUCCESS;
1203}
1204
1205//----------------------------------------------------------------------------------//
1207 const xAOD::EventInfo& eventInfo, const EventContext& ctx) const {
1208//----------------------------------------------------------------------------------//
1209
1210 auto IntLum = Monitored::Scalar<float>("IntLum", 0.0);
1211 auto LBvsLum = Monitored::Scalar<float>("LBvsLum", 0.0);
1212 auto LBvsTime_x = Monitored::Scalar<float>("LBvsTime_x", 0.0);
1213 auto LBvsTime_y = Monitored::Scalar<float>("LBvsTime_y", 0.0);
1214 auto IntLumWeight = Monitored::Scalar<float>("IntLumWeight", 0.0);
1215 auto LBvsLumWeight = Monitored::Scalar<float>("LBvsLumWeight", 0.0);
1216
1217 int lumiBlockNumber;
1218 int timeStamp;
1219 lumiBlockNumber = eventInfo.lumiBlock();
1220 timeStamp = eventInfo.timeStamp();
1221
1222 int runNumber;
1223 runNumber = eventInfo.runNumber();
1224 // get Online Luminosity
1225 double intLum = (lbDuration(ctx) * lbAverageLuminosity(ctx));
1226 IntLum = 0.5;
1227 IntLumWeight = intLum;
1228 fill("SmryHistograms", IntLumWeight, IntLum);
1229 LBvsLum = lumiBlockNumber;
1230 LBvsLumWeight = intLum;
1231 fill("SmryHistograms", LBvsLumWeight, LBvsLum);
1232 LBvsTime_x = lumiBlockNumber;
1233 LBvsTime_y = timeStamp;
1234 fill("SmryHistograms", LBvsTime_x, LBvsTime_y);
1235
1236 ATH_MSG_VERBOSE("Filling TRT Aging Histos");
1237
1238 auto Trackr_HT = Monitored::Scalar<float>("Trackr_HT", 0.0);
1239 auto Trackr_All = Monitored::Scalar<float>("Trackr_All", 0.0);
1240 auto Trackz_HT = Monitored::Scalar<float>("Trackz_HT", 0.0);
1241 auto Trackz_All = Monitored::Scalar<float>("Trackz_All", 0.0);
1242
1243 auto p_trk = trackCollection.begin();
1244 const Trk::Perigee *perigee = nullptr;
1245 const DataVector<const Trk::TrackParameters> *AllTrkPar(nullptr);
1247
1248 for (; p_trk != trackCollection.end(); ++p_trk) {
1249 AllTrkPar = ((*p_trk)->track())->trackParameters();
1250
1251 for (p_trkpariter = AllTrkPar->begin(); p_trkpariter != AllTrkPar->end(); ++p_trkpariter) {
1252 if ((perigee = dynamic_cast<const Trk::Perigee *>(*p_trkpariter))) break;
1253 }
1254
1255 // If you went through all of the track parameters and found no perigee mearsurement
1256 // then something is wrong with the track and so don't use the track.
1257 // i.e. continue to the next track.
1258 if (!perigee) {
1259 ATH_MSG_DEBUG("No perigee mearsurement found for the track. This entry will not be propogated to aging histograms.");
1260 continue;
1261 }
1262
1263 float track_eta = perigee->eta();
1264 float track_p = (perigee->parameters()[Trk::qOverP] != 0.) ? std::abs(1. / (perigee->parameters()[Trk::qOverP])) : 10e7;
1265 const Trk::TrackStates *trackStates = ((*p_trk)->track())->trackStateOnSurfaces();
1266
1267 if (trackStates == nullptr) continue;
1268
1269 Trk::TrackStates::const_iterator TSOSItBegin = trackStates->begin();
1270 Trk::TrackStates::const_iterator TSOSItEnd = trackStates->end();
1271
1272 uint8_t tempHitsVariable = 0;
1273 (*p_trk)->summaryValue(tempHitsVariable, xAOD::SummaryType::numberOfTRTHits);
1274 int trt_hits = unsigned(tempHitsVariable);
1275
1276 const bool passed_track_preselection = (static_cast<bool>(m_trackSelTool->accept(**p_trk)) || m_isCosmics) &&
1277 trt_hits >= 6. &&
1278 std::abs(track_p) >= 5000.;
1279
1280 if (!passed_track_preselection) continue;
1281
1282 // Now we have hit informations
1283 const Trk::TrackStates *track_states = ((*p_trk)->track())->trackStateOnSurfaces();
1284
1285 if (track_states) {
1286 ATH_MSG_DEBUG("This track has " << track_states->size() << " track states on surface.");
1287 } else {
1288 ATH_MSG_DEBUG("This track has null track states on surface.");
1289 continue;
1290 }
1291
1292 int barrel_ec_side = 0;
1293 int layer_or_wheel = 0;
1294 int phi_module = 0;
1295 int straw_layer = 0;
1296
1297 for (; TSOSItBegin != TSOSItEnd; ++TSOSItBegin) {
1298 if ((*TSOSItBegin) == nullptr) continue;
1299 if ( !((*TSOSItBegin)->type(Trk::TrackStateOnSurface::Measurement)) ) continue;
1300
1301 const InDet::TRT_DriftCircleOnTrack *trtCircle = dynamic_cast<const InDet::TRT_DriftCircleOnTrack *>((*TSOSItBegin)->measurementOnTrack());
1302 const Trk::TrackParameters *aTrackParam = dynamic_cast<const Trk::TrackParameters *>((*TSOSItBegin)->trackParameters());
1303
1304 if (!trtCircle) continue;
1305 if (!aTrackParam) continue;
1306
1307
1308 Identifier DCoTId = trtCircle->identify();
1309 barrel_ec_side = m_pTRTHelper->barrel_ec(DCoTId);
1310 layer_or_wheel = m_pTRTHelper->layer_or_wheel(DCoTId);
1311 phi_module = m_pTRTHelper->phi_module(DCoTId);
1312 straw_layer = m_pTRTHelper->straw_layer(DCoTId);
1313 int Ba_Ec = abs(barrel_ec_side) - 1; // Ba_Ec: 0 is barrel, 1 is Endcap
1314 int Side = barrel_ec_side > 0 ? 0 : 1; // Side : 0 is side_A, 1 is side_C
1315 double xPos = trtCircle->globalPosition().x(); // Global x coordinate
1316 double yPos = trtCircle->globalPosition().y(); // Global y coordinate
1317 double zPos = trtCircle->globalPosition().z(); // Global z coordinate
1318 double RPos = sqrt(xPos * xPos + yPos * yPos);
1319 Identifier surfaceID;
1320 surfaceID = trtCircle->identify();
1321 // Assume always Xe if m_ArgonXenonSplitter is not enabled, otherwise check the straw status (good is Xe, non-good is Ar)
1322 const InDet::TRT_DriftCircle *RawDriftCircle = dynamic_cast<const InDet::TRT_DriftCircle *>(trtCircle->prepRawData());
1323
1324 if (!RawDriftCircle) { //coverity 25097
1325 // This shouldn't happen in normal conditions because trtCircle is a TRT_DriftCircleOnTrack object
1326 ATH_MSG_WARNING("RawDriftCircle object returned null");
1327 continue;
1328 }
1329
1330 int middleHTbit = RawDriftCircle->getWord() & 0x00020000;
1331 //0x00020000 = 0000 0000 0000 0000 0000 0010 0000 0000 0000 0000
1332 bool is_middleHTbit_high = (middleHTbit != 0);
1333 // bool isHighLevel= RawDriftCircle->highLevel();
1334 bool isHighLevel = is_middleHTbit_high; // Hardcoded HT Middle Bit
1335 bool shortStraw = false;
1336 int InputBar = 0;
1337
1338 if (std::abs(track_eta) < 2. && Ba_Ec == 0.) {
1339 if ((layer_or_wheel == 0) && (phi_module < 4 || (phi_module > 7 && phi_module < 12) || (phi_module > 15 && phi_module < 20) || (phi_module > 23 && phi_module < 28))) InputBar = 1;
1340 else if ((runNumber >= 296939) && (layer_or_wheel == 0) && (phi_module > 27)) InputBar = 1;
1341 else if (layer_or_wheel == 0)
1342 InputBar = 0;
1343 else if ((layer_or_wheel == 1) && ((phi_module > 1 && phi_module < 6) || (phi_module > 9 && phi_module < 14) || (phi_module > 17 && phi_module < 22) || (phi_module > 25 && phi_module < 30)))
1344 InputBar = 1;
1345 else if (layer_or_wheel == 1)
1346 InputBar = 0;
1347 else if (layer_or_wheel == 2 && phi_module % 2 != 0)
1348 InputBar = 1;
1349 else if (layer_or_wheel == 2)
1350 InputBar = 0;
1351 else {
1352 ATH_MSG_WARNING("Should not pass here");
1353 continue;
1354 }
1355
1356 if ((layer_or_wheel == 0) && straw_layer < 9.)
1357 shortStraw = true;
1358 }
1359
1360 // Fill Barrel Plots
1361 if ((!shortStraw) && (Ba_Ec == 0)) {
1362 Trackz_All = zPos;
1363 fill("TRTAgingHistograms0"+std::to_string(layer_or_wheel)+std::to_string(InputBar), Trackz_All);
1364 if (isHighLevel) {
1365 Trackz_HT = zPos;
1366 fill("TRTAgingHistograms"+std::to_string(Ba_Ec)+std::to_string(layer_or_wheel)+std::to_string(InputBar), Trackz_HT);
1367 }
1368 }
1369
1370 if (shortStraw) {
1371 if (zPos > 0.) {
1372 Trackz_All = zPos;
1373 fill("TRTAgingHistograms"+std::to_string(Ba_Ec)+"3"+std::to_string(InputBar), Trackz_All);
1374 if (isHighLevel) {
1375 Trackz_HT = zPos;
1376 fill("TRTAgingHistograms"+std::to_string(Ba_Ec)+"3"+std::to_string(InputBar), Trackz_HT);
1377 }
1378 } else {
1379 Trackz_All = zPos;
1380 fill("TRTAgingHistograms"+std::to_string(Ba_Ec)+"4"+std::to_string(InputBar), Trackz_All);
1381
1382 if (isHighLevel) {
1383 Trackz_HT = zPos;
1384 fill("TRTAgingHistograms"+std::to_string(Ba_Ec)+"4"+std::to_string(InputBar), Trackz_HT);
1385 }
1386 }
1387 }
1388
1389 // End of Barrel plots, moving to Endcap plots
1390 int WType = -1;
1391
1392 if ((Ba_Ec == 1) && (layer_or_wheel < 6) &&
1393 ((straw_layer > 3 && straw_layer < 8) ||
1394 (straw_layer > 11))) {
1395 WType = 0;
1396 }
1397 if ((Ba_Ec == 1) && (layer_or_wheel >= 6) &&
1398 (straw_layer > 3)) {
1399 WType = 3;
1400 }
1401 if ((Ba_Ec == 1) && (layer_or_wheel < 6) &&
1402 ((straw_layer > -1 && straw_layer < 4) ||
1403 (straw_layer > 7 && straw_layer < 12))) {
1404 WType = 2;
1405 }
1406 if ((Ba_Ec == 1) && (layer_or_wheel >= 6) &&
1407 ((straw_layer > -1 && straw_layer < 4))) {
1408 WType = 1;
1409 }
1410
1411 if (WType < 0 && Ba_Ec == 1) { // Coverity CID 25096
1412 ATH_MSG_WARNING("The variable \"WType\" is less than zero!.");
1413 continue;
1414 }
1415
1416 if (Ba_Ec == 1) {
1417 Trackr_All = RPos;
1418 fill("TRTAgingHistograms"+std::to_string(Ba_Ec)+std::to_string(WType)+std::to_string(Side), Trackr_All);
1419 if (isHighLevel) {
1420 Trackr_HT = RPos;
1421 fill("TRTAgingHistograms"+std::to_string(Ba_Ec)+std::to_string(WType)+std::to_string(Side), Trackr_HT);
1422 }
1423 }
1424 }
1425 }
1426
1427 return StatusCode::SUCCESS;
1428}
1429
1430
1431StatusCode TRTMonitoringRun3ESD_Alg::fillHistograms( const EventContext& ctx ) const {
1432 using namespace Monitored;
1433 bool passEventBurst = true;
1434
1435 // Declare the quantities which should be monitored
1436
1437 // Set the values of the monitored variables for the event
1438
1439 ATH_MSG_VERBOSE("Monitoring Histograms being filled");
1440
1444
1445 if (!xAODEventInfo.isValid()) {
1446 ATH_MSG_ERROR("Could not find event info object " << m_xAODEventInfoKey.key() <<
1447 " in store");
1448 return StatusCode::FAILURE;
1449 }
1450
1451 if (m_doTracksMon) {
1452 if (!trackCollection.isValid()) {
1453 ATH_MSG_ERROR("Could not find track collection " << m_trackCollectionKey.key() <<
1454 " in store");
1455 return StatusCode::FAILURE;
1456 }
1457 const xAOD::TrigDecision* trigDecision = nullptr;
1458 if (! m_trigDecisionKey.empty()) {
1459 trigDecision = SG::get(m_trigDecisionKey, ctx);
1460 if (!trigDecision) {
1461 ATH_MSG_INFO("Could not find trigger decision object " << m_trigDecisionKey.key() <<
1462 " in store");
1463 }
1464 }
1465 const ComTime *comTimeObject=nullptr;
1466 if (!m_comTimeObjectKey.empty()) {
1467 SG::ReadHandle<ComTime> tmp_comTimeObject(m_comTimeObjectKey, ctx);
1468 if (!tmp_comTimeObject.isValid()) {
1469 // NOTE: failing to retrieve ComTime from store for some reason
1470 ATH_MSG_DEBUG("Could not find com time object " << m_comTimeObjectKey.key() << " in store" );
1471 }
1472 else {
1473 comTimeObject = tmp_comTimeObject.cptr();
1474 }
1475 }
1476 ATH_CHECK( fillTRTTracks(ctx, *trackCollection, trigDecision, comTimeObject, *xAODEventInfo) );
1477 }
1478
1479 if (!m_doTracksMon) {
1480 if (!trackCollection.isValid()) {
1481 ATH_MSG_ERROR("Could not find track collection " << m_trackCollectionKey.key() <<
1482 " in store");
1483 return StatusCode::FAILURE;
1484 }
1485 }
1486
1487 if (passEventBurst) { // ESD files does not have an RDO container to pass event burst, what to do?
1488 ATH_CHECK( fillTRTHighThreshold(*trackCollection, *xAODEventInfo, ctx) );
1489 }
1490
1491
1492
1493 return StatusCode::SUCCESS;
1494}
#define M_PI
Scalar theta() const
theta method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
This class provides an interface to generate or decode an identifier for the upper levels of the dete...
An STL vector of pointers that by default owns its pointed-to elements.
#define AmgSymMatrix(dim)
abstract interface to TRT calibration constants
Abstract interface to information on straws electronic grouping.
static Double_t t0
This is an Identifier helper class for the TRT subdetector.
InDetRawDataContainer< InDetRawDataCollection< TRT_RDORawData > > TRT_RDO_Container
const ServiceHandle< StoreGateSvc > & detStore() const
virtual StatusCode initialize() override
initialize
AthMonitorAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor.
double getTime() const
Definition ComTime.h:44
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.
size_type size() const noexcept
Returns the number of elements in the collection.
const_iterator end() const
return const_iterator for end of container
const_iterator begin() const
return const_iterator for first entry
This is a "hash" representation of an Identifier.
Virtual base class of TRT readout elements.
unsigned int nStraws() const
Number of straws in the element.
const Amg::Vector3D & strawCenter(int straw) const
Straw Surface: Local -> global transform of the straw via integer.
Represents 'corrected' measurements from the TRT (for example, corrected for wire sag).
virtual const InDetDD::TRT_BaseElement * detectorElement() const override final
returns the detector element, assoicated with the PRD of this class
virtual const TRT_DriftCircle * prepRawData() const override final
returns the PrepRawData - is a TRT_DriftCircle in this scope
virtual const Amg::Vector3D & globalPosition() const override final
return the global position of this RIO_OnTrack
unsigned int getWord() const
returns the TRT dataword
bool driftTimeValid() const
return true if the corrected drift time is OK
int driftTimeBin() const
returns the leading edge bin defined as in TRT_LoLumRawData to be the first 0-1 transition
bool lastBinHigh() const
returns true if the last bin is high
bool firstBinHigh() const
returns true if the first bin is high
double rawDriftTime() const
returns the raw driftTime
int trailingEdge() const
returns the trailing edge bin
Declare a monitored scalar variable.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type cptr()
Dereference the pointer.
Base class for rt-relations in the TRT.
Definition RtRelation.h:27
virtual float drifttime(float radius) const =0
drifttime for given radius
ToolHandle< InDet::IInDetTrackSelectionTool > m_trackSelTool
float radToDegrees(float radValue) const
SG::ReadHandleKey< xAOD::TrigDecision > m_trigDecisionKey
int strawNumberEndCap(int strawNumber, int strawLayerNumber, int LayerNumber, int phi_stack, int side) const
const AtlasDetectorID * m_idHelper
Gaudi::Property< int > m_min_trt_hits
Gaudi::Property< bool > m_doStraws
Gaudi::Property< float > m_min_pT
SG::ReadHandleKey< InDetTimeCollection > m_TRT_BCIDCollectionKey
virtual StatusCode fillHistograms(const EventContext &ctx) const override
adds event to the monitoring histograms
ToolHandle< Trk::ITrackSummaryTool > m_TrackSummaryTool
bool checkEventBurst(const TRT_RDO_Container &rdoContainer) const
std::vector< std::vector< unsigned char > > m_mat_chip_E
Gaudi::Property< int > m_minTRThits
Gaudi::Property< bool > m_doShift
StatusCode fillTRTHighThreshold(const xAOD::TrackParticleContainer &trackCollection, const xAOD::EventInfo &eventInfo, const EventContext &ctx) const
Gaudi::Property< float > m_DistToStraw
virtual StatusCode initialize() override
initialize
std::vector< std::vector< unsigned char > > m_mat_chip_B
ToolHandle< ITRT_CalDbTool > m_TRTCalDbTool
Gaudi::Property< float > m_minP
ServiceHandle< ITRT_StrawNeighbourSvc > m_TRTStrawNeighbourSvc
SG::ReadHandleKey< xAOD::EventInfo > m_xAODEventInfoKey
int strawNumber(int strawNumber, int strawlayerNumber, int LayerNumber) const
SG::ReadHandleKey< ComTime > m_comTimeObjectKey
Gaudi::Property< bool > m_doExpert
Gaudi::Property< bool > m_doTracksMon
ToolHandle< ITRT_DriftFunctionTool > m_drifttool
TRTMonitoringRun3ESD_Alg(const std::string &name, ISvcLocator *pSvcLocator)
Gaudi::Property< bool > m_doChips
const InDetDD::TRT_DetectorManager * m_mgr
int strawLayerNumber(int strawLayerNumber, int LayerNumber) const
GasType Straw_Gastype(int stat) const
StatusCode fillTRTTracks(const EventContext &ctx, const xAOD::TrackParticleContainer &trackCollection, const xAOD::TrigDecision *trigDecision, const ComTime *comTimeObject, const xAOD::EventInfo &eventInfo) const
ToolHandle< ITRT_StrawStatusSummaryTool > m_sumTool
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_trackCollectionKey
virtual bool highLevel() const override final
This class is the pure abstract base class for all fittable tracking measurements.
const LocalParameters & localParameters() const
Interface method to get the LocalParameters.
const Amg::MatrixX & localCovariance() const
Interface method to get the localError.
double eta() const
Access method for pseudorapidity - from momentum.
double pT() const
Access method for transverse momentum.
Identifier identify() const
return the identifier -extends MeasurementBase
@ Measurement
This is a measurement, and will at least contain a Trk::MeasurementBase.
@ Outlier
This TSoS contains an outlier, that is, it contains a MeasurementBase/RIO_OnTrack which was not used ...
@ Hole
A hole on the track - this is defined in the following way.
uint32_t lumiBlock() const
The current event's luminosity block number.
uint32_t timeStamp() const
POSIX time in seconds from 1970. January 1st.
uint32_t runNumber() const
The current event's run number.
const std::vector< uint32_t > & tav() const
Get the Trigger After Veto bits.
int trailingEdge(unsigned int m_word)
Definition driftCircle.h:64
int driftTimeBin(unsigned int m_word)
Definition driftCircle.h:50
bool lastBinHigh(unsigned int m_word)
bool firstBinHigh(unsigned int m_word)
virtual float lbAverageLuminosity(const EventContext &ctx=Gaudi::Hive::currentContext()) const
Calculate average luminosity (in ub-1 s-1 => 10^30 cm-2 s-1).
virtual double lbDuration(const EventContext &ctx=Gaudi::Hive::currentContext()) const
Calculate the duration of the luminosity block (in seconds)
double error(const Amg::MatrixX &mat, int index)
return diagonal error of the matrix caller should ensure the matrix is symmetric and the index is in ...
bool hasPositiveDiagElems(const AmgSymMatrix(N) &mat)
Returns true if all diagonal elements of the covariance matrix are finite aka sane in the above defin...
Generic monitoring tool for athena components.
const T * get(const ReadCondHandleKey< T > &key, const EventContext &ctx)
Convenience function to retrieve an object given a ReadCondHandleKey.
DataVector< const Trk::TrackStateOnSurface > TrackStates
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ driftRadius
trt, straws
Definition ParamDefs.h:53
@ locX
Definition ParamDefs.h:37
@ locR
Definition ParamDefs.h:44
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67
ParametersBase< TrackParametersDim, Charged > TrackParameters
Definition index.py:1
STL namespace.
EventInfo_v1 EventInfo
Definition of the latest event info version.
TrigDecision_v1 TrigDecision
Define the latest version of the trigger decision class.
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".
@ numberOfTRTHits
number of TRT hits [unit8_t].
void fill(H5::Group &out_file, size_t iterations)