ATLAS Offline Software
Loading...
Searching...
No Matches
TRT_SegmentToTrackTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
10
11#include <cmath>
12
13using Amg::Vector3D;
14using CLHEP::mm;
15
16namespace {
17 Amg::Vector2D project2D(const Amg::Vector3D &vec) {
18 Amg::Vector2D tmp{vec[0],vec[1] };
19 return tmp;
20 }
21}
22namespace InDet {
23
24 // Constructor with parameters:
26 const std::string &name,
27 const IInterface *parent) :
28 AthAlgTool(type,name,parent)
29 {
30 declareInterface<InDet::ITRT_SegmentToTrackTool>( this );
31 }
32
34 = default;
35
36
38
39 ATH_CHECK( AthAlgTool::initialize() );
40
41 ATH_CHECK(m_extrapolator.retrieve() );
42
43 ATH_CHECK(m_fitterTool.retrieve(DisableTool{ !m_doRefit }));
44 ATH_CHECK(m_assoTool.retrieve( DisableTool{m_assoTool.name().empty()} ));
45 ATH_CHECK(m_trackSummaryTool.retrieve( DisableTool{m_trackSummaryTool.name().empty()} ));
46
47 // Get the scoring tool
48 ATH_CHECK( m_scoringTool.retrieve() );
49
50 ATH_CHECK( detStore()->retrieve(m_trtId, "TRT_ID") );
54
55 // Get output print level
56 //
57 ATH_MSG_DEBUG( *this );
58
59 return StatusCode::SUCCESS;
60
61 }
62
64 StatusCode sc = AthAlgTool::finalize();
65 return sc;
66 }
67
68
70 // Dumps relevant information into the MsgStream
72 MsgStream& InDet::TRT_SegmentToTrackTool::dump( MsgStream& out ) const {
73 out<<std::endl;
74 return dumpconditions(out);
75 }
76
78 // Dumps conditions information into the MsgStream
80
81 MsgStream& InDet::TRT_SegmentToTrackTool::dumpconditions( MsgStream& out ) const
82 {
83 int n = 65-m_fitterTool.type().size();
84 std::string s1; for(int i=0; i<n; ++i) s1.append(" "); s1.append("|");
85 n = 65-m_assoTool.type().size();
86 std::string s2; for(int i=0; i<n; ++i) s2.append(" "); s2.append("|");
87 n = 65-m_scoringTool.type().size();
88 std::string s5; for(int i=0; i<n; ++i) s5.append(" "); s5.append("|");
89
90
91 out<<std::endl
92 <<"|----------------------------------------------------------------------"
93 <<"-------------------|"<<std::endl;
94 out<<"| Tool for final track refitting | "<<m_fitterTool.type() <<s1<<std::endl;
95 out<<"| Tool for track scoring | "<<m_scoringTool.type() <<s5<<std::endl;
96 out<<"| Association service | "<<m_assoTool.type() <<s2<<std::endl;
97 out<<"|----------------------------------------------------------------------"
98 <<"-------------------|";
99 return out;
100 }
101
102
104 // Dumps event information into the MsgStream
106 MsgStream& InDet::TRT_SegmentToTrackTool::dumpevent( MsgStream& out ) {
107 return out;
108 }
109
111 // Dumps relevant information into the ostream
113
114 std::ostream& InDet::TRT_SegmentToTrackTool::dump( std::ostream& out ) const
115 {
116 return out;
117 }
118
119
120 Trk::Track* TRT_SegmentToTrackTool::segToTrack(const EventContext& ctx, const Trk::TrackSegment& tS) const {
121
122 ATH_MSG_DEBUG("Transforming the TRT segment into a track...");
123
124 //
125 // Get the track segment fit quality. If not there drop segment
126 //
127 if (!tS.fitQuality()) {
128 ATH_MSG_DEBUG("Segment has no fit quality ! Discard...");
129 return nullptr;
130 }
131 auto fq = tS.fitQuality()->uniqueClone();
132
133 //
134 // Get the track segment information about the initial track parameters
135 //
136 const AmgVector(5)& par = tS.localParameters();
138 const Trk::TrackParameters* segPar =
141 par[Trk::loc2],
142 par[Trk::phi],
143 par[Trk::theta],
144 par[Trk::qOverP],
145 std::move(ep))
146 .release();
147
148 if (segPar) {
149 ATH_MSG_VERBOSE("Initial TRT Segment Parameters : " << (*segPar));
150 } else {
151 ATH_MSG_DEBUG("Could not get initial TRT segment parameters! ");
152 // clean up
153 return nullptr;
154 }
155
156 // --- create new track state on surface vector
157 auto ntsos = std::make_unique<Trk::TrackStates>();
158
159 //
160 // if no refit, make it a perigee
161 //
162 if (!m_doRefit) {
163
164 const Trk::TrackStateOnSurface* par_tsos = nullptr;
165
166 // --- create surface at perigee
167 Amg::Vector3D perigeePosition(0., 0., 0.);
168 Trk::PerigeeSurface perigeeSurface(perigeePosition);
169 // --- turn parameters into perigee...
170 std::unique_ptr<Trk::TrackParameters> tmp =
171 m_extrapolator->extrapolate(ctx, *segPar, perigeeSurface);
172 std::unique_ptr<Trk::Perigee> perParm = nullptr;
173 //pass ownership if of the right type
174 if (tmp && tmp->associatedSurface().type() == Trk::SurfaceType::Perigee) {
175 perParm.reset(static_cast<Trk::Perigee*>(tmp.release()));
176 }
177 if (perParm) {
178 ATH_MSG_VERBOSE("Perigee version of Parameters : " << (*segPar));
179 } else {
180 ATH_MSG_DEBUG("Failed to build perigee parameters.Discard...");
181 ntsos->clear();
182 delete segPar;
183 segPar = nullptr;
184 return nullptr;
185 }
186
187 // now create a perigee TSOS
188 std::bitset<Trk::TrackStateOnSurface::NumberOfTrackStateOnSurfaceTypes>
189 typePattern;
190 typePattern.set(Trk::TrackStateOnSurface::Perigee);
191 par_tsos = new Trk::TrackStateOnSurface(
192 nullptr, std::move(perParm), nullptr, typePattern);
193 // push new TSOS into the list
194 ntsos->push_back(par_tsos);
195 }
196
197 //
198 // now loop over the TSOS and copy them in
199 //
200
201 // psuedo measurement information
202 int npseudo = 0;
203 double pseudotheta = 0;
204 const Trk::MeasurementBase* pseudo = nullptr;
205 // other measurement information
206 const Trk::Surface *firstsurf = nullptr, *firstecsurf = nullptr,
207 *lastsurf = nullptr;
208 const Trk::MeasurementBase* firstmeas = nullptr;
209 // counters for barrel and endcap
210 int nbarrel = 0, nendcap = 0;
211 // some variables for endcaps
212 std::vector<std::pair<double, double>> points;
213 points.reserve(40);
214 double oldphi = 0;
215
216 const Trk::Surface *surf_zmin=nullptr;
217 const Trk::Surface *surf_zmax=nullptr;
218 int meas_zmin_i=int(tS.numberOfMeasurementBases());
219 int meas_zmax_i=int(tS.numberOfMeasurementBases());
220
221 // loop over the measurements in track segment (tS)
222 for (int it = 0; it < int(tS.numberOfMeasurementBases()); it++) {
223
224 // the track state on service we like to constuct ...
225 const Trk::TrackStateOnSurface* seg_tsos = nullptr;
226
227 // is this ROT a psuedo-measurement ?
228 if (dynamic_cast<const Trk::PseudoMeasurementOnTrack*>(
229 tS.measurement(it))) {
230 // did we have a speudo-measurement before ?
231 if (pseudo) {
232 // get theta from pseudo measurements
233 pseudotheta =
234 std::atan2(tS.measurement(it)->associatedSurface().center().perp(),
235 tS.measurement(it)->associatedSurface().center().z());
236 }
237 // keep this pseudo measurement
238 pseudo = tS.measurement(it);
239 // update counter
240 npseudo++;
241
242 if (m_doRefit) {
243 // refit means we can simply copy the state, otherwise we skip it
244 seg_tsos =
245 new Trk::TrackStateOnSurface(tS.measurement(it)->uniqueClone(), nullptr);
246 }
247
248 } else {
249 //
250 // normal measurement, not a pseudo measurement
251 //
252 // copy measurement
253 seg_tsos =
254 new Trk::TrackStateOnSurface(tS.measurement(it)->uniqueClone(), nullptr);
255 //
256 // --- following is for the hack below
257 //
258
259 // remember first measurement
260 if (!firstmeas)
261 firstmeas = tS.measurement(it);
262 if (!firstsurf)
263 firstsurf = &tS.measurement(it)->associatedSurface();
264 // it is always the last one
265 lastsurf = &tS.measurement(it)->associatedSurface();
266
267 // detect endcap elemeents
269 // increase counter and keep some information
270 nendcap++;
271 if (!firstecsurf)
272 firstecsurf = &tS.measurement(it)->associatedSurface();
273
274 double tmpphi =
275 tS.measurement(it)->associatedSurface().center().phi();
276 if (!points.empty() &&
277 std::abs(tmpphi - oldphi) > M_PI) { // correct for boundary at +/- pi
278 if (tmpphi < 0)
279 tmpphi += 2 * M_PI;
280 else
281 tmpphi -= 2 * M_PI;
282 }
283 // remember oldphi
284 oldphi = tmpphi;
285
286 if (!surf_zmax || std::abs(tS.measurement(it)->associatedSurface().center().z())>std::abs(surf_zmax->center().z())) {
287 surf_zmax = &tS.measurement(it)->associatedSurface();
288 meas_zmax_i = it;
289 }
290 if (!surf_zmin || std::abs(tS.measurement(it)->associatedSurface().center().z())<std::abs(surf_zmin->center().z())) {
291 surf_zmin = &tS.measurement(it)->associatedSurface();
292 meas_zmin_i = it;
293 }
294 // copy the points
295 points.emplace_back(
296 tS.measurement(it)->associatedSurface().center().z(), tmpphi);
297
298 } else
299 nbarrel++;
300
301 //
302 // --- end of hack stuff
303 //
304 }
305
306 // push new TSOS into the list
307 if (seg_tsos)
308 ntsos->push_back(seg_tsos);
309 }
310
311 // Construct the new track
312 Trk::TrackInfo info;
313 info.setPatternRecognitionInfo(Trk::TrackInfo::TRTStandalone);
314
315 // create new track candidate
316 if (!m_doRefit) {
317 return new Trk::Track(info, std::move(ntsos), std::move(fq));
318 } else {
319 //
320 // ----------------------------- this is a horrible hack to make the
321 // segments fittable
322 //
323
324 // in case of only 1 pseudo measurement, use the theta from it.
325 if (npseudo == 1)
326 pseudotheta = pseudo->localParameters()[Trk::theta];
327
328 // we need new perigee parameters
329 double myqoverp = 0, myphi = 0, myd0 = 0, myz0 = 0, mytheta = pseudotheta;
330
331 if (nendcap < 4) {
332 //
333 // --- are we in the barrel mostly
334 //
335
336 // momentum
337 myqoverp = par[4] * std::sin(pseudotheta) / std::sin(par[3]);
338
339 // --- create surface at perigee
340 Amg::Vector3D perigeePosition(0., 0., 0.);
341 Trk::PerigeeSurface perigeeSurface(perigeePosition);
342 // -- get perigee
343 std::unique_ptr<const Trk::TrackParameters> tmp =
344 m_extrapolator->extrapolate(ctx, *segPar, perigeeSurface);
345 std::unique_ptr<const Trk::Perigee> tempper = nullptr;
346 if (tmp && tmp->associatedSurface().type() == Trk::SurfaceType::Perigee) {
347 tempper.reset(static_cast<const Trk::Perigee*>(tmp.release()));
348 }
349 if (!tempper) {
350 ATH_MSG_DEBUG("Could not produce perigee");
351 delete segPar;
352 segPar = nullptr;
353 return nullptr;
354 }
355
356 // keep some values
357 myd0 = tempper->parameters()[Trk::d0];
358 myphi = tempper->parameters()[Trk::phi0];
359
360 } else {
361 //
362 // --- endcap or transition track
363 //
364
365 // get estimate of parameters
366 double sx = 0, sy = 0, sxx = 0, sxy = 0, d = 0;
367 float zmin = points.empty() ? 0 : points.begin()->first, zmax = 0;
368 // loop over all points
369
370 for (auto & point : points) {
371 sx += point.first;
372 sy += point.second;
373 sxy += point.first * point.second;
374 sxx += point.first * point.first;
375 if (std::abs(point.first) > std::abs(zmax)) {
376 zmax = point.first;
377 }
378 if (std::abs(point.first) < std::abs(zmin)) {
379 zmin = point.first;
380 }
381 }
382
383 if (std::abs(pseudotheta) < 1.e-6) {
384 ATH_MSG_DEBUG("pseudomeasurements missing on the segment?");
385
386 if ( meas_zmin_i < int(tS.numberOfMeasurementBases())
387 && meas_zmax_i < int(tS.numberOfMeasurementBases())) {
388 Amg::Vector3D meas_zmin = tS.measurement(meas_zmin_i)->globalPosition();
389 Amg::Vector3D meas_zmax = tS.measurement(meas_zmax_i)->globalPosition();
390 double diff_z (meas_zmax.z() - meas_zmin.z());
391 if (std::abs(diff_z)>1e-6) {
392 double diff_r( meas_zmax.perp() - meas_zmin.perp());
393 pseudotheta = std::atan2(diff_r, diff_z);
394 ATH_MSG_DEBUG("Initial pseudo theta is zero. Compute pseudo from inner- and outermost endcap measurements. delta R: " <<
395 meas_zmax.perp() << " - " << meas_zmin.perp() << " = " << diff_r
396 << " , diff Z: " << meas_zmax.z() << " - " << meas_zmin.z()
397 << " = " << diff_z
398 << " " << pseudotheta);
399 }
400 }
401
402 if (std::abs(pseudotheta) < 1.e-6) {
403 constexpr std::array<float,2> boundary_r {644., 1004.};
404 // true if track extends from endcap A to endcap C:
405 bool is_across_sides=std::signbit(zmin)!=std::signbit(zmax);
406 double r_path=is_across_sides ? 2 * boundary_r[1] : boundary_r[1]-boundary_r[0];
407 double r1 =0.;
408 double r2 =0.;
409 if (surf_zmin && surf_zmax) {
410 bool is_reverse=std::abs(surf_zmin->center().z())>std::abs(surf_zmax->center().z());
411 // use the innermost and outermost surfaces (defined by z-coordinate) to compute the
412 // travel distance in R-direction assume that the outermost surface is always crossed
413 // at the outermost radial position assume more over that the the innermost surface
414 // is either crossed at the innermost or outermost position what ever leads to the
415 // largest travel distance in r.
416
417 std::array<const Trk::Surface *,2> boundary_surface{
418 is_reverse ? surf_zmax : surf_zmin,
419 is_reverse ? surf_zmin : surf_zmax};
420 // true if track goes from the outside to the insde :
421 std::array<Amg::Vector2D,2> translation;
422 std::array<Amg::Vector2D,2> height;
423 for (unsigned int boundary_i=boundary_surface.size(); boundary_i-->0;) {
424 if ( boundary_surface[boundary_i]->type() == Trk::SurfaceType::Line
425 && boundary_surface[boundary_i]->bounds().type() == Trk::SurfaceBounds::Cylinder) {
426 const Trk::CylinderBounds &
427 cylinder_bounds = dynamic_cast<const Trk::CylinderBounds &>(boundary_surface[boundary_i]->bounds());
428 height[boundary_i]= project2D( boundary_surface[boundary_i]->transform().rotation()
429 * Amg::Vector3D(0,0, cylinder_bounds.halflengthZ()) );
430 translation[boundary_i]=project2D( boundary_surface[boundary_i]->transform().translation() );
431 }
432 }
433 r1 = (translation[1] + height[1] - ( translation[0] + height[0])).perp();
434 r2 = (translation[1] + height[1] - ( translation[0] - height[0])).perp();
435 r_path=std::max(r1,r2);
436
437 }
438 pseudotheta = std::atan2(r_path, zmax - zmin);
439 ATH_MSG_DEBUG("Initial pseudo theta is zero. Pseudo theta from inner- and outermost surfaces. Deleta r " << r2 << " - " << r1 << " -> " << r_path
440 << " , delta Z " << zmax << " - " << zmin << " -> " << (zmax-zmin)
441 << " " << pseudotheta);
442 }
443 }
444
445 // get q/p
446 d = (points.size() * sxx - sx * sx);
447 double dphidz = ((points.size() * sxy - sy * sx) / d);
448 myqoverp = (std::abs(pseudotheta) > 1e-6)
449 ? -dphidz / (0.6 * std::tan(pseudotheta))
450 : 1000.;
451
452 // some geometry stuff to estimate further paramters...
453 double halfz = 200.;
454 const Trk::CylinderBounds* cylb =
455 dynamic_cast<const Trk::CylinderBounds*>(&firstsurf->bounds());
456 if (!cylb)
457 ATH_MSG_ERROR("Cast of bounds failed, should never happen");
458 else
459 halfz = cylb->halflengthZ();
460 const Trk::CylinderBounds* cylb2 =
461 dynamic_cast<const Trk::CylinderBounds*>(&lastsurf->bounds());
462 double halfz2 = 200.;
463 if (!cylb2)
464 ATH_MSG_ERROR("Cast of bounds failed, should never happen");
465 else
466 halfz2 = cylb2->halflengthZ();
467 Amg::Vector3D strawdir1 = -firstsurf->transform().rotation().col(2);
468 Amg::Vector3D strawdir2 = -lastsurf->transform().rotation().col(2);
469 Amg::Vector3D pos1;
470 Amg::Vector3D pos2;
471
472 // ME: this is hardcoding, not nice and should be fixed
473 if (std::abs(lastsurf->center().z()) < 2650 * mm) {
474 pos2 = lastsurf->center() + halfz2 * strawdir2;
475 if (nbarrel == 0) {
476 double dr = std::abs(std::tan(pseudotheta) * (lastsurf->center().z() -
477 firstsurf->center().z()));
478 pos1 = firstsurf->center() + (halfz - dr) * strawdir1;
479 } else {
480 double dz = std::abs((pos2.perp() - firstsurf->center().perp()) /
481 std::tan(pseudotheta));
482 if (pos2.z() > 0)
483 dz = -dz;
484 double z1 = pos2.z() + dz;
485 pos1 = Amg::Vector3D(
486 firstsurf->center().x(), firstsurf->center().y(), z1);
487 }
488 } else {
489 double dr = std::abs(std::tan(pseudotheta) *
490 (lastsurf->center().z() - firstsurf->center().z()));
491 pos2 = lastsurf->center() + (dr - halfz2) * strawdir2;
492 pos1 = firstsurf->center() - halfz * strawdir1;
493 }
494
495 // ME: I don't understand this yet, why is this done only if barrel ==
496 // 0, while above this nendcap < 4 ?
497 if (nbarrel == 0 &&
498 std::abs(std::tan(pseudotheta) * (firstsurf->center().z() -
499 lastsurf->center().z())) < 250 * mm &&
500 std::abs(firstsurf->center().z()) > 1000 * mm) {
501
502 // ME: wow this is hacking the vector ...
503 const Trk::MeasurementBase* firstmeas =
504 (**ntsos->begin()).measurementOnTrack();
505 Amg::MatrixX newcov(2, 2);
506 newcov.setZero();
507 newcov(0, 0) = (firstmeas->localCovariance())(0, 0);
508 newcov(1, 1) = (myqoverp != 0) ? .0001 * myqoverp * myqoverp : 1.e-8;
509 Trk::LocalParameters newpar(std::make_pair(0, Trk::locZ),
510 std::make_pair(myqoverp, Trk::qOverP));
511 auto newpseudo =
512 std::make_unique<Trk::PseudoMeasurementOnTrack>(
513 std::move(newpar),
514 std::move(newcov),
515 firstmeas->associatedSurface());
516 // hack replace first measurement with pseudomeasurement
517 ntsos->erase(ntsos->begin());
518 ntsos->insert(ntsos->begin(),
519 new Trk::TrackStateOnSurface(std::move(newpseudo), nullptr));
520 }
521
522 Amg::Vector3D field1;
523
524 MagField::AtlasFieldCache fieldCache;
525
526 // Get field cache object
529 };
530 const AtlasFieldCacheCondObj* fieldCondObj{ *readHandle };
531 if (fieldCondObj == nullptr) {
533 "segToTrack: Failed to retrieve AtlasFieldCacheCondObj with key "
535 return nullptr;
536 }
537 fieldCondObj->getInitializedCache(fieldCache);
538
539 // MT version uses cache
540 fieldCache.getField(Amg::Vector3D(.5 * (pos1 + pos2)).data(),
541 field1.data());
542
543 field1 *= m_fieldUnitConversion; // field in Tesla
544
545 double phideflection =
546 -.3 * (pos2 - pos1).perp() * field1.z() * myqoverp / std::sin(pseudotheta);
547 double precisephi = (nbarrel == 0)
548 ? (pos2 - pos1).phi() - .5 * phideflection
549 : (pos2 - pos1).phi() + .5 * phideflection;
550 double radius = (myqoverp != 0. && field1.z() != 0.)
551 ? -std::sin(pseudotheta) / (.3 * field1.z() * myqoverp)
552 : 1.e6;
553 double precisetheta =
554 (myqoverp != 0.)
555 ? std::atan2(std::abs(radius * phideflection), pos2.z() - pos1.z())
556 : pseudotheta;
557 if (precisetheta < 0)
558 precisetheta += M_PI;
559 if (precisephi < -M_PI)
560 precisephi += 2 * M_PI;
561 if (precisephi > M_PI)
562 precisephi -= 2 * M_PI;
563
564 // now extrapolate backwards from the first surface to get starting
565 // parameters
566 const Trk::StraightLineSurface* surfforpar = nullptr;
567 if (nbarrel == 0)
568 surfforpar = dynamic_cast<const Trk::StraightLineSurface*>(firstsurf);
569 else
570 surfforpar = dynamic_cast<const Trk::StraightLineSurface*>(lastsurf);
571 if (!surfforpar) {
572 ATH_MSG_ERROR("Cast of surface failed, should never happen");
573 return nullptr;
574 }
575
576 Trk::AtaStraightLine ataline(((nbarrel == 0) ? pos1 : pos2),
577 precisephi,
578 precisetheta,
579 myqoverp,
580 *surfforpar);
581 Trk::PerigeeSurface persurf;
582 const Trk::TrackParameters* extrappar =
583 m_extrapolator->extrapolateDirectly(ctx, ataline, persurf).release();
584
585 // now get parameters
586 if (extrappar) {
587 if (nendcap >= 4) {
588 myphi = extrappar->parameters()[Trk::phi0];
589 myd0 = extrappar->parameters()[Trk::d0];
590 }
591
592 // construct theta again
593 double z0 = extrappar->parameters()[Trk::z0];
594 if (nbarrel == 0)
595 mytheta = std::atan(std::tan(extrappar->parameters()[Trk::theta]) *
596 std::abs((z0 - pos1.z()) / pos1.z()));
597 else
598 mytheta = std::atan(std::tan(extrappar->parameters()[Trk::theta]) *
599 std::abs((z0 - pos2.z()) / pos2.z()));
600
601 if (mytheta < 0)
602 mytheta += M_PI;
603
604 delete extrappar;
605 extrappar = nullptr;
606 }
607 }
608 while (myphi > M_PI)
609 myphi -= 2 * M_PI;
610 while (myphi < -M_PI)
611 myphi += 2 * M_PI;
612
613 double P[5] = { myd0, myz0, myphi, mytheta, myqoverp };
614
615 // create perigee TSOS and add as first (!) TSOS
616
617 auto per =
618 std::make_unique<Trk::Perigee>(P[0], P[1], P[2], P[3], P[4], Trk::PerigeeSurface());
619 std::bitset<Trk::TrackStateOnSurface::NumberOfTrackStateOnSurfaceTypes>
620 typePattern;
621 typePattern.set(Trk::TrackStateOnSurface::Perigee);
623 nullptr, std::move(per), nullptr, typePattern);
624 ntsos->insert(ntsos->begin(), seg_tsos);
625
626 //
627 // ------------------------------------------------------- now refit the
628 // track
629 //
630
631 Trk::Track newTrack (info, std::move(ntsos), std::move(fq));
632 Trk::Track* fitTrack =
633 m_fitterTool->fit(ctx,newTrack, true, Trk::nonInteracting).release();
634
635 // cleanup
636 if (segPar) {
637 delete segPar;
638 segPar = nullptr;
639 }
640
641 if (!fitTrack) {
642 ATH_MSG_DEBUG("Refit of TRT track segment failed!");
643 return nullptr;
644 }
645
646 //
647 // -------------------------------------- hack the covarinace back to
648 // something reasonable
649 //
650 const Trk::TrackParameters* firstmeaspar = nullptr;
652 fitTrack->trackParameters()->begin();
653 do {
654 // skip pesudo measurements on perigee
655 if ((*parit)->covariance() &&
656 ((*parit)->associatedSurface() == tS.associatedSurface()))
657 firstmeaspar = *parit;
658 ++parit;
659 } while (firstmeaspar == nullptr &&
660 parit != fitTrack->trackParameters()->end());
661
662 // Create new perigee starting from the modified first measurement that
663 // has a more reasonable covariance matrix
664 // const Trk::Perigee* perTrack=dynamic_cast<const
665 // Trk::Perigee*>(fitTrack->perigeeParameters());
666 const Trk::Perigee* perTrack = fitTrack->perigeeParameters();
667
668 if (!perTrack || !perTrack->covariance()) {
669 ATH_MSG_ERROR("Cast of perigee fails, should never happen !");
670 return nullptr;
671 } else {
672 ATH_MSG_VERBOSE ("Perigee after refit with fudges to make it converge : " << (*perTrack) );
673
674 if(firstmeaspar && firstmeaspar->position().perp()<2000*mm && std::abs(firstmeaspar->position().z())<3000*mm){
675
676 // Modify first measurement so that it has reasonable errors on z and theta
677 AmgSymMatrix(5) fcovmat = AmgSymMatrix(5)(*(firstmeaspar->covariance()));
678 // factors by which we like to scale the cov, this takes the original segment errors into account
679 double scaleZ = std::sqrt(tS.localCovariance()(1,1))/std::sqrt( (fcovmat)(1,1));
680 double scaleTheta = std::sqrt(tS.localCovariance()(3,3))/std::sqrt( (fcovmat)(3,3));
681 // now do it
682 fcovmat(1,0)=scaleZ*((fcovmat)(1,0));
683 fcovmat(0,1) = (fcovmat)(1,0);
684 fcovmat(1,1)=tS.localCovariance()(1,1);
685 fcovmat(2,1)=scaleZ*((fcovmat)(2,1));
686 fcovmat(1,2) = (fcovmat)(2,1);
687 fcovmat(3,1)=scaleZ*scaleTheta*((fcovmat)(3,1));
688 fcovmat(1,3) = (fcovmat)(3,1);
689 fcovmat(4,1)=scaleZ*((fcovmat)(4,1));
690 fcovmat(1,4) = (fcovmat)(4,1);
691 fcovmat(3,0)=scaleTheta*((fcovmat)(3,0));
692 fcovmat(0,3) = (fcovmat)(3,0);
693 fcovmat(3,2)=scaleTheta*((fcovmat)(3,2));
694 fcovmat(2,3) = (fcovmat)(3,2);
695 fcovmat(3,3)=tS.localCovariance()(3,3);
696 fcovmat(4,3)=scaleTheta*((fcovmat)(4,3));
697 fcovmat(3,4) = (fcovmat)(4,3);
698
699 // const Amg::VectorX& par = firstmeaspar->parameters();
700 const AmgVector(5)& par = firstmeaspar->parameters();
701 const Trk::TrackParameters* updatedPars =
703 par[Trk::loc1],
704 par[Trk::loc2],
705 par[Trk::phi],
706 par[Trk::theta],
707 par[Trk::qOverP],
708 std::move(fcovmat)).release();
709
710 // now take parameters at first measurement and exptrapolate to perigee
711 const Trk::TrackParameters* newperpar =
712 m_extrapolator->extrapolate(ctx,
713 *updatedPars,
714 perTrack->associatedSurface(),
716 false,
717 Trk::nonInteracting).release();
718 delete updatedPars; updatedPars = nullptr;
719
720 if (!newperpar || !newperpar->covariance()) {
721 ATH_MSG_WARNING ("Can not hack perigee parameters, extrapolation failed");
722 delete newperpar; newperpar = nullptr;
723 } else {
724 // this is a HACK !!!
725 // perTrack is owned by fitTrack which is not const here
726 // thus the const-ness is only removed from somthting which is not strictly const here.
727 AmgSymMatrix(5)& errmat ATLAS_THREAD_SAFE = const_cast<AmgSymMatrix(5)&>(*perTrack->covariance());
728 // overwrite cov in perTrack
729 errmat = *newperpar->covariance();
730 delete newperpar; newperpar = nullptr;
731 // check that new cov makes sense !
732 const AmgSymMatrix(5)& CM = *perTrack->covariance();
733 if( CM(1,1)==0.||CM(3,3)==0. ) {
734 ATH_MSG_DEBUG ("Hacked perigee covariance is CRAP, reject track");
735 delete fitTrack; return nullptr;
736 } else {
737 ATH_MSG_VERBOSE ("Perigee after fit with scaled covariance matrix : " << *perTrack);
738 }
739 }
740 }
741 }
742 // return fitted track
743 return fitTrack;
744 }
745 }
746
747
749 const Trk::PRDtoTrackMap *prd_to_track_map) const {
750
751 ATH_MSG_DEBUG ("Checking whether the TRT segment has already been used...");
752
753 // some counters to be handled
754 int nShared = 0; // Shared drift circles in segment
755 int nHits = 0; // Number of TRT measurements
756
757 if (m_assoTool.isEnabled()&& !prd_to_track_map) ATH_MSG_ERROR("PRDtoTrackMap to be used but not provided by the client");
758 // loop over the track states
759 for(int it=0; it<int(tS.numberOfMeasurementBases()); ++it){
760
761 // remove pseudo measurements
762 if ( dynamic_cast<const Trk::PseudoMeasurementOnTrack*>(tS.measurement(it)) )
763 continue;
764
765 // get the measurment
766 const InDet::TRT_DriftCircleOnTrack* trtcircle = dynamic_cast<const InDet::TRT_DriftCircleOnTrack*>(tS.measurement(it));
767 if (!trtcircle) continue;
768
769 // get PRD measurement
770 const InDet::TRT_DriftCircle* RawDataClus=dynamic_cast<const InDet::TRT_DriftCircle*>(trtcircle->prepRawData());
771 if(!RawDataClus) continue;
772
773 // count up number of hits
774 nHits++;
775
776 if(m_assoTool.isEnabled() && prd_to_track_map && prd_to_track_map->isUsed(*RawDataClus)) nShared++;
777 }
778
779 if(nShared >= int(m_sharedFrac * nHits)) {
780 ATH_MSG_DEBUG ("Too many shared hits.Will drop the TRT segment");
781 return true;
782 } else {
783 return false;
784 }
785
786
787 }
788
790
791 ATH_MSG_DEBUG ("Try to recover low TRT DC segments in crack...");
792
793 // counters
794 int nEC = 0; int nBRL = 0; int firstWheel = -999; int lastLayer = -999;
795
796 // loop over the track states
797 for(int it=0; it<int(tS.numberOfMeasurementBases()); ++it){
798
799 //test if it is a pseudo measurement
800 if ( dynamic_cast<const Trk::PseudoMeasurementOnTrack*>(tS.measurement(it)) )
801 continue;
802
803 // get measurement
804 const InDet::TRT_DriftCircleOnTrack* trtcircle = dynamic_cast<const InDet::TRT_DriftCircleOnTrack*>(tS.measurement(it));
805 if(!trtcircle) continue;
806
807 // get identifier
808 Identifier id = trtcircle->detectorElement()->identify();
809 // barrel or endcap
810 int isB = m_trtId->barrel_ec(id);
811 if (isB==2 || isB==-2) {
812 nEC++;
813 if(nEC == 1)
814 firstWheel = m_trtId->layer_or_wheel(id);
815 }
816 else if (isB==1 || isB==-1) {
817 nBRL++;
818 lastLayer = m_trtId->layer_or_wheel(id);
819 }
820
821 }
822
823 // now the logic
824 return (nEC>0 && nBRL>0) ||
825
826 (nEC==0 && nBRL>0 && lastLayer<2) ||
827
828 (nEC>0 && nBRL==0 && (firstWheel>10 || firstWheel<2));
829
830 }
831
832 void TRT_SegmentToTrackTool::addNewTrack(const EventContext& ctx, Trk::Track* trk, ITRT_SegmentToTrackTool::EventData &event_data) const {
833 // @TODO avoid non const member m_trackScoreTrackMap
834 ATH_MSG_DEBUG ("Add track to the scoring multimap...");
835 //Score the track under investigation
837 bool passBasicSelections = m_scoringTool->passBasicSelections(*trk);
838 if(passBasicSelections){
839 if (m_trackSummaryTool.isEnabled()) {
840 m_trackSummaryTool->computeAndReplaceTrackSummary(ctx, *trk,
842 }
843 score = m_scoringTool->score(*trk);
844 }
845 ATH_MSG_DEBUG ("TRT-only: score is " << score);
846
847 if (score==0) {
848 // statistics...
849 ATH_MSG_DEBUG ("Track score is zero, reject it");
851 // clean up
852 delete trk;
853 } else {
854 // add track to map, map is sorted small to big !
855 event_data.m_trackScores.emplace_back(-score, trk );
856 }
857
858
859 }
860
861
863 ITRT_SegmentToTrackTool::EventData &event_data) const {
864
865 ATH_MSG_DEBUG ("Resolving the TRT tracks in score map...");
866
867
868 std::unique_ptr<Trk::PRDtoTrackMap> prd_to_track_map;
869 if (m_assoTool.isEnabled()) {
870 prd_to_track_map=m_assoTool->createPRDtoTrackMap();
871 if (prd_to_track_map_in) {
872 *prd_to_track_map = *prd_to_track_map_in;
873 }
874 }
875
876 //final copy - ownership is passed out of algorithm
877 std::unique_ptr<TrackCollection> final_tracks = std::make_unique<TrackCollection>();
878 final_tracks->reserve( event_data.m_trackScores.size());
879 std::stable_sort(event_data.m_trackScores.begin(),
880 event_data.m_trackScores.end(),
881 []( const std::pair< Trk::TrackScore, Trk::Track* > &a,
882 const std::pair< Trk::TrackScore, Trk::Track* > &b)
883 { return a.first < b.first; });
884
885 for (std::pair< Trk::TrackScore, Trk::Track* > &track_score : event_data.m_trackScores) {
886
887 ATH_MSG_DEBUG ("--- Trying next track "<<track_score.second<<"\t with score "<<-track_score.first);
888
889 // some counters to be handled
890 int nShared = 0; // Shared drift circles in segment
891 int nHits = 0; // Number of TRT measurements
892
893 // get vector of TSOS
894 const Trk::TrackStates* tsos = (track_score.second)->trackStateOnSurfaces();
895
896 // loop over vector of TSOS
897 for ( const Trk::TrackStateOnSurface *a_tsos : *tsos) {
898
899 // get measurment from TSOS
900 const Trk::MeasurementBase* meas = a_tsos->measurementOnTrack();
901 if (!meas)
902 continue;
903
904 // make sure it is a TRT_DC and not a pseudo measurement
906 dynamic_cast<const InDet::TRT_DriftCircleOnTrack*>(meas);
907 if (!rot)
908 continue;
909
910 // get to the PRD object
911 const InDet::TRT_DriftCircle* RawDataClus =
912 dynamic_cast<const InDet::TRT_DriftCircle*>(rot->prepRawData());
913 if (!RawDataClus)
914 continue;
915
916 // count up number of hits
917 nHits++;
918
919 // count up number of shared hits
920 if (m_assoTool.isEnabled() && prd_to_track_map->isUsed(*RawDataClus))
921 nShared++;
922 }
923
924 ATH_MSG_DEBUG ("TRT-only has " << nHits << " hits and " << nShared << " of them are shared");
925
926 // cut on the number of shared hits with the max fraction
927 if(nShared >= int(m_sharedFrac * nHits)) {
928 // statistics
930 ATH_MSG_DEBUG ("Too many shared hits, remove it !");
931 delete track_score.second;
932 continue;
933 }
934
935 // ok, this seems like a useful track
936 final_tracks->push_back(track_score.second);
937
939 ATH_MSG_DEBUG ("TRT-only is accepted");
940
941 //Register the track with the association tool
942 if(m_assoTool.isEnabled()) {
943 if(m_assoTool->addPRDs(*prd_to_track_map,*(track_score.second)).isFailure()) {
944 ATH_MSG_WARNING ("addPRDs() failed!");
945 }
946 }
947
948 }
949
950 return final_tracks.release();
951
952 }
953
954}
#define M_PI
Scalar perp() const
perp method - perpendicular length
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
std::vector< size_t > vec
#define AmgSymMatrix(dim)
#define AmgVector(rows)
char data[hepevt_bytes_allocation_ATLAS]
Definition HepEvt.cxx:11
static Double_t a
static Double_t P(Double_t *tt, Double_t *par)
static Double_t sc
static const uint32_t nHits
This is an Identifier helper class for the TRT subdetector.
DataVector< Trk::Track > TrackCollection
This typedef represents a collection of Trk::Track objects.
#define ATLAS_THREAD_SAFE
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
const ServiceHandle< StoreGateSvc > & detStore() const
void getInitializedCache(MagField::AtlasFieldCache &cache) const
get B field cache for evaluation as a function of 2-d or 3-d position.
DataModel_detail::const_iterator< DataVector > const_iterator
Standard const_iterator.
Definition DataVector.h:838
virtual Identifier identify() const override final
identifier of this detector element:
Represents 'corrected' measurements from the TRT (for example, corrected for wire sag).
virtual const InDetDD::TRT_BaseElement * detectorElement() const override final
returns the detector element, assoicated with the PRD of this class
virtual const TRT_DriftCircle * prepRawData() const override final
returns the PrepRawData - is a TRT_DriftCircle in this scope
SG::ReadCondHandleKey< AtlasFieldCacheCondObj > m_fieldCacheCondObjInputKey
virtual Trk::Track * segToTrack(const EventContext &ctx, const Trk::TrackSegment &) const override
virtual StatusCode initialize() override
ToolHandle< Trk::IExtrapolator > m_extrapolator
ToolHandle< Trk::IExtendedTrackSummaryTool > m_trackSummaryTool
static constexpr double m_fieldUnitConversion
virtual TrackCollection * resolveTracks(const Trk::PRDtoTrackMap *prd_to_track_map, ITRT_SegmentToTrackTool::EventData &event_data) const override
Resolve the standalone TRT tracks based on the number of shared TRT hits.
virtual bool segIsUsed(const Trk::TrackSegment &, const Trk::PRDtoTrackMap *prd_to_track_map) const override
Check if the TRT segment has already been assigned a Si extension.
static MsgStream & dumpevent(MsgStream &out)
Tracks that will be passed out of AmbiProcessor.
TRT_SegmentToTrackTool(const std::string &type, const std::string &name, const IInterface *parent)
const TRT_ID * m_trtId
ID TRT helper.
virtual StatusCode finalize() override
ToolHandle< Trk::ITrackScoringTool > m_scoringTool
ToolHandle< Trk::IPRDtoTrackMapTool > m_assoTool
ToolHandle< Trk::ITrackFitter > m_fitterTool
virtual bool toLower(const Trk::TrackSegment &) const override
MsgStream & dumpconditions(MsgStream &out) const
virtual MsgStream & dump(MsgStream &out) const override
virtual void addNewTrack(const EventContext &, Trk::Track *, ITRT_SegmentToTrackTool::EventData &event_data) const override
Add track into the track-score multimap.
Local cache for magnetic field (based on MagFieldServices/AtlasFieldSvcTLS.h).
void getField(const double *ATH_RESTRICT xyz, double *ATH_RESTRICT bxyz, double *ATH_RESTRICT deriv=nullptr)
get B field value at given position xyz[3] is in mm, bxyz[3] is in kT if deriv[9] is given,...
Bounds for a cylindrical Surface.
double halflengthZ() const
This method returns the halflengthZ.
std::unique_ptr< FitQuality > uniqueClone() const
NVI uniqueClone.
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.
virtual const Amg::Vector3D & globalPosition() const =0
Interface method to get the global Position.
const Amg::MatrixX & localCovariance() const
Interface method to get the localError.
std::unique_ptr< MeasurementBase > uniqueClone() const
NVI Clone giving up unique pointer.
bool isUsed(const PrepRawData &prd) const
does this PRD belong to at least one track?
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.
virtual const S & associatedSurface() const override final
Access to the Surface method.
Class describing the Line to which the Perigee refers to.
Class to handle pseudo-measurements in fitters and on track objects.
const MeasurementBase * measurement(unsigned int) const
returns the Trk::MeasurementBase objects depending on the integer
const FitQuality * fitQuality() const
return the FitQuality object, returns NULL if no FitQuality is defined
unsigned int numberOfMeasurementBases() const
Return the number of contained Trk::MeasurementBase (s).
Class for a StraightLineSurface in the ATLAS detector to describe dirft tube and straw like detectors...
Abstract Base Class for tracking surfaces.
Definition Surface.h:79
virtual ChargedTrackParametersUniquePtr createUniqueTrackParameters(double l1, double l2, double phi, double theat, double qop, std::optional< AmgSymMatrix(5)> cov=std::nullopt) const =0
Use the Surface as a ParametersBase constructor, from local parameters - charged.
const Amg::Transform3D & transform() const
Returns HepGeom::Transform3D by reference.
Identifier associatedDetectorElementIdentifier() const
return Identifier of the associated Detector Element
virtual const SurfaceBounds & bounds() const =0
Surface Bounds method.
const Amg::Vector3D & center() const
Returns the center position of the Surface.
Contains information about the 'fitter' of this track.
Class for a generic track segment that holdes polymorphic Trk::MeasurementBase objects,...
const Surface & associatedSurface() const override final
returns the surface for the local to global transformation
represents the track state (measurement, material, fit parameters and quality) at a surface.
@ Perigee
This represents a perigee, and so will contain a Perigee object only.
const DataVector< const TrackParameters > * trackParameters() const
Return a pointer to a vector of TrackParameters.
const Perigee * perigeeParameters() const
return Perigee.
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Eigen::Matrix< double, 2, 1 > Vector2D
Eigen::Matrix< double, 3, 1 > Vector3D
Primary Vertex Finder.
@ anyDirection
float TrackScore
Definition TrackScore.h:10
DataVector< const Trk::TrackStateOnSurface > TrackStates
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
ParametersT< TrackParametersDim, Charged, StraightLineSurface > AtaStraightLine
@ phi0
Definition ParamDefs.h:65
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67
@ loc2
generic first and second local coordinate
Definition ParamDefs.h:35
@ phi
Definition ParamDefs.h:75
@ loc1
Definition ParamDefs.h:34
@ locZ
local cylindrical
Definition ParamDefs.h:42
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64
ParametersBase< TrackParametersDim, Charged > TrackParameters
void stable_sort(DataModel_detail::iterator< DVL > beg, DataModel_detail::iterator< DVL > end)
Specialization of stable_sort for DataVector/List.
@ knTRTTrk
Number of TRT-only tracks on output.
@ knTrkSegUsed
Number of excluded segments by other TRT segments.
@ knTrkScoreZero
Number of tracks rejected by score zero.
std::vector< std::pair< Trk::TrackScore, Trk::Track * > > m_trackScores