ATLAS Offline Software
Loading...
Searching...
No Matches
TRT_TrackSegmentsMaker_ATLxk.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6// Implementation file for class TRT_TrackSegmentsMaker_ATLxk
8// (c) ATLAS Detector software
10// AlgTool used for TRT_DriftCircleOnTrack object production
12// Version 1.0 21/04/2004 I.Gavrilenko
14
16
26#include <cmath>
28// Constructor
30
32(const std::string& t,const std::string& n,const IInterface* p)
33 : AthAlgTool(t,n,p)
34{
35 declareInterface<ITRT_TrackSegmentsMaker>(this);
36}
37
39// Destructor
42= default;
43
44
45
47// Initialisation
49
51{
52 StatusCode sc = AlgTool::initialize();
53
54 // Initiate magnetic field properties
55 //
57
58 // Initialize ReadHandle
59 //
60 ATH_CHECK(m_trtname.initialize());
61
62 // Get tool for track extension to TRT
63 //
64 if(m_extensionTool.retrieve().isFailure()) {
65 ATH_MSG_FATAL("Failed to retrieve tool "<< m_extensionTool);
66 return StatusCode::FAILURE;
67 }
68 ATH_MSG_DEBUG("Retrieved tool " << m_extensionTool);
69
70 // PRD-to-track association (optional)
71 ATH_CHECK( m_prdToTrackMap.initialize( !m_prdToTrackMap.key().empty()));
72 ATH_CHECK( m_condDataKey.initialize());
73
74 // TRT
75 if (detStore()->retrieve(m_trtid, "TRT_ID").isFailure()) {
76 ATH_MSG_FATAL("Could not get TRT ID helper");
77 return StatusCode::FAILURE;
78 }
79
80 //Define set of private variables
82 m_Psi = 2./float(m_nMom-1) ;
83 m_Psi128 = m_Psi*128. ;
84 m_Ns128 = m_nPhi*128 ;
85 m_A = float(m_nPhi)/(2.*M_PI) ;
86 m_Ai = 1./m_A ;
87 m_histsize = (m_nPhi*m_nMom*26/4) ;
88
89 // Get output print level
90 //
91 m_outputlevel = msg().level()-MSG::DEBUG;
92 return sc;
93}
94
96// Finalize
98
100{
101 StatusCode sc = AlgTool::finalize(); return sc;
102}
103
105// Initialize tool for new event
107
108std::unique_ptr<InDet::ITRT_TrackSegmentsMaker::IEventData>
110{
111
112 const float pi2 = 2.*M_PI;
113
115
116 // Get drift circles collection
117 //
119 if(not trtcontainer.isValid()) {
120 std::stringstream msg;
121 msg << name() << " Missing TRT_DriftCircleContainer " << m_trtname.key();
122 throw std::runtime_error(msg.str());
123 }
124
125 SG::ReadHandle<Trk::PRDtoTrackMap> prd_to_track_map;
126 const Trk::PRDtoTrackMap *prd_to_track_map_cptr = nullptr;
127 if (!m_prdToTrackMap.key().empty()) {
129 if (!prd_to_track_map.isValid()) {
130 ATH_MSG_ERROR("Failed to read PRD to track association map: " << m_prdToTrackMap.key());
131 }
132 prd_to_track_map_cptr = prd_to_track_map.cptr();
133 }
134
135
136 std::unique_ptr<TRT_TrackSegmentsMaker_ATLxk::EventData>
137 event_data = std::make_unique<TRT_TrackSegmentsMaker_ATLxk::EventData>(trtcontainer.cptr(), condData.m_cirsize);
138
139 // Initiate extension tool
140 //
141 event_data->m_extEventData = m_extensionTool->newEvent(ctx);
142
143 InDet::TRT_DriftCircleContainer::const_iterator
144 w = trtcontainer->begin(),we = trtcontainer->end();
145 int n = 0;
146 if(w!=we) {
147
148 eraseHistogramm(*event_data);
149
150 for(; w!=we; ++w) {
151
152 if(n >= condData.m_cirsize) break;
153
154 Identifier ID = (*w)->identify();
155 int be = m_trtid->barrel_ec (ID);
156 int lw = m_trtid->layer_or_wheel(ID);
157 int sl = m_trtid->straw_layer (ID);
158
159 int b; be < 0 ? b = be+2 : b = be+1;
160 int l = condData.m_flayers[b][lw]+sl;
161 unsigned int sb = condData.m_begin[b][l];
162 unsigned int se = condData.m_end [b][l];
163 unsigned int ad = 1000*b+l;
164
165 // Loop through all clusters from given detector element
166 //
167 InDet::TRT_DriftCircleCollection::const_iterator
168 c = (*w)->begin(), ce = (*w)->end();
169 for(; c!=ce; ++c) {
170
171 if(prd_to_track_map_cptr && prd_to_track_map_cptr->isUsed(*(*c))) continue;
172 if(m_removeNoise && (*c)->isNoise() ) continue;
173
174 if(n >= condData.m_cirsize) break;
175
176 int ns = m_trtid->straw((*c)->identify());
177 const Amg::Vector3D& sc = (*c)->detectorElement()->strawCenter(ns);
178 float Fs = std::atan2(sc.y(),sc.x()); if(Fs<0.) Fs+=pi2;
179 event_data->m_circles[n].set((*c),Fs,ad);
180
181 // Loop through all dz/dr for given cluster
182 //
183 for(unsigned int s=sb; s<=se; ++s) fillHistogramm (Fs,s,*event_data);
184 ++n;
185 }
186 }
187 }
188 event_data->m_clusters = n;
189 return std::unique_ptr<InDet::ITRT_TrackSegmentsMaker::IEventData>(event_data.release());
190}
191
193// Initialize tool for new region
195
196std::unique_ptr<InDet::ITRT_TrackSegmentsMaker::IEventData>
198(const EventContext& ctx, const std::vector<IdentifierHash>& vTRT) const
199{
200 const float pi2 = 2.*M_PI;
201
203
204 // Get drift cilrcles collection
205 //
207 if(not trtcontainer.isValid()) {
208 ATH_MSG_DEBUG("Could not get TRT_DriftCircleContainer");
209 }
210
211 SG::ReadHandle<Trk::PRDtoTrackMap> prd_to_track_map;
212 const Trk::PRDtoTrackMap *prd_to_track_map_cptr = nullptr;
213 if (!m_prdToTrackMap.key().empty()) {
215 if (!prd_to_track_map.isValid()) {
216 ATH_MSG_ERROR("Failed to read PRD to track association map: " << m_prdToTrackMap.key());
217 }
218 prd_to_track_map_cptr = prd_to_track_map.cptr();
219 }
220
221 std::unique_ptr<TRT_TrackSegmentsMaker_ATLxk::EventData>
222 event_data = std::make_unique<TRT_TrackSegmentsMaker_ATLxk::EventData>(trtcontainer.cptr(), condData.m_cirsize);
223
224 if(trtcontainer.isValid()) {
225 // Initiate extension tool
226 //
227 event_data->m_extEventData = m_extensionTool->newEvent(ctx);
228
229 eraseHistogramm(*event_data);
230
231 std::vector<IdentifierHash>::const_iterator d=vTRT.begin(),de=vTRT.end();
232 int n = 0;
233 for(; d!=de; ++d) {
234
235 const auto *w = trtcontainer->indexFindPtr((*d));
236
237 if(w!=nullptr) {
238
239 Identifier ID = w->identify();
240 int be = m_trtid->barrel_ec (ID);
241 int lw = m_trtid->layer_or_wheel(ID);
242 int sl = m_trtid->straw_layer (ID);
243
244 int b; be < 0 ? b = be+2 : b = be+1;
245 int l = condData.m_flayers[b][lw]+sl;
246 unsigned int sb = condData.m_begin[b][l];
247 unsigned int se = condData.m_end [b][l];
248 unsigned int ad = 1000*b+l;
249
250 InDet::TRT_DriftCircleCollection::const_iterator
251 c = w->begin(), ce = w->end();
252
253 for(; c!=ce; ++c) {
254
255 if(prd_to_track_map_cptr && prd_to_track_map_cptr->isUsed(*(*c))) continue;
256 if(m_removeNoise && (*c)->isNoise() ) continue;
257
258 if(n >= condData.m_cirsize) break;
259
260 int ns = m_trtid->straw((*c)->identify());
261 const Amg::Vector3D& sc = (*c)->detectorElement()->strawCenter(ns);
262 float Fs = std::atan2(sc.y(),sc.x()); if(Fs<0.) Fs+=pi2;
263 event_data->m_circles[n].set((*c),Fs,ad);
264
265 // Loop through all dz/dr for given cluster
266 //
267 for(unsigned int s=sb; s<=se; ++s) fillHistogramm (Fs,s,*event_data);
268 ++n;
269 }
270 }
271 }
272 event_data->m_clusters = n;
273 }
274 return std::unique_ptr<InDet::ITRT_TrackSegmentsMaker::IEventData>( event_data.release() );
275}
276
278// Inform tool about end of event or region investigation
280
282{
283 if (msgLvl(MSG::DEBUG)) {
286 dumpEvent(msg(MSG::DEBUG),event_data);
287 dumpConditions(msg(MSG::DEBUG));
288 msg(MSG::DEBUG) << endmsg;
289 }
290}
291
293// Methods for seeds production without vertex constraint
295void InDet::TRT_TrackSegmentsMaker_ATLxk::find(const EventContext &ctx,
298{
301
302 event_data.m_sizebin_iterator = event_data.m_sizebin.rbegin();
303
304 if(event_data.m_clusters<m_clustersCut) return;
305
307
308 event_data.m_clusterSegment.clear();
309 event_data.m_qualitySegment.clear();
310
311 unsigned int mc = event_data.m_clusters;
312
313 for(unsigned int n=0; n!=mc; ++n) {
314
315 unsigned int b = event_data.m_circles[n].buffer();
316 unsigned int l = event_data.m_circles[n].layer ();
317 unsigned int sb = condData.m_begin[b][l];
318 unsigned int se = condData.m_end [b][l];
319
320 // Loop through all dz/dr for given cluster
321 //
322 unsigned char max = 0;
323 unsigned int maxbin = 0;
324
325 for(unsigned int s=sb; s<=se; ++s) {
326 analyseHistogramm(max,maxbin,event_data.m_circles[n].phi(),s,event_data);
327 }
328
329 if(int(max) > m_clustersCut) {
330 event_data.m_bincluster.insert(std::make_pair(localMaximum(maxbin,event_data),n));
331 }
332 }
333
334 std::multimap<unsigned int,unsigned int>::iterator
335 bc,bce =event_data.m_bincluster.end();
336
337 unsigned int nbins = 0 ;
338 unsigned int fbin = 99999999;
339 for(bc = event_data.m_bincluster.begin(); bc!=bce; ++bc) {
340
341 if((*bc).first==fbin) ++nbins;
342 else {
343 if(fbin!=99999999 && nbins>=5) event_data.m_sizebin.insert(std::make_pair(nbins,fbin));
344 fbin=(*bc).first; nbins = 1;
345 }
346 }
347 if(fbin!=99999999 && nbins>=5) event_data.m_sizebin.insert(std::make_pair(nbins,fbin));
348 event_data.m_sizebin_iterator = event_data.m_sizebin.rbegin();
349
350 // @TODO add to event data
351 SG::ReadHandle<Trk::PRDtoTrackMap> prd_to_track_map;
352 const Trk::PRDtoTrackMap *prd_to_track_map_cptr = nullptr;
353 if (!m_prdToTrackMap.key().empty()) {
355 if (!prd_to_track_map.isValid()) {
356 ATH_MSG_ERROR("Failed to read PRD to track association map: " << m_prdToTrackMap.key());
357 }
358 prd_to_track_map_cptr = prd_to_track_map.cptr();
359 }
360
361 // Local reconstruction and track segments production
362 //
363 while(event_data.m_sizebin_iterator!=event_data.m_sizebin.rend()) {
364
365 unsigned int bin =(*event_data.m_sizebin_iterator++).second;
366 findLocaly(ctx, bin,prd_to_track_map_cptr, event_data, used);
367 }
368
369 // Final segments preparation
370 //
371 segmentsPreparation(event_data);
372
373 event_data.m_segiterator = event_data.m_segments.begin();
374}
375
377// Pseudo iterator
379
381{
384
385 if(event_data.m_segiterator!=event_data.m_segments.end()) return (*event_data.m_segiterator++);
386 return nullptr;
387}
388
390// Dumps relevant information into the MsgStream
392
393MsgStream& InDet::TRT_TrackSegmentsMaker_ATLxk::dump( MsgStream& out ) const
394{
395 out<<std::endl;
396 return dumpConditions(out);
397}
398
400// Dumps conditions information into the MsgStream
402
404{
405
407
408 std::string fieldmode[9] ={"NoField" ,"ConstantField","SolenoidalField",
409 "ToroidalField" ,"Grid3DField" ,"RealisticField" ,
410 "UndefinedField","AthenaField" , "?????" };
411
412 int mode = m_fieldprop.magneticFieldMode();
413 if(mode<0 || mode>8 ) mode = 8;
414
415 int n = 62-fieldmode[mode].size();
416 std::string s3; for(int i=0; i<n; ++i) s3.append(" "); s3.append("|");
417
418 n = 62-m_trtname.key().size();
419 std::string s4; for(int i=0; i<n; ++i) s4.append(" "); s4.append("|");
420
421 n = 62-m_extensionTool.type().size();
422 std::string s7; for(int i=0; i<n; ++i) s7.append(" "); s7.append("|");
423
424
425 out<<"|----------------------------------------------------------------------"
426 <<"-------------------|"
427 <<std::endl;
428 out<<"| Tool tracks extension | "<<m_extensionTool.type()<<s7<<std::endl;
429 out<<"| Magnetic field mode | "<<fieldmode[mode] <<s3<<std::endl;
430 out<<"| TRT container | "<<m_trtname.key().size() <<s4<<std::endl;
431 out<<"| Min. number straws | "
432 <<std::setw(12)<<m_clustersCut
433 <<" |"<<std::endl;
434 out<<"| Number neg. bar. layers | "
435 <<std::setw(12)<<condData.m_nlayers[1]
436 <<" |"<<std::endl;
437 out<<"| Number pos. bar. layers | "
438 <<std::setw(12)<<condData.m_nlayers[2]
439 <<" |"<<std::endl;
440 out<<"| Number neg. end. layers | "
441 <<std::setw(12)<<condData.m_nlayers[0]
442 <<" |"<<std::endl;
443 out<<"| Number pos. end. layers | "
444 <<std::setw(12)<<condData.m_nlayers[3]
445 <<" |"<<std::endl;
446 out<<"| Number neg. bar. straws | "
447 <<std::setw(12)<<condData.m_nstraws[1]
448 <<" |"<<std::endl;
449 out<<"| Number pos. bar. straws | "
450 <<std::setw(12)<<condData.m_nstraws[2]
451 <<" |"<<std::endl;
452 out<<"| Number neg. end. straws | "
453 <<std::setw(12)<<condData.m_nstraws[0]
454 <<" |"<<std::endl;
455 out<<"| Number pos. end. straws | "
456 <<std::setw(12)<<condData.m_nstraws[3]
457 <<" |"<<std::endl;
458 out<<"| Number azimut. channels | "
459 <<std::setw(12)<<m_nPhi
460 <<" |"<<std::endl;
461 out<<"| Number moment. channels | "
462 <<std::setw(12)<<m_nMom
463 <<" |"<<std::endl;
464 out<<"| Number histog. channels | "
465 <<std::setw(12)<<m_histsize
466 <<" |"<<std::endl;
467 out<<"| Number cluster links | "
468 <<std::setw(12)<<condData.m_cirsize
469 <<" |"<<std::endl;
470 out<<"| Use PRD-to-track assoc.?| "
471 <<std::setw(12)<< (!m_prdToTrackMap.key().empty() ? "yes" : "no ")
472 <<" |"<<std::endl;
473 out<<"| Remove noise ? | "
474 <<std::setw(12)<<m_removeNoise
475 <<" |"<<std::endl;
476 out<<"|----------------------------------------------------------------------"
477 <<"-------------------|"
478 <<std::endl;
479 return out;
480}
481
483// Dumps event information into the MsgStream
485
488{
491 out<<"|----------------------------------------------------------------------"
492 <<"-------------------|"
493 <<std::endl;
494 out<<"| Number drift circles | "
495 <<std::setw(12)<<event_data.m_clusters
496 <<" |"<<std::endl;
497 out<<"| Number local calls | "
498 <<std::setw(12)<<event_data.m_nlocal
499 <<" |"<<std::endl;
500 out<<"| Number found segments | "
501 <<std::setw(12)<<event_data.m_nsegments
502 <<" |"<<std::endl;
503 out<<"| Number save segments | "
504 <<std::setw(12)<<event_data.m_segments.size()
505 <<" |"<<std::endl;
506 out<<"|----------------------------------------------------------------------"
507 <<"-------------------|"
508 <<std::endl;
509 return out;
510}
511
513// Dumps relevant information into the ostream
515
516std::ostream& InDet::TRT_TrackSegmentsMaker_ATLxk::dump( std::ostream& out ) const
517{
518 return out;
519}
520
522// Erase histogramm
524
526{
527 for(int i=0; i!=m_histsize; i+=10) {
528
529 event_data.m_U.H4[i ]=0;
530
531 event_data.m_U.H4[i+1]=0;
532 event_data.m_U.H4[i+2]=0;
533 event_data.m_U.H4[i+3]=0;
534 event_data.m_U.H4[i+4]=0;
535 event_data.m_U.H4[i+5]=0;
536 event_data.m_U.H4[i+6]=0;
537 event_data.m_U.H4[i+7]=0;
538 event_data.m_U.H4[i+8]=0;
539 event_data.m_U.H4[i+9]=0;
540 }
541}
542
544// Fill histogramm
546
548(float Fs,int s, TRT_TrackSegmentsMaker_ATLxk::EventData &event_data) const
549{
551
552 int s0 = condData.m_ndzdr[s]*m_Ts;
553 int f = int((Fs*m_A-condData.m_slope[s])*128.); if(f<0) f+=m_Ns128;
554 int sf = condData.m_islope[s];
555
556 // Loop through all momentum slopes
557 //
558 for(int i=s0; i!=s0+m_Ts; i+=m_nPhi) {
559 int k =(f>>7); f+=sf; k<m_nPhi ? ++event_data.m_U.H[k+i] : ++event_data.m_U.H[k+i-m_nPhi];
560 }
561}
562
563
565// Analyse histogramm
567
569(unsigned char& max,unsigned int& maxbin,float Fs,int s, TRT_TrackSegmentsMaker_ATLxk::EventData &event_data) const
570{
572
573 int s0 = condData.m_ndzdr[s]*m_Ts;
574 int f = int((Fs*m_A-condData.m_slope[s])*128.); if(f<0) f+=m_Ns128;
575 int sf = condData.m_islope[s];
576
577 // Loop through all momentum slopes
578 //
579 for(int i=s0; i!=s0+m_Ts; i+=m_nPhi) {
580 int k =(f>>7); f+=sf; if(k>=m_nPhi) k-=m_nPhi;
581 if(event_data.m_U.H[k+i] > max) max = event_data.m_U.H[maxbin = k+i];
582 }
583}
584
586// TRT seeds production
588
590 unsigned int bin,
591 const Trk::PRDtoTrackMap *prd_to_track_map,
594{
596
597 const double pi=M_PI, pi2 = 2.*M_PI;
598
599 std::multimap<const InDet::TRT_DriftCircle*,Trk::TrackSegment*>::const_iterator
600 cse = event_data.m_clusterSegment.end();
601
602 std::multimap<unsigned int,unsigned int>::iterator
603 bc ,
604 bcb = event_data.m_bincluster.find(bin),
605 bce = event_data.m_bincluster.end() ;
606 int nfree = 0;
607 for(bc=bcb; bc!=bce; ++bc) {
608
609 if((*bc).first!=bin) break;
610 if(event_data.m_clusterSegment.find(event_data.m_circles[(*bc).second].circle())==cse) ++nfree;
611 }
612 if(nfree<5) return;
613 unsigned int ndzdr = bin/m_Ts;
614 unsigned int b = bin-ndzdr*m_Ts;
615 unsigned int m = b/m_nPhi;
616 float c0 = (1.-float(m)*m_Psi)*m_Ai;
617
618 bce=bc;
619 float fm = 0. ;
620 float Fo = 0. ;
621 bool first = false;
622 for(bc=bcb; bc!=bce; ++bc) {
623
624 unsigned int n = (*bc).second;
625
626 unsigned int b = event_data.m_circles[n].buffer();
627 unsigned int l = event_data.m_circles[n].layer ();
628 unsigned int s = condData.m_begin[b][l];
629 unsigned int se = condData.m_end [b][l];
630
631 for(; s<= se; ++s) {if(condData.m_ndzdr[s]==ndzdr) break;}
632 if(s>se) continue;
633 float F = event_data.m_circles[n].phi()-condData.m_slope[s]*c0;
634
635 if(!first) {
636 Fo = F; first = true;
637 }
638 else {
639 float df = F-Fo;
640 if (df > pi) df-=pi2;
641 else if(df <-pi) df+=pi2;
642 fm+=df;
643 }
644 }
645 fm = Fo+fm/float(nfree);
646
647 if (fm > pi) fm = fmod(fm+pi,pi2)-pi;
648 else if(fm <-pi) fm = fmod(fm-pi,pi2)+pi;
649
650 double pT = m_pTmin/(double(m)*m_Psi-1.);
651
652 double pin = 1./(pT*std::sqrt((1.+condData.m_dzdr[ndzdr]*condData.m_dzdr[ndzdr])));
653
654 Amg::Vector3D PSV(0.,0.,0.); Trk::PerigeeSurface PS(PSV);
655 auto Tp = PS.createUniqueTrackParameters(
656 0., 0., fm, std::atan2(1., condData.m_dzdr[ndzdr]), pin, std::nullopt);
657 ++event_data.m_nlocal;
658
659 Trk::TrackSegment* seg = m_extensionTool->findSegment(ctx, Tp.get(), *(event_data.m_extEventData),used);
660 if(!seg) return;
661
662 // Momentum cut
663 //
664 double T = seg->localParameters().get(Trk::theta );
665 double iP = seg->localParameters().get(Trk::qOverP);
666
667 // ME: let's not use a soft cut here
668 if(std::sin(T) < 0.9*m_pTmin*std::abs(iP)) {delete seg; return;}
669
670
671 ++event_data.m_nsegments;
672
673 bool isbarrel=false;
674 const Trk::MeasurementBase *lastmeas=seg->containedMeasurements().back();
675 if (std::abs(lastmeas->globalPosition().z())<750.) isbarrel=true;
676 // Number drift circles test
677 //
678 unsigned int size = isbarrel ? seg->numberOfMeasurementBases()-1 : seg->numberOfMeasurementBases()-2;
679
680 if(int(size) < m_clustersCut) {delete seg; return;}
681
682 std::vector<const Trk::MeasurementBase*>::const_iterator
683 s = seg->containedMeasurements().begin(), se = seg->containedMeasurements().end();
684
685 if(prd_to_track_map) {
686
687 // Test number unused drift circles
688 //
689 int nu = 0, ntot = 0;
690 for(++s; s!=se; ++s) {
691
692 const Trk::RIO_OnTrack* rio = dynamic_cast<const Trk::RIO_OnTrack*>(*s);
693 if(rio) {
694 ++ntot;
695 if (prd_to_track_map->isUsed(*rio->prepRawData())) ++nu;
696 }
697 }
698 // ME: use fraction cut
699 if( nu > int(ntot*m_sharedfrac)) {delete seg; return;}
700 }
701
702 // Save information about qiality and segment
703 //
704 int quality = size;
705 double Xi2 = seg->fitQuality()->chiSquared()/double(seg->fitQuality()->numberDoF());
706
707 if(Xi2 > .5) {
708 if (Xi2 > 6) quality-=9;
709 else if(Xi2 > 5.) quality-=7;
710 else if(Xi2 > 3.) quality-=5;
711 else if(Xi2 > 2.) quality-=3;
712 else if(Xi2 > 1.) quality-=1;
713 }
714
715 event_data.m_qualitySegment.insert(std::make_pair(quality,seg));
716
717 s = seg->containedMeasurements().begin();
718
719 // Save information about TRT clusters and segments
720 //
721 for(++s; s!=se; ++s) {
722 const Trk::RIO_OnTrack* rio = dynamic_cast<const Trk::RIO_OnTrack*>(*s);
723 if (!rio) continue;
724 const InDet::TRT_DriftCircleOnTrack* trt = static_cast<const InDet::TRT_DriftCircleOnTrack*>(*s);
725 const InDet::TRT_DriftCircle* dc = trt->prepRawData();
726 if(dc) event_data.m_clusterSegment.insert(std::make_pair(dc,seg));
727 }
728}
729
731// Remove fake TRT segments
733
735{
736 std::multimap<int,Trk::TrackSegment*>::reverse_iterator
737 qs = event_data.m_qualitySegment.rbegin();
738
739 for(; qs!=event_data.m_qualitySegment.rend(); ++qs) {
740
741 int nfree = 0;
742
743 std::vector<const Trk::MeasurementBase*>::const_iterator
744 s =(*qs).second->containedMeasurements().begin(),
745 se =(*qs).second->containedMeasurements().end ();
746
747 for(++s; s!=se; ++s) {
748 const Trk::RIO_OnTrack* rio = dynamic_cast<const Trk::RIO_OnTrack*>(*s);
749 if (!rio) continue;
750 const InDet::TRT_DriftCircleOnTrack* trt = static_cast<const InDet::TRT_DriftCircleOnTrack*>(*s);
751 const InDet::TRT_DriftCircle* dc = trt->prepRawData();
752 if(dc && event_data.m_clusterSegment.erase(dc)) ++nfree;
753 }
754 if(nfree >= 7) event_data.m_segments.push_back((*qs).second);
755 else delete (*qs).second;
756 }
757}
758
760// Local maximum
762
764(unsigned int bin, TRT_TrackSegmentsMaker_ATLxk::EventData &event_data) const
765{
766 int b = bin-(bin/m_Ts)*m_Ts ;
767 int m = b/m_nPhi ;
768 int f = b-m*m_nPhi ;
769 unsigned int maxb = bin ;
770 unsigned char max = event_data.m_U.H[bin];
771
772 int a1 = bin-1; if(f== 0 ) a1+=m_nPhi; if(event_data.m_U.H[a1]>max) {max=event_data.m_U.H[a1]; maxb=a1;}
773 int a2 = bin+1; if(f==m_nPhi-1) a2-=m_nPhi; if(event_data.m_U.H[a2]>max) {max=event_data.m_U.H[a2]; maxb=a2;}
774
775 if ( m < m_nMom-1) {
776
777 int a = bin+m_nPhi; if(event_data.m_U.H[a]>max) {max=event_data.m_U.H[a]; maxb=a;}
778 a = a1 +m_nPhi; if(event_data.m_U.H[a]>max) {max=event_data.m_U.H[a]; maxb=a;}
779 a = a2 +m_nPhi; if(event_data.m_U.H[a]>max) {max=event_data.m_U.H[a]; maxb=a;}
780
781 }
782 if ( m > 0 ) {
783
784 int a = bin-m_nPhi; if(event_data.m_U.H[a]>max) {max=event_data.m_U.H[a]; maxb=a;}
785 a = a1 -m_nPhi; if(event_data.m_U.H[a]>max) {max=event_data.m_U.H[a]; maxb=a;}
786 a = a2 -m_nPhi; if(event_data.m_U.H[a]>max) {max=event_data.m_U.H[a]; maxb=a;}
787 }
788 return maxb;
789}
790
#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_DEBUG(x)
std::vector< Identifier > ID
static Double_t s0
static Double_t a
static Double_t Tp(Double_t *t, Double_t *par)
static Double_t sc
#define F(x, y, z)
Definition MD5.cxx:112
Handle class for reading from StoreGate.
#define pi
static const Attributes_t empty
#define max(a, b)
Definition cfImp.cxx:41
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
const TRT_DriftCircle * circle() const
Represents 'corrected' measurements from the TRT (for example, corrected for wire sag).
virtual const TRT_DriftCircle * prepRawData() const override final
returns the PrepRawData - is a TRT_DriftCircle in this scope
std::unique_ptr< InDet::ITRT_TrackExtensionTool::IEventData > m_extEventData
std::list< Trk::TrackSegment * >::iterator m_segiterator
std::multimap< unsigned int, unsigned int > m_bincluster
std::multimap< unsigned int, unsigned int > m_sizebin
std::multimap< int,Trk::TrackSegment * > m_qualitySegment
std::multimap< const InDet::TRT_DriftCircle *, Trk::TrackSegment * > m_clusterSegment
union InDet::TRT_TrackSegmentsMaker_ATLxk::EventData::@317121305332267356337266061057375304276022267376 m_U
std::multimap< unsignedint, unsignedint >::reverse_iterator m_sizebin_iterator
MsgStream & dumpConditions(MsgStream &out) const
void findLocaly(const EventContext &ctx, unsigned int, const Trk::PRDtoTrackMap *prd_to_track_map, TRT_TrackSegmentsMaker_ATLxk::EventData &event_data, InDet::TRT_DetElementLink_xk::TRT_DetElemUsedMap &used) const
const TRT_TrackSegmentsToolCondData_xk * getConditionsData() const
void eraseHistogramm(TRT_TrackSegmentsMaker_ATLxk::EventData &event_data) const
ToolHandle< ITRT_TrackExtensionTool > m_extensionTool
virtual MsgStream & dump(MsgStream &out) const override
SG::ReadCondHandleKey< InDet::TRT_TrackSegmentsToolCondData_xk > m_condDataKey
SG::ReadHandleKey< InDet::TRT_DriftCircleContainer > m_trtname
void fillHistogramm(float, int, TRT_TrackSegmentsMaker_ATLxk::EventData &event_data) const
virtual Trk::TrackSegment * next(InDet::ITRT_TrackSegmentsMaker::IEventData &event_data) const override
void analyseHistogramm(unsigned char &, unsigned int &, float, int, TRT_TrackSegmentsMaker_ATLxk::EventData &event_data) const
virtual std::unique_ptr< InDet::ITRT_TrackSegmentsMaker::IEventData > newRegion(const EventContext &ctx, const std::vector< IdentifierHash > &) const override
unsigned int localMaximum(unsigned int, TRT_TrackSegmentsMaker_ATLxk::EventData &event_data) const
static MsgStream & dumpEvent(MsgStream &out, InDet::ITRT_TrackSegmentsMaker::IEventData &event_data)
TRT_TrackSegmentsMaker_ATLxk(const std::string &, const std::string &, const IInterface *)
virtual std::unique_ptr< InDet::ITRT_TrackSegmentsMaker::IEventData > newEvent(const EventContext &ctx) const override
void endEvent(InDet::ITRT_TrackSegmentsMaker::IEventData &event_data) const override
SG::ReadHandleKey< Trk::PRDtoTrackMap > m_prdToTrackMap
static void segmentsPreparation(TRT_TrackSegmentsMaker_ATLxk::EventData &event_data)
virtual void find(const EventContext &ctx, InDet::ITRT_TrackSegmentsMaker::IEventData &event_data, InDet::TRT_DetElementLink_xk::TRT_DetElemUsedMap &used) const override
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type cptr()
Dereference the pointer.
static EventData & getPrivateEventData(InDet::ITRT_TrackSegmentsMaker::IEventData &virt_event_data)
int numberDoF() const
returns the number of degrees of freedom of the overall track or vertex fit as integer
Definition FitQuality.h:60
double chiSquared() const
returns the of the overall track fit
Definition FitQuality.h:56
double get(ParamDefs par) const
Retrieve specified parameter (const version).
magnetic field properties to steer the behavior of the extrapolation
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 Amg::Vector3D & globalPosition() const =0
Interface method to get the global Position.
bool isUsed(const PrepRawData &prd) const
does this PRD belong to at least one track?
Class describing the Line to which the Perigee refers to.
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.
const FitQuality * fitQuality() const
return the FitQuality object, returns NULL if no FitQuality is defined
const std::vector< const Trk::MeasurementBase * > & containedMeasurements() const
returns the vector of Trk::MeasurementBase objects
unsigned int numberOfMeasurementBases() const
Return the number of contained Trk::MeasurementBase (s)
Class for a generic track segment that holdes polymorphic Trk::MeasurementBase objects,...
holding In fact this class is here in order to allow STL container for all features This class is sho...
Eigen::Matrix< double, 3, 1 > Vector3D
@ FastField
call the fast field access method of the FieldSvc
@ NoField
Field is set to 0., 0., 0.,.
@ FullField
Field is set to be realistic, but within a given Volume.
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67