ATLAS Offline Software
Loading...
Searching...
No Matches
AFPSiLayerAlgorithm.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*
4*
5* AFPSiLayerAlgorithm
6*
7*
8*/
9
13
14
15namespace {
16 constexpr int reorganizePlanes(const int station, const int layer)
17 {
18 bool reverse = station == 0 || station == 1;
19 return station * 4 + (reverse ? 3 - layer : layer);
20 }
21}
22
23
24AFPSiLayerAlgorithm::AFPSiLayerAlgorithm( const std::string& name, ISvcLocator* pSvcLocator )
25:AthMonitorAlgorithm(name,pSvcLocator)
26, m_afpHitContainerKey("AFPSiHitContainer")
27{
28 declareProperty("AFPSiHitContainer", m_afpHitContainerKey);
29}
30
31
33
34
36
37 using namespace Monitored;
38
41
42
43 // We must declare to the framework in initialize what SG objects we are going to use:
44 SG::ReadHandleKey<xAOD::AFPSiHitContainer> afpHitContainerKey("AFPSiHits");
45 ATH_CHECK(m_afpHitContainerKey.initialize());
46
47 ATH_MSG_INFO( "BunchCrossingKey initialization (SiT)" );
48 ATH_CHECK(m_bunchCrossingKey.initialize());
49 ATH_MSG_INFO( "initialization completed (SiT)" );
51}
52
53StatusCode AFPSiLayerAlgorithm::fillHistograms( const EventContext& ctx ) const {
54 using namespace Monitored;
55
56 auto bcidAll = Monitored::Scalar<int>("bcidAll", 0);
57 auto bcidFront = Monitored::Scalar<int>("bcidFront", 0);
58 auto bcidMiddle = Monitored::Scalar<int>("bcidMiddle", 0);
59 auto bcidEnd = Monitored::Scalar<int>("bcidEnd", 0);
60
61 auto numberOfEventsPerLumiblockFront = Monitored::Scalar<int>("numberOfEventsPerLumiblockFront", 0);
62 auto numberOfEventsPerLumiblockMiddle = Monitored::Scalar<int>("numberOfEventsPerLumiblockMiddle", 0);
63 auto numberOfEventsPerLumiblockEnd = Monitored::Scalar<int>("numberOfEventsPerLumiblockEnd", 0);
64
66 numberOfEventsPerLumiblockFront = eventInfo->lumiBlock();
67 numberOfEventsPerLumiblockMiddle = eventInfo->lumiBlock();
68 numberOfEventsPerLumiblockEnd = eventInfo->lumiBlock();
69
70 // BCX handler
71 const unsigned int bcid = eventInfo->bcid();
73 if (!bcidHdl.isValid()) {
74 ATH_MSG_ERROR( "Unable to retrieve BunchCrossing conditions object (SiT)" );
75 }
76 const BunchCrossingCondData* bcData{*bcidHdl};
77
78 // Classifying bunches by position in train (Front, Middle, End)
79 enum { FRONT, MIDDLE, END, NPOS } position = NPOS;
80 if(bcData->isFilled(bcid))
81 {
82 bcidAll = bcid;
83 fill("AFPSiLayerTool", bcidAll);
84 if(!bcData->isFilled(bcid-1))
85 {
86 position = FRONT;
87 bcidFront = bcid;
88 fill("AFPSiLayerTool", bcidFront);
89 fill("AFPSiLayerTool", numberOfEventsPerLumiblockFront);
90 }
91 else
92 {
93 if(bcData->isFilled(bcid+1))
94 {
95 position = MIDDLE;
96 bcidMiddle = bcid;
97 fill("AFPSiLayerTool", bcidMiddle);
98 fill("AFPSiLayerTool", numberOfEventsPerLumiblockMiddle);
99 }
100 else
101 {
102 position = END;
103 bcidEnd = bcid;
104 fill("AFPSiLayerTool", bcidEnd);
105 fill("AFPSiLayerTool", numberOfEventsPerLumiblockEnd);
106 }
107 }
108 }
109
110
111 // Declare the quantities which should be monitored:
112 auto lb = Monitored::Scalar<int>("lb", 0);
113 auto muPerBX = Monitored::Scalar<float>("muPerBX", 0.0);
114 //auto run = Monitored::Scalar<int>("run",0);
115
116 auto nSiHits = Monitored::Scalar<int>("nSiHits", 1);
117
118 auto pixelRowIDChip = Monitored::Scalar<int>("pixelRowIDChip", 0);
119 auto pixelColIDChip = Monitored::Scalar<int>("pixelColIDChip", 0);
120
121 auto timeOverThreshold = Monitored::Scalar<int>("timeOverThreshold", 0);
122
123 auto clusterX = Monitored::Scalar<float>("clusterX", 0.0);
124 auto clusterY = Monitored::Scalar<float>("clusterY", 0.0);
125 auto clustersInPlanes = Monitored::Scalar<int>("clustersInPlanes", 0);
126
127 auto trackX = Monitored::Scalar<float>("trackX", 0.0);
128 auto trackY = Monitored::Scalar<float>("trackY", 0.0);
129
130 auto planeHits = Monitored::Scalar<int>("planeHits", 0);
131
132 auto numberOfHitsPerStation = Monitored::Scalar<int>("numberOfHitsPerStation", 0);
133
134 auto lbEvents = Monitored::Scalar<int>("lbEvents", 0);
135 auto lbHits = Monitored::Scalar<int>("lbHits", 0);
136 auto lbEventsStations = Monitored::Scalar<int>("lbEventsStations", 0);
137 auto lbEventsStationsAll = Monitored::Scalar<int>("lbEventsStationsAll", 0);
138
139 auto planes = Monitored::Scalar<int>("planes", 0);
140
141 auto eventsPerStation = Monitored::Scalar<int>("eventsPerStation", 0);
142
143 auto clusterToT = Monitored::Scalar<int>("clusterToT", 0);
144
145 lb = eventInfo->lumiBlock();
146 lbEvents = eventInfo->lumiBlock();
147 //muPerBX = lbAverageInteractionsPerCrossing(ctx);
148 muPerBX = lbInteractionsPerCrossing(ctx);
149 if (muPerBX == 0.0) {
150 ATH_MSG_DEBUG("AverageInteractionsPerCrossing is 0, forcing to 1.0");
151 muPerBX=1.0;
152 }
153 fill("AFPSiLayerTool", lb, muPerBX);
154 fill("AFPSiLayerTool", lbEvents);
155
156
158 if(! afpHitContainer.isValid())
159 {
160 ATH_MSG_WARNING("evtStore() does not contain hits collection with name " << m_afpHitContainerKey);
161 return StatusCode::SUCCESS;
162 }
163
164 ATH_CHECK( afpHitContainer.initialize() );
165
166 nSiHits = afpHitContainer->size();
167 fill("AFPSiLayerTool", lb, nSiHits);
168
169 int eventsInStations[4] = {};
170 int numberOfHitsPerPlane[4][4] = {};
171
172 for(const xAOD::AFPSiHit *hitsItr: *afpHitContainer)
173 {
174 lb = eventInfo->lumiBlock();
175 lbHits = eventInfo->lumiBlock();
176 lbEventsStations = eventInfo->lumiBlock();
177 lbEventsStationsAll = eventInfo->lumiBlock();
178 pixelRowIDChip = hitsItr->pixelRowIDChip();
179 pixelColIDChip = hitsItr->pixelColIDChip();
180 timeOverThreshold = hitsItr->timeOverThreshold();
181
182
183 if (hitsItr->stationID()<4 && hitsItr->stationID()>=0 && hitsItr->pixelLayerID()<4 && hitsItr->pixelLayerID()>=0)
184 {
185 ++eventsInStations[hitsItr->stationID()];
186
187 fill(m_tools[m_StationPlaneGroup.at(m_stationnames.at(hitsItr->stationID())).at(m_pixlayers.at(hitsItr->pixelLayerID()))], pixelRowIDChip, pixelColIDChip);
188 fill(m_tools[m_StationPlaneGroup.at(m_stationnames.at(hitsItr->stationID())).at(m_pixlayers.at(hitsItr->pixelLayerID()))], pixelRowIDChip);
189 fill(m_tools[m_StationPlaneGroup.at(m_stationnames.at(hitsItr->stationID())).at(m_pixlayers.at(hitsItr->pixelLayerID()))], pixelColIDChip);
190 fill(m_tools[m_StationPlaneGroup.at(m_stationnames.at(hitsItr->stationID())).at(m_pixlayers.at(hitsItr->pixelLayerID()))], timeOverThreshold);
191
192 planeHits = hitsItr->pixelLayerID();
193 fill(m_tools[m_StationGroup.at(m_stationnames.at(hitsItr->stationID()))], planeHits);
194
195 ++numberOfHitsPerPlane[hitsItr->stationID()][hitsItr->pixelLayerID()];
196 numberOfHitsPerStation = hitsItr->stationID();
197 fill("AFPSiLayerTool", numberOfHitsPerStation);
198
199 fill("AFPSiLayerTool", lbHits);
200 }
201 else ATH_MSG_WARNING("Unrecognised station index: " << hitsItr->stationID());
202 }
203
204 auto hitsPerPlaneProfile = Monitored::Scalar<float>("hitsPerPlaneProfile", 0.0);
205 auto lbhitsPerPlaneProfile = Monitored::Scalar<int>("lbhitsPerPlaneProfile", 0);
206 auto hitsPerPlaneEventsMu = Monitored::Scalar<float>("hitsPerPlaneEventsMu", 0.0);
207 auto hitPerPlaneEventMuIndex = Monitored::Scalar<int>("hitPerPlaneEventMuIndex", 0);
208
209 lbhitsPerPlaneProfile = eventInfo->lumiBlock();
210 for(int i_station = 0; i_station < 4; i_station++)
211 for(int j_layer = 0; j_layer < 4; j_layer++)
212 {
213 hitsPerPlaneProfile = numberOfHitsPerPlane[i_station][j_layer]/muPerBX;
214 fill(m_tools[m_StationPlaneGroup.at(m_stationnames.at(i_station)).at(m_pixlayers.at(j_layer))], lbhitsPerPlaneProfile, hitsPerPlaneProfile);
215 if (muPerBX != 0.0) {
216 hitsPerPlaneEventsMu = numberOfHitsPerPlane[i_station][j_layer] / muPerBX;
217 }
218 hitPerPlaneEventMuIndex = reorganizePlanes(i_station, j_layer);
219 fill("AFPSiLayerTool", hitPerPlaneEventMuIndex, hitsPerPlaneEventsMu);
220 }
221
222 bool noEventsInStations = true;
223 for(int i=0; i<4; i++)
224 {
225 if(eventsInStations[i]>0) {
226 fill(m_tools[m_StationGroup.at(m_stationnames.at(i))], lbEventsStations);
227
228 eventsPerStation = i * 4;
229 fill("AFPSiLayerTool", eventsPerStation);
230 ++eventsPerStation;
231 fill("AFPSiLayerTool", eventsPerStation);
232 ++eventsPerStation;
233 fill("AFPSiLayerTool", eventsPerStation);
234 ++eventsPerStation;
235 fill("AFPSiLayerTool", eventsPerStation);
236
237 noEventsInStations = false;
238 }
239 }
240 if(!noEventsInStations)
241 {
242 fill("AFPSiLayerTool", lbEventsStationsAll);
243 }
244
245 // Filling of cluster and track 2D histograms
246 AFPMon::AFPFastReco fast(afpHitContainer.get());
247 fast.reco();
248
249 // Track histograms:
250 unsigned int totalTracksAll[4] = {};
251 unsigned int totalTracksFront[4] = {};
252 unsigned int totalTracksMiddle[4] = {};
253 unsigned int totalTracksEnd[4] = {};
254
255 for (const auto& track : fast.tracks())
256 {
257 trackX = track.x * 1.0;
258 trackY = track.y * 1.0;
259 fill(m_tools[m_StationGroup.at(m_stationnames.at(track.station))], trackY, trackX);
260
261 if (position == FRONT)
262 {
263 ++totalTracksFront[track.station];
264 ++totalTracksAll[track.station];
265 }
266 else if (position == MIDDLE)
267 {
268 ++totalTracksMiddle[track.station];
269 ++totalTracksAll[track.station];
270 }
271 else if (position == END)
272 {
273 ++totalTracksEnd[track.station];
274 ++totalTracksAll[track.station];
275 }
276 }
277
278 auto lbTracksAll = Monitored::Scalar<int>("lbTracksAll", 0);
279 auto lbTracksFront = Monitored::Scalar<int>("lbTracksFront", 0);
280 auto lbTracksMiddle = Monitored::Scalar<int>("lbTracksMiddle", 0);
281 auto lbTracksEnd = Monitored::Scalar<int>("lbTracksEnd", 0);
282
283 auto Total_tracks_All_profile = Monitored::Scalar<float>("Total_tracks_All_profile", 0.0);
284 auto Total_tracks_Front_profile = Monitored::Scalar<float>("Total_tracks_Front_profile", 0.0);
285 auto Total_tracks_Middle_profile = Monitored::Scalar<float>("Total_tracks_Middle_profile", 0.0);
286 auto Total_tracks_End_profile = Monitored::Scalar<float>("Total_tracks_End_profile", 0.0);
287
288 lbTracksAll = eventInfo->lumiBlock();
289 lbTracksFront = eventInfo->lumiBlock();
290 lbTracksMiddle = eventInfo->lumiBlock();
291 lbTracksEnd = eventInfo->lumiBlock();
292
293 for(int i = 0; i < 4; i++)
294 {
295 Total_tracks_All_profile = totalTracksAll[i] / muPerBX;
296 fill(m_tools[m_StationGroup.at(m_stationnames.at(i))], lbTracksAll, Total_tracks_All_profile);
297 totalTracksAll[i] = 0;
298
299 Total_tracks_Front_profile = totalTracksFront[i] / muPerBX;
300 if (position == FRONT)
301 fill(m_tools[m_StationGroup.at(m_stationnames.at(i))], lbTracksFront, Total_tracks_Front_profile);
302 totalTracksFront[i] = 0;
303
304 Total_tracks_Middle_profile = totalTracksMiddle[i] / muPerBX;
305 if (position == MIDDLE)
306 fill(m_tools[m_StationGroup.at(m_stationnames.at(i))], lbTracksMiddle, Total_tracks_Middle_profile);
307 totalTracksMiddle[i] = 0;
308
309 Total_tracks_End_profile = totalTracksEnd[i] / muPerBX;
310 if (position == END)
311 fill(m_tools[m_StationGroup.at(m_stationnames.at(i))], lbTracksEnd, Total_tracks_End_profile);
312 totalTracksEnd[i] = 0;
313 }
314
315 // Cluster histograms
316 unsigned int totalClustersAll[4][4] = {};
317 unsigned int totalClustersFront[4][4] = {};
318 unsigned int totalClustersMiddle[4][4] = {};
319 unsigned int totalClustersEnd[4][4] = {};
320
321 auto clustersPerPlaneAllPP = Monitored::Scalar<float>("clustersPerPlaneAllPP", 0.0);
322 auto clustersPerPlaneFrontPP = Monitored::Scalar<float>("clustersPerPlaneFrontPP", 0.0);
323 auto clustersPerPlaneMiddlePP = Monitored::Scalar<float>("clustersPerPlaneMiddlePP", 0.0);
324 auto clustersPerPlaneEndPP = Monitored::Scalar<float>("clustersPerPlaneEndPP", 0.0);
325
326 auto lbClustersPerPlanesAll = Monitored::Scalar<int>("lbClustersPerPlanesAll", 0);
327 auto lbClustersPerPlanesFront = Monitored::Scalar<int>("lbClustersPerPlanesFront", 0);
328 auto lbClustersPerPlanesMiddle = Monitored::Scalar<int>("lbClustersPerPlanesMiddle", 0);
329 auto lbClustersPerPlanesEnd = Monitored::Scalar<int>("lbClustersPerPlanesEnd", 0);
330
331 lbClustersPerPlanesAll = eventInfo->lumiBlock();
332 lbClustersPerPlanesFront = eventInfo->lumiBlock();
333 lbClustersPerPlanesMiddle = eventInfo->lumiBlock();
334 lbClustersPerPlanesEnd = eventInfo->lumiBlock();
335
336 for(const auto& cluster : fast.clusters())
337 {
338 clusterX = cluster.x * 1.0;
339 clusterY = cluster.y * 1.0;
340 fill(m_tools[m_StationPlaneGroup.at(m_stationnames.at(cluster.station)).at(m_pixlayers.at(cluster.layer))], clusterY, clusterX);
341 if (cluster.station == 0 || cluster.station == 1)
342 {
343 clustersInPlanes = reorganizePlanes(cluster.station, cluster.layer);
344 }
345 else
346 {
347 clustersInPlanes = (cluster.station*4)+cluster.layer;
348 }
349 fill("AFPSiLayerTool", clustersInPlanes);
350
351 clusterToT = cluster.sumToT;
352 fill(m_tools[m_StationPlaneGroup.at(m_stationnames.at(cluster.station)).at(m_pixlayers.at(cluster.layer))], clusterToT);
353
354 if (position == FRONT)
355 {
356 ++totalClustersFront[cluster.station][cluster.layer];
357 ++totalClustersAll[cluster.station][cluster.layer];
358 }
359 else if (position == MIDDLE)
360 {
361 ++totalClustersMiddle[cluster.station][cluster.layer];
362 ++totalClustersAll[cluster.station][cluster.layer];
363 }
364 else if (position == END)
365 {
366 ++totalClustersEnd[cluster.station][cluster.layer];
367 ++totalClustersAll[cluster.station][cluster.layer];
368 }
369 }
370
371 for(int i_station = 0; i_station < 4; i_station++)
372 for(int j_layer = 0; j_layer < 4; j_layer++)
373 {
374 clustersPerPlaneAllPP = totalClustersAll[i_station][j_layer] / muPerBX;
375 fill(m_tools[m_StationPlaneGroup.at(m_stationnames.at(i_station)).at(m_pixlayers.at(j_layer))], lbClustersPerPlanesAll, clustersPerPlaneAllPP);
376 totalClustersAll[i_station][j_layer] = 0;
377
378 clustersPerPlaneFrontPP = totalClustersFront[i_station][j_layer] / muPerBX;
379 if (position == FRONT)
380 fill(m_tools[m_StationPlaneGroup.at(m_stationnames.at(i_station)).at(m_pixlayers.at(j_layer))], lbClustersPerPlanesFront, clustersPerPlaneFrontPP);
381 totalClustersFront[i_station][j_layer] = 0;
382
383 clustersPerPlaneMiddlePP = totalClustersMiddle[i_station][j_layer] / muPerBX;
384 if (position == MIDDLE)
385 fill(m_tools[m_StationPlaneGroup.at(m_stationnames.at(i_station)).at(m_pixlayers.at(j_layer))], lbClustersPerPlanesMiddle, clustersPerPlaneMiddlePP);
386 totalClustersMiddle[i_station][j_layer] = 0;
387
388 clustersPerPlaneEndPP = totalClustersEnd[i_station][j_layer] / muPerBX;
389 if (position == END)
390 fill(m_tools[m_StationPlaneGroup.at(m_stationnames.at(i_station)).at(m_pixlayers.at(j_layer))], lbClustersPerPlanesEnd, clustersPerPlaneEndPP);
391 totalClustersEnd[i_station][j_layer] = 0;
392 }
393
394 return(fillHistogramsPlaneEff(*afpHitContainer));
395} // end of fillHistograms
396
398 using namespace Monitored;
399
400 // Define 2D histograms for tries and successes
401 Monitored::Scalar<bool> sit_plane_eff_passed("sit_plane_eff_passed", false);
402 Monitored::Scalar<float> sit_plane_eff_triedX("sit_plane_eff_triedX", 0.0);
403 Monitored::Scalar<float> sit_plane_eff_triedY("sit_plane_eff_triedY", 0.0);
404
405 auto triesX = Monitored::Scalar<float>("triesX", 0.0);
406 auto triesY = Monitored::Scalar<float>("triesY", 0.0);
407
408 auto successX = Monitored::Scalar<float>("successX", 0.0);
409 auto successY = Monitored::Scalar<float>("successY", 0.0);
410
411
412 auto clusterXLocal = Monitored::Scalar<float>("clusterXLocal", 0.0);
413 auto clusterYLocal = Monitored::Scalar<float>("clusterYLocal", 0.0);
414
415 auto clusterXTag = Monitored::Scalar<float>("clusterXTag", 0.0);
416 auto clusterYTag = Monitored::Scalar<float>("clusterYTag", 0.0);
417
418 AFPMon::AFPFastReco fast(&afpHitContainer);
419 fast.reco();
420
421 int min_hits[4] = {2, 3, 3, 3};
422 int numStations = 4;
423 int numPlanes = 4;
424
425 // Count clusters per station and plane
426 std::vector<std::vector<int>> nclusters_planes(numStations, std::vector<int>(numPlanes, 0));
427 for (const auto& cluster : fast.clusters()) {
428 ++nclusters_planes[cluster.station][cluster.layer];
429 }
430
431 // Count planes with clusters at each station
432 std::vector<std::vector<int>> nclusters_other_planes(numStations, std::vector<int>(numPlanes, 0));
433 for (int iStation = 0; iStation < numStations; ++iStation) {
434 for (int iPlane = 0; iPlane < numPlanes; ++iPlane) {
435 for (int ip = 0; ip < numPlanes; ++ip) {
436 if (iPlane != ip && nclusters_planes[iStation][ip] > 0) {
437 ++nclusters_other_planes[iStation][iPlane];
438 }
439 }
440 }
441 }
442
443 // Precomputed tag planes
444 std::array<std::set<int>, 4> precomputed_tag_planes = {
445 std::set<int>{1, 2, 3},
446 std::set<int>{0, 1, 2, 3},
447 std::set<int>{0, 1, 2, 3},
448 std::set<int>{0, 1, 2, 3}
449 };
450
451 // Precomputed clusters by Stations and Planes
452 using ClusterType = std::decay_t<decltype(*fast.clusters().begin())>;
453 std::vector<std::vector<std::vector<ClusterType>>> clusters_by_station_layer(numStations,
454 std::vector<std::vector<ClusterType>>(numPlanes));
455
456 for (const auto& cluster : fast.clusters()) {
457 clusters_by_station_layer[cluster.station][cluster.layer].push_back(cluster);
458 }
459
460 for (int iStation = 0; iStation < numStations; ++iStation) {
461 for (int iPlane = 0; iPlane < numPlanes; ++iPlane) {
462 // Skip if no enough hits in other planes
463 if (nclusters_other_planes[iStation][iPlane] < min_hits[iStation]) {continue;}
464
465 std::set<int> tag_planes = precomputed_tag_planes[iStation];
466 tag_planes.erase(iPlane);
467
468 std::vector<int> v_tag_planes(tag_planes.begin(), tag_planes.end());
469 int seed = v_tag_planes[0]; // Take first plane as a seed
470
471 // Iterate only stations coresponding to the given Station and Layer
472 for (const auto& cluster : clusters_by_station_layer[iStation][seed]) {
473 int found_in_other_tag_planes = 0;
474
475 // Check other tag planes using lookup
476 for (size_t i = 1; i < v_tag_planes.size(); ++i) {
477 int tag_plane = v_tag_planes[i];
478 bool found = false;
479
480 for (const auto& clusterTag : clusters_by_station_layer[iStation][tag_plane]) {
481 if (std::fabs(cluster.x - clusterTag.x) < 2.0 && std::fabs(cluster.y - clusterTag.y) < 2.0) {
482 found = true;
483 ++found_in_other_tag_planes;
484 break;
485 }
486 }
487 // if not found in current tag plane
488 if (!found) {break;}
489 }
490
491 // Check if valid tag
492 if (found_in_other_tag_planes == static_cast<int>(v_tag_planes.size()) - 1) {
493 sit_plane_eff_triedX = cluster.x;
494 sit_plane_eff_triedY = cluster.y;
495
496 // Check probe plane
497 for (const auto& clusterTag : clusters_by_station_layer[iStation][iPlane]) {
498 if (std::fabs(cluster.x - clusterTag.x) < 2.0 && std::fabs(cluster.y - clusterTag.y) < 2.0) {
499 sit_plane_eff_passed = true;
500 break;
501 }
502 }
503
504 fill(m_tools[m_StationPlaneGroup.at(m_stationnames.at(iStation)).at(m_pixlayers.at(iPlane))],
505 sit_plane_eff_passed, sit_plane_eff_triedY, sit_plane_eff_triedX);
506 sit_plane_eff_passed = false;
507 }
508 }
509 }
510 }
511 return StatusCode::SUCCESS;
512}
513
514
Definitions of AFP stations identification numbers.
#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)
Property holding a SG store/key/clid from which a ReadHandle is made.
SG::ReadCondHandleKey< BunchCrossingCondData > m_bunchCrossingKey
std::map< std::string, int > m_StationGroup
virtual StatusCode fillHistogramsPlaneEff(const xAOD::AFPSiHitContainer &) const
SG::ReadHandleKey< xAOD::AFPSiHitContainer > m_afpHitContainerKey
std::vector< std::string > m_stationnames
virtual StatusCode fillHistograms(const EventContext &ctx) const override
adds event to the monitoring histograms
std::vector< std::string > m_pixlayers
virtual StatusCode initialize() override
initialize
AFPSiLayerAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
std::map< std::string, std::map< std::string, int > > m_StationPlaneGroup
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
virtual StatusCode initialize() override
initialize
SG::ReadHandle< xAOD::EventInfo > GetEventInfo(const EventContext &) const
Return a ReadHandle for an EventInfo object (get run/event numbers, etc.)
AthMonitorAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor.
ToolHandleArray< GenericMonitoringTool > m_tools
Array of Generic Monitoring Tools.
bool isFilled(const bcid_type bcid) const
The simplest query: Is the bunch crossing filled or not?
Declare a monitored scalar variable.
Property holding a SG store/key/clid from which a ReadHandle is made.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type get() const
Dereference the pointer, but don't cache anything.
StatusCode initialize(bool used=true)
Verify that the handle has been configured properly.
double timeOverThreshold(unsigned int m_word)
int lb
Definition globals.cxx:23
void fill(const ToolHandle< GenericMonitoringTool > &groupHandle, std::vector< std::reference_wrapper< Monitored::IMonitoredVariable > > &&variables) const
Fills a vector of variables to a group by reference.
virtual float lbInteractionsPerCrossing(const EventContext &ctx=Gaudi::Hive::currentContext()) const
Calculate instantaneous number of interactions, i.e.
Generic monitoring tool for athena components.
std::vector< V > buildToolMap(const ToolHandleArray< GenericMonitoringTool > &tools, const std::string &baseName, int nHist)
Builds an array of indices (base case)
@ layer
Definition HitInfo.h:79
AFPSiHit_v2 AFPSiHit
Definition AFPSiHit.h:12
AFPSiHitContainer_v2 AFPSiHitContainer