ATLAS Offline Software
Loading...
Searching...
No Matches
InDet::PixelClusterOnTrackTool Class Referencefinal

#include <PixelClusterOnTrackTool.h>

Inheritance diagram for InDet::PixelClusterOnTrackTool:
Collaboration diagram for InDet::PixelClusterOnTrackTool:

Public Member Functions

 PixelClusterOnTrackTool (const std::string &, const std::string &, const IInterface *)
 AlgTool constructor.
virtual ~PixelClusterOnTrackTool ()
virtual StatusCode initialize () override
 AlgTool initialisation.
virtual StatusCode finalize () override
 AlgTool termination.
virtual InDet::PixelClusterOnTrackcorrect (const Trk::PrepRawData &, const Trk::TrackParameters &, const EventContext &ctx=Gaudi::Hive::currentContext()) const override
 produces a PixelClusterOnTrack (object factory!).
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysInitialize () override
 Perform system initialization for an algorithm.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const

Static Public Member Functions

static const InterfaceID & interfaceID ()
 The AlgTool InterfaceID.

Protected Member Functions

void correctBow (const Identifier &, Amg::Vector2D &locpos, const double tanphi, const double taneta, const EventContext &ctx) const
InDet::PixelClusterOnTrackcorrectDefault (const Trk::PrepRawData &, const Trk::TrackParameters &, const EventContext &ctx) const
 The correct method produces a PixelClusterOnTrack using the measured PixelCluster and the track prediction.
InDet::PixelClusterOnTrackcorrectNN (const Trk::PrepRawData &, const Trk::TrackParameters &, const EventContext &ctx) const
bool getErrorsDefaultAmbi (const InDet::PixelCluster *, const Trk::TrackParameters &, Amg::Vector2D &, Amg::MatrixX &) const
bool getErrorsTIDE_Ambi (const InDet::PixelCluster *, const Trk::TrackParameters &, Amg::Vector2D &, Amg::MatrixX &) const
InDet::PixelClusterOnTrackcorrect (const Trk::PrepRawData &, const Trk::TrackParameters &, const InDet::PixelClusterStrategy) const
const Trk::ClusterSplitProbabilityContainer::ProbabilityInfogetClusterSplittingProbability (const InDet::PixelCluster *pix) const
void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

SG::ReadCondHandleKey< PixelDistortionDatam_distortionKey {this, "PixelDistortionData", "PixelDistortionData", "Output readout distortion data"}
ToolHandle< ISiLorentzAngleToolm_lorentzAngleTool {this, "LorentzAngleTool", "SiLorentzAngleTool", "Tool to retrieve Lorentz angle"}
SG::ReadCondHandleKey< PixelCalib::PixelOfflineCalibDatam_clusterErrorKey {this, "PixelOfflineCalibData", "PixelOfflineCalibData", "Output key of pixel cluster"}
SG::ReadCondHandleKey< RIO_OnTrackErrorScalingm_pixelErrorScalingKey {this,"PixelErrorScalingKey", "/Indet/TrkErrorScalingPixel", "Key for pixel error scaling conditions data."}
BooleanProperty m_disableDistortions
 toolhandle for central error scaling flag storing if errors need scaling or should be kept nominal
IntegerProperty m_positionStrategy
std::atomic_int m_errorStrategy {2}
IntegerProperty m_errorStrategyProperty {this, "ErrorStrategy", 2, "Which calibration of cluster position errors"}
const PixelIDm_pixelid = nullptr
 Flag controlling how module distortions are taken into account:
BooleanProperty m_applyNNcorrectionProperty {this, "applyNNcorrection", false}
 Enable NN based calibration (do only if NN calibration is applied)
bool m_applyNNcorrection = false
BooleanProperty m_NNIBLcorrection {this, "NNIBLcorrection", false}
bool m_IBLAbsent = true
ToolHandle< NnClusterizationFactorym_NnClusterizationFactory
 NN clusterizationi factory for NN based positions and errors.
ServiceHandle< IIBLParameterSvcm_IBLParameterSvc {this, "IBLParameterSvc", "IBLParameterSvc"}
BooleanProperty m_doNotRecalibrateNN {this, "doNotRecalibrateNN", false}
BooleanProperty m_noNNandBroadErrors {this, "noNNandBroadErrors", false}
BooleanProperty m_usingTIDE_Ambi {this, "RunningTIDE_Ambi", false}
 Enable different treatment of cluster errors based on NN information (do only if TIDE ambi is run)
