ATLAS Offline Software
Loading...
Searching...
No Matches
SimpleTRT_SeededSpacePointFinder_ATL.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
7// Implementation file for class TRT_SeededSpacePointFinder_ATL
9// (c) ATLAS Detector software
11// Version 1.0 04/2007 Martin Siebel
13
14#include <iostream>
15#include <iomanip>
16#include <set>
17#include "GaudiKernel/MsgStream.h"
18#include "CLHEP/Vector/ThreeVector.h"
19#include "CLHEP/Units/SystemOfUnits.h"
24#include "TrkSurfaces/Surface.h"
26
28// Constructor
30
32(const std::string& t,const std::string& n,const IInterface* p)
33 : AthAlgTool(t,n,p)
34{
35 declareInterface<ITRT_SeededSpacePointFinder>(this);
36}
37
39// Destructor
41
43= default;
44
46// Initialisation
48
50{
51
52 // msg().setLevel(outputLevel());
53
54 // PRD-to-track association (optional)
55 ATH_CHECK( m_prdToTrackMap.initialize( !m_prdToTrackMap.key().empty()));
56
57 // get region selector
58 StatusCode sc = m_pRegionSelector.retrieve();
59 if( sc.isFailure() )
60 {
61 msg(MSG::FATAL) << "Failed to retrieve RegionSelector Service";
62 return sc;
63 }
64
65 sc = detStore()->retrieve(m_sctId, "SCT_ID");
66 if (sc.isFailure())
67 {
68 msg(MSG::ERROR) << "Could not get SCT_ID helper !" << endmsg;
69 return StatusCode::FAILURE;
70 }
71
72 sc = detStore()->retrieve(m_trtId, "TRT_ID");
73 if (sc.isFailure())
74 {
75 msg(MSG::ERROR) << "Could not get TRT_ID helper !" << endmsg;
76 return StatusCode::FAILURE;
77 }
78
80
81 ATH_CHECK(m_spacepointsSCTname.initialize());
83
84 return sc;
85}
86
88// Finalize
90
92{
93 StatusCode sc = AthAlgTool::finalize();
94 return sc;
95}
96//============================================================================================
97
99// Initialize tool for new region
101
102std::unique_ptr<InDet::ITRT_SeededSpacePointFinder::IEventData> InDet::SimpleTRT_SeededSpacePointFinder_ATL::newRegion
103(const std::vector<IdentifierHash>& /*vPixel*/, const std::vector<IdentifierHash>& /*vSCT*/) const
104{
105 return {};
106}
107
108std::list<std::pair<const Trk::SpacePoint*, const Trk::SpacePoint*> >
110 const Trk::TrackParameters& directionTRT,
112{
131
132
133 msg(MSG::VERBOSE) << "Enter getListOfSpacePointPairs, TrackParameter given is: " << endmsg;
134 msg(MSG::VERBOSE) << directionTRT << endmsg;
135
136 // clear output buffer
137
138 std::list<std::pair<const Trk::SpacePoint*, const Trk::SpacePoint*> > listOfSpacePointPairsBuffer;
139
140 // IdHash in ROI
141 std::set<IdentifierHash> setOfSCT_Hashes;
142
143 if (m_useROI)
144 {
145 // fill IdHashes in ROI
146 getHashesInROI(ctx,directionTRT,setOfSCT_Hashes);
147 msg(MSG::VERBOSE) << "Retrieved " << setOfSCT_Hashes.size() << " potentially interesting detector elements." << endmsg;
148 }
149
150 if ( !m_useROI || !setOfSCT_Hashes.empty())
151 {
152 // map of < SCT layer number, SP* >
153 std::multimap<int,const Trk::SpacePoint*> relevantSpacePoints;
154
155 // fill the map of relevant SP as defined by hashes from ROI
156 int modulTRT = TRT_Module(directionTRT);
157 getSpacePointsInROI(ctx, setOfSCT_Hashes, modulTRT, relevantSpacePoints);
158
159 msg(MSG::VERBOSE) << "Retrieved " << relevantSpacePoints.size() << " potentially interesting SpacePoints" << endmsg;
160
161 // build pairs of the relevant SP according to the look-up table
162 combineSpacePoints(relevantSpacePoints, directionTRT, modulTRT, listOfSpacePointPairsBuffer);
163
164 /* output for debug purposes, deactivated now. Once development is finished, it will be removed.
165 */
166 Amg::Vector3D r0 = directionTRT.position();
167 const Amg::Vector3D& v0(directionTRT.momentum());
168 msg(MSG::VERBOSE) << "------------------------------------------------------------------------------------------" << endmsg;
169 msg(MSG::VERBOSE) << "Final SpacePoint pairs: " << listOfSpacePointPairsBuffer.size() << endmsg;
170 msg(MSG::VERBOSE) << " Position of initial vector: ( " << r0.x() << ", " << r0.y() << ", "<< r0.z() << " ) " << endmsg;
171 msg(MSG::VERBOSE) << " Direction of initial vector: ( " << v0.unit().x() << ", " << v0.unit().y() << ", "<< v0.unit().z() << " ) , phi = "
172 << v0.phi() << " theta = " << v0.theta() << " eta = "<< v0.eta() << endmsg;
173 msg(MSG::VERBOSE) << "------------------------------------------------------------------------------------------" << endmsg;
174 msg(MSG::VERBOSE) << " Direction of space Point vectors: "<< endmsg;
175 msg(MSG::VERBOSE) << ". . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . " << endmsg;
176 for (auto & it : listOfSpacePointPairsBuffer)
177 {
178 Amg::Vector3D s1 = it.first->globalPosition();
179 Amg::Vector3D s2 = it.second->globalPosition();
180 Amg::Vector3D s1s2 = s2-s1; // vector from s1 to s2
181 msg(MSG::VERBOSE) << " Positions: ( " << s1.x() << " , "<< s1.y() << " , "<< s1.z() << " ) " << endmsg;
182 msg(MSG::VERBOSE) << " ( " << s2.x() << " , "<< s2.y() << " , "<< s2.z() << " ) " << endmsg;
183 msg(MSG::VERBOSE) << " direction: ( "
184 << s1s2.unit().x() << ", " << s1s2.unit().y() << ", " << s1s2.unit().z() <<" ) , phi = "
185 << s1s2.phi() << " theta = " << s1s2.theta() << " eta = "<< s1s2.eta() << endmsg;
186 msg(MSG::VERBOSE) << ". . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . " << endmsg;
187 }
188 msg(MSG::VERBOSE) << "------------------------------------------------------------------------------------------" << endmsg;
189 }
190
191 return listOfSpacePointPairsBuffer;
192}
193
194//=====================================================================================================
195
196void InDet::SimpleTRT_SeededSpacePointFinder_ATL::getHashesInROI(const EventContext& ctx, const Trk::TrackParameters& directionTRT, std::set<IdentifierHash>& setOfSCT_Hashes) const
197{
202
203
204 //double eta = directionTRT.eta();
205 //double phi = directionTRT.position().phi();
206 double phi = directionTRT.parameters()[Trk::phi];
207 double eta = directionTRT.position().eta();
208 double deltaPhi = 0.;
209 double deltaEta = 0.;
210
211 // define tolerance for ROI
212 //getSearchRange(directionTRT, deltaPhi, deltaEta);
213 getSearchRange(deltaPhi, deltaEta);
214
215 // retrieve SCT hashes in Region of interest
216 std::vector<IdentifierHash> listOfSCT_Hashes;
217
219 RoiDescriptor roi( eta-deltaEta, eta+deltaEta, phi-deltaPhi, phi+deltaPhi);
220
221 m_pRegionSelector->lookup(ctx)->HashIDList( roi, listOfSCT_Hashes );
222
223 // copy Hashes into Set to be able to search them
224 for (auto listOfSCT_Hashe : listOfSCT_Hashes)
225 setOfSCT_Hashes.insert(listOfSCT_Hashe);
226}
227
228//=====================================================================================================
229
231 std::set<IdentifierHash>& setOfSCT_Hashes,
232 int modulTRT, std::multimap<int,
233 const Trk::SpacePoint*>& relevantSpacePoints) const
234{
250
251
252 const std::set<IdentifierHash>::const_iterator endSCT_Hashes = setOfSCT_Hashes.end();
253
254 SG::ReadHandle<Trk::PRDtoTrackMap> prd_to_track_map;
255 const Trk::PRDtoTrackMap *prd_to_track_map_cptr = nullptr;
256 if (!m_prdToTrackMap.key().empty()) {
258 if (!prd_to_track_map.isValid()) {
259 ATH_MSG_ERROR("Failed to read PRD to track association map: " << m_prdToTrackMap.key());
260 }
261 prd_to_track_map_cptr = prd_to_track_map.cptr();
262 }
263
264 // retrieve SP Container
266 if(spacepointsSCT.isValid())
267 {
268 // loop over SP collections in SP container
269 SpacePointContainer::const_iterator itCont = spacepointsSCT->begin();
270 SpacePointContainer::const_iterator endCont = spacepointsSCT->end ();
271 for(; itCont != endCont; ++itCont)
272 {
273 bool acceptCollection = true;
274
275 if (m_useROI)
276 {
277 // Check if IdHash of the Collection is in the list of ROI hashes
278 IdentifierHash idHash = (*itCont)->identifyHash();
279 acceptCollection = ( idHash.is_valid() && setOfSCT_Hashes.find( idHash ) != endSCT_Hashes );
280 }
281
282 int SCT_LayerNumber = -99;
283 if (acceptCollection)
284 {
285 // check if the module is on a relevant layer
286 Identifier id = (*itCont)->identify();
287 int detRegionIndex = m_sctId->barrel_ec(id) ;
288 if ( detRegionIndex == -2 || detRegionIndex == 0 || detRegionIndex == 2 )
289 {
290 // conversion factors for layer numbering
291 int detLayerOffset = 1;
292 int detLayerSign = 1;
293 if (detRegionIndex == 0) detLayerOffset = 10;
294 if (detRegionIndex == -2) detLayerSign = -1;
295
296 // retrieve layer number from IdHelper
297 int detLayerIndex = m_sctId->layer_disk(id);
298 SCT_LayerNumber = detLayerSign*(detLayerIndex+detLayerOffset);
299 }
300 else
301 {
302 msg(MSG::VERBOSE) << "detRegionIndex " << detRegionIndex << " recieved. " << endmsg;
303 acceptCollection = false;
304 }
305 }
306
307 if (acceptCollection)
308 {
309 if (modulTRT > -14 && modulTRT < 15 )
310 {
311 std::set<int>::const_iterator pos = m_modulLookupTableIndex[modulTRT+SIMPLE_TRT_INDEX_OFFSET].find(SCT_LayerNumber);
312 if ( pos == m_modulLookupTableIndex[modulTRT+SIMPLE_TRT_INDEX_OFFSET].end() ) acceptCollection = false;
313 }
314 else
315 {
316 msg(MSG::WARNING) << "Received TRTmodul number " << modulTRT << endmsg;
317 acceptCollection = false;
318 }
319 }
320
321 if (acceptCollection)
322 {
323 // Loop over SP Collection and add SP if they are not yet used (or if this does not matter).
324 SpacePointCollection::const_iterator itColl = (*itCont)->begin();
325 SpacePointCollection::const_iterator endColl = (*itCont)->end ();
326 for(; itColl != endColl; ++itColl)
327 if ( !prd_to_track_map_cptr
328 || !( prd_to_track_map_cptr->isUsed(*((*itColl)->clusterList().first))
329 || prd_to_track_map_cptr->isUsed(*((*itColl)->clusterList().second)) ) )
330 {
331 relevantSpacePoints.insert(std::make_pair( SCT_LayerNumber ,*itColl));
332 msg(MSG::VERBOSE) << "Added SpacePoint for layer " << SCT_LayerNumber << " at ( "
333 << (*itColl)->globalPosition().x() << " , "
334 << (*itColl)->globalPosition().y() << " , "
335 << (*itColl)->globalPosition().z() << " ) " << endmsg;
336 }
337 }
338 }
339 }
340
341 // retrieve the overlap collection
343 if(spacepointsOverlap.isValid())
344 {
345
346 // Loop over Overlap SP
347 SpacePointOverlapCollection::const_iterator itColl = spacepointsOverlap->begin();
348 SpacePointOverlapCollection::const_iterator endColl = spacepointsOverlap->end ();
349 for (; itColl != endColl; ++itColl)
350 {
351 // check if SP is in ROI
352 std::pair<IdentifierHash, IdentifierHash> idHashPair = (*itColl)->elementIdList();
353 if ( !m_useROI
354 || (idHashPair.first.is_valid() && setOfSCT_Hashes.find( idHashPair.first ) != endSCT_Hashes)
355 || (idHashPair.second.is_valid() && setOfSCT_Hashes.find( idHashPair.second ) != endSCT_Hashes) )
356 {
357
358 // find out if one of the Clusters has already been used, if relevant
359 if(prd_to_track_map_cptr)
360 {
361 bool u1=false;
362 bool u2=false;
363 const Trk::PrepRawData* p1=(*itColl)->clusterList().first;
364 u1=prd_to_track_map->isUsed(*p1);
365 const Trk::PrepRawData* p2=(*itColl)->clusterList().second;
366 u2=prd_to_track_map->isUsed(*p2);
367 if(u1 || u2) continue;
368 }
369
370 // retrieve identifier and fill SP in SP buffer.
371 Identifier id = (*itColl)->associatedSurface().associatedDetectorElementIdentifier();
372 if ( id.is_valid() )
373 {
374 int detRegionIndex = m_sctId->barrel_ec(id) ;
375 if ( detRegionIndex == -2 || detRegionIndex == 0 || detRegionIndex == 2 )
376 {
377 int detLayerOffset = 1;
378 int detLayerSign = 1;
379 if (detRegionIndex == 0) detLayerOffset = 10;
380 if (detRegionIndex == -2) detLayerSign = -1;
381 int detLayerIndex = m_sctId->layer_disk(id);
382
383 relevantSpacePoints.insert(std::make_pair(detLayerSign*(detLayerIndex+detLayerOffset),*itColl));
384 msg(MSG::VERBOSE) << "Added OverlapSpacePoint for layer " << detLayerSign*(detLayerIndex+detLayerOffset) << " at ( "
385 << (*itColl)->globalPosition().x() << " , "
386 << (*itColl)->globalPosition().y() << " , "
387 << (*itColl)->globalPosition().z() << " ) " << endmsg;
388 }
389 else
390 msg(MSG::VERBOSE) << "detRegionIndex " << detRegionIndex << " recieved. " << endmsg;
391 }
392 else
393 msg(MSG::DEBUG) << "Invalid Id from OverlapCollection recieved. " << endmsg;
394 }
395 }
396 }
397}
398
399
400//=====================================================================================================
401
402//void InDet::SimpleTRT_SeededSpacePointFinder_ATL::getSearchRange(const Trk::TrackParameters& directionTRT, double& deltaPhi, double& deltaEta)
404{
409
410 // first guess: the ROI tolerance is constant
411 deltaPhi = 0.05;
412 deltaEta = .2;
413}
414
415//=====================================================================================================
416
417void InDet::SimpleTRT_SeededSpacePointFinder_ATL::combineSpacePoints(const std::multimap<int,const Trk::SpacePoint*>& relevantSpacePoints,
418 const Trk::TrackParameters& directionTRT, int modulNumber,
419 std::list<std::pair<const Trk::SpacePoint*, const Trk::SpacePoint*> > &listOfSpacePointPairsBuffer) const
420{
432
433 if ( modulNumber+SIMPLE_TRT_INDEX_OFFSET < 29 && -1 < modulNumber+SIMPLE_TRT_INDEX_OFFSET)
434 {
435 modulLookupTable::const_iterator itTable = m_modulLookupTable[modulNumber+SIMPLE_TRT_INDEX_OFFSET].begin();
436 modulLookupTable::const_iterator endTable = m_modulLookupTable[modulNumber+SIMPLE_TRT_INDEX_OFFSET].end();
437 for ( ; itTable != endTable ; ++itTable)
438 {
439 msg(MSG::VERBOSE) << " Combining Space Poinst from modules " << itTable->first << " and " << itTable->second << endmsg;
440 // SPs from SCT layer 1
441 std::pair< std::multimap<int,const Trk::SpacePoint*>::const_iterator,
442 std::multimap<int,const Trk::SpacePoint*>::const_iterator >
443 range1 = relevantSpacePoints.equal_range(itTable->first);
444 // if (range1.second != relevantSpacePoints.end() ) ++(range1.second); // transforming the 2nd pointer to an end of range
445
446 //SPs from SCT layer 2
447 std::pair< std::multimap<int,const Trk::SpacePoint*>::const_iterator,
448 std::multimap<int,const Trk::SpacePoint*>::const_iterator >
449 range2 = relevantSpacePoints.equal_range(itTable->second);
450 // if (range2.second != relevantSpacePoints.end() ) ++(range2.second); // transforming the 2nd pointer to an end of range
451
452 // add the SP pairs
453 for ( std::multimap<int,const Trk::SpacePoint*>::const_iterator it1 = range1.first;
454 it1 != range1.second ; ++it1 )
455 for ( std::multimap<int,const Trk::SpacePoint*>::const_iterator it2 = range2.first;
456 it2 != range2.second ; ++it2 )
457 if ( pairIsOk(it1->second,it2->second,directionTRT) )
458 listOfSpacePointPairsBuffer.emplace_back(it2->second,it1->second);
459
460 }
461 }
462 else
463 msg(MSG::WARNING) << "TRT module not in look-up table --> no SP pairs formed" << endmsg;
464
465}
466//=====================================================================================================
467
469{
477
478
479 Amg::Vector3D s1 = sp1->globalPosition();
480 Amg::Vector3D s2 = sp2->globalPosition();
481 //Amg::Vector3D r0 = directionTRT.position();
482 const Amg::Vector3D& v0(directionTRT.momentum());
483 Amg::Vector3D s1s2 = s2-s1; // vector from s1 to s2
484
485 msg(MSG::VERBOSE) << "Checking Space Point Pair at ( "<< s1.x() << ", " << s1.y() << ", " << s1.z()
486 << " ) and ( "<< s2.x() << ", " << s2.y() << ", " << s2.z() << " )" << endmsg;
487
488 // Cut on closest approach to z-axis
489 double t = ( s1.x()-s2.x() )*( s1.x()-s2.x() ) + ( s1.y()-s2.y() )*( s1.y()-s2.y() ) ;
490 if (t)
491 {
492 t = ( s1.x()*(s2.x()-s1.x()) + s1.y()*(s2.y()-s1.y()) ) / t;
493 Amg::Vector3D perigee = s1 + t*(s1-s2);
494 msg(MSG::VERBOSE) << " closest approach to beam-pipe at ( "<< perigee.x() << ", " << perigee.y() << ", " << perigee.z() << " ) transversal: " << perigee.perp()<<endmsg;
495 if (perigee.perp() > m_perigeeCut ) return false;
496 }
497 msg(MSG::VERBOSE) << " Passed cut on r-phi impact parameter" <<endmsg;
498
499
500 // Cut on angle between s1-s2 and parameter direction
501 double diffPhi= std::abs(s1s2.phi()-v0.phi());
502
503 if (diffPhi > CLHEP::pi ) diffPhi = 2.*CLHEP::pi - diffPhi;
504 if (diffPhi > CLHEP::pi/2. ) diffPhi = CLHEP::pi - diffPhi;
505 msg(MSG::VERBOSE) << "Phi directions differ by "<< diffPhi << endmsg;
506 if ( diffPhi > m_directionPhiCut ) return false;
507
508 msg(MSG::VERBOSE) << " Passed cut on direction phi deviation" <<endmsg;
509
510 return true;
511}
512//=====================================================================================================
514{
549
550 m_modulLookupTable[0+SIMPLE_TRT_INDEX_OFFSET].emplace_back(12,13);
551 if (m_maxHoles > 0 )
552 {
553 m_modulLookupTable[0+SIMPLE_TRT_INDEX_OFFSET].emplace_back(11,12);
554 m_modulLookupTable[0+SIMPLE_TRT_INDEX_OFFSET].emplace_back(11,13);
555 m_modulLookupTable[0+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,12);
556 }
557 if (m_maxHoles > 1 )
558 {
559 m_modulLookupTable[0+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,13);
560 m_modulLookupTable[0+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,11);
561 }
562
563 m_modulLookupTable[1+SIMPLE_TRT_INDEX_OFFSET].emplace_back(12,13);
564 m_modulLookupTable[1+SIMPLE_TRT_INDEX_OFFSET].emplace_back(11,12);
565 m_modulLookupTable[1+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,11);
566 m_modulLookupTable[1+SIMPLE_TRT_INDEX_OFFSET].emplace_back(13,1);
567 m_modulLookupTable[1+SIMPLE_TRT_INDEX_OFFSET].emplace_back(12,1);
568 m_modulLookupTable[1+SIMPLE_TRT_INDEX_OFFSET].emplace_back(11,1);
569 m_modulLookupTable[1+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,1);
570 m_modulLookupTable[1+SIMPLE_TRT_INDEX_OFFSET].emplace_back(1,2);
571 if (m_maxHoles > 0 )
572 {
573 m_modulLookupTable[1+SIMPLE_TRT_INDEX_OFFSET].emplace_back(11,13);
574 m_modulLookupTable[1+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,12);
575 m_modulLookupTable[1+SIMPLE_TRT_INDEX_OFFSET].emplace_back(13,2);
576 m_modulLookupTable[1+SIMPLE_TRT_INDEX_OFFSET].emplace_back(12,2);
577 m_modulLookupTable[1+SIMPLE_TRT_INDEX_OFFSET].emplace_back(11,2);
578 m_modulLookupTable[1+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,2);
579 }
580 if (m_maxHoles > 1 )
581 {
582 m_modulLookupTable[1+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,13);
583 }
584
585 // TRT wheels 2 and 3
586 for (int i = 2; i<4 ; ++ i)
587 {
588 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(12,13);
589 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(13,1);
590 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(12,1);
591 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(11,1);
592 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,1);
593 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(1,2);
594 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(2,3);
595 if (m_maxHoles > 0 )
596 {
597 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(11,12);
598 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(13,2);
599 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(12,2);
600 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(11,2);
601 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,2);
602 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(1,3);
603 }
604 if (m_maxHoles > 1 )
605 {
606 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,11);
607 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,12);
608 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(13,3);
609 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(12,3);
610 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(11,3);
611 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,3);
612 }
613 }
614
615 // TRT wheels 4 - 6
616 for (int i = 4 ; i < 7 ; ++i )
617 {
618 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(12,13);
619 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(13,1);
620 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(12,1);
621 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(11,1);
622 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,1);
623 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(1,2);
624 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(2,3);
625 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(3,4);
626 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(4,5);
627 if (m_maxHoles > 0 )
628 {
629 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(11,13);
630 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(11,12);
631 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(13,2);
632 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(12,2);
633 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(11,2);
634 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,2);
635 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(1,3);
636 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(2,4);
637 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(3,5);
638 }
639 if (m_maxHoles > 1 )
640 {
641 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(13,3);
642 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(12,3);
643 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(11,3);
644 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(10,3);
645 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(1,4);
646 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(2,5);
647 }
648 }
649
650
651 // TRT wheels 7 - 9
652 for ( int i = 7 ; i < 10 ; ++i)
653 {
654 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(2,3);
655 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(3,4);
656 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(4,5);
657 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(5,6);
658 if (m_maxHoles > 0 )
659 {
660 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(1,2);
661 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(1,3);
662 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(2,4);
663 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(3,5);
664 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(4,6);
665 }
666 if (m_maxHoles > 1 )
667 {
668 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(13,3);
669 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(12,3);
670 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(1,4);
671 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(2,5);
672 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(3,6);
673 }
674 }
675
676
677 // TRT wheels 10 - 12
678 for (int i = 10 ; i < 13 ; ++i )
679 {
680 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(4,5);
681 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(5,6);
682 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(6,7);
683 if (m_maxHoles > 0 )
684 {
685 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(3,4);
686 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(3,5);
687 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(4,6);
688 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(5,7);
689 }
690 if (m_maxHoles > 1 )
691 {
692 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(2,3);
693 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(2,4);
694 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(2,5);
695 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(3,6);
696 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(4,7);
697 }
698 }
699
700
701 // TRT - wheels 13 and 14
702 for (int i = 13 ; i < 15 ; ++i )
703 {
704 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(5,6);
705 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(6,7);
706 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(7,8);
707 if (m_maxHoles > 0 )
708 {
709 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(4,5);
710 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(4,6);
711 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(5,7);
712 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(6,8);
713 }
714 if (m_maxHoles > 0 )
715 {
716 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(3,4);
717 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(3,5);
718 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(3,6);
719 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(4,7);
720 m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].emplace_back(5,8);
721 }
722 }
723
724
725 // mirror the endcap information
726 for ( int iTable = 1 ; iTable < 15 ; ++iTable )
727 {
728 std::list<std::pair<int,int> >::const_iterator itList = m_modulLookupTable[iTable+SIMPLE_TRT_INDEX_OFFSET].begin();
729 std::list<std::pair<int,int> >::const_iterator endList = m_modulLookupTable[iTable+SIMPLE_TRT_INDEX_OFFSET].end();
730 for ( ; itList != endList ; ++itList )
731 {
732 int i = itList->first;
733 int j = itList->second;
734 if (i<10) i*=-1; // don't change the sign of the barrel modules...
735 if (j<10) j*=-1;
736 m_modulLookupTable[-1*iTable+SIMPLE_TRT_INDEX_OFFSET].emplace_back(i,j);
737 }
738 }
739
740 // fill the short index for the lookup table
741 for ( int iTable = -14 ; iTable < 15 ; ++iTable )
742 {
743 std::list<std::pair<int,int> >::const_iterator itList = m_modulLookupTable[iTable+SIMPLE_TRT_INDEX_OFFSET].begin();
744 std::list<std::pair<int,int> >::const_iterator endList = m_modulLookupTable[iTable+SIMPLE_TRT_INDEX_OFFSET].end();
745 for ( ; itList != endList ; ++itList )
746 {
747 m_modulLookupTableIndex[iTable+SIMPLE_TRT_INDEX_OFFSET].insert(itList->first);
748 m_modulLookupTableIndex[iTable+SIMPLE_TRT_INDEX_OFFSET].insert(itList->second);
749 }
750 }
751
753
754}
755
756//=====================================================================================================
757
759{
760 out << "to be implemented soon..." << std::endl;
761
762 return out;
763}
764//=====================================================================================================
765
766std::ostream& InDet::SimpleTRT_SeededSpacePointFinder_ATL::dump( std::ostream& out ) const
767{
768 out << "to be implemented soon..." << std::endl;
769 return out;
770}
771
772
773//=====================================================================================================
774
776{
778
779
780 msg(MSG::VERBOSE) << "=====================================================================================================================" << endmsg;
781 msg(MSG::VERBOSE) << " Module Lookup table: " << endmsg;
782 msg(MSG::VERBOSE) << "---------------------------------------------------------------------------------------------------------------------" << endmsg;
783 for (int i = -14; i<15 ; ++i)
784 {
785 msg(MSG::VERBOSE) << "Module "<< i << " : " << endmsg;
786 std::list<std::pair<int,int> >::const_iterator itList = m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].begin();
787 std::list<std::pair<int,int> >::const_iterator endList = m_modulLookupTable[i+SIMPLE_TRT_INDEX_OFFSET].end();
788 for ( ; itList != endList ; ++itList )
789 msg(MSG::VERBOSE) << " (" << itList->first << ", " << itList->second << ") " ;
790 msg(MSG::VERBOSE) << " " << endmsg;
791 std::set<int>::const_iterator itIndex = m_modulLookupTableIndex[i+SIMPLE_TRT_INDEX_OFFSET].begin();
792 std::set<int>::const_iterator endIndex = m_modulLookupTableIndex[i+SIMPLE_TRT_INDEX_OFFSET].end();
793 msg(MSG::VERBOSE) << "Needs Layers: { " ;
794 for ( ; itIndex != endIndex ; ++ itIndex )
795 msg(MSG::VERBOSE) << *itIndex << ", " ;
796 msg(MSG::VERBOSE) << " } " << endmsg;
797 msg(MSG::VERBOSE) << "...................................................................................................................." << endmsg;
798
799 }
800 msg(MSG::VERBOSE) << "=====================================================================================================================" << endmsg;
801
802}
803
804//=====================================================================================================
805
807{
811
812 // find out from which SCT part the segment is built
813 Identifier id ;
814
816
817 if (!id.is_valid())
818 {
819 msg(MSG::WARNING) << " Id not valid "<<endmsg;
820 return -999;
821 }
822 int modulNumber;
823 int detRegionIndex = m_trtId->barrel_ec(id) ;
824 if ( detRegionIndex == -1 || detRegionIndex == 1 )
825 modulNumber = 0; // central TRT
826 else if ( detRegionIndex == -2 || detRegionIndex == 2 )
827 {
828 modulNumber = m_trtId->layer_or_wheel(id) + 1 ; // forward/bachward wheels
829 modulNumber *= (detRegionIndex/2) ; // add sign to distinguish forward and backward
830 }
831 else
832 {
833 msg(MSG::WARNING) << "TRT barrel-endcap id not in {-2;-1;1;2}" << endmsg;
834 return -999;
835 }
836
837 return modulNumber;
838
839}
Scalar eta() const
pseudorapidity method
Scalar deltaPhi(const MatrixBase< Derived > &vec) const
Scalar phi() const
phi method
#define endmsg
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
static Double_t sc
This is an Identifier helper class for the SCT subdetector.
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
MsgStream & msg() const
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
This is a "hash" representation of an Identifier.
bool is_valid() const
Check if id is in a valid state.
void combineSpacePoints(const std::multimap< int, const Trk::SpacePoint * > &relevantSpacePoints, const Trk::TrackParameters &directionTRT, int modulTRT, std::list< std::pair< const Trk::SpacePoint *, const Trk::SpacePoint * > > &listOfSpacePointPairsBuffer) const
builds pairs of SP according to the kook-up table
SimpleTRT_SeededSpacePointFinder_ATL(const std::string &, const std::string &, const IInterface *)
ToolHandle< IRegSelTool > m_pRegionSelector
Region Selector.
static void getSearchRange(double &deltaPhi, double &deltaEta)
List with SP pairs as seed for the Si part of the back-track.
IntegerProperty m_maxHoles
controls how many not considered SCT layers are allowed between two SP in order to form a seed pair
void getSpacePointsInROI(const EventContext &ctx, std::set< IdentifierHash > &setOfSCT_Hashes, int modulTRT, std::multimap< int, const Trk::SpacePoint * > &relevantSpacePoints) const
retrieves SP Collections of modules in the ROI and sorts them by SCT layer
bool pairIsOk(const Trk::SpacePoint *sp1, const Trk::SpacePoint *sp2, const Trk::TrackParameters &directionTRT) const
applies rough cuts on the quality of a SP pair
BooleanProperty m_useROI
Controls, if SP have to be checked with the AssociationTool of the forward tracking and to avoid doub...
std::set< int > m_modulLookupTableIndex[2 *SIMPLE_TRT_INDEX_OFFSET+1]
std::list< std::pair< const Trk::SpacePoint *, const Trk::SpacePoint * > > find2Sp(const EventContext &ctx, const Trk::TrackParameters &, InDet::ITRT_SeededSpacePointFinder::IEventData &event_data) const override
main method, calls the private methods and returns a pointer to the list of SpacePointpairs.
SG::ReadHandleKey< SpacePointContainer > m_spacepointsSCTname
virtual std::unique_ptr< InDet::ITRT_SeededSpacePointFinder::IEventData > newRegion(const std::vector< IdentifierHash > &, const std::vector< IdentifierHash > &) const override
SG::ReadHandleKey< SpacePointOverlapCollection > m_spacepointsOverlapname
DoubleProperty m_perigeeCut
rough cuts on the quality of the suggested SP pair
void getHashesInROI(const EventContext &ctx, const Trk::TrackParameters &directionTRT, std::set< IdentifierHash > &setOfSCT_Hashes) const
obtains the hashes of modules in the ROI
int TRT_Module(const Trk::TrackParameters &directionTRT) const
returns the number of the TRT wheel or barrel where the TP belongs to
modulLookupTable m_modulLookupTable[2 *SIMPLE_TRT_INDEX_OFFSET+1]
Lookup table that contains the SCT Layers to be considered to provide SP for the pairing in dependenc...
Describes the Region of Ineterest geometry It has basically 9 parameters.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type cptr()
Dereference the pointer.
bool isUsed(const PrepRawData &prd) const
does this PRD belong to at least one track?
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.
virtual const Amg::Vector3D & globalPosition() const override final
Interface method to get the global Position.
Identifier associatedDetectorElementIdentifier() const
return Identifier of the associated Detector Element
Eigen::Matrix< double, 3, 1 > Vector3D
@ phi
Definition ParamDefs.h:75
ParametersBase< TrackParametersDim, Charged > TrackParameters