ATLAS Offline Software
Loading...
Searching...
No Matches
InDetPhysHitDecoratorAlg.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
9
11#include "safeDecorator.h"
12#include "TrkParameters/TrackParameters.h" // Contains typedef to Trk::CurvilinearParameters
16// for the identifiers
22//
23#include <tuple>
24#include <limits>
25
26
27InDetPhysHitDecoratorAlg::InDetPhysHitDecoratorAlg(const std::string& name, ISvcLocator* pSvcLocator) :
28 AthReentrantAlgorithm(name,pSvcLocator) {}
29
32
33StatusCode
35 ATH_CHECK(m_holeSearchTool.retrieve());
36 if (not (m_updatorHandle.empty())) {
37 ATH_CHECK(m_updatorHandle.retrieve());
38 }
39 ATH_CHECK(m_lorentzAngleTool.retrieve());
40
41 std::vector<std::string> float_decor_names(kNFloatDecorators);
42 std::vector<std::string> int_decor_names(kNIntDecorators);
43 std::vector<std::string> uint64_decor_names(kNUInt64Decorators);
44
45 int_decor_names[kDecorRegion]="measurement_region";
46 int_decor_names[kDecorDet]="measurement_det";
47 int_decor_names[kDecorILayer]="measurement_iLayer";
48 int_decor_names[kDecorType]="measurement_type";
49 int_decor_names[kDecorPhiWidth]="hitResiduals_phiWidth";
50 int_decor_names[kDecorEtaWidth]="hitResiduals_etaWidth";
51
52 float_decor_names[kDecorResidualLocX]="hitResiduals_residualLocX";
53 float_decor_names[kDecorPullLocX]="hitResiduals_pullLocX";
54 float_decor_names[kDecorMeasLocX]="measurementLocX";
55 float_decor_names[kDecorTrkParamLocX]="trackParamLocX";
56 float_decor_names[kDecorMeasLocCovX]="measurementLocCovX";
57
58 float_decor_names[kDecorResidualLocY]="hitResiduals_residualLocY";
59 float_decor_names[kDecorPullLocY]="hitResiduals_pullLocY";
60 float_decor_names[kDecorMeasLocY]="measurementLocY";
61 float_decor_names[kDecorTrkParamLocY]="trackParamLocY";
62 float_decor_names[kDecorMeasLocCovY]="measurementLocCovY";
63
64 float_decor_names[kDecorAngle]="angle";
65 float_decor_names[kDecorEtaLoc]="etaloc";
66
67 uint64_decor_names[kDecorID]="surfaceID";
68
69 ATH_CHECK( m_trkParticleName.initialize() );
73 assert( m_intDecor.size() == kNIntDecorators);
74 assert( m_uint64Decor.size() == kNUInt64Decorators);
75 assert( m_floatDecor.size() == kNFloatDecorators);
76
77 // Get the dictionary manager from the detector store
78 ATH_CHECK(detStore()->retrieve(m_idHelper, "AtlasID"));
79 ATH_CHECK(detStore()->retrieve(m_pixelID, "PixelID"));
80 ATH_CHECK(detStore()->retrieve(m_sctID, "SCT_ID"));
81 ATH_CHECK(detStore()->retrieve(m_trtID, "TRT_ID"));
82
83 if (m_residualPullCalculator.empty()) {
84 ATH_MSG_INFO("No residual/pull calculator for general hit residuals configured.");
85 ATH_MSG_INFO("It is recommended to give R/P calculators to the det-specific tool handle lists then.");
86 } else if (m_residualPullCalculator.retrieve().isFailure()) {
87 ATH_MSG_FATAL( "Could not retrieve " << m_residualPullCalculator << " (to calculate residuals and pulls) " );
88 } else {
89 ATH_MSG_INFO("Generic hit residuals & pulls will be calculated in one or both available local coordinates");
90 }
91 return StatusCode::SUCCESS;
92}
93
94StatusCode
96 return StatusCode::SUCCESS;
97}
98
99// to migrate to AthReentrantAlgorithm later
100StatusCode
101InDetPhysHitDecoratorAlg::execute(const EventContext &ctx) const {
103 if ((not ptracks.isValid())) {
104 ATH_MSG_ERROR("Cannot get ReadHandle " << m_trkParticleName);
105 return StatusCode::FAILURE;
106 }
107
108 std::vector< SG::WriteDecorHandle<xAOD::TrackParticleContainer,std::vector<float> > >
109 float_decor( IDPVM::createDecorators<xAOD::TrackParticleContainer, std::vector<float> >(m_floatDecor, ctx) );
110 std::vector< SG::WriteDecorHandle<xAOD::TrackParticleContainer,std::vector<int> > >
111 int_decor( IDPVM::createDecorators<xAOD::TrackParticleContainer,std::vector<int> >(m_intDecor, ctx) );
112 std::vector< SG::WriteDecorHandle<xAOD::TrackParticleContainer,std::vector<uint64_t> > >
113 uint64_decor( IDPVM::createDecorators<xAOD::TrackParticleContainer,std::vector<uint64_t> >(m_uint64Decor, ctx) );
114
115 for (const xAOD::TrackParticle *trk_particle : *ptracks) {
116 if (not decorateTrack(*trk_particle, float_decor, int_decor, uint64_decor) ) {
117 ATH_MSG_ERROR("Could not decorate track");
118 return StatusCode::FAILURE;
119 }
120 }
121 return StatusCode::SUCCESS;
122}
123
124bool
126 std::vector< SG::WriteDecorHandle<xAOD::TrackParticleContainer,std::vector<float> > > &float_decor,
127 std::vector< SG::WriteDecorHandle<xAOD::TrackParticleContainer,std::vector<int> > > &int_decor,
128 std::vector< SG::WriteDecorHandle<xAOD::TrackParticleContainer,std::vector<uint64_t> > > &uint64_decor) const
129{
130 int trackNumber(0);
131
132 using SingleResult_t = std::tuple<int, int, int, // detector, region, layer index,
133 float, float, // residual local position X, pull local position X
134 float, float, // residual local position Y, pull local position Y
135 int, int, int, // cluster dimensions in phi and eta directions, measurement type
136 float, float, // measurement local position X, measurement local position Y
137 float, float,// track parameter local X, track parameter local Y
138 float, float, // track angle in local module frame and associated eta
139 float, float, uint64_t>; // measurement covariance local X, measurement covariance local Y, surface identifier
140 using TrackResult_t = std::vector<SingleResult_t>;
141 constexpr float invalidFloat(-1.);
142 constexpr int invalidInteger(-1);
143 constexpr uint64_t invalidID(0);
144 const SingleResult_t invalidResult = std::make_tuple(invalidInteger, invalidInteger, invalidInteger,
145 invalidFloat, invalidFloat, invalidFloat, invalidFloat,
146 invalidInteger, invalidInteger, invalidInteger,
147 invalidFloat, invalidFloat, invalidFloat, invalidFloat,
148 invalidFloat, invalidFloat, invalidFloat, invalidFloat,
149 invalidID);
150 bool isUnbiased(true);
151 // get element link to the original track
152 const ElementLink< TrackCollection >& trackLink = particle.trackLink();
153 if (trackLink.isValid()) {
154 ATH_MSG_VERBOSE("Track link found ");
155 std::unique_ptr<const Trk::Track> trackWithHoles(m_holeSearchTool->getTrackWithHoles(**trackLink));
156 const auto& allTrackStates = *(trackWithHoles->trackStateOnSurfaces());
157 const int numberOfHits(allTrackStates.size());
158 unsigned int trackParametersCounter(numberOfHits);
159 TrackResult_t result;
160 result.reserve(numberOfHits);
161 if (!m_updatorHandle.empty()) {
162 isUnbiased = true;
163 } else {
164 ATH_MSG_WARNING("The updater handle is empty, now using biased estimators");
165 isUnbiased = false;
166 }
167 ATH_MSG_DEBUG("Num. track states in track " << ++trackNumber << ": " << allTrackStates.size());
168
169 for (const auto *const thisTrackState: allTrackStates) {
170 // Copy logic from InDetRttPerformance to get hits/outliers/holes
171 // Variable specifying measurement type filled
172 SingleResult_t thisResult(invalidResult);
173 if (not thisTrackState) { // is this check needed?
174 msg(MSG::ERROR) << "TSOS is NULL" << (thisTrackState) << endmsg;
175 continue;
176 }
177 Identifier surfaceID;
178 const Trk::MeasurementBase* mesb = (thisTrackState)->measurementOnTrack();
179 const Trk::RIO_OnTrack* hit = mesb ? dynamic_cast<const Trk::RIO_OnTrack*>(mesb) : nullptr;
180 if (mesb && !hit) {
181 continue; // skip pseudomeasurements
182 }
183 // Get surfaceID, different for measuremnt hits & outliers, and holes
184 if (mesb && mesb->associatedSurface().associatedDetectorElement()) {
186 } // Holes
187 else {
188 if (not (thisTrackState)->trackParameters()) {
189 msg(MSG::INFO) << "TSOS surface is NULL" << endmsg;
190 continue;
191 }
192 surfaceID = (thisTrackState)->trackParameters()->associatedSurface().associatedDetectorElementIdentifier();
193 }
194 bool isMesb = (thisTrackState)->type(Trk::TrackStateOnSurface::Measurement);
195 bool isOutl = (thisTrackState)->type(Trk::TrackStateOnSurface::Outlier);
196 bool isHole = (thisTrackState)->type(Trk::TrackStateOnSurface::Hole);
197
198 int measureType = -1;
199 if (isMesb) {
200 measureType = 0;
201 }
202 if (isOutl) {
203 measureType = 1;
204 }
205 if (isHole) {
206 measureType = 2;
207 }
208
209 bool anyHit = isMesb || isOutl || isHole;
210 if (!anyHit) {
211 continue;
212 }
215 int iLayer(invalidInteger);
216 const bool successfulIdentification = decideDetectorRegion(surfaceID, det, r, iLayer);
217 if (not successfulIdentification) {
218 ATH_MSG_DEBUG("Could not identify surface");
219 continue;
220 }
221 uint64_t Surface_ID = surfaceID.get_compact();
222 // defining invalid values, they will be filled when and where needed
223 float residualLocX(invalidFloat), pullLocX(invalidFloat);
224 float residualLocY(invalidFloat), pullLocY(invalidFloat);
225 float measurementLocX(invalidFloat), trackParamLocX(invalidFloat), measurementLocCovX(invalidFloat);
226 float measurementLocY(invalidFloat), trackParamLocY(invalidFloat), measurementLocCovY(invalidFloat);
227 float angle(0.), etaloc(0.);
228 int phiWidth(invalidInteger), etaWidth(invalidInteger);
229 std::optional<Trk::ResidualPull> residualPull(std::nullopt);
230 const Trk::TrackParameters* biasedTrackParameters = thisTrackState->trackParameters();
231 if (biasedTrackParameters) {
232 ATH_MSG_VERBOSE("biased track parameters ok");
233 }
234 ATH_MSG_VERBOSE("checking mesb and track parameters");
235 if (mesb && biasedTrackParameters) {
236 ATH_MSG_DEBUG("mesb and biased track parameters are ok");
237 // for outliers, the measurement is not part of the fit, so track parameters are already unbiased
238 std::unique_ptr<const Trk::TrackParameters> cleanup_trackparam;
239 const Trk::TrackParameters* trackParameters =
241 biasedTrackParameters, mesb, isUnbiased) : biasedTrackParameters;
242
243 if (trackParameters != biasedTrackParameters) {
244 cleanup_trackparam.reset(trackParameters);
245 }
246 if (not trackParameters) {
247 ATH_MSG_DEBUG("unbiased track parameters pointer is NULL");
248 }
249
252
253 residualPull= m_residualPullCalculator->residualPull(hit, trackParameters, resType);
254 ATH_MSG_VERBOSE("checking residual pull");
255 if (not residualPull) {
256 ATH_MSG_DEBUG("residualPull is NULL");
257 continue;
258 }
259
260 ATH_MSG_DEBUG("residualPull is OK");
261
262 residualLocX = 1000. * residualPull->residual()[Trk::loc1]; // residuals in microns
263 pullLocX = residualPull->pull()[Trk::loc1];
264 measurementLocX = hit->localParameters()[Trk::loc1];
265 trackParamLocX = trackParameters->parameters()[Trk::loc1];
266 measurementLocCovX = hit->localCovariance()(Trk::loc1,Trk::loc1);
267
268 if (residualPull->dimension() > 1) {
269 residualLocY = 1000. * residualPull->residual()[Trk::loc2];
270 pullLocY = residualPull->pull()[Trk::loc2];
271 measurementLocY = hit->localParameters()[Trk::loc2];
272 trackParamLocY = trackParameters->parameters()[Trk::loc2];
273 measurementLocCovY = hit->localCovariance()(Trk::loc2,Trk::loc2);
274 }
275
276 // Unbiased residuals
277 measureType = 4;
278
279 if (hit && isUnbiased) {
280 // Cluster width determination
281 if ((det == L0PIXBARR)or(det == PIXEL) or(det == SCT)) {
282 assert(dynamic_cast <const InDet::SiCluster*>(hit->prepRawData())!=nullptr);
283 const InDet::SiCluster* pCluster = static_cast <const InDet::SiCluster*>(hit->prepRawData());
284 InDet::SiWidth width = pCluster->width();
285 phiWidth = int(width.colRow().x());
286 etaWidth = int(width.colRow().y());
287
288 // get candidate track angle in module local frame
289 Amg::Vector3D my_track = trackParameters->momentum();
290 const InDetDD::SiDetectorElement* element = pCluster->detectorElement();
291 Amg::Vector3D my_normal = element->normal();
292 Amg::Vector3D my_phiax = element->phiAxis();
293 Amg::Vector3D my_etaax = element->etaAxis();
294 double trkphicomp = my_track.dot(my_phiax);
295 double trketacomp = my_track.dot(my_etaax);
296 double trknormcomp = my_track.dot(my_normal);
297 double bowphi = std::atan2(trkphicomp,trknormcomp);
298 double boweta = std::atan2(trketacomp,trknormcomp);
299
300 double tanl = m_lorentzAngleTool->getTanLorentzAngle(element->identifyHash(), Gaudi::Hive::currentContext());
301 int readoutside = element->design().readoutSide();
302
303 // map the angles of inward-going tracks onto [-PI/2, PI/2]
304 if(bowphi > M_PI/2) bowphi -= M_PI;
305 if(bowphi < -M_PI/2) bowphi += M_PI;
306
307 // finally, subtract the Lorentz angle effect
308 // the readoutside term is needed because of a bug in old
309 // geometry versions (CSC-01-* and CSC-02-*)
310 angle = std::atan(std::tan(bowphi)-readoutside*tanl);
311
312 double thetaloc=-999.;
313 if(boweta > -0.5*M_PI && boweta < M_PI/2.) {
314 thetaloc = M_PI/2.-boweta;
315 } else if(boweta > M_PI/2. && boweta < M_PI) {
316 thetaloc = 1.5*M_PI-boweta;
317 } else { // 3rd quadrant
318 thetaloc = -0.5*M_PI-boweta;
319 }
320 etaloc = -1*log(tan(thetaloc/2.));
321 }
322 ATH_MSG_VERBOSE("hit and isUnbiased ok");
323 }
324 } else {
325 if (not mesb) {
326 ATH_MSG_VERBOSE("mesb not ok");
327 }
328 if (not biasedTrackParameters) {
329 ATH_MSG_VERBOSE("biasedTrackParameters were not found");
330 }
331 --trackParametersCounter;
332 }
333 thisResult = std::make_tuple(det, r, iLayer,
334 residualLocX, pullLocX, residualLocY, pullLocY,
335 phiWidth, etaWidth, measureType,
336 measurementLocX, measurementLocY,
337 trackParamLocX, trackParamLocY,
338 angle, etaloc,
339 measurementLocCovX, measurementLocCovY,
340 Surface_ID);
341 result.push_back(thisResult);
342 } // end of for loop*/
343 ATH_MSG_DEBUG( "Out of " << numberOfHits << " hits, " << trackParametersCounter
344 << " had track params, and " << result.size() << " had residuals.");
345 if (not result.empty()) {
346 const unsigned int arraySize = result.size();
347 std::vector<int> result_det;
348 result_det.reserve(arraySize);
349 std::vector<int> result_r;
350 result_r.reserve(arraySize);
351 std::vector<int> result_iLayer;
352 result_iLayer.reserve(arraySize);
353 std::vector<float> result_residualLocX;
354 result_residualLocX.reserve(arraySize);
355 std::vector<float> result_pullLocX;
356 result_pullLocX.reserve(arraySize);
357 std::vector<float> result_residualLocY;
358 result_residualLocY.reserve(arraySize);
359 std::vector<float> result_pullLocY;
360 result_pullLocY.reserve(arraySize);
361 std::vector<int> result_phiWidth;
362 result_phiWidth.reserve(arraySize);
363 std::vector<int> result_etaWidth;
364 result_etaWidth.reserve(arraySize);
365 std::vector<int> result_measureType;
366 result_measureType.reserve(arraySize);
367 std::vector<float> result_measurementLocX;
368 result_measurementLocX.reserve(arraySize);
369 std::vector<float> result_measurementLocY;
370 result_measurementLocY.reserve(arraySize);
371 std::vector<float> result_trackParamLocX;
372 result_trackParamLocX.reserve(arraySize);
373 std::vector<float> result_trackParamLocY;
374 result_trackParamLocY.reserve(arraySize);
375 std::vector<float> result_angle;
376 result_angle.reserve(arraySize);
377 std::vector<float> result_etaloc;
378 result_etaloc.reserve(arraySize);
379 std::vector<float> result_measurementLocCovX;
380 result_measurementLocCovX.reserve(arraySize);
381 std::vector<float> result_measurementLocCovY;
382 result_measurementLocCovY.reserve(arraySize);
383 std::vector<uint64_t> result_surfaceID;
384 result_surfaceID.reserve(arraySize);
385
386 for (const SingleResult_t& single_result : result) {
387 result_det.push_back(std::get<0>(single_result));
388 result_r.push_back(std::get<1>(single_result));
389 result_iLayer.push_back(std::get<2>(single_result));
390 result_residualLocX.push_back(std::get<3>(single_result));
391 result_pullLocX.push_back(std::get<4>(single_result));
392 result_residualLocY.push_back(std::get<5>(single_result));
393 result_pullLocY.push_back(std::get<6>(single_result));
394 result_phiWidth.push_back(std::get<7>(single_result));
395 result_etaWidth.push_back(std::get<8>(single_result));
396 result_measureType.push_back(std::get<9>(single_result));
397 result_measurementLocX.push_back(std::get<10>(single_result));
398 result_measurementLocY.push_back(std::get<11>(single_result));
399 result_trackParamLocX.push_back(std::get<12>(single_result));
400 result_trackParamLocY.push_back(std::get<13>(single_result));
401 result_angle.push_back(std::get<14>(single_result));
402 result_etaloc.push_back(std::get<15>(single_result));
403 result_measurementLocCovX.push_back(std::get<16>(single_result));
404 result_measurementLocCovY.push_back(std::get<17>(single_result));
405 result_surfaceID.push_back(std::get<18>(single_result));
406 }
407
408 int_decor[kDecorRegion](particle) = result_r;
409 int_decor[kDecorDet](particle) = result_det;
410 int_decor[kDecorILayer](particle) = result_iLayer;
411 int_decor[kDecorPhiWidth](particle) = result_phiWidth;
412 int_decor[kDecorEtaWidth](particle) = result_etaWidth;
413 int_decor[kDecorType](particle) = result_measureType;
414
415 float_decor[kDecorResidualLocX](particle) = result_residualLocX;
416 float_decor[kDecorPullLocX](particle) = result_pullLocX;
417 float_decor[kDecorMeasLocX](particle) = result_measurementLocX;
418 float_decor[kDecorTrkParamLocX](particle) = result_trackParamLocX;
419 float_decor[kDecorMeasLocCovX](particle) = result_measurementLocCovX;
420
421 float_decor[kDecorResidualLocY](particle) = result_residualLocY;
422 float_decor[kDecorPullLocY](particle) = result_pullLocY;
423 float_decor[kDecorMeasLocY](particle) = result_measurementLocY;
424 float_decor[kDecorTrkParamLocY](particle) = result_trackParamLocY;
425 float_decor[kDecorMeasLocCovY](particle) = result_measurementLocCovY;
426
427 float_decor[kDecorAngle](particle) = result_angle;
428 float_decor[kDecorEtaLoc](particle) = result_etaloc;
429
430 uint64_decor[kDecorID](particle) = result_surfaceID;
431
432 return true;
433 } else {
434 // particle below pt threshold for decoration. Since this is not an error now "true" is returned.
435 // If "false" is returned the job would be aborted.
436 return true;
437 }
438 } else {
439 ATH_MSG_ERROR("No valid track link found ");
440 }
441 return false;
442}
443
445 Region& region, int& layer) const {
446 bool success(false);
447 const int pixelSctBarrelIndex(0);
448 const int trtBarrelIndex(1);
449
450 detector = INVALID_DETECTOR;
451 region = INVALID_REGION;
452
453 if (m_idHelper->is_pixel(id)) {
454 detector = PIXEL;
455 region = (m_pixelID->barrel_ec(id) == pixelSctBarrelIndex) ? (BARREL) : (ENDCAP);
456 layer = m_pixelID->layer_disk(id);
457 if (BARREL == region and layer == 0) {
458 detector = L0PIXBARR;
459 }
460 }
461 else if (m_idHelper->is_sct(id)) {
462 detector = SCT;
463 region = (m_sctID->barrel_ec(id) == pixelSctBarrelIndex) ? (BARREL) : (ENDCAP);
464 layer = m_sctID->layer_disk(id);
465 }
466 else if (m_idHelper->is_trt(id)) {
467 detector = TRT;
468 region = (std::abs(m_trtID->barrel_ec(id)) == trtBarrelIndex) ? (BARREL) : (ENDCAP);
469 layer = m_trtID->layer_or_wheel(id);
470 }
471 success = (detector != INVALID_DETECTOR) and (region != INVALID_REGION);
472
473 return success;
474}
475
478 const Trk::MeasurementBase* measurement,
479 bool &isUnbiased) const {
480 const Trk::TrackParameters* unbiasedTrkParameters(trkParameters);
481
482 if (!m_updatorHandle.empty() && (isUnbiased)) {
483 if (trkParameters->covariance()) {
484 // Get unbiased state
485 unbiasedTrkParameters = m_updatorHandle->removeFromState(*trkParameters,
486 measurement->localParameters(),
487 measurement->localCovariance()).release();
488 if (!unbiasedTrkParameters) {
489 ATH_MSG_INFO( "Could not get unbiased track parameters, use normal parameters" );
490 isUnbiased = false;
491 }
492 } else if (not m_alreadyWarned) {
493 // warn only once!
494 ATH_MSG_WARNING("TrackParameters contain no covariance, unbiased track states can not be calculated "
495 "(ie. pulls and residuals will be too small)" );
496 m_alreadyWarned = true;
497 isUnbiased = false;
498 } else {
499 isUnbiased = false;
500 } // end if no measured track parameter
501 }
502 return unbiasedTrkParameters;
503}
#define M_PI
#define endmsg
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#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...
header file for class of same name
This is an Identifier helper class for the Pixel subdetector.
This is an Identifier helper class for the SCT subdetector.
double angle(const GeoTrf::Vector2D &a, const GeoTrf::Vector2D &b)
This is an Identifier helper class for the TRT subdetector.
const double width
const ServiceHandle< StoreGateSvc > & detStore() const
An algorithm that can be simultaneously executed in multiple threads.
value_type get_compact() const
Get the compact id.
int readoutSide() const
ReadoutSide.
Class to hold geometrical description of a silicon detector element.
virtual const SiDetectorDesign & design() const override final
access to the local description (inline):
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)
bool decideDetectorRegion(const Identifier &id, Subdetector &det, Region &r, int &layer) const
bool decorateTrack(const xAOD::TrackParticle &particle, std::vector< SG::WriteDecorHandle< xAOD::TrackParticleContainer, std::vector< float > > > &float_decor, std::vector< SG::WriteDecorHandle< xAOD::TrackParticleContainer, std::vector< int > > > &int_decor, std::vector< SG::WriteDecorHandle< xAOD::TrackParticleContainer, std::vector< uint64_t > > > &uint64_decor) const
InDetPhysHitDecoratorAlg(const std::string &name, ISvcLocator *pSvcLocator)
Gaudi::Property< std::string > m_prefix
ToolHandle< Trk::ITrackHoleSearchTool > m_holeSearchTool
const Trk::TrackParameters * getUnbiasedTrackParameters(const Trk::TrackParameters *trkParameters, const Trk::MeasurementBase *measurement, bool &isUnbiased) const
virtual StatusCode initialize() override
const AtlasDetectorID * m_idHelper
ToolHandle< Trk::IResidualPullCalculator > m_residualPullCalculator
std::vector< SG::WriteDecorHandleKey< xAOD::TrackParticleContainer > > m_intDecor
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_trkParticleName
ToolHandle< ISiLorentzAngleTool > m_lorentzAngleTool
ToolHandle< Trk::IUpdator > m_updatorHandle
std::vector< SG::WriteDecorHandleKey< xAOD::TrackParticleContainer > > m_floatDecor
std::vector< SG::WriteDecorHandleKey< xAOD::TrackParticleContainer > > m_uint64Decor
virtual StatusCode execute(const EventContext &ctx) const override
virtual StatusCode finalize() override
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...
virtual bool isValid() override final
Can the handle be successfully dereferenced?
Handle class for adding a decoration to an object.
This class is the pure abstract base class for all fittable tracking measurements.
const LocalParameters & localParameters() const
Interface method to get the LocalParameters.
virtual const Surface & associatedSurface() const =0
Interface method to get the associated Surface.
const Amg::MatrixX & localCovariance() const
Interface method to get the localError.
const Amg::Vector3D & momentum() const
Access method for the momentum.
Class to handle RIO On Tracks ROT) for InDet and Muons, it inherits from the common MeasurementBase.
Definition RIO_OnTrack.h:70
virtual const Trk::PrepRawData * prepRawData() const =0
returns the PrepRawData (also known as RIO) object to which this RIO_OnTrack is associated.
@ Biased
RP with track state including the hit.
@ Unbiased
RP with track state that has measurement not included.
const TrkDetElementBase * associatedDetectorElement() const
return associated Detector Element
@ Measurement
This is a measurement, and will at least contain a Trk::MeasurementBase.
@ Outlier
This TSoS contains an outlier, that is, it contains a MeasurementBase/RIO_OnTrack which was not used ...
@ Hole
A hole on the track - this is defined in the following way.
virtual Identifier identify() const =0
Identifier.
int r
Definition globals.cxx:22
Eigen::Matrix< double, 3, 1 > Vector3D
void createDecoratorKeys(T_Parent &parent, const SG::ReadHandleKey< T_Cont > &container_key, const std::string &prefix, const std::vector< std::string > &decor_names, std::vector< SG::WriteDecorHandleKey< T_Cont > > &decor_out)
std::vector< SG::WriteDecorHandle< T_Cont, T > > createDecorators(const std::vector< SG::WriteDecorHandleKey< T_Cont > > &keys, const EventContext &ctx)
Definition HitInfo.h:33
@ loc2
generic first and second local coordinate
Definition ParamDefs.h:35
@ loc1
Definition ParamDefs.h:34
ParametersBase< TrackParametersDim, Charged > TrackParameters
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".
implementation file for function of same name