ATLAS Offline Software
Loading...
Searching...
No Matches
ITkPixelClusterOnTrackTool.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// Implementation file for class ITk::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
21
22#include <cmath>
24
25//clustermap is most likely to be removed at later date
26#define __clustermap
27
28
30// Constructor
32
33namespace
34{
35 // using x*x might be quicker than pow(x,2), depends on compiler optimisation
36 inline double
37 square(const double x) {
38 return x * x;
39 }
40
41 double
42 distance(const std::vector<Amg::Vector2D> &vectorOfPositions,
43 const std::vector<Amg::Vector2D> &allLocalPositions,
44 const std::vector<Amg::MatrixX> &allErrorMatrix,
45 const int i, const int j, const int k) {
46 return
47 square(vectorOfPositions[i][0] - allLocalPositions[0][0]) / allErrorMatrix[0](0, 0) +
48 square(vectorOfPositions[j][0] - allLocalPositions[1][0]) / allErrorMatrix[1](0, 0) +
49 square(vectorOfPositions[k][0] - allLocalPositions[2][0]) / allErrorMatrix[2](0, 0) +
50 square(vectorOfPositions[i][1] - allLocalPositions[0][1]) / allErrorMatrix[0](1, 1) +
51 square(vectorOfPositions[j][1] - allLocalPositions[1][1]) / allErrorMatrix[1](1, 1) +
52 square(vectorOfPositions[k][1] - allLocalPositions[2][1]) / allErrorMatrix[2](1, 1);
53 }
54}
55
56namespace ITk
57{
58
60 (const std::string &t, const std::string &n, const IInterface *p) :
61 ::AthAlgTool(t, n, p)
62{
63 declareInterface<IRIO_OnTrackCreator>(this);
64}
65
66
68// Initialisation
70
71StatusCode
73
74 ATH_MSG_DEBUG(name() << " initialize()");
75
77 ATH_MSG_DEBUG("Error strategy is" << m_errorStrategy);
78
80
81 ATH_CHECK(m_clusterITkErrorKey.initialize());
82
84
85 // get the error scaling tool
87 if (!m_pixelErrorScalingKey.key().empty()) ATH_MSG_DEBUG("Detected need for scaling Pixel errors.");
88
89 ATH_CHECK (detStore()->retrieve(m_pixelid, "PixelID"));
90
93 ATH_CHECK(m_NnClusterizationFactory.retrieve( DisableTool{!m_applyNNcorrection} ));
94
95 ATH_CHECK(m_lorentzAngleTool.retrieve());
96 return StatusCode::SUCCESS;
97}
98
99
101// Trk::SiClusterOnTrack production
103
104
107 (const Trk::PrepRawData &rio, const Trk::TrackParameters &trackPar, const EventContext& /*ctx*/) const {
108
109 if (not m_applyNNcorrection){
110 return correctDefault(rio, trackPar);
111 }else {
112 if (m_errorStrategy == 0 || m_errorStrategy == 1) {
113 // version from Giacinto
115 return nullptr;
116 }
117 // if we try broad errors, get Pixel Cluster to test if it is split
118 const InDet::PixelCluster *pix = nullptr;
120 pix = static_cast<const InDet::PixelCluster *>(&rio);
121 }
122 if (!pix) {
123 return nullptr;
124 }
126 if (splitProb.isSplit()) {
127 return correctNN(rio, trackPar);
128 } else {
129 return correctDefault(rio, trackPar);
130 }
131 } else {
132 return correctNN(rio, trackPar);
133 }
134 }
135}
136
137
138
144 (const Trk::PrepRawData &rio, const Trk::TrackParameters &trackPar) const {
145 using CLHEP::micrometer;
146
147
148 const double TOPHAT_SIGMA = 1. / std::sqrt(12.);
149 const InDet::PixelCluster *pix = nullptr;
151 pix = static_cast<const InDet::PixelCluster *>(&rio);
152 }
153 else{
154 return nullptr;
155 }
156
157 ATH_MSG_VERBOSE("Correct called with Error strategy " << m_errorStrategy);
158
159 // PixelClusterOnTrack production
160 //
162 Amg::Vector3D glob(pix->globalPosition());
163
164
165 // Get pointer to detector element
166 const InDetDD::SiDetectorElement *element = pix->detectorElement();
167 if (!element) {
168 return nullptr;
169 }
170 IdentifierHash idHash = element->identifyHash();
171
172 double errphi = -1;
173 double erreta = -1;
174
175 if (pix->rdoList().empty()) {
176 ATH_MSG_WARNING("Pixel RDO-list size is 0, check integrity of pixel clusters! stop ROT creation.");
177 return nullptr;
178 } else {
179 const InDetDD::PixelModuleDesign *design =
180 dynamic_cast<const InDetDD::PixelModuleDesign *>(&element->design());
181
182 // get candidate track angle in module local frame
183 const Amg::Vector3D& my_track = trackPar.momentum();
184 const Amg::Vector3D& my_normal = element->normal();
185 const Amg::Vector3D& my_phiax = element->phiAxis();
186 const Amg::Vector3D& my_etaax = element->etaAxis();
187 float trkphicomp = my_track.dot(my_phiax);
188 float trketacomp = my_track.dot(my_etaax);
189 float trknormcomp = my_track.dot(my_normal);
190 double bowphi = std::atan2(trkphicomp, trknormcomp);
191 double boweta = std::atan2(trketacomp, trknormcomp);
192
193 float tanl = m_lorentzAngleTool->getTanLorentzAngle(idHash, Gaudi::Hive::currentContext());
194 int readoutside = element->design().readoutSide();
195
196 // map the angles of inward-going tracks onto [-PI/2, PI/2]
197 if (bowphi > M_PI *0.5) {
198 bowphi -= M_PI;
199 }
200 if (bowphi < -M_PI *0.5) {
201 bowphi += M_PI;
202 }
203 // finally, subtract the Lorentz angle effect
204 // the readoutside term is needed because of a bug in old
205 // geometry versions (CSC-01-* and CSC-02-*)
206 double angle = std::atan(std::tan(bowphi) - readoutside * tanl);
207
208 // settle the sign/pi periodicity issues
209 double thetaloc = -999.;
210 if (boweta > -0.5 * M_PI && boweta < M_PI / 2.) { //M_PI_2 in cmath
211 thetaloc = M_PI_2 - boweta;
212 }else if (boweta > M_PI_2 && boweta < M_PI) {
213 thetaloc = 1.5 * M_PI - boweta;
214 } else { // 3rd quadrant
215 thetaloc = -M_PI_2 - boweta;
216 }
217 double etaloc = -1 * log(tan(thetaloc * 0.5));
218
219 // try to understand...
220 const Identifier element_id = element->identify();
221 int PixEtaModule = m_pixelid->eta_module(element_id);
222 int PixPhiModule = m_pixelid->phi_module(element_id);
223 double PixTrkPt = trackPar.pT();
224 double PixTrkEta = trackPar.eta();
225 ATH_MSG_VERBOSE("tanl = " << tanl << " readout side is " << readoutside <<
226 " module " << PixEtaModule << " " << PixPhiModule <<
227 " track pt, eta = " << PixTrkPt << " " << PixTrkEta <<
228 " track momentum phi, norm = " << trkphicomp << " " <<
229 trknormcomp << " bowphi = " << bowphi << " angle = " << angle);
230
231 float omegaphi = pix->omegax();
232 float omegaeta = pix->omegay();
233 double localphi = -9999.;
234 double localeta = -9999.;
235
236 const std::vector<Identifier> & rdos = pix->rdoList();
237 InDetDD::SiLocalPosition meanpos(0, 0, 0);
238 int rowmin = 9999;
239 int rowmax = -9999;
240 int colmin = 9999;
241 int colmax = -9999;
242 for (const auto & rId:rdos) {
243 const int row = m_pixelid->phi_index(rId);
244 const int col = m_pixelid->eta_index(rId);
245 rowmin = std::min(rowmin, row);
246 rowmax = std::max(rowmax,row);
247 colmin = std::min(colmin, col);
248 colmax = std::max(colmax, col);
249 meanpos += design->positionFromColumnRow(col, row);
250 }
251 meanpos = meanpos / rdos.size();
253 design->positionFromColumnRow(colmin, rowmin);
255 design->positionFromColumnRow(colmax, rowmin);
257 design->positionFromColumnRow(colmin, rowmax);
259 design->positionFromColumnRow(colmax, rowmax);
260
261 InDetDD::SiLocalPosition centroid = 0.25 * (pos1 + pos2 + pos3 + pos4);
262 double shift = m_lorentzAngleTool->getLorentzShift(idHash, Gaudi::Hive::currentContext());
263 int nrows = rowmax - rowmin + 1;
264 int ncol = colmax - colmin + 1;
265
266 // TOT interpolation for collision data
268
269 if (m_positionStrategy > 0 && omegaphi > -0.5 && omegaeta > -0.5) {
270 localphi = centroid.xPhi() + shift;
271 localeta = centroid.xEta();
272
273 std::pair<double,double> delta = offlineITkCalibDataHandle->getClusterErrorData()->getDelta(idHash,nrows,angle,ncol,etaloc);
274 double delta_phi = nrows != 1 ? delta.first : 0.;
275 double delta_eta = ncol != 1 ? delta.second : 0.;
276 localphi += delta_phi*(omegaphi-0.5);
277 localeta += delta_eta*(omegaeta-0.5);
278 }
279 // digital
280 else {
281 localphi = meanpos.xPhi() + shift;
282 localeta = meanpos.xEta();
283 }
284
285 const InDet::SiWidth& width = pix->width();
286
287 // Error strategies
288
289 // For very shallow tracks the cluster can easily break as
290 // the average charge per pixel is of the order of the threshold
291 // In this case, an error equal to the geometrical projection
292 // of the track path in silicon onto the module surface seems
293 // appropriate
294 if (std::abs(angle) > 1) {
295 errphi = 250 * micrometer * std::tan(std::abs(angle)) * TOPHAT_SIGMA;
296 erreta = width.z() > 250 * micrometer * std::tan(std::abs(boweta)) ?
297 width.z() * TOPHAT_SIGMA : 250 * micrometer * std::tan(std::abs(boweta)) * TOPHAT_SIGMA;
298 ATH_MSG_VERBOSE("Shallow track with tanl = " << tanl << " bowphi = " <<
299 bowphi << " angle = " << angle << " width.z = " << width.z() <<
300 " errphi = " << errphi << " erreta = " << erreta);
301 }else if (m_errorStrategy == 0) {
302 errphi = width.phiR() * TOPHAT_SIGMA;
303 erreta = width.z() * TOPHAT_SIGMA;
304 }else if (m_errorStrategy == 1) {
305 errphi = (width.phiR() / nrows) * TOPHAT_SIGMA;
306 erreta = (width.z() / ncol) * TOPHAT_SIGMA;
307 }else if (m_errorStrategy == 2) {
308 std::pair<double,double> delta_err = offlineITkCalibDataHandle->getClusterErrorData()->getDeltaError(idHash);
309 errphi = nrows != 1 ? delta_err.first : width.phiR()*TOPHAT_SIGMA;
310 erreta = ncol != 1 ? delta_err.second : width.z()*TOPHAT_SIGMA;
311 }
312
313 Amg::Vector2D locpos = Amg::Vector2D(localphi, localeta);
314 locpar = Trk::LocalParameters(locpos);
315 centroid = InDetDD::SiLocalPosition(localeta, localphi, 0.);
316 glob = element->globalPosition(centroid);
317 }
318
319 // Error matrix production
320
321 Amg::MatrixX cov = pix->localCovariance();
322
323 // corrected phi error
324 if (errphi > 0) {
325 cov(0, 0) = errphi * errphi;
326 }
327 if (erreta > 0) {
328 cov(1, 1) = erreta * erreta;
329 }
330
331 ATH_MSG_VERBOSE(" errphi = " << errphi << " erreta = " << erreta);
332
333 // create new copy of error matrix
334 if (!m_pixelErrorScalingKey.key().empty()) {
337 ->getScaledCovariance(std::move(cov), *m_pixelid,
338 element->identify());
339 }
340 bool isbroad = m_errorStrategy == 0;
341 return new InDet::PixelClusterOnTrack(pix, std::move(locpar),
342 std::move(cov),
343 idHash, glob, pix->gangedPixel(), isbroad);
344}
345
346
349 (const Trk::PrepRawData &rio, const Trk::TrackParameters &trackPar,
350 const ITk::PixelClusterStrategy strategy) const {
351 int initial_errorStrategy;
353
354 switch (strategy) {
355 case PixelClusterStrategy::OUTLIER: // if cluster is outlier, increase errors
357 initial_errorStrategy = m_errorStrategy;
358 m_errorStrategy = 0; // error as size of cluster /sqrt(12)
359 newROT = correct(rio, trackPar);
360 m_errorStrategy = initial_errorStrategy;
361 return newROT;
362
363 default:
364 return correct(rio, trackPar);
365 }
366}
367
368// GP: NEW correct() method in case of NN based calibration */
371 (const Trk::PrepRawData &rio,
372 const Trk::TrackParameters &trackPar) const {
373
374 const InDet::PixelCluster *pixelPrepCluster = nullptr;
376 pixelPrepCluster = static_cast<const InDet::PixelCluster *>(&rio);
377 }
378
379 if (pixelPrepCluster == nullptr) {
380 ATH_MSG_WARNING("This is not a pixel cluster, return 0.");
381 return nullptr;
382 }
383
384 const InDetDD::SiDetectorElement *element = pixelPrepCluster->detectorElement();
385 if (!element) {
386 ATH_MSG_WARNING("Cannot access detector element. Aborting cluster correction...");
387 return nullptr;
388 }
389
390 IdentifierHash idHash = element->identifyHash();
391
393 Amg::Vector3D glob(pixelPrepCluster->globalPosition());
394
395 Amg::Vector2D locpos = pixelPrepCluster->localPosition();
397 Amg::MatrixX cov = pixelPrepCluster->localCovariance();
398
399 return new InDet::PixelClusterOnTrack(pixelPrepCluster, std::move(locpar),
400 std::move(cov), idHash, glob,
401 pixelPrepCluster->gangedPixel(), false);
402 }
403
404
405
406 Amg::Vector2D finalposition;
407 Amg::MatrixX finalerrormatrix;
408
409 if (m_usingTIDE_Ambi) {
410 if (!getErrorsTIDE_Ambi(pixelPrepCluster, trackPar, finalposition, finalerrormatrix)) {
411 return correctDefault(rio, trackPar);
412 }
413 }else {
414 if (!getErrorsDefaultAmbi(pixelPrepCluster, trackPar, finalposition, finalerrormatrix)) {
415 return correctDefault(rio, trackPar);
416 }
417 }
418
419 ATH_MSG_DEBUG( " Old position x: " << pixelPrepCluster->localPosition()[0]
420 << " +/- " << std::sqrt(pixelPrepCluster->localCovariance()(0, 0))
421 << " y: " << pixelPrepCluster->localPosition()[1]
422 << " +/- " << std::sqrt(pixelPrepCluster->localCovariance()(1, 1)) <<"\n"
423 << " Final position x: " << finalposition[0]
424 << " +/- " << std::sqrt(finalerrormatrix(0, 0))
425 << " y: " << finalposition[1] << " +/- "
426 <<std::sqrt(finalerrormatrix(1, 1)) );
427
428 Amg::MatrixX cov = finalerrormatrix;
429 // create new copy of error matrix
430 if (!m_pixelErrorScalingKey.key().empty()) {
433 ->getScaledCovariance(std::move(cov), *m_pixelid,
434 element->identify());
435 }
436
437 InDetDD::SiLocalPosition centroid = InDetDD::SiLocalPosition(finalposition[1],
438 finalposition[0],
439 0);
440 Trk::LocalParameters locpar = Trk::LocalParameters(finalposition);
441
442 const Amg::Vector3D &glob = element->globalPosition(centroid);
443
444
445 return new InDet::PixelClusterOnTrack(pixelPrepCluster,
446 std::move(locpar),
447 std::move(cov), idHash,
448 glob,
449 pixelPrepCluster->gangedPixel(),
450 false);
451}
452
453bool
455 const Trk::TrackParameters &trackPar,
456 Amg::Vector2D &finalposition,
457 Amg::MatrixX &finalerrormatrix) const {
458 std::vector<Amg::Vector2D> vectorOfPositions;
459 int numberOfSubclusters = 1;
460 vectorOfPositions.push_back(pixelPrepCluster->localPosition());
461
464 InDet::PixelGangedClusterAmbiguities::const_iterator mapBegin = splitClusterMap->begin();
465 InDet::PixelGangedClusterAmbiguities::const_iterator mapEnd = splitClusterMap->end();
466 for (InDet::PixelGangedClusterAmbiguities::const_iterator mapIter = mapBegin; mapIter != mapEnd; ++mapIter) {
467 const InDet::SiCluster *first = (*mapIter).first;
468 const InDet::SiCluster *second = (*mapIter).second;
469 if (first == pixelPrepCluster && second != pixelPrepCluster) {
470 ATH_MSG_DEBUG("Found additional split cluster in ambiguity map (+=1).");
471 numberOfSubclusters += 1;
472 const InDet::SiCluster *otherOne = second;
473 const InDet::PixelCluster *pixelAddCluster = nullptr;
475 pixelAddCluster = static_cast<const InDet::PixelCluster *>(otherOne);
476 }
477 if (pixelAddCluster == nullptr) {
478 ATH_MSG_WARNING("Pixel ambiguity map has empty pixel cluster. Please DEBUG!");
479 continue;
480 }
481 vectorOfPositions.push_back(pixelAddCluster->localPosition());
482
483 ATH_MSG_DEBUG( "Found one more pixel cluster. Position x: "
484 << pixelAddCluster->localPosition()[0] << "y: " << pixelAddCluster->localPosition()[1]);
485 }// find relevant element of map
486 }// iterate over map
487 }
488
489 // now you have numberOfSubclusters and the vectorOfPositions (Amg::Vector2D)
490
491 if (trackPar.surfaceType() != Trk::SurfaceType::Plane ||
492 trackPar.type() != Trk::AtaSurface) {
494 "Parameters are not at a plane ! Aborting cluster correction... ");
495 return false;
496 }
497
498 std::vector<Amg::Vector2D> allLocalPositions;
499 std::vector<Amg::MatrixX> allErrorMatrix;
500 allLocalPositions =
501 m_NnClusterizationFactory->estimatePositions(*pixelPrepCluster,
502 trackPar.associatedSurface(),
503 trackPar,
504 allErrorMatrix,
505 numberOfSubclusters);
506
507 if (allLocalPositions.empty()) {
508 ATH_MSG_DEBUG( " Cluster cannot be treated by NN. Giving back to default clusterization " );
509
510 return false;
511 }
512
513 if (allLocalPositions.size() != size_t(numberOfSubclusters)) {
514 ATH_MSG_WARNING( "Returned position vector size " << allLocalPositions.size() <<
515 " not according to expected number of subclusters: " << numberOfSubclusters << ". Abort cluster correction..." );
516 return false;
517 }
518
519
520 // GP: now the not so nice part of matching the new result with the old one...
521 // Takes the error into account to improve the matching
522
523 if (numberOfSubclusters == 1) {
524 finalposition = allLocalPositions[0];
525 finalerrormatrix = allErrorMatrix[0];
526 }
527
528 else if (numberOfSubclusters == 2) {
529 double distancesq1 =
530 square(vectorOfPositions[0][0] - allLocalPositions[0][0]) / allErrorMatrix[0](0, 0) +
531 square(vectorOfPositions[1][0] - allLocalPositions[1][0]) / allErrorMatrix[1](0, 0) +
532 square(vectorOfPositions[0][1] - allLocalPositions[0][1]) / allErrorMatrix[0](1, 1) +
533 square(vectorOfPositions[1][1] - allLocalPositions[1][1]) / allErrorMatrix[1](1, 1);
534
535 double distancesq2 =
536 square(vectorOfPositions[1][0] - allLocalPositions[0][0]) / allErrorMatrix[0](0, 0) +
537 square(vectorOfPositions[0][0] - allLocalPositions[1][0]) / allErrorMatrix[1](0, 0) +
538 square(vectorOfPositions[1][1] - allLocalPositions[0][1]) / allErrorMatrix[0](1, 1) +
539 square(vectorOfPositions[0][1] - allLocalPositions[1][1]) / allErrorMatrix[1](1, 1);
540
542 " Old pix (1) x: " << vectorOfPositions[0][0] << " y: " << vectorOfPositions[0][1] << "\n"
543 << " Old pix (2) x: " << vectorOfPositions[1][0] << " y: " << vectorOfPositions[1][1] << "\n"
544 << " Pix (1) x: " << allLocalPositions[0][0] << " +/- " << std::sqrt(allErrorMatrix[0](0, 0))
545 << " y: " << allLocalPositions[0][1] << " +/- " << std::sqrt(allErrorMatrix[0](1, 1)) <<"\n"
546 << " Pix (2) x: " << allLocalPositions[1][0] << " +/- " << std::sqrt(allErrorMatrix[1](0, 0))
547 << " y: " << allLocalPositions[1][1] << " +/- " << std::sqrt(allErrorMatrix[1](1, 1)) << "\n"
548 << " Old (1) new (1) dist: " << std::sqrt(distancesq1) << " Old (1) new (2) " << std::sqrt(distancesq2) );
549
550
551 if (distancesq1 < distancesq2) {
552 finalposition = allLocalPositions[0];
553 finalerrormatrix = allErrorMatrix[0];
554 }else {
555 finalposition = allLocalPositions[1];
556 finalerrormatrix = allErrorMatrix[1];
557 }
558 }
559
560
561 else if (numberOfSubclusters == 3) {
562 double distances[6];
563
564 distances[0] = distance(vectorOfPositions, allLocalPositions, allErrorMatrix, 0, 1, 2);
565 distances[1] = distance(vectorOfPositions, allLocalPositions, allErrorMatrix, 0, 2, 1);
566 distances[2] = distance(vectorOfPositions, allLocalPositions, allErrorMatrix, 1, 0, 2);
567 distances[3] = distance(vectorOfPositions, allLocalPositions, allErrorMatrix, 1, 2, 0);
568 distances[4] = distance(vectorOfPositions, allLocalPositions, allErrorMatrix, 2, 0, 1);
569 distances[5] = distance(vectorOfPositions, allLocalPositions, allErrorMatrix, 2, 1, 0);
570
571 int smallestDistanceIndex = -10;
572 double minDistance = 1e10;
573
574 for (int i = 0; i < 6; i++) {
575 ATH_MSG_VERBOSE(" distance n.: " << i << " distance is: " << distances[i]);
576
577 if (distances[i] < minDistance) {
578 minDistance = distances[i];
579 smallestDistanceIndex = i;
580 }
581 }
582
583 ATH_MSG_DEBUG(" The minimum distance is : " << minDistance << " for index: " << smallestDistanceIndex);
584
585 if (smallestDistanceIndex == 0 || smallestDistanceIndex == 1) {
586 finalposition = allLocalPositions[0];
587 finalerrormatrix = allErrorMatrix[0];
588 }
589 if (smallestDistanceIndex == 2 || smallestDistanceIndex == 4) {
590 finalposition = allLocalPositions[1];
591 finalerrormatrix = allErrorMatrix[1];
592 }
593 if (smallestDistanceIndex == 3 || smallestDistanceIndex == 5) {
594 finalposition = allLocalPositions[2];
595 finalerrormatrix = allErrorMatrix[2];
596 }
597 }
598 return true;
599}
600
601bool
603 const Trk::TrackParameters &trackPar,
604 Amg::Vector2D &finalposition,
605 Amg::MatrixX &finalerrormatrix) const {
607 std::vector<Amg::Vector2D> vectorOfPositions;
608 int numberOfSubclusters = 1;
611 numberOfSubclusters = 1 + splitClusterMap->count(pixelPrepCluster);
612
613 if (splitClusterMap->count(pixelPrepCluster) == 0 && splitProb.isSplit()) {
614 numberOfSubclusters = 2;
615 }
616 if (splitClusterMap->count(pixelPrepCluster) != 0 && !splitProb.isSplit()) {
617 numberOfSubclusters = 1;
618 }
619 }
620
621 // now you have numberOfSubclusters and the vectorOfPositions (Amg::Vector2D)
622 if (trackPar.surfaceType() != Trk::SurfaceType::Plane ||
623 trackPar.type() != Trk::AtaSurface) {
624 ATH_MSG_WARNING("Parameters are not at a plane surface ! Aborting cluster "
625 "correction... ");
626 return false;
627 }
628
629 std::vector<Amg::Vector2D> allLocalPositions;
630 std::vector<Amg::MatrixX> allErrorMatrix;
631 allLocalPositions = m_NnClusterizationFactory->estimatePositions(
632 *pixelPrepCluster,
633 trackPar.associatedSurface(),
634 trackPar,
635 allErrorMatrix,
636 numberOfSubclusters);
637
638 if (allLocalPositions.empty()) {
640 " Cluster cannot be treated by NN. Giving back to default clusterization, too big: " <<
641 splitProb.isTooBigToBeSplit());
642 return false;
643 }
644
645 if (allLocalPositions.size() != size_t(numberOfSubclusters)) {
647 "Returned position vector size " << allLocalPositions.size() << " not according to expected number of subclusters: " << numberOfSubclusters <<
648 ". Abort cluster correction...");
649 return false;
650 }
651
652 // AKM: now the not so nice part find the best match position option
653 // Takes the error into account to scale the importance of the measurement
654
655 if (numberOfSubclusters == 1) {
656 finalposition = allLocalPositions[0];
657 finalerrormatrix = allErrorMatrix[0];
658 return true;
659 }
660
661 // Get the track parameters local position
662 const Amg::Vector2D localpos = trackPar.localPosition();
663 // Use the track parameters cov to weight distcances
664 Amg::Vector2D localerr(0.01, 0.05);
665 if (trackPar.covariance()) {
666 localerr = Amg::Vector2D(std::sqrt((*trackPar.covariance())(0, 0)), std::sqrt((*trackPar.covariance())(1, 1)));
667 }
668
669 double minDistance(1e300);
670 int index(0);
671
672 for (unsigned int i(0); i < allLocalPositions.size(); ++i) {
673 double distance =
674 square(localpos[0] - allLocalPositions[i][0]) / localerr[0]
675 + square(localpos[1] - allLocalPositions[i][1]) / localerr[1];
676
677 if (distance < minDistance) {
678 index = i;
679 minDistance = distance;
680 }
681 }
682
683 finalposition = allLocalPositions[index];
684 finalerrormatrix = allErrorMatrix[index];
685 return true;
686}
687
688} // namespace ITk
#define M_PI
#define ATH_CHECK
Evaluate an expression and check for errors.
#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
bool getErrorsDefaultAmbi(const InDet::PixelCluster *, const Trk::TrackParameters &, Amg::Vector2D &, Amg::MatrixX &) const
SG::ReadHandleKey< InDet::PixelGangedClusterAmbiguities > m_splitClusterMapKey
SG::ReadHandleKey< Trk::ClusterSplitProbabilityContainer > m_clusterSplitProbContainer
const Trk::ClusterSplitProbabilityContainer::ProbabilityInfo & getClusterSplittingProbability(const InDet::PixelCluster *pix) const
IntegerProperty m_positionStrategy
toolhandle for central error scaling flag storing if errors need scaling or should be kept nominal
BooleanProperty m_usingTIDE_Ambi
Enable different treatment of cluster errors based on NN information (do only if TIDE ambi is run)
InDet::PixelClusterOnTrack * correctNN(const Trk::PrepRawData &, const Trk::TrackParameters &) const
bool m_applyNNcorrection
Enable NN based calibration (do only if NN calibration is applied)
virtual StatusCode initialize() override
AlgTool initialisation.
InDet::PixelClusterOnTrack * correctDefault(const Trk::PrepRawData &, const Trk::TrackParameters &) const
The correct method produces a PixelClusterOnTrack using the measured PixelCluster and the track predi...
bool getErrorsTIDE_Ambi(const InDet::PixelCluster *, const Trk::TrackParameters &, Amg::Vector2D &, Amg::MatrixX &) const
const PixelID * m_pixelid
Flag controlling how module distortions are taken into account:
SG::ReadCondHandleKey< ITk::PixelOfflineCalibData > m_clusterITkErrorKey
ToolHandle< ISiLorentzAngleTool > m_lorentzAngleTool
ToolHandle< InDet::NnClusterizationFactory > m_NnClusterizationFactory
NN clusterizationi factory for NN based positions and errors.
SG::ReadCondHandleKey< RIO_OnTrackErrorScaling > m_pixelErrorScalingKey
PixelClusterOnTrackTool(const std::string &, const std::string &, const IInterface *)
AlgTool constructor.
virtual InDet::PixelClusterOnTrack * correct(const Trk::PrepRawData &, const Trk::TrackParameters &, const EventContext &ctx=Gaudi::Hive::currentContext()) const override
produces a PixelClusterOnTrack (object factory!).
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)
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.
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
PixelClusterStrategy
creates PixelClusterOnTrack objects allowing to calibrate cluster position and error using a given tr...
ParametersBase< TrackParametersDim, Charged > TrackParameters
const T_res * ErrorScalingCast(const T_src *src)
Definition index.py:1