ATLAS Offline Software
Loading...
Searching...
No Matches
PixelClusterOnTrackTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
6// Implementation file for class PixelClusterOnTrackTool
8// (c) ATLAS Detector software
10// AlgTool used for PixelClusterOnTrack object production
12// started 1.0 21/04/2004 I.Gavrilenko - see ChangeLog
14
22
23#include <cmath>
25
26//clustermap is most likely to be removed at later date
27#define __clustermap
28
29
31// Constructor
33
34namespace
35{
36 // using x*x might be quicker than pow(x,2), depends on compiler optimisation
37 inline double
38 square(const double x) {
39 return x * x;
40 }
41
42 double
43 distance(const std::vector<Amg::Vector2D> &vectorOfPositions,
44 const std::vector<Amg::Vector2D> &allLocalPositions,
45 const std::vector<Amg::MatrixX> &allErrorMatrix,
46 const int i, const int j, const int k) {
47 return
48 square(vectorOfPositions[i][0] - allLocalPositions[0][0]) / allErrorMatrix[0](0, 0) +
49 square(vectorOfPositions[j][0] - allLocalPositions[1][0]) / allErrorMatrix[1](0, 0) +
50 square(vectorOfPositions[k][0] - allLocalPositions[2][0]) / allErrorMatrix[2](0, 0) +
51 square(vectorOfPositions[i][1] - allLocalPositions[0][1]) / allErrorMatrix[0](1, 1) +
52 square(vectorOfPositions[j][1] - allLocalPositions[1][1]) / allErrorMatrix[1](1, 1) +
53 square(vectorOfPositions[k][1] - allLocalPositions[2][1]) / allErrorMatrix[2](1, 1);
54 }
55}
56
58 (const std::string &t, const std::string &n, const IInterface *p) :
59 ::AthAlgTool(t, n, p) {
60 declareInterface<IRIO_OnTrackCreator>(this);
61}
62
64// Destructor
66
68
70// Initialisation
72
73StatusCode
75
76 ATH_MSG_DEBUG(name() << " initialize()");
77
79 ATH_MSG_DEBUG("Error strategy is" << m_errorStrategy);
80
81 if (m_IBLParameterSvc.retrieve().isFailure()) {
82 ATH_MSG_WARNING("Could not retrieve IBLParameterSvc");
83 } else {
84 m_IBLParameterSvc->setBoolParameters(m_applyNNcorrectionProperty.value(), "doPixelClusterSplitting");
85 m_IBLParameterSvc->setBoolParameters(m_IBLAbsent, "IBLAbsent");
86 }
88
89 ATH_CHECK(m_clusterErrorKey.initialize());
91
92 // get the error scaling tool
94 if (!m_pixelErrorScalingKey.key().empty()) ATH_MSG_DEBUG("Detected need for scaling Pixel errors.");
95
96 // get the module distortions tool
98 if (m_disableDistortions) ATH_MSG_DEBUG("No PixelDistortions will be simulated.");
99
100 ATH_CHECK (detStore()->retrieve(m_pixelid, "PixelID"));
101
104 ATH_CHECK(m_NnClusterizationFactory.retrieve( DisableTool{!m_applyNNcorrection} ));
105
106 // Include IBL calibration constants
107 //Moved to initialize to remove statics and prevent repitition
108
109 constexpr double phimin=-0.27, phimax=0.27;
110 for (int i=0; i<=s_nbinphi; i++) m_phix[i]=phimin+i*(phimax-phimin)/s_nbinphi;
111 constexpr double etacen[s_nbineta]={-0.,1.,1.55,1.9,2.15,2.35};
112 m_etax[0]=0.; m_etax[s_nbineta]=2.7;
113 for (int i=0; i<s_nbineta-1; i++) m_etax[i+1]=(etacen[i]+etacen[i+1])/2.;
114
116#include "IBL_calibration.h"
118
119 ATH_CHECK(m_lorentzAngleTool.retrieve());
120 return StatusCode::SUCCESS;
121}
122
123
124
126// Finalize
128
129StatusCode
131return StatusCode::SUCCESS;
132}
133
135// Trk::SiClusterOnTrack production
137
138
141 (const Trk::PrepRawData &rio, const Trk::TrackParameters &trackPar, const EventContext& ctx) const {
142
143 const auto *element= dynamic_cast<const InDetDD::SiDetectorElement *>(rio.detectorElement());
144 if ((not m_applyNNcorrection) or (element and element->isBlayer() and (not m_NNIBLcorrection) and (not m_IBLAbsent))){
145 return correctDefault(rio, trackPar, ctx);
146 }else {
147 if (m_errorStrategy == 0 || m_errorStrategy == 1) {
148 // version from Giacinto
150 return nullptr;
151 }
152 // if we try broad errors, get Pixel Cluster to test if it is split
153 const InDet::PixelCluster *pix = nullptr;
155 pix = static_cast<const InDet::PixelCluster *>(&rio);
156 }
157 if (!pix) {
158 return nullptr;
159 }
161 if (splitProb.isSplit()) {
162 return correctNN(rio, trackPar, ctx);
163 } else {
164 return correctDefault(rio, trackPar, ctx);
165 }
166 } else {
167 return correctNN(rio, trackPar, ctx);
168 }
169 }
170}
171
172
173
179 (const Trk::PrepRawData &rio, const Trk::TrackParameters &trackPar, const EventContext& ctx) const {
180 using CLHEP::micrometer;
181
182 const double TOPHAT_SIGMA = 1. / std::sqrt(12.);
183
184 const InDet::PixelCluster *pix = nullptr;
186 pix = static_cast<const InDet::PixelCluster *>(&rio);
187 }
188 else{
189 return nullptr;
190 }
191
192 ATH_MSG_VERBOSE("Correct called with Error strategy " << m_errorStrategy);
193
194 // PixelClusterOnTrack production
195 //
197 Amg::Vector3D glob(pix->globalPosition());
198
199
200 // Get pointer to detector element
201 const InDetDD::SiDetectorElement *element = pix->detectorElement();
202 if (!element) {
203 return nullptr;
204 }
205 bool blayer = element->isBlayer();
206 IdentifierHash iH = element->identifyHash();
207
208 double errphi = -1;
209 double erreta = -1;
210
211 if (pix->rdoList().empty()) {
212 ATH_MSG_WARNING("Pixel RDO-list size is 0, check integrity of pixel clusters! stop ROT creation.");
213 return nullptr;
214 } else {
215 const InDetDD::PixelModuleDesign *design =
216 dynamic_cast<const InDetDD::PixelModuleDesign *>(&element->design());
217
218 // get candidate track angle in module local frame
219 const Amg::Vector3D& my_track = trackPar.momentum();
220 const Amg::Vector3D& my_normal = element->normal();
221 const Amg::Vector3D& my_phiax = element->phiAxis();
222 const Amg::Vector3D& my_etaax = element->etaAxis();
223 float trkphicomp = my_track.dot(my_phiax);
224 float trketacomp = my_track.dot(my_etaax);
225 float trknormcomp = my_track.dot(my_normal);
226 double bowphi = std::atan2(trkphicomp, trknormcomp);
227 double boweta = std::atan2(trketacomp, trknormcomp);
228 float etatrack = trackPar.eta();
229
230 float tanl = m_lorentzAngleTool->getTanLorentzAngle(iH, ctx);
231 int readoutside = element->design().readoutSide();
232
233 // map the angles of inward-going tracks onto [-PI/2, PI/2]
234 if (bowphi > M_PI *0.5) {
235 bowphi -= M_PI;
236 }
237 if (bowphi < -M_PI *0.5) {
238 bowphi += M_PI;
239 }
240 // finally, subtract the Lorentz angle effect
241 // the readoutside term is needed because of a bug in old
242 // geometry versions (CSC-01-* and CSC-02-*)
243 double angle = std::atan(std::tan(bowphi) - readoutside * tanl);
244
245 // try to understand...
246 const Identifier element_id = element->identify();
247 int PixEtaModule = m_pixelid->eta_module(element_id);
248 int PixPhiModule = m_pixelid->phi_module(element_id);
249 double PixTrkPt = trackPar.pT();
250 double PixTrkEta = trackPar.eta();
251 ATH_MSG_VERBOSE("tanl = " << tanl << " readout side is " << readoutside <<
252 " module " << PixEtaModule << " " << PixPhiModule <<
253 " track pt, eta = " << PixTrkPt << " " << PixTrkEta <<
254 " track momentum phi, norm = " << trkphicomp << " " <<
255 trknormcomp << " bowphi = " << bowphi << " angle = " << angle);
256
257 float omegaphi = pix->omegax();
258 float omegaeta = pix->omegay();
259 double localphi = -9999.;
260 double localeta = -9999.;
261
262 const std::vector<Identifier> & rdos = pix->rdoList();
263 InDetDD::SiLocalPosition meanpos(0, 0, 0);
264 int rowmin = 9999;
265 int rowmax = -9999;
266 int colmin = 9999;
267 int colmax = -9999;
268 for (const auto & rId:rdos) {
269 const int row = m_pixelid->phi_index(rId);
270 const int col = m_pixelid->eta_index(rId);
271 rowmin = std::min(rowmin, row);
272 rowmax = std::max(rowmax,row);
273 colmin = std::min(colmin, col);
274 colmax = std::max(colmax, col);
275 meanpos += design->positionFromColumnRow(col, row);
276 }
277 meanpos = meanpos / rdos.size();
279 design->positionFromColumnRow(colmin, rowmin);
281 design->positionFromColumnRow(colmax, rowmin);
283 design->positionFromColumnRow(colmin, rowmax);
285 design->positionFromColumnRow(colmax, rowmax);
286
287 InDetDD::SiLocalPosition centroid = 0.25 * (pos1 + pos2 + pos3 + pos4);
288 double shift = m_lorentzAngleTool->getLorentzShift(iH, ctx);
289 int nrows = rowmax - rowmin + 1;
290 int ncol = colmax - colmin + 1;
291 double ang = 999.;
292
293 // TOT interpolation for collision data
294 // Force IBL to use digital clustering and broad errors.
296 if (m_positionStrategy > 0 && omegaphi > -0.5 && omegaeta > -0.5) {
297 localphi = centroid.xPhi() + shift;
298 localeta = centroid.xEta();
299 // barrel
300 if (element->isBarrel()) {
301 ang = 180 * angle * M_1_PI; //M_1_PI in cmath, = 1/pi
302 double delta = 0.;
303 if (m_IBLAbsent || !blayer) {
304 delta = offlineCalibData->getPixelChargeInterpolationParameters()->getDeltaXbarrel(nrows, ang, 1);
305 } else { // special calibration for IBL
306 if (angle < m_phix[0] || angle > m_phix[s_nbinphi] || nrows != 2) {
307 delta = 0.;
308 }else {
309 int bin = -1;
310 while (angle > m_phix[bin + 1]) {
311 bin++;
312 }
313 if ((bin >= 0)and(bin < s_nbinphi)) {
314 delta = m_calphi[bin];
315 } else {
316 ATH_MSG_ERROR("bin out of range in line " << __LINE__ << " of PixelClusterOnTrackTool.cxx.");
317 }
318 }
319
320 if (offlineCalibData->getPixelChargeInterpolationParameters()->getVersion()<-1) {
321 delta = offlineCalibData->getPixelChargeInterpolationParameters()->getDeltaXbarrel(nrows, ang, 0);
322 }
323 }
324 localphi += delta * (omegaphi - 0.5);
325 // settle the sign/pi periodicity issues
326 double thetaloc = -999.;
327 if (boweta > -0.5 * M_PI && boweta < M_PI / 2.) { //M_PI_2 in cmath
328 thetaloc = M_PI_2 - boweta;
329 }else if (boweta > M_PI_2 && boweta < M_PI) {
330 thetaloc = 1.5 * M_PI - boweta;
331 } else { // 3rd quadrant
332 thetaloc = -M_PI_2 - boweta;
333 }
334 double etaloc = -1 * log(tan(thetaloc * 0.5));
335 if (m_IBLAbsent || !blayer) {
336 delta = offlineCalibData->getPixelChargeInterpolationParameters()->getDeltaYbarrel(ncol, etaloc, 1);
337 } else { // special calibration for IBL
338 etaloc = std::abs(etaloc);
339 if (etaloc < m_etax[0] || etaloc > m_etax[s_nbineta]) {
340 delta = 0.;
341 } else {
342 int bin = -1;
343 while (etaloc > m_etax[bin + 1]) {
344 bin++;
345 }
346 if ((bin >= 0)and(bin < s_nbineta)) {
347 if (ncol == bin) {
348 delta = m_caleta[bin][0];
349 } else if (ncol == bin + 1) {
350 delta = m_caleta[bin][1];
351 } else if (ncol == bin + 2) {
352 delta = m_caleta[bin][2];
353 } else {
354 delta = 0.;
355 }
356 } else {// bin out of range of array indices
357 ATH_MSG_ERROR("bin out of range in line " << __LINE__ << " of PixelClusterOnTrackTool.cxx.");
358 }
359 }
360 if (offlineCalibData->getPixelChargeInterpolationParameters()->getVersion()<-1) {
361 delta = offlineCalibData->getPixelChargeInterpolationParameters()->getDeltaYbarrel(ncol, std::abs(etaloc), 0);
362 }
363 }
364 localeta += delta * (omegaeta - 0.5);
365 }else {
366 // collision endcap data
367 if (m_positionStrategy == 1) {
370 localphi += deltax * (omegaphi - 0.5);
371 localeta += deltay * (omegaeta - 0.5);
372 }
373 // SR1 cosmics endcap data
374 // some parametrization dependent on track angle
375 // would be better here
376 else if (m_positionStrategy == 20) {
377 double deltax = 35 * micrometer;
378 double deltay = 35 * micrometer;
379 localphi += deltax * (omegaphi - 0.5);
380 localeta += deltay * (omegaeta - 0.5);
381 }
382 }
383 }
384// digital
385 else {
386 localphi = meanpos.xPhi() + shift;
387 localeta = meanpos.xEta();
388 }
389
390 const InDet::SiWidth& width = pix->width();
391
392 // Error strategies
393
394 // For very shallow tracks the cluster can easily break as
395 // the average charge per pixel is of the order of the threshold
396 // In this case, an error equal to the geometrical projection
397 // of the track path in silicon onto the module surface seems
398 // appropriate
399 if (std::abs(angle) > 1) {
400 errphi = 250 * micrometer * std::tan(std::abs(angle)) * TOPHAT_SIGMA;
401 erreta = width.z() > 250 * micrometer * std::tan(std::abs(boweta)) ?
402 width.z() * TOPHAT_SIGMA : 250 * micrometer * std::tan(std::abs(boweta)) * TOPHAT_SIGMA;
403 ATH_MSG_VERBOSE("Shallow track with tanl = " << tanl << " bowphi = " <<
404 bowphi << " angle = " << angle << " width.z = " << width.z() <<
405 " errphi = " << errphi << " erreta = " << erreta);
406 }else if (m_errorStrategy == 0) {
407 errphi = width.phiR() * TOPHAT_SIGMA;
408 erreta = width.z() * TOPHAT_SIGMA;
409 }else if (m_errorStrategy == 1) {
410 errphi = (width.phiR() / nrows) * TOPHAT_SIGMA;
411 erreta = (width.z() / ncol) * TOPHAT_SIGMA;
412 }else if (m_errorStrategy == 2) {
413 if (element->isBarrel()) {
414 if (m_IBLAbsent || !blayer) {
415 int ibin = offlineCalibData->getPixelClusterOnTrackErrorData()->getBarrelBinPhi(ang, nrows);
416 errphi = offlineCalibData->getPixelClusterOnTrackErrorData()->getPixelBarrelPhiError(ibin);
417 } else { // special calibration for IBL
418 if (angle < m_phix[0] || angle > m_phix[s_nbinphi]) {
419 errphi = width.phiR() * TOPHAT_SIGMA;
420 } else {
421 int bin = -1;// cannot be used as array index, which will happen if angle<m_phix[bin+1]
422 while (angle > m_phix[bin + 1]) {
423 bin++;
424 }
425 if ((bin >= 0)and(bin < s_nbinphi)) {
426 if (nrows == 1) {
427 errphi = m_calerrphi[bin][0];
428 } else if (nrows == 2) {
429 errphi = m_calerrphi[bin][1];
430 } else {
431 errphi = m_calerrphi[bin][2];
432 }
433 } else {
434 ATH_MSG_ERROR("bin out of range in line " << __LINE__ << " of PixelClusterOnTrackTool.cxx.");
435 }
436 }
437 }
438
439 if (m_IBLAbsent || !blayer) {
440 int ibin = offlineCalibData->getPixelClusterOnTrackErrorData()->getBarrelBinEta(std::abs(etatrack), ncol, nrows);
441 erreta = offlineCalibData->getPixelClusterOnTrackErrorData()->getPixelBarrelEtaError(ibin);
442 } else { // special calibration for IBL
443 double etaloc = std::abs(etatrack);
444 if (etaloc < m_etax[0] || etaloc > m_etax[s_nbineta]) {
445 erreta = width.z() * TOPHAT_SIGMA;
446 } else {
447 int bin = 0;
448 while (bin < s_nbineta && etaloc > m_etax[bin + 1]) {
449 ++bin;
450 }
451 if (bin >= s_nbineta) {
452 ATH_MSG_ERROR("bin out of range in line " << __LINE__ << " of PixelClusterOnTrackTool.cxx.");
453 } else {
454 if (ncol == bin) {
455 erreta = m_calerreta[bin][0];
456 } else if (ncol == bin + 1) {
457 erreta = m_calerreta[bin][1];
458 } else if (ncol == bin + 2) {
459 erreta = m_calerreta[bin][2];
460 } else {
461 erreta = width.z() * TOPHAT_SIGMA;
462 }
463 }
464 }
465 }
466 }else {
467 int ibin = offlineCalibData->getPixelClusterErrorData()->getEndcapBin(ncol, nrows);
468 errphi = offlineCalibData->getPixelClusterErrorData()->getPixelEndcapPhiError(ibin);
469 erreta = offlineCalibData->getPixelClusterErrorData()->getPixelEndcapRError(ibin);
470 }
471 if (errphi > erreta) {
472 erreta = width.z() * TOPHAT_SIGMA;
473 }
474 }
475
476 Amg::Vector2D locpos = Amg::Vector2D(localphi, localeta);
477 if (element->isBarrel() && !m_disableDistortions) {
478 correctBow(element->identify(), locpos, bowphi, boweta,ctx);
479 }
480
481
482 locpar = Trk::LocalParameters(locpos);
483 centroid = InDetDD::SiLocalPosition(localeta, localphi, 0.);
484 glob = element->globalPosition(centroid);
485 }
486
487 // Error matrix production
488
489 Amg::MatrixX cov = pix->localCovariance();
490
491 // corrected phi error
492 if (errphi > 0) {
493 cov(0, 0) = errphi * errphi;
494 }
495 if (erreta > 0) {
496 cov(1, 1) = erreta * erreta;
497 }
498
499 ATH_MSG_VERBOSE(" errphi = " << errphi << " erreta = " << erreta);
500
501 // create new copy of error matrix
502 if (!m_pixelErrorScalingKey.key().empty()) {
505 ->getScaledCovariance(std::move(cov), *m_pixelid,
506 element->identify());
507 }
508 bool isbroad = m_errorStrategy == 0;
510 std::move(locpar),
511 std::move(cov), iH,
512 glob, pix->gangedPixel(), isbroad);
513}
514
515void
517 const double theta, const EventContext& ctx) const {
518 Amg::Vector3D dir(tan(phi), tan(theta), 1.);
519 Amg::Vector2D newpos = SG::ReadCondHandle<PixelDistortionData>(m_distortionKey, ctx)->correctReconstruction(m_pixelid->wafer_hash(id), localpos, dir);
520
521 localpos = newpos;
522}
523
526 (const Trk::PrepRawData &rio, const Trk::TrackParameters &trackPar,
527 const InDet::PixelClusterStrategy strategy) const {
528 int initial_errorStrategy;
530
531 switch (strategy) {
532 case InDet::PIXELCLUSTER_OUTLIER: // if cluster is outlier, increase errors
533 case InDet::PIXELCLUSTER_SHARED:
534 initial_errorStrategy = m_errorStrategy;
535 m_errorStrategy = 0; // error as size of cluster /sqrt(12)
536 newROT = correct(rio, trackPar);
537 m_errorStrategy = initial_errorStrategy;
538 return newROT;
539
540 default:
541 return correct(rio, trackPar);
542 }
543}
544
545// GP: NEW correct() method in case of NN based calibration */
548 (const Trk::PrepRawData &rio,
549 const Trk::TrackParameters &trackPar,
550 const EventContext& ctx) const {
551
552 const InDet::PixelCluster *pixelPrepCluster = nullptr;
554 pixelPrepCluster = static_cast<const InDet::PixelCluster *>(&rio);
555 }
556
557 if (pixelPrepCluster == nullptr) {
558 ATH_MSG_WARNING("This is not a pixel cluster, return 0.");
559 return nullptr;
560 }
561
562 const InDetDD::SiDetectorElement *element = pixelPrepCluster->detectorElement();
563 if (!element) {
564 ATH_MSG_WARNING("Cannot access detector element. Aborting cluster correction...");
565 return nullptr;
566 }
567
568 IdentifierHash iH = element->identifyHash();
569
571 Amg::Vector3D glob(pixelPrepCluster->globalPosition());
572
573 const Amg::Vector3D& my_track = trackPar.momentum();
574 const Amg::Vector3D& my_normal = element->normal();
575 const Amg::Vector3D& my_phiax = element->phiAxis();
576 const Amg::Vector3D& my_etaax = element->etaAxis();
577 float trkphicomp = my_track.dot(my_phiax);
578 float trketacomp = my_track.dot(my_etaax);
579 float trknormcomp = my_track.dot(my_normal);
580 double bowphi = std::atan2(trkphicomp, trknormcomp);
581 double boweta = std::atan2(trketacomp, trknormcomp);
582
583 Amg::Vector2D locpos = pixelPrepCluster->localPosition();
584 if (element->isBarrel() && !m_disableDistortions) {
585 correctBow(element->identify(), locpos, bowphi, boweta,ctx);
586 }
587
589 Amg::MatrixX cov = pixelPrepCluster->localCovariance();
590
591 return new InDet::PixelClusterOnTrack(pixelPrepCluster, std::move(locpar), std::move(cov), iH, glob,
592 pixelPrepCluster->gangedPixel(), false);
593 }
594
595
596
597 Amg::Vector2D finalposition;
598 Amg::MatrixX finalerrormatrix;
599
600 if (m_usingTIDE_Ambi) {
601 if (!getErrorsTIDE_Ambi(pixelPrepCluster, trackPar, finalposition, finalerrormatrix)) {
602 return correctDefault(rio, trackPar, ctx);
603 }
604 }else {
605 if (!getErrorsDefaultAmbi(pixelPrepCluster, trackPar, finalposition, finalerrormatrix)) {
606 return correctDefault(rio, trackPar, ctx);
607 }
608 }
609
610 ATH_MSG_DEBUG( " Old position x: " << pixelPrepCluster->localPosition()[0]
611 << " +/- " << std::sqrt(pixelPrepCluster->localCovariance()(0, 0))
612 << " y: " << pixelPrepCluster->localPosition()[1]
613 << " +/- " << std::sqrt(pixelPrepCluster->localCovariance()(1, 1)) <<"\n"
614 << " Final position x: " << finalposition[0]
615 << " +/- " << std::sqrt(finalerrormatrix(0, 0))
616 << " y: " << finalposition[1] << " +/- "
617 <<std::sqrt(finalerrormatrix(1, 1)) );
618
619 const Amg::Vector3D& my_track = trackPar.momentum();
620 Amg::Vector3D my_normal = element->normal();
621 Amg::Vector3D my_phiax = element->phiAxis();
622 Amg::Vector3D my_etaax = element->etaAxis();
623 float trkphicomp = my_track.dot(my_phiax);
624 float trketacomp = my_track.dot(my_etaax);
625 float trknormcomp = my_track.dot(my_normal);
626 double bowphi = std::atan2(trkphicomp, trknormcomp);
627 double boweta = std::atan2(trketacomp, trknormcomp);
628
629 if (element->isBarrel() && !m_disableDistortions) {
630 correctBow(element->identify(), finalposition, bowphi, boweta,ctx);
631 }
632
633 Amg::MatrixX cov = finalerrormatrix;
634 // create new copy of error matrix
635 if (!m_pixelErrorScalingKey.key().empty()) {
639 ->getScaledCovariance(std::move(cov), *m_pixelid,
640 element->identify());
641 }
642
643 InDetDD::SiLocalPosition centroid = InDetDD::SiLocalPosition(finalposition[1],
644 finalposition[0],
645 0);
646 Trk::LocalParameters locpar = Trk::LocalParameters(finalposition);
647
648 const Amg::Vector3D &glob = element->globalPosition(centroid);
649
650
651 return new InDet::PixelClusterOnTrack(pixelPrepCluster, std::move(locpar),
652 std::move(cov), iH,
653 glob,
654 pixelPrepCluster->gangedPixel(),
655 false);
656}
657
658bool
660 const Trk::TrackParameters &trackPar,
661 Amg::Vector2D &finalposition,
662 Amg::MatrixX &finalerrormatrix) const {
663 std::vector<Amg::Vector2D> vectorOfPositions;
664 int numberOfSubclusters = 1;
665 vectorOfPositions.push_back(pixelPrepCluster->localPosition());
666
669 InDet::PixelGangedClusterAmbiguities::const_iterator mapBegin = splitClusterMap->begin();
670 InDet::PixelGangedClusterAmbiguities::const_iterator mapEnd = splitClusterMap->end();
671 for (InDet::PixelGangedClusterAmbiguities::const_iterator mapIter = mapBegin; mapIter != mapEnd; ++mapIter) {
672 const SiCluster *first = (*mapIter).first;
673 const SiCluster *second = (*mapIter).second;
674 if (first == pixelPrepCluster && second != pixelPrepCluster) {
675 ATH_MSG_DEBUG("Found additional split cluster in ambiguity map (+=1).");
676 numberOfSubclusters += 1;
677 const SiCluster *otherOne = second;
678 const InDet::PixelCluster *pixelAddCluster = nullptr;
680 pixelAddCluster = static_cast<const InDet::PixelCluster *>(otherOne);
681 }
682 if (pixelAddCluster == nullptr) {
683 ATH_MSG_WARNING("Pixel ambiguity map has empty pixel cluster. Please DEBUG!");
684 continue;
685 }
686 vectorOfPositions.push_back(pixelAddCluster->localPosition());
687
688 ATH_MSG_DEBUG( "Found one more pixel cluster. Position x: "
689 << pixelAddCluster->localPosition()[0] << "y: " << pixelAddCluster->localPosition()[1]);
690 }// find relevant element of map
691 }// iterate over map
692 }
693
694 // now you have numberOfSubclusters and the vectorOfPositions (Amg::Vector2D)
695
696 if (trackPar.surfaceType() != Trk::SurfaceType::Plane ||
697 trackPar.type() != Trk::AtaSurface) {
699 "Parameters are not at a plane ! Aborting cluster correction... ");
700 return false;
701 }
702
703 std::vector<Amg::Vector2D> allLocalPositions;
704 std::vector<Amg::MatrixX> allErrorMatrix;
705 allLocalPositions =
706 m_NnClusterizationFactory->estimatePositions(*pixelPrepCluster,
707 trackPar.associatedSurface(),
708 trackPar,
709 allErrorMatrix,
710 numberOfSubclusters);
711
712 if (allLocalPositions.empty()) {
713 ATH_MSG_DEBUG( " Cluster cannot be treated by NN. Giving back to default clusterization " );
714
715 return false;
716 }
717
718 if (allLocalPositions.size() != size_t(numberOfSubclusters)) {
719 ATH_MSG_WARNING( "Returned position vector size " << allLocalPositions.size() <<
720 " not according to expected number of subclusters: " << numberOfSubclusters << ". Abort cluster correction..." );
721 return false;
722 }
723
724
725 // GP: now the not so nice part of matching the new result with the old one...
726 // Takes the error into account to improve the matching
727
728 if (numberOfSubclusters == 1) {
729 finalposition = allLocalPositions[0];
730 finalerrormatrix = allErrorMatrix[0];
731 }
732
733 else if (numberOfSubclusters == 2) {
734 double distancesq1 =
735 square(vectorOfPositions[0][0] - allLocalPositions[0][0]) / allErrorMatrix[0](0, 0) +
736 square(vectorOfPositions[1][0] - allLocalPositions[1][0]) / allErrorMatrix[1](0, 0) +
737 square(vectorOfPositions[0][1] - allLocalPositions[0][1]) / allErrorMatrix[0](1, 1) +
738 square(vectorOfPositions[1][1] - allLocalPositions[1][1]) / allErrorMatrix[1](1, 1);
739
740 double distancesq2 =
741 square(vectorOfPositions[1][0] - allLocalPositions[0][0]) / allErrorMatrix[0](0, 0) +
742 square(vectorOfPositions[0][0] - allLocalPositions[1][0]) / allErrorMatrix[1](0, 0) +
743 square(vectorOfPositions[1][1] - allLocalPositions[0][1]) / allErrorMatrix[0](1, 1) +
744 square(vectorOfPositions[0][1] - allLocalPositions[1][1]) / allErrorMatrix[1](1, 1);
745
747 " Old pix (1) x: " << vectorOfPositions[0][0] << " y: " << vectorOfPositions[0][1] << "\n"
748 << " Old pix (2) x: " << vectorOfPositions[1][0] << " y: " << vectorOfPositions[1][1] << "\n"
749 << " Pix (1) x: " << allLocalPositions[0][0] << " +/- " << std::sqrt(allErrorMatrix[0](0, 0))
750 << " y: " << allLocalPositions[0][1] << " +/- " << std::sqrt(allErrorMatrix[0](1, 1)) <<"\n"
751 << " Pix (2) x: " << allLocalPositions[1][0] << " +/- " << std::sqrt(allErrorMatrix[1](0, 0))
752 << " y: " << allLocalPositions[1][1] << " +/- " << std::sqrt(allErrorMatrix[1](1, 1)) << "\n"
753 << " Old (1) new (1) dist: " << std::sqrt(distancesq1) << " Old (1) new (2) " << std::sqrt(distancesq2) );
754
755
756 if (distancesq1 < distancesq2) {
757 finalposition = allLocalPositions[0];
758 finalerrormatrix = allErrorMatrix[0];
759 }else {
760 finalposition = allLocalPositions[1];
761 finalerrormatrix = allErrorMatrix[1];
762 }
763 }
764
765
766 else if (numberOfSubclusters == 3) {
767 double distances[6];
768
769 distances[0] = distance(vectorOfPositions, allLocalPositions, allErrorMatrix, 0, 1, 2);
770 distances[1] = distance(vectorOfPositions, allLocalPositions, allErrorMatrix, 0, 2, 1);
771 distances[2] = distance(vectorOfPositions, allLocalPositions, allErrorMatrix, 1, 0, 2);
772 distances[3] = distance(vectorOfPositions, allLocalPositions, allErrorMatrix, 1, 2, 0);
773 distances[4] = distance(vectorOfPositions, allLocalPositions, allErrorMatrix, 2, 0, 1);
774 distances[5] = distance(vectorOfPositions, allLocalPositions, allErrorMatrix, 2, 1, 0);
775
776 int smallestDistanceIndex = -10;
777 double minDistance = 1e10;
778
779 for (int i = 0; i < 6; i++) {
780 ATH_MSG_VERBOSE(" distance n.: " << i << " distance is: " << distances[i]);
781
782 if (distances[i] < minDistance) {
783 minDistance = distances[i];
784 smallestDistanceIndex = i;
785 }
786 }
787
788 ATH_MSG_DEBUG(" The minimum distance is : " << minDistance << " for index: " << smallestDistanceIndex);
789
790 if (smallestDistanceIndex == 0 || smallestDistanceIndex == 1) {
791 finalposition = allLocalPositions[0];
792 finalerrormatrix = allErrorMatrix[0];
793 }
794 if (smallestDistanceIndex == 2 || smallestDistanceIndex == 4) {
795 finalposition = allLocalPositions[1];
796 finalerrormatrix = allErrorMatrix[1];
797 }
798 if (smallestDistanceIndex == 3 || smallestDistanceIndex == 5) {
799 finalposition = allLocalPositions[2];
800 finalerrormatrix = allErrorMatrix[2];
801 }
802 }
803 return true;
804}
805
806bool
808 const Trk::TrackParameters &trackPar,
809 Amg::Vector2D &finalposition,
810 Amg::MatrixX &finalerrormatrix) const {
812 std::vector<Amg::Vector2D> vectorOfPositions;
813 int numberOfSubclusters = 1;
816 numberOfSubclusters = 1 + splitClusterMap->count(pixelPrepCluster);
817
818 if (splitClusterMap->count(pixelPrepCluster) == 0 && splitProb.isSplit()) {
819 numberOfSubclusters = 2;
820 }
821 if (splitClusterMap->count(pixelPrepCluster) != 0 && !splitProb.isSplit()) {
822 numberOfSubclusters = 1;
823 }
824 }
825
826 // Position NN expects 3 clusters at most
827 if(numberOfSubclusters>3) numberOfSubclusters = 3;
828
829 // now you have numberOfSubclusters and the vectorOfPositions (Amg::Vector2D)
830 if (trackPar.surfaceType() != Trk::SurfaceType::Plane ||
831 trackPar.type() != Trk::AtaSurface) {
832 ATH_MSG_WARNING("Parameters are not at a plane surface ! Aborting cluster "
833 "correction... ");
834 return false;
835 }
836
837 std::vector<Amg::Vector2D> allLocalPositions;
838 std::vector<Amg::MatrixX> allErrorMatrix;
839 allLocalPositions = m_NnClusterizationFactory->estimatePositions(
840 *pixelPrepCluster,
841 trackPar.associatedSurface(),
842 trackPar,
843 allErrorMatrix,
844 numberOfSubclusters);
845
846 if (allLocalPositions.empty()) {
848 " Cluster cannot be treated by NN. Giving back to default clusterization, too big: " <<
849 splitProb.isTooBigToBeSplit());
850 return false;
851 }
852
853 if (allLocalPositions.size() != size_t(numberOfSubclusters)) {
855 "Returned position vector size " << allLocalPositions.size() << " not according to expected number of subclusters: " << numberOfSubclusters <<
856 ". Abort cluster correction...");
857 return false;
858 }
859
860 // AKM: now the not so nice part find the best match position option
861 // Takes the error into account to scale the importance of the measurement
862
863 if (numberOfSubclusters == 1) {
864 finalposition = allLocalPositions[0];
865 finalerrormatrix = allErrorMatrix[0];
866 return true;
867 }
868
869 // Get the track parameters local position
870 const Amg::Vector2D localpos = trackPar.localPosition();
871 // Use the track parameters cov to weight distcances
872 Amg::Vector2D localerr(0.01, 0.05);
873 if (trackPar.covariance()) {
874 localerr = Amg::Vector2D(std::sqrt((*trackPar.covariance())(0, 0)), std::sqrt((*trackPar.covariance())(1, 1)));
875 }
876
877 double minDistance(1e300);
878 int index(0);
879
880 for (unsigned int i(0); i < allLocalPositions.size(); ++i) {
881 double distance =
882 square(localpos[0] - allLocalPositions[i][0]) / localerr[0]
883 + square(localpos[1] - allLocalPositions[i][1]) / localerr[1];
884
885 if (distance < minDistance) {
886 index = i;
887 minDistance = distance;
888 }
889 }
890
891 finalposition = allLocalPositions[index];
892 finalerrormatrix = allErrorMatrix[index];
893 return true;
894}
#define M_PI
Scalar phi() const
phi method
Scalar theta() const
theta method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
This is an Identifier helper class for the Pixel subdetector.
double angle(const GeoTrf::Vector2D &a, const GeoTrf::Vector2D &b)
const double width
#define x
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
const ServiceHandle< StoreGateSvc > & detStore() const
This is a "hash" representation of an Identifier.
int readoutSide() const
ReadoutSide.
Class used to describe the design of a module (diode segmentation and readout scheme)
SiLocalPosition positionFromColumnRow(const int column, const int row) const
Given row and column index of a diode, return position of diode center ALTERNATIVE/PREFERED way is to...
Class to hold geometrical description of a silicon detector element.
virtual const SiDetectorDesign & design() const override final
access to the local description (inline):
Class to represent a position in the natural frame of a silicon sensor, for Pixel and SCT For Pixel: ...
double xPhi() const
position along phi direction:
double xEta() const
position along eta direction:
virtual const Amg::Vector3D & normal() const override final
Get reconstruction local normal axes in global frame.
virtual IdentifierHash identifyHash() const override final
identifier hash (inline)
HepGeom::Point3D< double > globalPosition(const HepGeom::Point3D< double > &localPos) const
transform a reconstruction local position into a global position (inline):
virtual Identifier identify() const override final
identifier of this detector element (inline)
const PixelID * m_pixelid
Flag controlling how module distortions are taken into account:
SG::ReadCondHandleKey< RIO_OnTrackErrorScaling > m_pixelErrorScalingKey
BooleanProperty m_disableDistortions
toolhandle for central error scaling flag storing if errors need scaling or should be kept nominal
InDet::PixelClusterOnTrack * correctDefault(const Trk::PrepRawData &, const Trk::TrackParameters &, const EventContext &ctx) const
The correct method produces a PixelClusterOnTrack using the measured PixelCluster and the track predi...
SG::ReadHandleKey< InDet::PixelGangedClusterAmbiguities > m_splitClusterMapKey
bool getErrorsDefaultAmbi(const InDet::PixelCluster *, const Trk::TrackParameters &, Amg::Vector2D &, Amg::MatrixX &) const
virtual InDet::PixelClusterOnTrack * correct(const Trk::PrepRawData &, const Trk::TrackParameters &, const EventContext &ctx=Gaudi::Hive::currentContext()) const override
produces a PixelClusterOnTrack (object factory!).
SG::ReadCondHandleKey< PixelCalib::PixelOfflineCalibData > m_clusterErrorKey
virtual StatusCode finalize() override
AlgTool termination.
InDet::PixelClusterOnTrack * correctNN(const Trk::PrepRawData &, const Trk::TrackParameters &, const EventContext &ctx) const
virtual StatusCode initialize() override
AlgTool initialisation.
SG::ReadCondHandleKey< PixelDistortionData > m_distortionKey
const Trk::ClusterSplitProbabilityContainer::ProbabilityInfo & getClusterSplittingProbability(const InDet::PixelCluster *pix) const
ToolHandle< ISiLorentzAngleTool > m_lorentzAngleTool
PixelClusterOnTrackTool(const std::string &, const std::string &, const IInterface *)
AlgTool constructor.
BooleanProperty m_usingTIDE_Ambi
Enable different treatment of cluster errors based on NN information (do only if TIDE ambi is run)
BooleanProperty m_applyNNcorrectionProperty
Enable NN based calibration (do only if NN calibration is applied)
SG::ReadHandleKey< Trk::ClusterSplitProbabilityContainer > m_clusterSplitProbContainer
bool getErrorsTIDE_Ambi(const InDet::PixelCluster *, const Trk::TrackParameters &, Amg::Vector2D &, Amg::MatrixX &) const
ServiceHandle< IIBLParameterSvc > m_IBLParameterSvc
ToolHandle< NnClusterizationFactory > m_NnClusterizationFactory
NN clusterizationi factory for NN based positions and errors.
void correctBow(const Identifier &, Amg::Vector2D &locpos, const double tanphi, const double taneta, const EventContext &ctx) const
Specific class to represent the pixel measurements.
const Amg::Vector3D & globalPosition() const
return global position reference
bool gangedPixel() const
return the flag of this cluster containing a gangedPixel
virtual const InDetDD::SiDetectorElement * detectorElement() const override final
return the detector element corresponding to this PRD The pointer will be zero if the det el is not d...
double eta() const
Access method for pseudorapidity - from momentum.
const Amg::Vector3D & momentum() const
Access method for the momentum.
virtual constexpr SurfaceType surfaceType() const override=0
Returns the Surface Type enum for the surface used to define the derived class.
virtual constexpr ParametersType type() const override=0
Return the ParametersType enum.
virtual const Surface & associatedSurface() const override=0
Access to the Surface associated to the Parameters.
double pT() const
Access method for transverse momentum.
Amg::Vector2D localPosition() const
Access method for the local coordinates, local parameter definitions differ for each surface type.
virtual const TrkDetElementBase * detectorElement() const =0
return the detector element corresponding to this PRD The pointer will be zero if the det el is not d...
const Amg::Vector2D & localPosition() const
return the local position reference
virtual bool type(PrepRawDataType type) const
Interface method checking the type.
const Amg::MatrixX & localCovariance() const
return const ref to the error matrix
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Eigen::Matrix< double, 2, 1 > Vector2D
Eigen::Matrix< double, 3, 1 > Vector3D
ParametersBase< TrackParametersDim, Charged > TrackParameters
const T_res * ErrorScalingCast(const T_src *src)
Definition index.py:1