ATLAS Offline Software
Loading...
Searching...
No Matches
PixelAthMVAMonAlg.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
9#include "PixelAthMVAMonAlg.h"
16#include "TTree.h"
17#include "TFile.h"
18
19PixelAthMVAMonAlg::PixelAthMVAMonAlg( const std::string& name, ISvcLocator* pSvcLocator ) :
20 AthMonitorAlgorithm(name, pSvcLocator),
21 m_holeSearchTool("InDet::InDetTrackHoleSearchTool/InDetHoleSearchTool", this),
22 m_trackSelTool("InDet::InDetTrackSelectionTool/TrackSelectionTool", this),
23 m_trkextrapolator("Trk::Extrapolator/InDetExtrapolator", this),
24 m_atlasid(nullptr)
25{
26 declareProperty("HoleSearchTool", m_holeSearchTool);
27 declareProperty("TrackSelectionTool", m_trackSelTool);
28 declareProperty("Extrapolator", m_trkextrapolator);
29 declareProperty("calibFolder", m_calibFolder = "20220503");
30 declareProperty("dumpTree", m_dumpTree = false);
31}
32
34
37 ATH_CHECK( detStore()->retrieve(m_atlasid, "AtlasID") );
38 ATH_CHECK( m_pixelRDOName.initialize() );
39 if ( !m_holeSearchTool.empty() ) ATH_CHECK( m_holeSearchTool.retrieve() );
40 if ( !m_trackSelTool.empty() ) ATH_CHECK( m_trackSelTool.retrieve() );
41 if ( !m_trkextrapolator.empty() ) ATH_CHECK( m_trkextrapolator.retrieve() );
42 ATH_CHECK( m_trackParticlesKey.initialize() );
43 ATH_CHECK( m_clustersKey.initialize() );
44
45 for (int ii = 0; ii < PixLayers::NBASELAYERS; ++ii) {
46
47 std::string partitionLabel("Disks");
48 if (ii>=PixLayers::kBLayer) partitionLabel = pixBaseLayersLabel[ii];
49
50 // std::string fullPathToTrainingFile = partitionLabel + "_training.root"; //TEST
51 std::string fullPathToTrainingFile = PathResolverFindCalibFile("PixelDQMonitoring/" + m_calibFolder + "/" + partitionLabel + "_training.root");
52
53 std::unique_ptr<TFile> calibFile(TFile::Open(fullPathToTrainingFile.c_str(), "READ"));
54 if ( !calibFile ) {
55 ATH_MSG_ERROR("Can not retrieve PixelMonitoring calibration file: " << fullPathToTrainingFile);
56 return StatusCode::FAILURE;
57 }
58 std::unique_ptr<TTree> trainingWeights( (TTree*)calibFile->Get("lgbm") );
59 m_classBDT.emplace( std::make_pair(partitionLabel, std::make_unique<MVAUtils::BDT>( trainingWeights.get())) );
60
64 }
65 return StatusCode::SUCCESS;
66}
67
68
69StatusCode PixelAthMVAMonAlg::fillHistograms( const EventContext& ctx ) const {
70
71 ATH_MSG_DEBUG("Pixel MVAMon: Filling MVA Monitoring Histograms");
72
73 using namespace Monitored;
74 int lb = GetEventInfo(ctx)->lumiBlock();
75 auto lbval = Monitored::Scalar<int>( "pixmvamontool_lb", lb );
76
77 const int MAXHASH = m_pixelid->wafer_hash_max(); // 2048
78
79 //
80 // fill status vector from conditions summary
81 //
84 std::vector<float> status(MAXHASH*2);
85
86 PixelID::const_id_iterator idIt = m_pixelid->wafer_begin();
87 PixelID::const_id_iterator idItEnd = m_pixelid->wafer_end();
88 for (; idIt != idItEnd; ++idIt) {
89
90 Identifier waferID = *idIt;
91 IdentifierHash modHash = m_pixelid->wafer_hash(waferID);
92 int pixlayer = getPixLayersID(m_pixelid->barrel_ec(waferID), m_pixelid->layer_disk(waferID) );
93 if (pixlayer == 99) continue;
94
95 if (pixlayer == PixLayers::kIBL)
96 {
97 // per FE for IBL
98 //
99 int nFE = getNumberOfFEs(pixlayer, m_pixelid->eta_module(waferID));
100 for (int iFE=0; iFE<nFE; iFE++) {
101 Identifier pixID = m_pixelReadout->getPixelIdfromHash(modHash, iFE, 1, 1);
102 if (pixID.is_valid()) {
103 auto [is_active,is_good] = isChipGood( !m_pixelDetElStatusActiveOnly.empty() ? pixel_active.cptr() : nullptr,
104 !m_pixelDetElStatus.empty() ? pixel_status.cptr() : nullptr,
105 modHash, iFE);
106 if (is_active)
107 {
108 if (is_good) status[modHash+MAXHASH*iFE] = 0;
109 else status[modHash+MAXHASH*iFE] = 1;
110 }
111 else status[modHash+MAXHASH*iFE] = 2;
112 }
113 }
114 }
115 else
116 {
117 // per module for the old pixel layers
118 //
119 bool is_active = isActive( !m_pixelDetElStatusActiveOnly.empty() ? pixel_active.cptr() : nullptr, modHash);
120 bool is_good = isGood( !m_pixelDetElStatus.empty() ? pixel_status.cptr() : nullptr, modHash);
121 if (is_active)
122 {
123 if (is_good) status[modHash] = 0;
124 else status[modHash] = 1;
125 }
126 else status[modHash] = 2;
127 }
128 }
129 //
130 // input data from tracks
131 //
132
133 auto trackParticles = SG::makeHandle(m_trackParticlesKey, ctx);
134 if ( !(trackParticles.isValid()) ) {
135 ATH_MSG_ERROR("Pixel MVAMon: TrackParticle container "<< m_trackParticlesKey.key() << " could not be found.");
136 return StatusCode::RECOVERABLE;
137 } else {
138 ATH_MSG_DEBUG("Pixel MVAMon: TrackParticle container "<< trackParticles.name() <<" is found.");
139 }
140
141 bool eventHasPixHits(false);
142 std::vector<std::pair<Identifier, double> > ClusterIDs;
143
144 std::vector<float> holes(MAXHASH*2);
145 std::vector<float> outliers(MAXHASH*2);
146 std::vector<float> measurements(MAXHASH*2);
147 std::vector<float> trkalpha(MAXHASH*2);
148 std::vector<float> trkchi2byndf(MAXHASH*2);
149 std::vector<float> trknpixdead(MAXHASH*2);
150 std::vector<float> trknblayerhits(MAXHASH*2);
151
152 uint8_t iSummaryValue(0); // Dummy counter to retrieve summary values
153
154 for (auto trackPart: *trackParticles) {
155 const Trk::Track * track = trackPart->track();
156 int numberOfPixelHits = trackPart->summaryValue(iSummaryValue, xAOD::numberOfPixelHits) ? iSummaryValue : 0;
157
158 if (track == nullptr || track->perigeeParameters() == nullptr || track->trackSummary() == nullptr || numberOfPixelHits == 0) {
159 ATH_MSG_DEBUG("Pixel MVAMon: Track either invalid or it does not contain pixel hits, continuing...");
160 continue;
161 }
162
163 int nPixelHits = 0;
164 bool hasIBLTSOS(false);
165
166 bool passJOTrkCut = static_cast<bool>( m_trackSelTool->accept(*track) );
167 if (!passJOTrkCut) continue;
168 float trkfitndf = track->fitQuality()->numberDoF();
169 float trkfitchi2byndf = track->fitQuality()->chiSquared();
170 if (trkfitndf != 0) trkfitchi2byndf /= trkfitndf;
171 else continue;
172 float npixdead = trackPart->summaryValue(iSummaryValue, xAOD::numberOfPixelDeadSensors) ? iSummaryValue : 0;
173 float nblayerhits(0);
174 int expBLhit = trackPart->summaryValue(iSummaryValue, xAOD::expectNextToInnermostPixelLayerHit) ? iSummaryValue : 0;
175 if (expBLhit==1) {
176 nblayerhits = trackPart->summaryValue(iSummaryValue, xAOD::numberOfNextToInnermostPixelLayerHits) ? iSummaryValue : 0;
177 }
178
179 const Trk::Track* trackWithHoles( track );
180 std::unique_ptr<const Trk::Track> trackWithHolesUnique = nullptr;
181 int nPixHoles = trackPart->summaryValue(iSummaryValue, xAOD::numberOfPixelHoles) ? iSummaryValue : 0;
182 if ( nPixHoles > 0 ) {
183 trackWithHolesUnique.reset( m_holeSearchTool->getTrackWithHoles(*track) );
184 trackWithHoles = trackWithHolesUnique.get();
185 }
186 const Trk::TrackStates *trackStates = trackWithHoles->trackStateOnSurfaces();
187 for (auto trackStateOnSurface: *trackStates) {
188
189 const Trk::MeasurementBase* mesBase = trackStateOnSurface->measurementOnTrack();
190
191 // skip pseudomeasurements (only hits, holes, outliers are considered)
192 //
193 const Trk::RIO_OnTrack* RIOOnTrack = nullptr;
194 if ( mesBase && mesBase->type(Trk::MeasurementBaseType::RIO_OnTrack) ) {
195 RIOOnTrack = static_cast<const Trk::RIO_OnTrack*>(mesBase);
196 }
197 if (mesBase && !RIOOnTrack) continue;
198
199 // obtaining surfaceID
200 //
201 Identifier surfaceID;
202 const Trk::TrackParameters* trkParameters = trackStateOnSurface->trackParameters();
203 if (mesBase && mesBase->associatedSurface().associatedDetectorElement()) {
204 surfaceID = mesBase->associatedSurface().associatedDetectorElement()->identify();
205 } else { // holes, perigee
206 if (trkParameters) {
207 surfaceID = trkParameters->associatedSurface().associatedDetectorElementIdentifier();
208 } else {
209 ATH_MSG_DEBUG("Pixel MVAMon: can't obtain track parameters for TSOS.");
210 continue;
211 }
212 }
213
214 // selecting only pixel part
215 //
216 if (!m_atlasid->is_pixel(surfaceID)) continue;
217 int pixlayer = getPixLayersID(m_pixelid->barrel_ec(surfaceID), m_pixelid->layer_disk(surfaceID) );
218 if (pixlayer == 99) continue;
219 if (pixlayer == PixLayers::kIBL) hasIBLTSOS = true;
220
221 int indexModule = static_cast<int>( m_pixelid->wafer_hash(surfaceID) ); // [0,2047]
222
223 const InDetDD::SiDetectorElement *sde = dynamic_cast<const InDetDD::SiDetectorElement *>(trkParameters->associatedSurface().associatedDetectorElement());
224 const Trk::AtaPlane *trackAtPlane = dynamic_cast<const Trk::AtaPlane *>(trkParameters);
225 const InDetDD::SiLocalPosition& trkLocalPos = trkParameters->localPosition();
226 Identifier locPosID;
227 int iFE(0);
228 bool goodPixelMeasurement(false);
229
230 if (trackStateOnSurface->type(Trk::TrackStateOnSurface::Outlier))
231 {
232 if ( isIBL2D(indexModule) )
233 {
234 const InDet::PixelClusterOnTrack *pixClus = dynamic_cast<const InDet::PixelClusterOnTrack *>(mesBase);
235 if ( mesBase && pixClus) {
236 locPosID = pixClus->identify();
237 if ( !(locPosID.is_valid()) ) {
238 ATH_MSG_WARNING("Pixel MVAMon: got invalid track local position on surface for an outlier.");
239 continue;
240 }
241 iFE = m_pixelReadout->getFE(locPosID, locPosID);
242 outliers[indexModule+MAXHASH*iFE] += 1;
243 }
244 }
245 else
246 {
247 outliers[indexModule] += 1;
248 }
249 }
250 else if (trackStateOnSurface->type(Trk::TrackStateOnSurface::Hole))
251 {
252 if ( isIBL2D(indexModule) )
253 {
254 locPosID = sde->identifierOfPosition(trkLocalPos);
255 if ( !(locPosID.is_valid()) ) {
256 ATH_MSG_WARNING("Pixel MVAMon: got invalid track local position on surface for a hole.");
257 continue;
258 }
259 iFE = m_pixelReadout->getFE(locPosID, locPosID);
260 holes[indexModule+MAXHASH*iFE] += 1;
261 }
262 else
263 {
264 holes[indexModule] += 1;
265 }
266 }
267 else if (trackStateOnSurface->type(Trk::TrackStateOnSurface::Measurement))
268 {
269 // making sure we get raw pixel clusters data
270 //
271 if (not mesBase) continue;
272 sde = dynamic_cast<const InDetDD::SiDetectorElement *>(mesBase->associatedSurface().associatedDetectorElement());
273 const InDet::PixelClusterOnTrack *pixClus = dynamic_cast<const InDet::PixelClusterOnTrack *>(mesBase);
274 if (!sde || !pixClus) continue;
275 const InDet::PixelCluster *rawPixClus = pixClus->prepRawData();
276 if (!rawPixClus) continue;
277
278 nPixelHits++;
279 goodPixelMeasurement = true;
280 locPosID = pixClus->identify();
281 if ( !(locPosID.is_valid()) ) {
282 ATH_MSG_WARNING("Pixel MVAMon: got invalid cluster on track ID.");
283 continue;
284 }
285 if ( isIBL2D(indexModule) ) iFE = m_pixelReadout->getFE(locPosID, locPosID);
286 measurements[indexModule+MAXHASH*iFE] += 1;
287 }
288 else continue;
289
290 if (trackAtPlane) {
291 Amg::Vector3D normal = sde->normal();
292 Amg::Vector3D trackp = trackAtPlane->momentum();
293 double trackpnormcomp = trackp.dot(normal);
294 double trackp_mag = trackp.mag();
295 double alpha = 0.5 * M_PI;
296 if (trackp_mag != 0) alpha = std::acos(std::abs(trackpnormcomp / trackp_mag));
297
298 trkalpha[indexModule+MAXHASH*iFE] += alpha;
299
300 if (goodPixelMeasurement) ClusterIDs.emplace_back(locPosID, alpha);
301 }
302 trkchi2byndf[indexModule+MAXHASH*iFE] += trkfitchi2byndf;
303 trknpixdead[indexModule+MAXHASH*iFE] += npixdead;
304 trknblayerhits[indexModule+MAXHASH*iFE] += nblayerhits;
305 } // end of TSOS loop
306 eventHasPixHits = eventHasPixHits || (nPixelHits > 0);
307
308 int expIBLhit = trackPart->summaryValue(iSummaryValue, xAOD::expectInnermostPixelLayerHit) ? iSummaryValue : 0;
309 if (!hasIBLTSOS && expIBLhit==1 && nPixHoles==0 )
310 {
311 const Trk::Perigee *perigee = track->perigeeParameters();
312 Amg::Transform3D transSurf;
313 transSurf.setIdentity();
314 const Trk::CylinderSurface BiggerThanIBLSurface(transSurf, 40.0, 400.0);
315 std::vector<std::unique_ptr<Trk::TrackParameters> > iblTrkParameters =
316 m_trkextrapolator->extrapolateStepwise(ctx, *perigee, BiggerThanIBLSurface, Trk::alongMomentum, false);
317 Identifier sensorPosID, chipPosID;
318 double alpha = 0.5 * M_PI;
319 bool foundIBLElem(false); //more than one IBL element possible, only the last one is considered
320
321 for (auto& tp: iblTrkParameters) {
322 const InDetDD::SiDetectorElement *sde = dynamic_cast<const InDetDD::SiDetectorElement *>(tp->associatedSurface().associatedDetectorElement());
323 if (!(sde != nullptr && sde->identify() != 0)) {
324 continue;
325 }
326 sensorPosID = sde->identify();
327 if (!m_atlasid->is_pixel(sensorPosID) || !m_pixelid->is_barrel(sensorPosID) || m_pixelid->layer_disk(sensorPosID)!= 0) {
328 ATH_MSG_DEBUG("Pixel MVAMon: found non-IBL element: "<< sensorPosID);
329 continue;
330 }
331 const InDetDD::SiLocalPosition& trkLocalPos = tp->localPosition();
332 chipPosID = sde->identifierOfPosition(trkLocalPos);
333 if ( !(chipPosID.is_valid()) ) {
334 foundIBLElem=false;
335 ATH_MSG_DEBUG("Pixel MVAMon: got invalid track local position on surface for an expected IBL hit.");
336 continue;
337 } else foundIBLElem=true;
338 double trackpnormcomp = tp->momentum().dot(tp->associatedSurface().normal());
339 double trackp_mag = tp->momentum().mag();
340 if (trackp_mag != 0) alpha = std::acos(std::abs(trackpnormcomp / trackp_mag));
341 }
342 if (!foundIBLElem) {
343 ATH_MSG_DEBUG("Pixel MVAMon: couldn't find IBL pos ID: "<< sensorPosID);
344 continue;
345 }
346 int indexModule = static_cast<int>( m_pixelid->wafer_hash(sensorPosID) );
347 int iFE = m_pixelReadout->getFE(chipPosID, chipPosID);
348 holes[indexModule+MAXHASH*iFE] += 1;
349 trkalpha[indexModule+MAXHASH*iFE] += alpha;
350 trkchi2byndf[indexModule+MAXHASH*iFE] += trkfitchi2byndf;
351 trknpixdead[indexModule+MAXHASH*iFE] += npixdead;
352 trknblayerhits[indexModule+MAXHASH*iFE] += nblayerhits;
353 }
354 } // end of track loop
355
356 if (!eventHasPixHits) {
357 ATH_MSG_DEBUG("Pixel MVAMon: event doesn't contain any pixel hits on tracks, skipping evaluation.");
358 return StatusCode::SUCCESS;
359 }
360
361
362 sort(ClusterIDs.begin(), ClusterIDs.end(), [](const std::pair<Identifier, double> &left, const std::pair<Identifier, double> &right) {
363 return left.first < right.first;
364 });
365
366 std::vector<float> misshitsf(MAXHASH*2);
367 float alltrackinfo;
368 for (int ih=0; ih<MAXHASH; ++ih) {
369 alltrackinfo = holes[ih] + outliers[ih] + measurements[ih];
370 if (alltrackinfo) {
371 misshitsf[ih] = (holes[ih]+outliers[ih])/alltrackinfo;
372 trkalpha[ih] = trkalpha[ih]/alltrackinfo;
373 trkchi2byndf[ih] = trkchi2byndf[ih]/alltrackinfo;
374 trknpixdead[ih] = trknpixdead[ih]/alltrackinfo;
375 trknblayerhits[ih] = trknblayerhits[ih]/alltrackinfo;
376 }
377 if ( isIBL2D(ih) ) {
378 alltrackinfo = holes[ih+MAXHASH]+outliers[ih+MAXHASH]+measurements[ih+MAXHASH];
379 if (alltrackinfo) {
380 misshitsf[ih+MAXHASH] = (holes[ih+MAXHASH]+outliers[ih+MAXHASH])/alltrackinfo;
381 trkalpha[ih+MAXHASH] = trkalpha[ih+MAXHASH]/alltrackinfo;
382 trkchi2byndf[ih+MAXHASH] = trkchi2byndf[ih+MAXHASH]/alltrackinfo;
383 trknpixdead[ih+MAXHASH] = trknpixdead[ih+MAXHASH]/alltrackinfo;
384 trknblayerhits[ih+MAXHASH] = trknblayerhits[ih+MAXHASH]/alltrackinfo;
385 }
386 }
387 }
388
389
390 //
391 // get CLUSTER information from container
392 //
393 std::vector<float> clsall(MAXHASH*2);
394 std::vector<float> clsontrkf(MAXHASH*2);
395 std::vector<float> clsontrk(MAXHASH*2);
396 std::vector<float> clsofftrk(MAXHASH*2);
397
398
399 std::vector<float> clsontrksize(MAXHASH*2);
400 std::vector<float> clsontrkrowsize(MAXHASH*2);
401 std::vector<float> clsontrkcolsize(MAXHASH*2);
402 std::vector<float> clsofftrksize(MAXHASH*2);
403 std::vector<float> clsofftrkrowsize(MAXHASH*2);
404 std::vector<float> clsofftrkcolsize(MAXHASH*2);
405
406 std::vector<float> clsontrktot(MAXHASH*2);
407 std::vector<float> clsofftrktot(MAXHASH*2);
408
409 auto pixel_clcontainer = SG::makeHandle(m_clustersKey, ctx);
410
411 if ( !(pixel_clcontainer.isValid()) ) {
412 ATH_MSG_ERROR("Pixel MVAMon: Pixel Cluster container "<< m_clustersKey.key() << " could not be found.");
413 return StatusCode::RECOVERABLE;
414 } else {
415 ATH_MSG_DEBUG("Pixel MVAMon: Pixel Cluster container "<< pixel_clcontainer.name() <<" is found.");
416 }
417
418 Identifier clusID;
419
420 for (auto colNext: *pixel_clcontainer) {
421
422 const InDet::PixelClusterCollection* ClusterCollection(colNext);
423 if (!ClusterCollection) {
424 ATH_MSG_DEBUG("Pixel MVAMon: Pixel Cluster container is empty.");
425 continue;
426 }
427
428 IdentifierHash clusIDHash = ClusterCollection->identifyHash();
429 int indexModule = static_cast<int>( clusIDHash );
430
431 for (auto p_clus: *ClusterCollection) {
432
433 clusID = p_clus->identify();
434
435 int pixlayer = getPixLayersID(m_pixelid->barrel_ec(clusID), m_pixelid->layer_disk(clusID) );
436 if (pixlayer == 99) continue;
437
438 int iFE=0;
439 if ( isIBL2D(indexModule) ) {
440 iFE = m_pixelReadout->getFE(clusID, clusID);
441 }
442 int idxCluster = indexModule+MAXHASH*iFE;
443
444 clsall[idxCluster] += 1;
445
446 const InDet::PixelCluster& cluster = *p_clus;
447 double alpha = 0.5 * M_PI;
448 if ( isClusterOnTrack(clusID, ClusterIDs, alpha) ) {
449 clsontrk[idxCluster] += 1;
450 clsontrksize[idxCluster] += cluster.rdoList().size();
451 clsontrkrowsize[idxCluster] += cluster.width().colRow().x();
452 clsontrkcolsize[idxCluster] += cluster.width().colRow().y();
453 clsontrktot[idxCluster] += cluster.totalToT();
454 } else {
455 clsofftrk[idxCluster] += 1;
456 clsofftrksize[idxCluster] += cluster.rdoList().size();
457 clsofftrkrowsize[idxCluster]+= cluster.width().colRow().x();
458 clsofftrkcolsize[idxCluster]+= cluster.width().colRow().y();
459 clsofftrktot[idxCluster] += cluster.totalToT();
460 }
461 }
462 }
463
464 ClusterIDs.clear();
465
466 //
467 // normalize cluster information
468 //
469 for (int ih=0; ih<MAXHASH; ++ih) {
470 if (clsontrk[ih]) {
471 clsontrksize[ih] /= clsontrk[ih];
472 clsontrkrowsize[ih]/= clsontrk[ih];
473 clsontrkcolsize[ih]/= clsontrk[ih];
474 clsontrktot[ih] /= clsontrk[ih];
475 }
476 if (clsofftrk[ih]) {
477 clsofftrksize[ih] /= clsofftrk[ih];
478 clsofftrkrowsize[ih]/= clsofftrk[ih];
479 clsofftrkcolsize[ih]/= clsofftrk[ih];
480 clsofftrktot[ih] /= clsofftrk[ih];
481 }
482 if (clsall[ih]) {
483 clsontrkf[ih] = clsontrk[ih]/clsall[ih];
484 }
485
486 if ( isIBL2D(ih) ) {
487 if (clsontrk[ih+MAXHASH]) {
488 clsontrksize[ih+MAXHASH] /= clsontrk[ih+MAXHASH];
489 clsontrkrowsize[ih+MAXHASH]/= clsontrk[ih+MAXHASH];
490 clsontrkcolsize[ih+MAXHASH]/= clsontrk[ih+MAXHASH];
491 clsontrktot[ih+MAXHASH] /= clsontrk[ih+MAXHASH];
492 }
493 if (clsofftrk[ih+MAXHASH]) {
494 clsofftrksize[ih+MAXHASH] /= clsofftrk[ih+MAXHASH];
495 clsofftrkrowsize[ih+MAXHASH]/= clsofftrk[ih+MAXHASH];
496 clsofftrkcolsize[ih+MAXHASH]/= clsofftrk[ih+MAXHASH];
497 clsofftrktot[ih+MAXHASH] /= clsofftrk[ih+MAXHASH];
498 }
499 if (clsall[ih+MAXHASH]) {
500 clsontrkf[ih+MAXHASH] = clsontrk[ih+MAXHASH]/clsall[ih+MAXHASH];
501 }
502 }
503 }
504
505 //
506 // compute BDT weights
507 //
508 std::vector<float> bdtweights(MAXHASH*2);
509 VecAccumulator2DMap BDT_Weights(*this, "BDTWeights");
510
511 for (int ih=12; ih<MAXHASH-12; ++ih)
512 {
513 int nFE(1);
514 int module(0);
515 float mod_eta(0);
516 float el_eta(0);
517
518 if ( ih>=156 && ih<=435 ) //IBL
519 {
520 module = (ih-156) % 20;
521 nFE = 2;
522 if (module<4 || module>15) continue; //3D's are out of acceptance
523 else mod_eta = module-9.5;
524 }
525 else if (ih>=436 && ih<=721 ) //BLayer
526 {
527 module = (ih-436) % 13;
528 el_eta = std::abs(module-6);
529 }
530 else if (ih>=722 && ih<=1215) //Layer1
531 {
532 module = (ih-722) % 13;
533 el_eta = std::abs(module-6);
534 }
535 else if (ih>=1216 && ih<=1891) //Layer2
536 {
537 module = (ih-1216) % 13;
538 el_eta = std::abs(module-6);
539 }
540 else if ( (ih>=12 && ih<=155) || (ih>=1892 && ih<=2035) ) //Disks
541 {
542 //int idxdisks = ih - 12;
543 //if (ih>=1892) idxdisks-=1880;
544 el_eta = 0; // idxdisks / 48;
545 }
546 else continue;
547
548 Identifier waferID = m_pixelid->wafer_id(ih);
549 int pixlayer = getPixLayersID(m_pixelid->barrel_ec(waferID), m_pixelid->layer_disk(waferID) );
550 if (pixlayer == 99) continue;
551 std::string partitionLabel("Disks");
552 if (pixlayer>=PixLayers::kBLayer) partitionLabel = pixBaseLayersLabel[pixlayer];
553
554 for (int iFE=0; iFE<nFE; iFE++) {
555 Identifier pixID = waferID;
556 if (pixlayer == PixLayers::kIBL) {
557 pixID = m_pixelReadout->getPixelIdfromHash(ih, iFE, 1, 1);
558 if (!pixID.is_valid()) {
559 ATH_MSG_ERROR("Pixel MVAMon: got invalid pixID " << pixID);
560 continue;
561 }
562 el_eta = mod_eta*2 + iFE;
563 if (mod_eta > 0) el_eta-=1.0;
564 el_eta = std::abs(el_eta);
565 }
566
567 int idx = ih+MAXHASH*iFE;
568 if ( status[idx]==2 ) BDT_Weights.add(pixlayer, pixID, 0);
569 if ( status[idx]!=0 || (measurements[idx]+holes[idx]+outliers[idx]==0) ) continue;
570
571 std::vector<float> bdtVars;
572 if (pixlayer == PixLayers::kIBL) {
573 bdtVars = { el_eta, misshitsf[idx], clsontrkf[idx], clsontrkrowsize[idx],
574 clsontrkcolsize[idx], clsofftrkrowsize[idx], clsofftrkcolsize[idx],
575 clsontrktot[idx], clsofftrktot[idx], trkalpha[idx],
576 trkchi2byndf[idx], trknpixdead[idx], trknblayerhits[idx] };
577 } else {
578 bdtVars = { el_eta, misshitsf[idx], clsontrkf[idx], clsontrkrowsize[idx],
579 clsontrkcolsize[idx], clsofftrkrowsize[idx], clsofftrkcolsize[idx],
580 clsontrktot[idx], clsofftrktot[idx], trkalpha[idx],
581 trkchi2byndf[idx], trknpixdead[idx] };
582 }
583 bdtweights[idx] = m_classBDT.at(partitionLabel)->GetClassification(bdtVars);
584 BDT_Weights.add(pixlayer, pixID, bdtweights[idx]);
585 }
586 }
587 fill2DProfLayerAccum( BDT_Weights );
588 fill1DModProfAccum( BDT_Weights, lb );
589
590 if (m_dumpTree) {
591 auto mvaGroup = getGroup("MVA");
592 auto mon_status_vec = Monitored::Collection("status_vec", status);
593 auto mon_h_vec = Monitored::Collection("holes_vec", holes);
594 auto mon_o_vec = Monitored::Collection("outliers_vec", outliers);
595 auto mon_m_vec = Monitored::Collection("meas_vec", measurements);
596
597 auto mon_clsontrkf_vec = Monitored::Collection("clsontrkf_vec", clsontrkf);
598 auto mon_clsontrk_vec = Monitored::Collection("clsontrk_vec", clsontrk);
599 auto mon_clsofftrk_vec = Monitored::Collection("clsofftrk_vec", clsofftrk);
600 auto mon_clsall_vec = Monitored::Collection("clsall_vec", clsall);
601
602 auto mon_clsontrksize_vec = Monitored::Collection("clsontrksize_vec", clsontrksize);
603 auto mon_clsontrkrowsize_vec = Monitored::Collection("clsontrkrowsize_vec", clsontrkrowsize);
604 auto mon_clsontrkcolsize_vec = Monitored::Collection("clsontrkcolsize_vec", clsontrkcolsize);
605 auto mon_clsofftrksize_vec = Monitored::Collection("clsofftrksize_vec", clsofftrksize);
606 auto mon_clsofftrkrowsize_vec = Monitored::Collection("clsofftrkrowsize_vec", clsofftrkrowsize);
607 auto mon_clsofftrkcolsize_vec = Monitored::Collection("clsofftrkcolsize_vec", clsofftrkcolsize);
608
609 auto mon_clsontrktot_vec = Monitored::Collection("clsontrktot_vec", clsontrktot);
610 auto mon_clsofftrktot_vec = Monitored::Collection("clsofftrktot_vec", clsofftrktot);
611 auto mon_trkalpha_vec = Monitored::Collection("trkalpha_vec", trkalpha);
612 auto mon_trkchi2byndf_vec = Monitored::Collection("trkchi2byndf_vec", trkchi2byndf);
613 auto mon_trknpixdead_vec = Monitored::Collection("trknpixdead_vec", trknpixdead);
614 auto mon_trknblayerhits_vec = Monitored::Collection("trknblayerhits_vec", trknblayerhits);
615
616 auto mon_mva_vec = Monitored::Collection("mva_vec", bdtweights);
617
618 fill( mvaGroup, lbval, mon_status_vec, mon_h_vec, mon_o_vec, mon_m_vec,
619 mon_clsontrkf_vec, mon_clsontrk_vec, mon_clsofftrk_vec, mon_clsall_vec,
620 mon_clsontrksize_vec, mon_clsontrkrowsize_vec, mon_clsontrkcolsize_vec,
621 mon_clsofftrksize_vec, mon_clsofftrkrowsize_vec, mon_clsofftrkcolsize_vec,
622 mon_clsontrktot_vec, mon_clsofftrktot_vec, mon_trkalpha_vec, mon_trkchi2byndf_vec,
623 mon_trknpixdead_vec, mon_trknblayerhits_vec, mon_mva_vec );
624 }
625
626 return StatusCode::SUCCESS;
627}
628
#define M_PI
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(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...
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
std::string PathResolverFindCalibFile(const std::string &logical_file_name)
const std::string pixBaseLayersLabel[PixLayers::NBASELAYERS]
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
const ServiceHandle< StoreGateSvc > & detStore() const
const ToolHandle< GenericMonitoringTool > & getGroup(const std::string &name) const
Get a specific monitoring tool from the tool handle array.
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.
This is a "hash" representation of an Identifier.
bool is_valid() const
Check if id is in a valid state.
Class to hold geometrical description of a silicon detector element.
Class to represent a position in the natural frame of a silicon sensor, for Pixel and SCT For Pixel: ...
virtual const Amg::Vector3D & normal() const override final
Get reconstruction local normal axes in global frame.
virtual Identifier identify() const override final
identifier of this detector element (inline)
Identifier identifierOfPosition(const Amg::Vector2D &localPos) const
Full identifier of the cell for a given position: assumes a raw local position (no Lorentz shift)
Specific class to represent the pixel measurements.
virtual const PixelCluster * prepRawData() const override final
returns the PrepRawData - is a SiCluster in this scope
const InDet::SiWidth & width() const
return width class reference
const Amg::Vector2D & colRow() const
Definition SiWidth.h:115
Declare a monitored scalar variable.
ToolHandle< Trk::IExtrapolator > m_trkextrapolator
const AtlasDetectorID * m_atlasid
ToolHandle< Trk::ITrackHoleSearchTool > m_holeSearchTool
virtual StatusCode initialize() override
initialize
SG::ReadHandleKey< PixelRDO_Container > m_pixelRDOName
std::map< std::string, std::unique_ptr< MVAUtils::BDT > > m_classBDT
virtual StatusCode fillHistograms(const EventContext &ctx) const override
adds event to the monitoring histograms
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_trackParticlesKey
std::string m_calibFolder
PixelAthMVAMonAlg(const std::string &name, ISvcLocator *pSvcLocator)
SG::ReadHandleKey< InDet::PixelClusterContainer > m_clustersKey
ToolHandle< InDet::IInDetTrackSelectionTool > m_trackSelTool
void fill2DProfLayerAccum(const VecAccumulator2DMap &accumulator) const
take VecAccumulator2DMap and fill the corresponding group
SG::ReadHandleKey< InDet::SiDetectorElementStatus > m_pixelDetElStatusActiveOnly
Optional read handle to get status data to test whether a pixel detector element is active.
int getNumberOfFEs(int pixlayer, int etaMod) const
helper function to get number of FEs per module
bool isClusterOnTrack(Identifier id, std::vector< std::pair< Identifier, double > > const &ClusterIDs) const
checks if cluster is on track
void fill1DModProfAccum(const VecAccumulator2DMap &accumulator, int lumiblock) const
take VecAccumulator2DMap and fill 3D arrays [layer, pm, em] with its values and lumiblock
ServiceHandle< InDetDD::IPixelReadoutManager > m_pixelReadout
bool isGood(const InDet::SiDetectorElementStatus *element_status, const IdentifierHash &module_hash) const
int getPixLayersID(int ec, int ld) const
helper function to get layers ID
std::vector< int > m_modData[PixLayers::NBASELAYERS]
bool isIBL2D(int hashID) const
helper function to check if module is IBL planar based on pixel hash ID
std::tuple< bool, bool > isChipGood(const IdentifierHash &module_hash, unsigned int chip_i) const
SG::ReadHandleKey< InDet::SiDetectorElementStatus > m_pixelDetElStatus
Optional read handle to get status data to test whether a pixel detector element is good.
bool isActive(const InDet::SiDetectorElementStatus *element_status, const IdentifierHash &module_hash) const
SG::ReadHandle< InDet::SiDetectorElementStatus > getPixelDetElStatus(const SG::ReadHandleKey< InDet::SiDetectorElementStatus > &key, const EventContext &ctx) const
virtual StatusCode initialize() override
initialize
std::vector< Identifier >::const_iterator const_id_iterator
Definition PixelID.h:72
const_pointer_type cptr()
Dereference the pointer.
Class for a CylinderSurface in the ATLAS detector.
This class is the pure abstract base class for all fittable tracking measurements.
virtual const Surface & associatedSurface() const =0
Interface method to get the associated Surface.
virtual bool type(MeasurementBaseType::Type type) const =0
Interface method checking the type.
const Amg::Vector3D & momentum() const
Access method for the momentum.
virtual const Surface & associatedSurface() const override=0
Access to the Surface associated to the Parameters.
Amg::Vector2D localPosition() const
Access method for the local coordinates, local parameter definitions differ for each surface type.
const std::vector< Identifier > & rdoList() const
return the List of rdo identifiers (pointers)
Class to handle RIO On Tracks ROT) for InDet and Muons, it inherits from the common MeasurementBase.
Definition RIO_OnTrack.h:70
Identifier identify() const
return the identifier -extends MeasurementBase
const TrkDetElementBase * associatedDetectorElement() const
return associated Detector Element
Identifier associatedDetectorElementIdentifier() const
return Identifier of the associated Detector Element
@ 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.
const Trk::TrackStates * trackStateOnSurfaces() const
return a pointer to a const DataVector of const TrackStateOnSurfaces.
virtual Identifier identify() const =0
Identifier.
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.
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 3, 1 > Vector3D
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)
ValuesCollection< T > Collection(std::string name, const T &collection)
Declare a monitored (double-convertible) collection.
const int pixPhiSteps[PixLayers::NBASELAYERS]
const int pixEtaSteps[PixLayers::NBASELAYERS]
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
@ alongMomentum
DataVector< const Trk::TrackStateOnSurface > TrackStates
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
ParametersBase< TrackParametersDim, Charged > TrackParameters
ParametersT< TrackParametersDim, Charged, PlaneSurface > AtaPlane
@ expectInnermostPixelLayerHit
Do we expect a 0th-layer barrel hit for this track?
@ numberOfPixelHoles
number of pixel layers on track with absence of hits [unit8_t].
@ numberOfNextToInnermostPixelLayerHits
these are the hits in the 1st pixel barrel layer
@ expectNextToInnermostPixelLayerHit
Do we expect a 1st-layer barrel hit for this track?
@ numberOfPixelHits
these are the pixel hits, including the b-layer [unit8_t].
@ numberOfPixelDeadSensors
number of dead pixel sensors crossed [unit8_t].