ATLAS Offline Software
Loading...
Searching...
No Matches
TRTTrackHoleSearchTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5// TRTTrackHoleSearchTool.cxx
6// author: Ryan D. Reece <ryan.reece@cern.ch>
7// created: Nov 2009
8
10
11// STD
12#include <bitset>
13#include <cstdlib> // abs
14#include <cmath> // fabs
15#include <fstream>
16
17// athena
18#include "TrkTrack/Track.h"
20
23#include "TrkSurfaces/Surface.h"
24#include "Identifier/Identifier.h"
30#include "CLHEP/Units/SystemOfUnits.h"
31#include "CLHEP/Geometry/Transform3D.h"
36
37
38//____________________________________________________________________________
39TRTTrackHoleSearchTool::TRTTrackHoleSearchTool(const std::string& type, const std::string& name, const IInterface* parent)
40 : AthAlgTool(type, name, parent)
41{
42 declareInterface<ITrackHoleSearchTool>(this);
43}
44
45
46//____________________________________________________________________________
48 ATH_MSG_DEBUG( "TRTTrackHoleSearchTool::initialize()" );
49
50 // retrieve extrapolator
51 ATH_CHECK(m_extrapolator.retrieve());
52
53 // retrieve TRT_ID
54 ATH_CHECK(detStore()->retrieve(m_TRT_ID, "TRT_ID"));
55
56 // retrieve ConditionsSummarySvc
57 ATH_CHECK(m_conditions_svc.retrieve());
58
59 m_trt_outer_surf = new Trk::CylinderSurface(Amg::Transform3D(Amg::Transform3D::Identity()),
61 // note: HepGeom::Translate3D is deleted by Trk::Surface destructor
62
63 return StatusCode::SUCCESS;
64}
65
66
67//____________________________________________________________________________
69 ATH_MSG_DEBUG( "TRTTrackHoleSearchTool::finalize()" );
70 delete m_trt_outer_surf;
71 return StatusCode::SUCCESS;
72}
73
74
75//____________________________________________________________________________
77 std::vector<int>& information ,
78 const Trk::ParticleHypothesis partHyp) const {
79 ATH_MSG_DEBUG( "TRTTrackHoleSearchTool::countHoles" );
80 const Trk::TrackStates* holes = getHolesOnTrack(track, partHyp);
81 if (holes) {
82 information[Trk::numberOfTRTHoles] = holes->size();
83 delete holes;
84 } else {
85 information[Trk::numberOfTRTHoles] = -1;
86 }
87}
88
89
90//____________________________________________________________________________
92 const Trk::Track& track,
93 const Trk::ParticleHypothesis partHyp) const {
94 ATH_MSG_DEBUG( "TRTTrackHoleSearchTool::getHolesOnTrack" );
95
96 // write out list of bad straws for debugging purposes.
97 // only dump log on first call of this function.
100 }
101 m_has_been_called = true;
102
103 if (track.perigeeParameters()) {
104 ATH_MSG_VERBOSE( " This track has perigee: \n"
105 << " pT = " << track.perigeeParameters()->pT()/CLHEP::GeV << " CLHEP::GeV\n"
106 << " eta = " << track.perigeeParameters()->eta() << "\n"
107 << " phi = " << track.perigeeParameters()->parameters()[Trk::phi0] << "\n"
108 << " d0 = " << track.perigeeParameters()->parameters()[Trk::d0] << "\n"
109 << " z0 = " << track.perigeeParameters()->parameters()[Trk::z0] );
110 }
111
112 // get TrackStateOnSurfaces
113 const Trk::TrackStates* track_states = track.trackStateOnSurfaces();
114 if (track_states) {
115 ATH_MSG_DEBUG( " This track has " << track_states->size() << " track states on surface." );
116 } else {
117 ATH_MSG_WARNING( " This track has null track states on surface. Returning 0." );
118 return nullptr;
119 }
120
121 if (track_states->size() < 2) {
122 ATH_MSG_WARNING( " Fewer than 2 TrackStatesOnSurface. Returning 0." );
123 return nullptr;
124 }
125
126 // set beginning point of extrapolation
127 Trk::TrackStates::const_iterator beginning_track_state;
129 // begin at first TRT hit
130 beginning_track_state = find_first_trt_hit(*track_states);
131 } else {
132 // begin at last Si hit
133 beginning_track_state = find_last_hit_before_trt(*track_states);
134 if (beginning_track_state == track_states->end()) {
135 ATH_MSG_WARNING( " beginning_track_state == track_states->end(). There must be no Si hits.\n"
136 << " Will try to begin at the first TRT hit." );
137 beginning_track_state = find_first_trt_hit(*track_states);
138 }
139 }
140
141 if (beginning_track_state == track_states->end()) {
142 ATH_MSG_WARNING( " beginning_track_state == track_states->end(). No where to extrapolate to. Returning 0." );
143 return nullptr;
144 }
145
146 // to be returned:
148
149 Trk::TrackStates::const_iterator track_state = beginning_track_state;
150 const Trk::TrackParameters* start_parameters = (*track_state)->trackParameters();
151 ++track_state;
152 // loop over TrackStateOnSurfaces (destination surfaces for extrapolation)
153 for (; track_state != track_states->end(); ++track_state) {
154 // skip track states that are not measurements
155 if (!(*track_state)->type(Trk::TrackStateOnSurface::Measurement)) {
156 ATH_MSG_VERBOSE( " TrackStateOnSurface is not of type Trk::TrackStateOnSurface::Measurement." );
157 continue;
158 }
159
160 const Trk::TrackParameters* end_parameters = (*track_state)->trackParameters();
161 if(end_parameters) {
162 ATH_MSG_VERBOSE( " TrackStateOnSurface: (x, y, z) = ("
163 << end_parameters->position().x() << ", "
164 << end_parameters->position().y() << ", "
165 << end_parameters->position().z() << "); (rho, eta, phi) = ("
166 << end_parameters->position().perp() << ", "
167 << end_parameters->position().eta() << ", "
168 << end_parameters->position().phi() << ")");
169 } else {
170 ATH_MSG_WARNING( " TrackStateOnSurface has no TrackParameters." );
171 continue;
172 }
173
174 const Trk::Surface& end_surf = end_parameters->associatedSurface();
175 if (!start_parameters) {
176 start_parameters = end_parameters;
177 } else {
178 extrapolateBetweenHits(start_parameters, end_surf, holes, partHyp);
179 start_parameters = end_parameters;
180 }
181 } // end loop over TrackStateOnSurfaces
182
183 if( !m_end_at_last_trt_hit ) {
184 // final extrapolation to the edge of TRT to check for holes in outermost straws
185 int trailing_hole_count = extrapolateBetweenHits(start_parameters, *m_trt_outer_surf, holes, partHyp);
186
187 // remove trailing holes
188 if(trailing_hole_count > m_max_trailing_holes) {
189 ATH_MSG_DEBUG("There are " << trailing_hole_count << " trailing holes removed.");
190 for(int i=0; i < trailing_hole_count; i++) {
191 holes->pop_back();
192 /*
193 Note: DataVector::pop_back() deletes pointers if it owns them.
194 */
195 }
196 }
197 }
198
199 return holes;
200}
201
202
203//____________________________________________________________________________
205 const Trk::ParticleHypothesis partHyp) const {
206 ATH_MSG_DEBUG( "TRTTrackHoleSearchTool::getTrackWithHoles" );
207 const Trk::TrackStates* holes = getHolesOnTrack(track, partHyp);
208 const Trk::Track* new_track = addHolesToTrack(track, holes);
209 delete holes;
210 return new_track;
211}
212
213
214//____________________________________________________________________________
216 const Trk::ParticleHypothesis partHyp) const {
217 ATH_MSG_DEBUG( "TRTTrackHoleSearchTool::getTrackWithHolesAndOutliers" );
218 return getTrackWithHoles(track, partHyp);
219}
220
221//----------------------------------------------------------------------------
222// Private Methods
223//----------------------------------------------------------------------------
224
225//____________________________________________________________________________
227 const Trk::Surface& end_surf,
228 Trk::TrackStates* holes,
229 const Trk::ParticleHypothesis partHyp/*=Trk::pion*/) const {
230 int hole_count = 0;
231 // initialize previous id
232 const Trk::Surface& start_surf = start_parameters->associatedSurface();
233
234 Identifier previous_id = start_surf.associatedDetectorElementIdentifier();
235
236 // get end id
237 const Identifier end_id = end_surf.associatedDetectorElementIdentifier();
238
239 // look for holes
240 std::vector<std::unique_ptr<Trk::TrackParameters> > steps =
241 m_extrapolator->extrapolateStepwise(Gaudi::Hive::currentContext(),
242 *start_parameters,
243 end_surf,
245 partHyp);
246
247 if(steps.empty()) {
248 ATH_MSG_DEBUG("extrapolateBetweenHits: extrapolateStepwise returned null");
249 } else {
250 // loop over parameters from extrapolation
251 // note: the last element in the vector is always the track parameters at the destination surface
252 for(std::vector<std::unique_ptr<Trk::TrackParameters> >::const_iterator step = steps.begin(); step != steps.end()-1; ++step) {
253 // check for surface
254 const Trk::Surface& surf = (*step)->associatedSurface();
255
256 // get id
258
259 // check that the surface is a TRT straw
260 if(!m_TRT_ID->is_trt(id)) {
261 ATH_MSG_DEBUG("extrapolateBetweenHits: surf is not a TRT straw. Skipping.");
262 continue;
263 }
264
265 // check that we are not still on the previous id
266 if(id == previous_id) {
267 ATH_MSG_DEBUG("extrapolateBetweenHits: id == previous_id");
268 continue;
269 }
270
271 // check that we haven't hit the destination surface
272 if(id == end_id) {
273 ATH_MSG_DEBUG("extrapolateBetweenHits: id == end_id");
274 continue;
275 }
276
277 // ignore id 0xffffffff from m_trt_outer_surf
278 if(id == 0xffffffff) {
279 ATH_MSG_DEBUG("extrapolateBetweenHits: id == 0xffffffff. Skipping.");
280 continue;
281 }
282
283 // don't count bad straws as holes
284 if( (m_use_conditions_svc)&&(!m_conditions_svc->isGood(id)) ) {
285 ATH_MSG_DEBUG("extrapolateBetweenHits: ConditionsSvc says this straw is bad. Skipping.");
286 continue;
287 }
288
289 // check fiducial |locR| by sigma
290 const float locR = (*step)->parameters()[Trk::locR];
291 if(m_locR_sigma_cut > 0.0) {
292 //const Trk::MeasuredTrackParameters* meas = dynamic_cast< const Trk::MeasuredTrackParameters* >(*step);
293 //if(meas)
294 //{
295
296 const AmgSymMatrix(5)* merr = (*step)->covariance();
297
298 if(merr){
299 //const float locR_error = (*merr)(Trk::locR,Trk::locR);
300
301 const float locR_error = Amg::error(*merr,Trk::locR);
302
303 // hole must be within the straw by some tolerance
304 if( (2.0*CLHEP::mm - fabs(locR))/locR_error < m_locR_sigma_cut ) {
305 continue;
306 }
307 } else {
308 ATH_MSG_WARNING("extrapolateBetweenHits: Track parameters failed to dynamic_cast< const Trk::MeasuredTrackParameters* >.");
309 continue;
310 }
311 }
312 // check fiducial |locR| by value
313 if(m_locR_cut > 0.0) {
314 if( fabs(locR) > m_locR_cut )
315 continue;
316 }
317 // check fiducial |locZ| by value
318 if(m_locZ_cut > 0.0) {
319 const float locZ = (*step)->parameters()[Trk::locZ];
320 const Trk::CylinderBounds* cylb = dynamic_cast<const Trk::CylinderBounds *>(&surf.bounds());
321 if(cylb) {
322 const float halfz=cylb->halflengthZ();
323 if( fabs(locZ) > halfz - m_locZ_cut ) continue;
324 } else {
325 ATH_MSG_WARNING("extrapolateBetweenHits: Surface failed to dynamic_cast to Trk::CylinderBounds. Skipping.");
326 continue;
327 }
328 }
329
330 // if we've gotten here, it is a legitimate hole
331 ATH_MSG_DEBUG(" HOLE Found! Identifier = " << id.getString() << " indicates:\n"
332 << " is_pixel = " << m_TRT_ID->is_pixel(id) << "\n"
333 << " is_sct = " << m_TRT_ID->is_sct(id) << "\n"
334 << " is_trt = " << m_TRT_ID->is_trt(id) << "\n"
335 << " barrel_ec = " << m_TRT_ID->barrel_ec(id) << "\n"
336 << " phi_module = " << m_TRT_ID->phi_module(id) << "\n"
337 << " layer_or_wheel = " << m_TRT_ID->layer_or_wheel(id) << "\n"
338 << " straw_layer = " << m_TRT_ID->straw_layer(id) << "\n"
339 << " straw = " << m_TRT_ID->straw(id) << "\n"
340 << " rho = " << (*step)->position().perp() );
341
342 std::bitset<Trk::TrackStateOnSurface::NumberOfTrackStateOnSurfaceTypes> typePattern;
343 typePattern.set(Trk::TrackStateOnSurface::Hole);
344 holes->push_back( new Trk::TrackStateOnSurface(nullptr, (*step)->uniqueClone(), nullptr, typePattern) );
345 hole_count++;
346 previous_id = id;
347 } // end loop over parameters from extrapolation
348 }
349 return hole_count;
350}
351
352
353//____________________________________________________________________________
355 ATH_MSG_DEBUG( "TRTTrackHoleSearchTool::dump_bad_straw_log" );
356 std::ofstream out("TRT_ConditionsSummarySvc_bad_straws.log");
357 out << "# id barrel_ec phi_module layer_or_wheel straw_layer straw" << std::endl;
358 for(std::vector<Identifier>::const_iterator it = m_TRT_ID->straw_layer_begin(); it != m_TRT_ID->straw_layer_end(); ++it) {
359 for(int i=0; i<= m_TRT_ID->straw_max(*it); ++i) {
360 Identifier id = m_TRT_ID->straw_id(*it, i);
361 if(!m_conditions_svc->isGood(id)) {
362 out << id.getString()
363 << std::setw(3) << m_TRT_ID->barrel_ec(id)
364 << std::setw(4) << m_TRT_ID->phi_module(id)
365 << std::setw(4) << m_TRT_ID->layer_or_wheel(id)
366 << std::setw(4) << m_TRT_ID->straw_layer(id)
367 << std::setw(4) << m_TRT_ID->straw(id)
368 << std::endl;
369 }
370 }
371 }
372 out.close();
373}
374
375
376//____________________________________________________________________________
379 ATH_MSG_DEBUG( "TRTTrackHoleSearchTool::find_first_trt_hit" );
380 // finds the first TRT hit or returns end
381 Trk::TrackStates::const_iterator track_state = track_states.begin();
382 for(; track_state != track_states.end(); ++track_state) {
383 // skip track states that are not measurements
384 if(!(*track_state)->type(Trk::TrackStateOnSurface::Measurement)) {
385 ATH_MSG_VERBOSE( " TrackStateOnSurface is not of type Trk::TrackStateOnSurface::Measurement." );
386 continue;
387 }
388 if( dynamic_cast<const InDet::TRT_DriftCircleOnTrack*>( (*track_state)->measurementOnTrack() ) ) {
389 ATH_MSG_DEBUG(" Found first TRT hit.");
390 break;
391 }
392 }
393 return track_state;
394}
395
396
397//____________________________________________________________________________
400 ATH_MSG_DEBUG( "TRTTrackHoleSearchTool::find_last_hit_before_trt" );
401 // start at first TRT hit
402 Trk::TrackStates::const_iterator track_state = find_first_trt_hit(track_states);
403 // if there is no room to step backwards in the hits, return end
404 if(track_states.size() < 2 || track_state == track_states.begin()) {
405 return track_states.end();
406 }
407 --track_state; // step back and look for last measurement before the TRT hit
408 for(; track_state != track_states.begin(); --track_state) {
409 if((*track_state)->type(Trk::TrackStateOnSurface::Measurement)) {
410 break;
411 }
412 }
413 // if the loop ended without finding a measurement, return end
414 if(!(*track_state)->type(Trk::TrackStateOnSurface::Measurement)) {
415 return track_states.end();
416 }
417 return track_state;
418}
419
420//____________________________________________________________________________
421const Trk::Track*
423 const Trk::Track& track,
424 const Trk::TrackStates* holes) const
425{
426 ATH_MSG_DEBUG("TRTTrackHoleSearchTool::addHolesToTrack");
427 /*
428 This method was basically coppied from here:
429 http://alxr.usatlas.bnl.gov/lxr-stb4/source/atlas/InnerDetector/InDetRecTools/InDetTrackHoleSearch/src/InDetTrackHoleSearchTool.cxx#931
430 */
431
432 // get states from track
433 auto tsos = std::make_unique<Trk::TrackStates>();
434 for (const auto *it : *track.trackStateOnSurfaces()) {
435 // veto old holes
436 if (!it->type(Trk::TrackStateOnSurface::Hole)) {
437 tsos->push_back(new Trk::TrackStateOnSurface(*it));
438 }
439 }
440
441 // if we have no holes on the old track and no holes found by search, then we
442 // just copy the track
443 if (track.trackStateOnSurfaces()->size() == tsos->size() && holes->empty()) {
444 // create copy of track
445 const Trk::Track* new_track = new Trk::Track(
446 track.info(),
447 std::move(tsos),
448 track.fitQuality() ? track.fitQuality()->uniqueClone() : nullptr);
449 return new_track;
450 }
451
452 // add new holes
453 tsos->insert(tsos->end(), holes->begin(), holes->end());
454
455 // sort
456 const Trk::TrackParameters* perigee = track.perigeeParameters();
457 if (!perigee)
458 perigee = (*(track.trackStateOnSurfaces()->begin()))->trackParameters();
459
460 if (perigee) {
462
463 if (fabs(perigee->parameters()[Trk::qOverP]) > 0.002) {
464 /* invest n*(logN)**2 sorting time for lowPt, coping with a possibly
465 not 100% transitive comparison functor. */
466 if (msgLvl(MSG::DEBUG)) {
467 msg() << "sorting vector with stable_sort" << endmsg;
468 }
469 std::stable_sort(tsos->begin(), tsos->end(), CompFunc);
470 } else {
471 tsos->sort(CompFunc); // respects DV object ownership
472 }
473 }
474
475 // create copy of track
476 const Trk::Track* new_track =
477 new Trk::Track(track.info(),
478 std::move(tsos),
479 track.fitQuality() ? track.fitQuality()->uniqueClone() : nullptr);
480
481 return new_track;
482}
483
484// EOF
#define endmsg
#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)
#define AmgSymMatrix(dim)
This is an Identifier helper class for the TRT subdetector.
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
const ServiceHandle< StoreGateSvc > & detStore() const
bool msgLvl(const MSG::Level lvl) const
MsgStream & msg() const
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
size_type size() const noexcept
Returns the number of elements in the collection.
Represents 'corrected' measurements from the TRT (for example, corrected for wire sag).
void countHoles(const Trk::Track &track, std::vector< int > &information, const Trk::ParticleHypothesis partHyp=Trk::pion) const
Input : track, partHyp Output: Changes in information This method fills the fields relevant to the ho...
Trk::TrackStates::const_iterator find_last_hit_before_trt(const Trk::TrackStates &track_states) const
BooleanProperty m_begin_at_first_trt_hit
const Trk::Track * getTrackWithHolesAndOutliers(const Trk::Track &track, const Trk::ParticleHypothesis partHyp=Trk::pion) const
Input : track, parthyp Return: A pointer to a new Trk::Track which containes the information of the i...
Trk::TrackStates::const_iterator find_first_trt_hit(const Trk::TrackStates &track_states) const
TRTTrackHoleSearchTool(const std::string &type, const std::string &name, const IInterface *parent)
ServiceHandle< IInDetConditionsSvc > m_conditions_svc
const Trk::Track * addHolesToTrack(const Trk::Track &track, const Trk::TrackStates *holes) const
BooleanProperty m_do_dump_bad_straw_log
BooleanProperty m_end_at_last_trt_hit
ToolHandle< Trk::IExtrapolator > m_extrapolator
const Trk::TrackStates * getHolesOnTrack(const Trk::Track &track, const Trk::ParticleHypothesis partHyp=Trk::pion) const
Input : track, parthyp Return: A DataVector containing pointers to TrackStateOnSurfaces which each re...
const Trk::Track * getTrackWithHoles(const Trk::Track &track, const Trk::ParticleHypothesis partHyp=Trk::pion) const
Input : track, parthyp Return: A pointer to a new Trk::Track which containes the information of the i...
int extrapolateBetweenHits(const Trk::TrackParameters *start_parameters, const Trk::Surface &end_surf, Trk::TrackStates *holes, const Trk::ParticleHypothesis partHyp=Trk::pion) const
Trk::CylinderSurface * m_trt_outer_surf
Bounds for a cylindrical Surface.
double halflengthZ() const
This method returns the halflengthZ.
Class for a CylinderSurface in the ATLAS detector.
const Amg::Vector3D & momentum() const
Access method for the momentum.
const Amg::Vector3D & position() const
Access method for the position.
virtual const Surface & associatedSurface() const override=0
Access to the Surface associated to the Parameters.
Abstract Base Class for tracking surfaces.
Identifier associatedDetectorElementIdentifier() const
return Identifier of the associated Detector Element
virtual const SurfaceBounds & bounds() const =0
Surface Bounds method.
Class providing comparison function, or relational definition, for sorting MeasurementBase objects.
represents the track state (measurement, material, fit parameters and quality) at a surface.
@ Measurement
This is a measurement, and will at least contain a Trk::MeasurementBase.
@ Hole
A hole on the track - this is defined in the following way.
double error(const Amg::MatrixX &mat, int index)
return diagonal error of the matrix caller should ensure the matrix is symmetric and the index is in ...
Eigen::Affine3d Transform3D
@ alongMomentum
DataVector< const Trk::TrackStateOnSurface > TrackStates
@ locR
Definition ParamDefs.h:44
@ phi0
Definition ParamDefs.h:65
@ qOverP
perigee
Definition ParamDefs.h:67
@ locZ
local cylindrical
Definition ParamDefs.h:42
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64
ParticleHypothesis
Enumeration for Particle hypothesis respecting the interaction with material.
ParametersBase< TrackParametersDim, Charged > TrackParameters
@ numberOfTRTHoles
number of TRT hits which pass the high threshold (only xenon counted) total number of TRT hits which ...
void stable_sort(DataModel_detail::iterator< DVL > beg, DataModel_detail::iterator< DVL > end)
Specialization of stable_sort for DataVector/List.