SG::ReadHandleKey< InDet::PixelGangedClusterAmbiguities > m_splitClusterMapKey {this, "SplitClusterAmbiguityMap", ""}
SG::ReadHandleKey< Trk::ClusterSplitProbabilityContainerm_clusterSplitProbContainer {this, "ClusterSplitProbabilityName", "",""}
double m_calphi [s_nbinphi] {}
double m_caleta [s_nbineta][3] {}
double m_calerrphi [s_nbinphi][3] {}
double m_calerreta [s_nbineta][3] {}
double m_phix [s_nbinphi+1] {}
double m_etax [s_nbineta+1] {}
StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default)
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Static Private Attributes

static constexpr int s_nbinphi =9
static constexpr int s_nbineta =6

Detailed Description

Definition at line 61 of file PixelClusterOnTrackTool.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< AlgTool > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ PixelClusterOnTrackTool()

InDet::PixelClusterOnTrackTool::PixelClusterOnTrackTool ( const std::string & t,
const std::string & n,
const IInterface * p )

AlgTool constructor.

Definition at line 57 of file PixelClusterOnTrackTool.cxx.

58 :
59 ::AthAlgTool(t, n, p) {
60 declareInterface<IRIO_OnTrackCreator>(this);
61}
AthAlgTool()
Default constructor:

◆ ~PixelClusterOnTrackTool()

InDet::PixelClusterOnTrackTool::~PixelClusterOnTrackTool ( )
virtualdefault

Member Function Documentation

◆ correct() [1/2]

InDet::PixelClusterOnTrack * InDet::PixelClusterOnTrackTool::correct ( const Trk::PrepRawData & rio,
const Trk::TrackParameters & trackPar,
const EventContext & ctx = Gaudi::Hive::currentContext() ) const
overridevirtual

produces a PixelClusterOnTrack (object factory!).

Depending on job options it changes the pixel cluster position and error according to the parameters (in particular, the angle) of the intersecting track.

Implements Trk::IRIO_OnTrackCreator.

Definition at line 140 of file PixelClusterOnTrackTool.cxx.

141 {
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 }
160 const Trk::ClusterSplitProbabilityContainer::ProbabilityInfo &splitProb = getClusterSplittingProbability(pix);
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}
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...
InDet::PixelClusterOnTrack * correctNN(const Trk::PrepRawData &, const Trk::TrackParameters &, const EventContext &ctx) const
const Trk::ClusterSplitProbabilityContainer::ProbabilityInfo & getClusterSplittingProbability(const InDet::PixelCluster *pix) const
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...
virtual bool type(PrepRawDataType type) const
Interface method checking the type.

◆ correct() [2/2]

InDet::PixelClusterOnTrack * InDet::PixelClusterOnTrackTool::correct ( const Trk::PrepRawData & rio,
const Trk::TrackParameters & trackPar,
const InDet::PixelClusterStrategy strategy ) const
protected

Definition at line 525 of file PixelClusterOnTrackTool.cxx.

