ATLAS Offline Software
Loading...
Searching...
No Matches
ClusterMakerTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5//***************************************************************************
6//
7// Implementation for ClusterMaker
8//
9//****************************************************************************
10
11#include "CLHEP/Units/SystemOfUnits.h"
12#include "CLHEP/Matrix/SymMatrix.h"
13#include "GaudiKernel/ToolHandle.h"
14#include "GaudiKernel/ServiceHandle.h"
22
25
27
29
30#include <memory>
31
33using CLHEP::micrometer;
34
35namespace {
36
37inline double square(const double x){
38 return x*x;
39}
40constexpr double ONE_TWELFTH = 1./12.;
41
42// Some methods below can be parameterized on the pixel cluster type,
43// The following functions allow using a function parameter for common
44// operations.
45// [ Omegax = TOT1/(TOT1+TOT2), where TOT1 and TOT2 are the sum of the
46// charges of the first and last row of the cluster respectively
47// Omegay: similar definition with columns rather than rows ]
48InDet::PixelCluster newInDetpixelCluster(const Identifier& RDOId,
49 const Amg::Vector2D& locpos,
50 const Amg::Vector3D& globpos,
51 std::vector<Identifier>&& rdoList,
52 const int lvl1a,
53 std::vector<int>&& totList,
54 std::vector<float>&& chargeList,
55 const InDet::SiWidth& width,
56 const InDetDD::SiDetectorElement* detEl,
57 Amg::MatrixX&& locErrMat,
58 const float omegax,
59 const float omegay,
60 bool split,
61 float splitProb1,
62 float splitProb2)
63{
64 return InDet::PixelCluster(RDOId,
65 locpos,
66 globpos,
67 std::move(rdoList),
68 lvl1a,
69 std::move(totList),
70 std::move(chargeList),
71 width,
72 detEl,
73 std::move(locErrMat),
74 omegax,
75 omegay,
76 split,
77 splitProb1,
78 splitProb2);
79}
80
81// Function-like class to add an xAOD::PixelCluster to an
82// xAOD::PixelClusterContainer. This is needed because the
83// PixelCluster object needs an aux store for the setMeasurement call
84// to not crash
85class AddNewxAODpixelCluster {
86public:
87 explicit AddNewxAODpixelCluster(xAOD::PixelCluster& cluster)
88 : m_cluster(&cluster) {}
89
90 xAOD::PixelCluster* operator()(const Identifier& /*RDOId*/,
91 const Amg::Vector2D& locpos,
92 const Amg::Vector3D& globpos,
93 const std::vector<Identifier>& rdoList,
94 const int lvl1a,
95 const std::vector<int>& totList,
96 const std::vector<float>& chargeList,
97 const InDet::SiWidth& width,
98 const InDetDD::SiDetectorElement* detEl,
99 const Amg::MatrixX& locErrMat,
100 bool split,
101 float splitProb1,
102 float splitProb2) {
103 IdentifierHash idHash = detEl->identifyHash();
104
105 Eigen::Matrix<float,2,1> localPosition(locpos.x(), locpos.y());
106 Eigen::Matrix<float,2,2> localCovariance = Eigen::Matrix<float,2,2>::Zero();
107 localCovariance(0, 0) = locErrMat(0, 0);
108 localCovariance(1, 1) = locErrMat(1, 1);
109
110 m_cluster->setMeasurement<2>(idHash, localPosition, localCovariance);
111 m_cluster->setRDOlist(rdoList);
112 m_cluster->globalPosition() = globpos.cast<float>();
113 m_cluster->setToTlist(totList);
114 m_cluster->setTotalToT( xAOD::xAODInDetMeasurement::Utilities::computeTotalToT(totList) );
115 m_cluster->setChargelist(chargeList);
116 m_cluster->setTotalCharge( xAOD::xAODInDetMeasurement::Utilities::computeTotalCharge(chargeList) );
117 m_cluster->setLVL1A(lvl1a);
118 m_cluster->setChannelsInPhiEta(width.colRow()[0], width.colRow()[1]);
119 m_cluster->setWidthInEta(static_cast<float>(width.widthPhiRZ()[1]));
120 m_cluster->setIsSplit(split);
121 m_cluster->setSplitProbabilities(splitProb1, splitProb2);
122
123 return m_cluster;
124 }
125
126private:
127 xAOD::PixelCluster *m_cluster;
128};
129
130
131}
132
133
134namespace InDet {
135
136// using namespace Trk;
137
138// Constructor with parameters:
140 const std::string& n,
141 const IInterface* p) :
142 AthAlgTool(t,n,p)
143{
144 declareInterface<ClusterMakerTool>(this);
145}
146
147//================ Initialisation =============================================
148
150 // Code entered here will be executed once at program start.
151
152 ATH_MSG_DEBUG ( name() << " initialize()" );
153
154 if (not m_pixelLorentzAngleTool.empty()) {
156 } else {
157 m_pixelLorentzAngleTool.disable();
158 }
159 if (not m_sctLorentzAngleTool.empty()) {
161 } else {
162 m_sctLorentzAngleTool.disable();
163 }
164
165
166 return StatusCode::SUCCESS;
167}
168
169
170// Compute the pixel cluster global position, and the error associated
171// to the position.
172// Called by the pixel clustering tools
173//
174// Input parameters
175// - the cluster Identifier
176// - the position in local reference frame
177// - the list of identifiers of the Raw Data Objects belonging to the cluster
178// - the width of the cluster
179// - the module the cluster belongs to
180// - wheter the cluster contains ganged pixels
181// - the error strategy, currently
182// 0: cluster width/sqrt(12.)
183// 1: pixel pitch/sqrt(12.)
184// 2: parametrized as a function ofpseudorapidity and cluster size
185// (default)
186// 10: CTB parametrization (as a function of module and cluster size)
187// no magnetic field
188// - const reference to a PixelID helper class
189
190template <typename ClusterType, typename IdentifierVec, typename ToTList>
192 const Identifier& clusterID,
193 const Amg::Vector2D& localPos,
194 IdentifierVec&& rdoList,
195 const int lvl1a,
196 ToTList&& totList,
197 const SiWidth& width,
198 const InDetDD::SiDetectorElement* element,
199 bool ganged,
200 int errorStrategy,
201 const PixelID& pixelID,
202 bool split,
203 double splitProb1,
204 double splitProb2,
205 const PixelChargeCalibCondData *calibData,
206 const PixelOfflineCalibData *offlineCalibData,
207 const EventContext& ctx,
208 xAOD::PixelCluster* cluster) const{
209
210 ATH_MSG_VERBOSE("ClusterMakerTool called, number ");
211 if ( errorStrategy==2 && m_forceErrorStrategy1B ) errorStrategy=1;
212
213 // Fill vector of charges and compute charge balance
214 const InDetDD::PixelModuleDesign* design = (dynamic_cast<const InDetDD::PixelModuleDesign*>(&element->design()));
215 if (not design){
216 throw std::runtime_error( "Dynamic cast failed for design in ClusterMakerTool.cxx");
217 }
218 int rowMin = design->rows();
219 int rowMax = 0;
220 int colMin = design->columns();
221 int colMax = 0;
222 float qRowMin = 0; float qRowMax = 0;
223 float qColMin = 0; float qColMax = 0;
224 std::vector<float> chargeList;
225 int nRDO=rdoList.size();
226 if (calibData) {
227 chargeList.reserve(nRDO);
228 IdentifierHash moduleHash = element->identifyHash(); // wafer hash
229 for (int i=0; i<nRDO; i++) {
230 Identifier pixid=rdoList[i];
231 int ToT=totList[i];
232
233 float charge = ToT;
234 assert( element->identifyHash() == pixelID.wafer_hash(pixelID.wafer_id(pixid)));
235 std::array<InDetDD::PixelDiodeTree::CellIndexType,2> diode_idx
237 pixelID.eta_index(pixid));
238 InDetDD::PixelDiodeTree::DiodeProxy si_param ( design->diodeProxyFromIdx(diode_idx));
239 std::uint32_t feValue = design->getFE(si_param);
240 auto diode_type = design->getDiodeType(si_param);
242 && design->numberOfConnectedCells( design->readoutIdOfCell(InDetDD::SiCellId(diode_idx[0],diode_idx[1])))>1) {
244 }
245
246 charge = calibData->getCharge(diode_type, moduleHash, feValue, ToT);
247 if (design->getReadoutTechnology() != InDetDD::PixelReadoutTechnology::RD53 && (moduleHash<12 || moduleHash>2035)) {
248 charge = ToT/8.0*(8000.0-1200.0)+1200.0;
249 }
250 chargeList.push_back(charge);
251 }
252 }
253
254 for (int i=0; i<nRDO; i++) {
255 Identifier pixid=rdoList[i];
256 int ToT=totList[i];
257
258 float charge = ToT;
259 if (calibData) { charge=chargeList[i]; }
260
261 // std::cout << "tot, charge = " << ToT << " " << charge << std::endl;
262 int row = pixelID.phi_index(pixid);
263 int col = pixelID.eta_index(pixid);
264 if (row == rowMin) qRowMin += charge;
265 if (row < rowMin){
266 rowMin = row;
267 qRowMin = charge;
268 }
269
270 if (row == rowMax) qRowMax += charge;
271 if (row > rowMax){
272 rowMax = row;
273 qRowMax = charge;
274 }
275 if (col == colMin) qColMin += charge;
276 if (col < colMin){
277 colMin = col;
278 qColMin = charge;
279 }
280
281 if (col == colMax) qColMax += charge;
282 if (col > colMax){
283 colMax = col;
284 qColMax = charge;
285 }
286 }
287
288 Identifier newClusterID = pixelID.pixel_id(pixelID.wafer_id(clusterID),rowMin,colMin);
289 // Compute omega for charge interpolation correction (if required)
290 // Two pixels may have charge=0 (very rarely, hopefully)
291 float omegax = -1;
292 float omegay = -1;
293 if(qRowMin+qRowMax > 0) omegax = qRowMax/float(qRowMin+qRowMax);
294 if(qColMin+qColMax > 0) omegay = qColMax/float(qColMin+qColMax);
295
296 ATH_MSG_VERBOSE("omega = " << omegax << " " << omegay);
297
298// ask for Lorentz correction, get global position
299 double shift = m_pixelLorentzAngleTool->getLorentzShift(element->identifyHash(), ctx);
300 Amg::Vector2D locpos(localPos[Trk::locX]+shift, localPos[Trk::locY]);
301// find global position of element
302 const Amg::Transform3D& T = element->surface().transform();
303 double Ax[3] = {T(0,0),T(1,0),T(2,0)};
304 double Ay[3] = {T(0,1),T(1,1),T(2,1)};
305 double R [3] = {T(0,3),T(1,3),T(2,3)};
306
307 const Amg::Vector2D& M = locpos;
308 Amg::Vector3D globalPos(M[0]*Ax[0]+M[1]*Ay[0]+R[0],M[0]*Ax[1]+M[1]*Ay[1]+R[1],M[0]*Ax[2]+M[1]*Ay[2]+R[2]);
309
310 // error matrix
311 const Amg::Vector2D& colRow = width.colRow();// made ref to avoid
312 // unnecessary copy EJWM
313 auto errorMatrix = Amg::MatrixX(2,2);
314 errorMatrix.setIdentity();
315
316 // switches are more readable **OPT**
317 // actually they're slower as well (so I'm told) so perhaps
318 // this should be re-written at some point EJWM
319 double eta = std::abs(globalPos.eta());
320 double zPitch = width.z()/colRow.y();
321
322 const AtlasDetectorID* aid = element->getIdHelper();
323
325 throw std::runtime_error( "Wrong helper type in ClusterMakerTool.cxx.");
326 }
327 const PixelID* pid = static_cast<const PixelID*>(aid);
328 int layer = pid->layer_disk(clusterID);
329 int phimod = pid->phi_module(clusterID);
330 switch (errorStrategy){
331 case 0:
332 errorMatrix.fillSymmetric(0,0,square(width.phiR())*ONE_TWELFTH);
333 errorMatrix.fillSymmetric(1,1,square(width.z())*ONE_TWELFTH);
334 break;
335 case 1:
336 errorMatrix.fillSymmetric(0,0,square(width.phiR()/colRow.x())*ONE_TWELFTH);
337 errorMatrix.fillSymmetric(1,1,square(width.z()/colRow.y())*ONE_TWELFTH);
338 break;
339 case 2:
340 // use parameterization only if the cluster does not
341 // contain long pixels or ganged pixels
342 // Also require calibration service is available....
343 if (!ganged && zPitch>399*micrometer && zPitch<401*micrometer) {
344 if (offlineCalibData) {
345 if (element->isBarrel()) {
346 int ibin = offlineCalibData->getPixelClusterErrorData()->getBarrelBin(eta,int(colRow.y()),int(colRow.x()));
347 double phiError = offlineCalibData->getPixelClusterErrorData()->getPixelBarrelPhiError(ibin);
348 double etaError = offlineCalibData->getPixelClusterErrorData()->getPixelBarrelEtaError(ibin);
349 errorMatrix.fillSymmetric(0,0,pow(phiError,2));
350 errorMatrix.fillSymmetric(1,1,pow(etaError,2));
351 }
352 else {
353 int ibin = offlineCalibData->getPixelClusterErrorData()->getEndcapBin(int(colRow.y()),int(colRow.x()));
354 double phiError = offlineCalibData->getPixelClusterErrorData()->getPixelEndcapPhiError(ibin);
355 double etaError = offlineCalibData->getPixelClusterErrorData()->getPixelEndcapRError(ibin);
356 errorMatrix.fillSymmetric(0,0,square(phiError));
357 errorMatrix.fillSymmetric(1,1,square(etaError));
358 }
359 }
360 }else{// cluster with ganged and/or long pixels
361 errorMatrix.fillSymmetric(0,0,square(width.phiR()/colRow.x())*ONE_TWELFTH);
362 errorMatrix.fillSymmetric(1,1,square(zPitch)*ONE_TWELFTH);
363 }
364 break;
365
366 case 10:
367 errorMatrix.fillSymmetric(0,0,square( getPixelCTBPhiError(layer,phimod,int(colRow.x()))));
368 errorMatrix.fillSymmetric(1,1,square(width.z()/colRow.y())*ONE_TWELFTH);
369 break;
370
371 default:
372 errorMatrix.fillSymmetric(0,0,square(width.phiR()/colRow.x())*ONE_TWELFTH);
373 errorMatrix.fillSymmetric(1,1,square(width.z()/colRow.y())*ONE_TWELFTH);
374 break;
375 }
376
377 //1) We want to move always for the Trk::PixelCluster
378 //2) We forward the rdoList and chargeList that could have be passed
379 // different ways.
380 //
381 static_assert(std::is_same_v<ClusterType, PixelCluster> ||
382 std::is_same_v<std::remove_pointer_t<ClusterType>,
384 "Not an InDet::PixelCluster or xAOD::PixelCluster");
385 if constexpr (std::is_same<ClusterType, InDet::PixelCluster>::value) {
386 return newInDetpixelCluster(newClusterID, locpos,
387 globalPos,
388 std::forward<IdentifierVec>(rdoList),
389 lvl1a,
390 std::forward<ToTList>(totList),
391 std::move(chargeList),
392 width,
393 element,
394 std::move(errorMatrix),
395 omegax,
396 omegay,
397 split,
398 splitProb1,
399 splitProb2);
400 } else{
401 return AddNewxAODpixelCluster(*cluster)(newClusterID,
402 locpos,
403 globalPos,
404 std::forward<IdentifierVec>(rdoList),
405 lvl1a,
406 std::forward<ToTList>(totList),
407 chargeList,
408 width,
409 element,
410 errorMatrix,
411 split,
412 splitProb1,
413 splitProb2);
414 }
415}
416
418 const Identifier& clusterID,
419 const Amg::Vector2D& localPos,
420 std::vector<Identifier>&& rdoList,
421 const int lvl1a,
422 std::vector<int>&& totList,
423 const SiWidth& width,
424 const InDetDD::SiDetectorElement* element,
425 bool ganged,
426 int errorStrategy,
427 const PixelID& pixelID,
428 bool split,
429 double splitProb1,
430 double splitProb2,
431 const PixelChargeCalibCondData *calibData,
432 const PixelOfflineCalibData *offlineCalibData,
433 const EventContext& ctx) const
434{
436 clusterID,
437 localPos,
438 std::move(rdoList),
439 lvl1a,
440 std::move(totList),
441 width,
442 element,
443 ganged,
444 errorStrategy,
445 pixelID,
446 split,
447 splitProb1,
448 splitProb2,
449 calibData,
450 offlineCalibData,
451 ctx);
452}
453
455 xAOD::PixelCluster& cluster,
456 const Amg::Vector2D& localPos,
457 const std::vector<Identifier>& rdoList,
458 const int lvl1a,
459 const std::vector<int>& totList,
460 const SiWidth& width,
461 const InDetDD::SiDetectorElement* element,
462 bool ganged,
463 int errorStrategy,
464 const PixelID& pixelID,
465 bool split,
466 double splitProb1,
467 double splitProb2,
468 const PixelChargeCalibCondData *calibData,
469 const PixelOfflineCalibData *offlineCalibData,
470 const EventContext& ctx) const
471{
473 Identifier(),
474 localPos,
475 rdoList,
476 lvl1a,
477 totList,
478 width,
479 element,
480 ganged,
481 errorStrategy,
482 pixelID,
483 split,
484 splitProb1,
485 splitProb2,
486 calibData,
487 offlineCalibData,
488 ctx,
489 &cluster);
490}
491
492// Computes global position and errors for SCT cluster.
493// Called by SCT Clustering tools
494//
495// Input parameters
496// - the cluster Identifier
497// - the position in local reference frame
498// - the list of identifiers of the Raw Data Objects belonging to the cluster
499// - the width of the cluster
500// - the module the cluster belongs to
501// - the error strategy, currently
502// 0: Cluster Width/sqrt(12.)
503// 1: Set to a different values for one and two-strip clusters (def.)
504// The scale factors were derived by the study reported on 25th September 2006.
505// https://indico.cern.ch/event/430391/contributions/1066157/attachments/929942/1317007/SCTSoft_25Sept06_clusters.pdf
506
509 const Amg::Vector2D& localPos,
510 std::vector<Identifier>&& rdoList,
511 const SiWidth& width,
512 const InDetDD::SiDetectorElement* element,
513 int errorStrategy) const
514{
515
516 double shift =
517 m_sctLorentzAngleTool->getLorentzShift(element->identifyHash(), Gaudi::Hive::currentContext());
518 Amg::Vector2D locpos(localPos[Trk::locX] + shift, localPos[Trk::locY]);
519
520 // error matrix
521 const Amg::Vector2D& colRow = width.colRow(); // made ref to avoid
522 // unnecessary copy EJWM
523
524 auto errorMatrix = Amg::MatrixX(2,2);
525 errorMatrix.setIdentity();
526
527 // switches are more readable **OPT**
528 // actually they're slower as well (so I'm told) so perhaps
529 // this should be re-written at some point EJWM
530
531 switch (errorStrategy){
532 case 0:
533 errorMatrix.fillSymmetric(0,0,square(width.phiR())*ONE_TWELFTH);
534 errorMatrix.fillSymmetric(1,1,square(width.z())*ONE_TWELFTH);
535 break;
536 case 1:
537 // mat(1,1) = pow(width.phiR()/colRow.x(),2)/12;
538 // single strip - resolution close to pitch/sqrt(12)
539 // two-strip hits: better resolution, approx. 40% lower
540 if(colRow.x() == 1){
541 errorMatrix.fillSymmetric(0,0,square(1.05*width.phiR())*ONE_TWELFTH);
542 }
543 else if(colRow.x() == 2){
544 errorMatrix.fillSymmetric(0,0,square(0.27*width.phiR())*ONE_TWELFTH);
545 }
546 else{
547 errorMatrix.fillSymmetric(0,0,square(width.phiR())*ONE_TWELFTH);
548 }
549 errorMatrix.fillSymmetric(1,1,square(width.z()/colRow.y())*ONE_TWELFTH);
550 break;
551 default:
552 // single strip - resolution close to pitch/sqrt(12)
553 // two-strip hits: better resolution, approx. 40% lower
554 if(colRow.x() == 1){
555 errorMatrix.fillSymmetric(0,0,square(width.phiR())*ONE_TWELFTH);
556 }
557 else if(colRow.x() == 2){
558 errorMatrix.fillSymmetric(0,0,square(0.27*width.phiR())*ONE_TWELFTH);
559 }
560 else{
561 errorMatrix.fillSymmetric(0,0,square(width.phiR())*ONE_TWELFTH);
562 }
563 errorMatrix.fillSymmetric(1,1,square(width.z()/colRow.y())*ONE_TWELFTH);
564 break;
565 }
566
567 auto designShape = element->design().shape();
568 // rotation for endcap SCT
569 if(designShape == InDetDD::Trapezoid || designShape == InDetDD::Annulus) {
570 double sn = element->sinStereoLocal(localPos);
571 double sn2 = sn*sn;
572 double cs2 = 1.-sn2;
573 double w = element->phiPitch(localPos)/element->phiPitch();
574 double v0 = (errorMatrix)(0,0)*w*w;
575 double v1 = (errorMatrix)(1,1);
576 errorMatrix.fillSymmetric(0,0,cs2*v0+sn2*v1);
577 errorMatrix.fillSymmetric(0,1,sn*sqrt(cs2)*(v0-v1));
578 errorMatrix.fillSymmetric(1,1,sn2*v0+cs2*v1);
579 } //else if (designShape == InDetDD::PolarAnnulus) {// Polar rotation for endcap}
580
581 return SCT_Cluster(clusterID, locpos, std::move(rdoList), width, element, std::move(errorMatrix));
582
583}
584
585//---------------------------------------------------------------------------
586// CTB parameterization, B field off
588 int phiClusterSize) const{
589
590 double sigmaL0Phi1[3] = { 8.2*micrometer, 9.7*micrometer, 14.6*micrometer};
591 double sigmaL1Phi1[3] = {14.6*micrometer, 9.3*micrometer, 14.6*micrometer};
592 double sigmaL2Phi1[3] = {14.6*micrometer, 8.6*micrometer, 14.6*micrometer};
593 double sigmaL0Phi0[3] = {14.6*micrometer, 13.4*micrometer, 13.0*micrometer};
594 double sigmaL1Phi0[3] = {14.6*micrometer, 8.5*micrometer, 11.0*micrometer};
595 double sigmaL2Phi0[3] = {14.6*micrometer, 11.6*micrometer, 9.3*micrometer};
596
597 if(phiClusterSize > 3) return 14.6*micrometer;
598
599 if(layer == 0 && phi == 0) return sigmaL0Phi0[phiClusterSize-1];
600 if(layer == 1 && phi == 0) return sigmaL1Phi0[phiClusterSize-1];
601 if(layer == 2 && phi == 0) return sigmaL2Phi0[phiClusterSize-1];
602 if(layer == 0 && phi == 1) return sigmaL0Phi1[phiClusterSize-1];
603 if(layer == 1 && phi == 1) return sigmaL1Phi1[phiClusterSize-1];
604 if(layer == 2 && phi == 1) return sigmaL2Phi1[phiClusterSize-1];
605
606 // shouldn't really happen...
607 ATH_MSG_WARNING("Unexpected layer and phi numbers: layer = "
608 << layer << " and phi = " << phi);
609 return 14.6*micrometer;
610
611}
612
613}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#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 class provides an interface to generate or decode an identifier for the upper levels of the dete...
double charge(const T &p)
Definition AtlasPID.h:997
This is an Identifier helper class for the Pixel subdetector.
const double width
#define x
constexpr int pow(int base, int exp) noexcept
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
This class provides an interface to generate or decode an identifier for the upper levels of the dete...
virtual HelperType helper() const
Type of helper, defaulted to 'Unimplemented'.
This is a "hash" representation of an Identifier.
virtual DetectorShape shape() const
Shape of element.
static constexpr std::array< PixelDiodeTree::CellIndexType, 2 > makeCellIndex(T local_x_idx, T local_y_idx)
Create a 2D cell index from the indices in local-x (phi, row) and local-y (eta, column) direction.
Class used to describe the design of a module (diode segmentation and readout scheme)
virtual int numberOfConnectedCells(const SiReadoutCellId &readoutId) const
readout id -> id of connected diodes
PixelReadoutTechnology getReadoutTechnology() const
int columns() const
Number of cell columns per module:
PixelDiodeTree::DiodeProxy diodeProxyFromIdx(const std::array< PixelDiodeTree::IndexType, 2 > &idx) const
int rows() const
Number of cell rows per module:
virtual SiReadoutCellId readoutIdOfCell(const SiCellId &cellId) const
diode id -> readout id
static InDetDD::PixelDiodeType getDiodeType(const PixelDiodeTree::DiodeProxy &diode_proxy)
static unsigned int getFE(const PixelDiodeTree::DiodeProxy &diode_proxy)
Identifier for the strip or pixel cell.
Definition SiCellId.h:29
Class to hold geometrical description of a silicon detector element.
virtual const SiDetectorDesign & design() const override final
access to the local description (inline):
double phiPitch() const
Pitch (inline methods)
double sinStereoLocal(const Amg::Vector2D &localPos) const
Angle of strip in local frame with respect to the etaAxis.
virtual IdentifierHash identifyHash() const override final
identifier hash (inline)
const AtlasDetectorID * getIdHelper() const
Returns the id helper (inline)
Trk::Surface & surface()
Element Surface.
PixelCluster pixelCluster(const Identifier &clusterID, const Amg::Vector2D &localPos, std::vector< Identifier > &&rdoList, const int lvl1a, std::vector< int > &&totList, const SiWidth &width, const InDetDD::SiDetectorElement *element, bool ganged, int errorStrategy, const PixelID &pixelID, bool split, double splitProb1, double splitProb2, const PixelChargeCalibCondData *calibData, const PixelCalib::PixelOfflineCalibData *offlineCalibData, const EventContext &ctx) const
ClusterType makePixelCluster(const Identifier &clusterID, const Amg::Vector2D &localPos, IdentifierVec &&rdoList, const int lvl1a, ToTList &&totList, const SiWidth &width, const InDetDD::SiDetectorElement *element, bool ganged, int errorStrategy, const PixelID &pixelID, bool split, double splitProb1, double splitProb2, const PixelChargeCalibCondData *calibData, const PixelCalib::PixelOfflineCalibData *offlineCalibData, const EventContext &ctx, xAOD::PixelCluster *cluster=nullptr) const
ToolHandle< ISiLorentzAngleTool > m_sctLorentzAngleTool
ClusterMakerTool(const std::string &type, const std::string &name, const IInterface *parent)
double getPixelCTBPhiError(int layer, int phi, int PhiClusterSize) const
xAOD::PixelCluster * xAODpixelCluster(xAOD::PixelCluster &cluster, const Amg::Vector2D &localPos, const std::vector< Identifier > &rdoList, const int lvl1a, const std::vector< int > &totList, const SiWidth &width, const InDetDD::SiDetectorElement *element, bool ganged, int errorStrategy, const PixelID &pixelID, bool split, double splitProb1, double splitProb2, const PixelChargeCalibCondData *calibData, const PixelCalib::PixelOfflineCalibData *offlineCalibData, const EventContext &ctx) const
ToolHandle< ISiLorentzAngleTool > m_pixelLorentzAngleTool
SCT_Cluster sctCluster(const Identifier &clusterID, const Amg::Vector2D &localPos, std::vector< Identifier > &&rdoList, const SiWidth &width, const InDetDD::SiDetectorElement *element, int errorStrategy) const
int getBarrelBin(double eta, int etaClusterSize, int phiClusterSize) const
int getEndcapBin(int etaClusterSize, int phiClusterSize) const
PixelClusterErrorData * getPixelClusterErrorData()
float getCharge(InDetDD::PixelDiodeType type, unsigned int moduleHash, unsigned int FE, float ToT) const
This is an Identifier helper class for the Pixel subdetector.
Definition PixelID.h:67
int eta_index(const Identifier &id) const
Definition PixelID.h:645
Identifier wafer_id(int barrel_ec, int layer_disk, int phi_module, int eta_module) const
For a single crystal.
Definition PixelID.h:360
Identifier pixel_id(int barrel_ec, int layer_disk, int phi_module, int eta_module, int phi_index, int eta_index) const
For an individual pixel.
Definition PixelID.h:428
IdentifierHash wafer_hash(Identifier wafer_id) const
wafer hash from id
Definition PixelID.h:383
int phi_index(const Identifier &id) const
Definition PixelID.h:639
const Amg::Transform3D & transform() const
Returns HepGeom::Transform3D by reference.
std::vector< std::string > split(const std::string &s, const std::string &t=":")
Definition hcg.cxx:177
constexpr double ONE_TWELFTH
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 2, 1 > Vector2D
Eigen::Matrix< double, 3, 1 > Vector3D
Primary Vertex Finder.
@ locY
local cartesian
Definition ParamDefs.h:38
@ locX
Definition ParamDefs.h:37
PixelCluster_v1 PixelCluster
Define the version of the pixel cluster class.
Helper class to access parameters of a diode.