527 {
528 int initial_errorStrategy;
529 InDet::PixelClusterOnTrack *newROT;
530
531 switch (strategy) {
532 case InDet::PIXELCLUSTER_OUTLIER: // if cluster is outlier, increase errors
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}
virtual InDet::PixelClusterOnTrack * correct(const Trk::PrepRawData &, const Trk::TrackParameters &, const EventContext &ctx=Gaudi::Hive::currentContext()) const override
produces a PixelClusterOnTrack (object factory!).

◆ correctBow()

void InDet::PixelClusterOnTrackTool::correctBow ( const Identifier & id,
Amg::Vector2D & locpos,
const double tanphi,
const double taneta,
const EventContext & ctx ) const
protected

Definition at line 516 of file PixelClusterOnTrackTool.cxx.

517 {
519 Amg::Vector2D newpos = SG::ReadCondHandle<PixelDistortionData>(m_distortionKey, ctx)->correctReconstruction(m_pixelid->wafer_hash(id), localpos, dir);
520
521 localpos = newpos;
522}
Scalar phi() const
phi method
Scalar theta() const
theta method
const PixelID * m_pixelid
Flag controlling how module distortions are taken into account:
SG::ReadCondHandleKey< PixelDistortionData > m_distortionKey
Eigen::Matrix< double, 2, 1 > Vector2D
Eigen::Matrix< double, 3, 1 > Vector3D

◆ correctDefault()

InDet::PixelClusterOnTrack * InDet::PixelClusterOnTrackTool::correctDefault ( const Trk::PrepRawData & rio,
const Trk::TrackParameters & trackPar,
const EventContext & ctx ) const
protected

The correct method produces a PixelClusterOnTrack using the measured PixelCluster and the track prediction.

Definition at line 178 of file PixelClusterOnTrackTool.cxx.

179 {
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 //
196 Trk::LocalParameters locpar;
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();
278 InDetDD::SiLocalPosition pos1 =
279 design->positionFromColumnRow(colmin, rowmin);
280 InDetDD::SiLocalPosition pos2 =
281 design->positionFromColumnRow(colmax, rowmin);
282 InDetDD::SiLocalPosition pos3 =
283 design->positionFromColumnRow(colmin, rowmax);
284 InDetDD::SiLocalPosition pos4 =
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.
295 SG::ReadCondHandle<PixelCalib::PixelOfflineCalibData> offlineCalibData(m_clusterErrorKey, ctx);
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
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()) {
503 SG::ReadCondHandle<RIO_OnTrackErrorScaling> error_scaling( m_pixelErrorScalingKey , ctx);
505 ->getScaledCovariance(std::move(cov), *m_pixelid,
506 element->identify());
507 }
508 bool isbroad = m_errorStrategy == 0;
509 return new InDet::PixelClusterOnTrack(pix,
510 std::move(locpar),
511 std::move(cov), iH,
512 glob, pix->gangedPixel(), isbroad);
513}
#define M_PI
#define ATH_MSG_ERROR(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
double angle(const GeoTrf::Vector2D &a, const GeoTrf::Vector2D &b)
const double width
int readoutSide() const
ReadoutSide.
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...
virtual const SiDetectorDesign & design() const override final
access to the local description (inline):
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)
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
SG::ReadCondHandleKey< PixelCalib::PixelOfflineCalibData > m_clusterErrorKey
ToolHandle< ISiLorentzAngleTool > m_lorentzAngleTool
void correctBow(const Identifier &, Amg::Vector2D &locpos, const double tanphi, const double taneta, const EventContext &ctx) const
const Amg::Vector3D & globalPosition() const
return global position reference
bool gangedPixel() const
return the flag of this cluster containing a gangedPixel
const InDet::SiWidth & width() const
return width class reference
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.
double pT() const
Access method for transverse momentum.
const Amg::MatrixX & localCovariance() const
return const ref to the error matrix
const std::vector< Identifier > & rdoList() const
return the List of rdo identifiers (pointers)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
row
Appending html table to final .html summary file.
const T_res * ErrorScalingCast(const T_src *src)

◆ correctNN()

InDet::PixelClusterOnTrack * InDet::PixelClusterOnTrackTool::correctNN ( const Trk::PrepRawData & rio,
const Trk::TrackParameters & trackPar,
const EventContext & ctx ) const
protected

Definition at line 547 of file PixelClusterOnTrackTool.cxx.

550 {
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
588 Trk::LocalParameters locpar = Trk::LocalParameters(locpos);
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()) {
636 SG::ReadCondHandle<RIO_OnTrackErrorScaling> error_scaling(
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}
#define ATH_MSG_DEBUG(x)
bool getErrorsDefaultAmbi(const InDet::PixelCluster *, const Trk::TrackParameters &, Amg::Vector2D &, Amg::MatrixX &) const
BooleanProperty m_usingTIDE_Ambi
Enable different treatment of cluster errors based on NN information (do only if TIDE ambi is run)
bool getErrorsTIDE_Ambi(const InDet::PixelCluster *, const Trk::TrackParameters &, Amg::Vector2D &, Amg::MatrixX &) const
const Amg::Vector2D & localPosition() const
return the local position reference

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ finalize()

StatusCode InDet::PixelClusterOnTrackTool::finalize ( )
overridevirtual

AlgTool termination.

Definition at line 130 of file PixelClusterOnTrackTool.cxx.

130 {
131return StatusCode::SUCCESS;
132}

◆ getClusterSplittingProbability()

const Trk::ClusterSplitProbabilityContainer::ProbabilityInfo & InDet::PixelClusterOnTrackTool::getClusterSplittingProbability ( const InDet::PixelCluster * pix) const
inlineprotected

Definition at line 113 of file PixelClusterOnTrackTool.h.

113 {
115
116 SG::ReadHandle<Trk::ClusterSplitProbabilityContainer> splitProbContainer(m_clusterSplitProbContainer);
117 if (!splitProbContainer.isValid()) {
118 ATH_MSG_FATAL("Failed to get cluster splitting probability container " << m_clusterSplitProbContainer);
119 }
120 return splitProbContainer->splitProbability(pix);
121 }
#define ATH_MSG_FATAL(x)
SG::ReadHandleKey< Trk::ClusterSplitProbabilityContainer > m_clusterSplitProbContainer
static const ProbabilityInfo & getNoSplitProbability()

◆ getErrorsDefaultAmbi()

bool InDet::PixelClusterOnTrackTool::getErrorsDefaultAmbi ( const InDet::PixelCluster * pixelPrepCluster,
const Trk::TrackParameters & trackPar,
Amg::Vector2D & finalposition,
Amg::MatrixX & finalerrormatrix ) const
protected

Definition at line 659 of file PixelClusterOnTrackTool.cxx.

662 {
663 std::vector<Amg::Vector2D> vectorOfPositions;
664 int numberOfSubclusters = 1;
665 vectorOfPositions.push_back(pixelPrepCluster->localPosition());
666
668 SG::ReadHandle<InDet::PixelGangedClusterAmbiguities> splitClusterMap(m_splitClusterMapKey);
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;
679 if (otherOne->type(Trk::PrepRawDataType::PixelCluster)) {
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}
SG::ReadHandleKey< InDet::PixelGangedClusterAmbiguities > m_splitClusterMapKey
ToolHandle< NnClusterizationFactory > m_NnClusterizationFactory
NN clusterizationi factory for NN based positions and errors.
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.
float distance(const Amg::Vector3D &p1, const Amg::Vector3D &p2)
calculates the distance between two point in 3D space
bool first
Definition DeMoScan.py:534

◆ getErrorsTIDE_Ambi()

bool InDet::PixelClusterOnTrackTool::getErrorsTIDE_Ambi ( const InDet::PixelCluster * pixelPrepCluster,
const Trk::TrackParameters & trackPar,
Amg::Vector2D & finalposition,
Amg::MatrixX & finalerrormatrix ) const
protected

Definition at line 807 of file PixelClusterOnTrackTool.cxx.

810 {
811 const Trk::ClusterSplitProbabilityContainer::ProbabilityInfo &splitProb = getClusterSplittingProbability(pixelPrepCluster);
812 std::vector<Amg::Vector2D> vectorOfPositions;
813 int numberOfSubclusters = 1;
815 SG::ReadHandle<InDet::PixelGangedClusterAmbiguities> splitClusterMap(m_splitClusterMapKey);
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}
Amg::Vector2D localPosition() const
Access method for the local coordinates, local parameter definitions differ for each surface type.
str index
Definition DeMoScan.py:362

◆ initialize()

StatusCode InDet::PixelClusterOnTrackTool::initialize ( )
overridevirtual

AlgTool initialisation.

UGLY!

Definition at line 74 of file PixelClusterOnTrackTool.cxx.

74 {
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}
#define ATH_CHECK
Evaluate an expression and check for errors.
const ServiceHandle< StoreGateSvc > & detStore() const
BooleanProperty m_applyNNcorrectionProperty
Enable NN based calibration (do only if NN calibration is applied)
ServiceHandle< IIBLParameterSvc > m_IBLParameterSvc
retrieve(aClass, aKey=None)
Definition PyKernel.py:110

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ interfaceID()

const InterfaceID & Trk::IRIO_OnTrackCreator::interfaceID ( )
inlinestaticinherited

The AlgTool InterfaceID.

Definition at line 42 of file IRIO_OnTrackCreator.h.

◆ msg()

MsgStream & AthCommonMsg< AlgTool >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msgLvl()

bool AthCommonMsg< AlgTool >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< AlgTool > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ sysInitialize()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysInitialize ( )
overridevirtualinherited

Perform system initialization for an algorithm.

We override this to declare all the elements of handle key arrays at the end of initialization. See comments on updateVHKA.

Reimplemented in asg::AsgMetadataTool, AthCheckedComponent< AthAlgTool >, AthCheckedComponent<::AthAlgTool >, and DerivationFramework::CfAthAlgTool.

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }
std::vector< SG::VarHandleKeyArray * > m_vhka

Member Data Documentation

◆ m_applyNNcorrection

bool InDet::PixelClusterOnTrackTool::m_applyNNcorrection = false
private

Definition at line 168 of file PixelClusterOnTrackTool.h.

◆ m_applyNNcorrectionProperty

BooleanProperty InDet::PixelClusterOnTrackTool::m_applyNNcorrectionProperty {this, "applyNNcorrection", false}
private

Enable NN based calibration (do only if NN calibration is applied)

Definition at line 167 of file PixelClusterOnTrackTool.h.

167{this, "applyNNcorrection", false};

◆ m_calerreta

double InDet::PixelClusterOnTrackTool::m_calerreta[s_nbineta][3] {}
private

Definition at line 195 of file PixelClusterOnTrackTool.h.

195{};

◆ m_calerrphi

double InDet::PixelClusterOnTrackTool::m_calerrphi[s_nbinphi][3] {}
private

Definition at line 194 of file PixelClusterOnTrackTool.h.

194{};

◆ m_caleta

double InDet::PixelClusterOnTrackTool::m_caleta[s_nbineta][3] {}
private

Definition at line 193 of file PixelClusterOnTrackTool.h.

193{};

◆ m_calphi

double InDet::PixelClusterOnTrackTool::m_calphi[s_nbinphi] {}
private

Definition at line 192 of file PixelClusterOnTrackTool.h.

192{};

◆ m_clusterErrorKey

SG::ReadCondHandleKey<PixelCalib::PixelOfflineCalibData> InDet::PixelClusterOnTrackTool::m_clusterErrorKey {this, "PixelOfflineCalibData", "PixelOfflineCalibData", "Output key of pixel cluster"}
private

Definition at line 135 of file PixelClusterOnTrackTool.h.

136{this, "PixelOfflineCalibData", "PixelOfflineCalibData", "Output key of pixel cluster"};

◆ m_clusterSplitProbContainer

SG::ReadHandleKey<Trk::ClusterSplitProbabilityContainer> InDet::PixelClusterOnTrackTool::m_clusterSplitProbContainer {this, "ClusterSplitProbabilityName", "",""}
private

Definition at line 186 of file PixelClusterOnTrackTool.h.

187{this, "ClusterSplitProbabilityName", "",""};

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_disableDistortions

BooleanProperty InDet::PixelClusterOnTrackTool::m_disableDistortions
private
Initial value:
{this, "DisableDistortions", false,
"Disable simulation of module distortions"}

toolhandle for central error scaling flag storing if errors need scaling or should be kept nominal

Definition at line 142 of file PixelClusterOnTrackTool.h.

142 {this, "DisableDistortions", false,
143 "Disable simulation of module distortions"};

◆ m_distortionKey

SG::ReadCondHandleKey<PixelDistortionData> InDet::PixelClusterOnTrackTool::m_distortionKey {this, "PixelDistortionData", "PixelDistortionData", "Output readout distortion data"}
private

Definition at line 129 of file PixelClusterOnTrackTool.h.

130{this, "PixelDistortionData", "PixelDistortionData", "Output readout distortion data"};

◆ m_doNotRecalibrateNN

BooleanProperty InDet::PixelClusterOnTrackTool::m_doNotRecalibrateNN {this, "doNotRecalibrateNN", false}
private

Definition at line 179 of file PixelClusterOnTrackTool.h.

179{this, "doNotRecalibrateNN", false};

◆ m_errorStrategy

std::atomic_int InDet::PixelClusterOnTrackTool::m_errorStrategy {2}
mutableprivate

Definition at line 146 of file PixelClusterOnTrackTool.h.

146{2};

◆ m_errorStrategyProperty

IntegerProperty InDet::PixelClusterOnTrackTool::m_errorStrategyProperty {this, "ErrorStrategy", 2, "Which calibration of cluster position errors"}
private

Definition at line 147 of file PixelClusterOnTrackTool.h.

148{this, "ErrorStrategy", 2, "Which calibration of cluster position errors"};

◆ m_etax

double InDet::PixelClusterOnTrackTool::m_etax[s_nbineta+1] {}
private

Definition at line 197 of file PixelClusterOnTrackTool.h.

197{};

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_IBLAbsent

bool InDet::PixelClusterOnTrackTool::m_IBLAbsent = true
private

Definition at line 170 of file PixelClusterOnTrackTool.h.

◆ m_IBLParameterSvc

ServiceHandle<IIBLParameterSvc> InDet::PixelClusterOnTrackTool::m_IBLParameterSvc {this, "IBLParameterSvc", "IBLParameterSvc"}
private

Definition at line 176 of file PixelClusterOnTrackTool.h.

177{this, "IBLParameterSvc", "IBLParameterSvc"};

◆ m_lorentzAngleTool

ToolHandle<ISiLorentzAngleTool> InDet::PixelClusterOnTrackTool::m_lorentzAngleTool {this, "LorentzAngleTool", "SiLorentzAngleTool", "Tool to retrieve Lorentz angle"}
private

Definition at line 132 of file PixelClusterOnTrackTool.h.

133{this, "LorentzAngleTool", "SiLorentzAngleTool", "Tool to retrieve Lorentz angle"};

◆ m_NnClusterizationFactory

ToolHandle<NnClusterizationFactory> InDet::PixelClusterOnTrackTool::m_NnClusterizationFactory
private
Initial value:
{this, "NnClusterizationFactory",
"InDet::NnClusterizationFactory/NnClusterizationFactory"}

NN clusterizationi factory for NN based positions and errors.

Definition at line 173 of file PixelClusterOnTrackTool.h.

174 {this, "NnClusterizationFactory",
175 "InDet::NnClusterizationFactory/NnClusterizationFactory"};

◆ m_NNIBLcorrection

BooleanProperty InDet::PixelClusterOnTrackTool::m_NNIBLcorrection {this, "NNIBLcorrection", false}
private

Definition at line 169 of file PixelClusterOnTrackTool.h.

169{this, "NNIBLcorrection", false};

◆ m_noNNandBroadErrors

BooleanProperty InDet::PixelClusterOnTrackTool::m_noNNandBroadErrors {this, "noNNandBroadErrors", false}
private

Definition at line 180 of file PixelClusterOnTrackTool.h.

180{this, "noNNandBroadErrors", false};

◆ m_phix

double InDet::PixelClusterOnTrackTool::m_phix[s_nbinphi+1] {}
private

Definition at line 196 of file PixelClusterOnTrackTool.h.

196{};

◆ m_pixelErrorScalingKey

SG::ReadCondHandleKey<RIO_OnTrackErrorScaling> InDet::PixelClusterOnTrackTool::m_pixelErrorScalingKey {this,"PixelErrorScalingKey", "/Indet/TrkErrorScalingPixel", "Key for pixel error scaling conditions data."}
private

Definition at line 137 of file PixelClusterOnTrackTool.h.

138{this,"PixelErrorScalingKey", "/Indet/TrkErrorScalingPixel", "Key for pixel error scaling conditions data."};

◆ m_pixelid

const PixelID* InDet::PixelClusterOnTrackTool::m_pixelid = nullptr
private

Flag controlling how module distortions are taken into account:

case 0 -----> No distorsions implemented;

case 1 -----> Set curvature (in 1/meter) and twist (in radiant) equal for all modules;

case 2 -----> Read curvatures and twists from textfile containing Survey data;

case 3 -----> Set curvature and twist from Gaussian random generator with mean and RMS coming from Survey data;

case 4 -----> Read curvatures and twists from database (not ready yet); identifier-helper

Definition at line 164 of file PixelClusterOnTrackTool.h.

◆ m_positionStrategy

IntegerProperty InDet::PixelClusterOnTrackTool::m_positionStrategy
private
Initial value:
{this, "PositionStrategy", 1,
"Which calibration of cluster positions"}

Definition at line 144 of file PixelClusterOnTrackTool.h.

144 {this, "PositionStrategy", 1,
145 "Which calibration of cluster positions"};

◆ m_splitClusterMapKey

SG::ReadHandleKey<InDet::PixelGangedClusterAmbiguities> InDet::PixelClusterOnTrackTool::m_splitClusterMapKey {this, "SplitClusterAmbiguityMap", ""}
private

Definition at line 183 of file PixelClusterOnTrackTool.h.

184{this, "SplitClusterAmbiguityMap", ""};

◆ m_usingTIDE_Ambi

BooleanProperty InDet::PixelClusterOnTrackTool::m_usingTIDE_Ambi {this, "RunningTIDE_Ambi", false}
private

Enable different treatment of cluster errors based on NN information (do only if TIDE ambi is run)

Definition at line 182 of file PixelClusterOnTrackTool.h.

182{this, "RunningTIDE_Ambi", false};

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< AlgTool > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< AlgTool > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.

◆ s_nbineta

int InDet::PixelClusterOnTrackTool::s_nbineta =6
staticconstexprprivate

Definition at line 191 of file PixelClusterOnTrackTool.h.

◆ s_nbinphi

int InDet::PixelClusterOnTrackTool::s_nbinphi =9
staticconstexprprivate

Definition at line 190 of file PixelClusterOnTrackTool.h.


The documentation for this class was generated from the following files: