ATLAS Offline Software
Loading...
Searching...
No Matches
InDet::SiTrackMaker_xk Class Referencefinal

InDet::SiTrackMaker_xk is algorithm which produce Trk::Track started from 3 space points information of SCT and Pixels in the road of InDetDD::SiDetectorElement* sorted in propagation order. More...

#include <SiTrackMaker_xk.h>

Inheritance diagram for InDet::SiTrackMaker_xk:

Public Member Functions

Standard tool methods
 SiTrackMaker_xk (const std::string &, const std::string &, const IInterface *)
virtual ~SiTrackMaker_xk ()=default
virtual StatusCode initialize () override
virtual StatusCode finalize () override
Main methods for local track finding
virtual std::list< Trk::Track * > getTracks (const EventContext &ctx, SiTrackMakerEventData_xk &data, const std::vector< const Trk::SpacePoint * > &Sp) const override
virtual std::list< Trk::Track * > getTracks (const EventContext &ctx, SiTrackMakerEventData_xk &data, const Trk::TrackParameters &Tp, const std::vector< Amg::Vector3D > &Gp) const override
virtual void newEvent (const EventContext &ctx, SiTrackMakerEventData_xk &data, bool PIX, bool SCT) const override
virtual void newTrigEvent (const EventContext &ctx, SiTrackMakerEventData_xk &data, bool PIX, bool SCT) const override
virtual void endEvent (SiTrackMakerEventData_xk &data) const override
Print internal tool parameters and status
MsgStream & dump (SiTrackMakerEventData_xk &data, MsgStream &out) const override

Private Member Functions

Disallow default constructor, copy constructor and assignment operator
 SiTrackMaker_xk ()=delete
 SiTrackMaker_xk (const SiTrackMaker_xk &)=delete
SiTrackMaker_xkoperator= (const SiTrackMaker_xk &)=delete

Private Attributes

Tool handles
ToolHandle< InDet::ISiDetElementsRoadMakerm_roadmaker {this, "RoadTool", "InDet::SiDetElementsRoadMaker_xk"}
ToolHandle< InDet::ISiCombinatorialTrackFinderm_tracksfinder {this, "CombinatorialTrackFinder", "InDet::SiCombinatorialTrackFinder_xk"}
ToolHandle< ITrigInDetTrackFollowingToolm_trigInDetTrackFollowingTool {this, "TrigTrackFollowingTool", "TrigInDetTrackFollowingTool"}
ToolHandle< ITrigInDetRoadPredictorToolm_trigInDetRoadPredictorTool {this, "TrigInDetRoadPredictorTool", "TrigInDetRoadPredictorTool_FTF"}
ToolHandle< InDet::ISeedToTrackConversionToolm_seedtrack {this, "SeedToTrackConversion", "InDet::SeedToTrackConversionTool"}
Data handles
SG::ReadCondHandleKey< InDet::BeamSpotDatam_beamSpotKey {this, "BeamSpotKey", "BeamSpotData", "SG key for beam spot"}
SG::ReadCondHandleKey< AtlasFieldCacheCondObjm_fieldCondObjInputKey {this, "AtlasFieldCacheCondObj", "fieldCondObj", "Name of the Magnetic Field conditions object key"}
SG::ReadHandleKey< ROIPhiRZContainerm_caloCluster {this, "EMROIPhiRZContainer", ""}
SG::ReadHandleKey< ROIPhiRZContainerm_caloHad {this, "HadROIPhiRZContainer", ""}
Properties
BooleanProperty m_seedsfilter {this, "UseSeedFilter", true, "Use seed filter"}
StringProperty m_fieldmode {this, "MagneticFieldMode", "MapSolenoid", "Mode of magnetic field"}
StringProperty m_patternName {this, "TrackPatternRecoInfo", "SiSPSeededFinder", "Name of the pattern recognition"}
BooleanProperty m_usePix {this, "usePixel", true, "flags to set whether to use pixel/sct cluster, irrespective of what is in event"}
BooleanProperty m_useSct {this, "useSCT", true}
BooleanProperty m_useassoTool {this, "UseAssociationTool", false, "Use prd-track association tool"}
BooleanProperty m_cosmicTrack {this, "CosmicTrack", false, "Is it cosmic track"}
BooleanProperty m_multitracks {this, "doMultiTracksProd", false}
BooleanProperty m_useBremModel {this, "useBremModel", false}
BooleanProperty m_useCaloSeeds {this, "doCaloSeededBrem", false}
BooleanProperty m_useSSSfilter {this, "useSSSseedsFilter", true}
BooleanProperty m_useHClusSeed {this, "doHadCaloSeedSSS", false, "Hadronic Calorimeter Seeds"}
BooleanProperty m_ITKGeometry {this, "ITKGeometry", false, "ITK geometry"}
BooleanProperty m_seedsegmentsWrite {this, "SeedSegmentsWrite", false, "Call seed to track conversion"}
BooleanProperty m_useTrigTrackFollowingTool {this, "useTrigTrackFollowingTool", false, "Option to use TrigInDetTrackFollowingTool instead of SiCombinatorialTrackFinder_xk"}
BooleanProperty m_useTrigInDetRoadPredictorTool {this, "useTrigInDetRoadPredictorTool", false, "Option to use TrigInDetRoadPredictorTool instead of ISiDetElementsRoadMaker"}
BooleanProperty m_LRTmode {this, "LRTMode", false}
DoubleProperty m_xi2max {this, "Xi2max", 15., "max Xi2 for updators"}
DoubleProperty m_xi2maxNoAdd {this, "Xi2maxNoAdd", 35., "max Xi2 for clusters"}
DoubleProperty m_xi2maxlink {this, "Xi2maxlink", 200., "max Xi2 for clusters"}
DoubleProperty m_pTmin {this, "pTmin", 500., "min pT"}
DoubleProperty m_pTminBrem {this, "pTminBrem", 1000., "min pT for Brem mode"}
DoubleProperty m_distmax {this, "MaxDistanceForSCTsp", 5.}
DoubleProperty m_xi2multitracks {this, "Xi2maxMultiTracks", 3., "max Xi2 for multi tracks"}
IntegerProperty m_nholesmax {this, "nHolesMax", 2, "Max number holes"}
IntegerProperty m_dholesmax {this, "nHolesGapMax", 2, "Max holes gap"}
IntegerProperty m_nclusmin {this, "nClustersMin", 6, "Min number clusters"}
IntegerProperty m_nwclusmin {this, "nWeightedClustersMin", 6, "Min umber weighted clusters(pix=2 sct=1)"}
IntegerProperty m_trackletPoints {this, "trackletPoints", 1, "Select which tracklet points to use"}
DoubleProperty m_phiWidth {this, "phiWidth", 0.3}
DoubleProperty m_etaWidth {this, "etaWidth", 0.3}
DoubleArrayProperty m_etabins {this, "etaBins", {}, "eta bins"}
DoubleArrayProperty m_ptbins {this, "pTBins", {}, "pT bins"}

Data members, which are updated only in initialize method

enum  statAllTypes {
  kTotalInputSeeds , kTotalUsedSeeds , kTotalNoTrackPar , kTotalBremSeeds ,
  kTwoClusters , kWrongInit , kWrongRoad , kNoTrack ,
  kNotNewTrack , kBremAttempt , kOutputTracks , kExtraTracks ,
  kBremTracks , kDESize , kSeedsWithTracks
}
enum  kNStatEtaTypes { kUsedSeedsEta , kSeedsWithTracksEta }
Trk::TrackInfo m_trackinfo
bool m_heavyion {false}
Trk::MagneticFieldMode m_fieldModeEnum {Trk::FullField}
std::mutex m_counterMutex
std::array< std::atomic< int >, SiCombinatorialTrackFinderData_xk::kNSeedTypes > m_totalInputSeeds ATLAS_THREAD_SAFE {}
std::array< std::atomic< int >, SiCombinatorialTrackFinderData_xk::kNSeedTypes > m_totalUsedSeeds ATLAS_THREAD_SAFE {}
std::array< std::atomic< int >, SiCombinatorialTrackFinderData_xk::kNSeedTypes > m_totalNoTrackPar ATLAS_THREAD_SAFE {}
std::array< std::atomic< int >, SiCombinatorialTrackFinderData_xk::kNSeedTypes > m_totalBremSeeds ATLAS_THREAD_SAFE {}
std::array< std::atomic< int >, SiCombinatorialTrackFinderData_xk::kNSeedTypes > m_twoClusters ATLAS_THREAD_SAFE {}
std::array< std::atomic< int >, SiCombinatorialTrackFinderData_xk::kNSeedTypes > m_wrongRoad ATLAS_THREAD_SAFE {}
std::array< std::atomic< int >, SiCombinatorialTrackFinderData_xk::kNSeedTypes > m_wrongInit ATLAS_THREAD_SAFE {}
std::array< std::atomic< int >, SiCombinatorialTrackFinderData_xk::kNSeedTypes > m_noTrack ATLAS_THREAD_SAFE {}
std::array< std::atomic< int >, SiCombinatorialTrackFinderData_xk::kNSeedTypes > m_notNewTrack ATLAS_THREAD_SAFE {}
std::array< std::atomic< int >, SiCombinatorialTrackFinderData_xk::kNSeedTypes > m_bremAttempt ATLAS_THREAD_SAFE {}
std::array< std::atomic< int >, SiCombinatorialTrackFinderData_xk::kNSeedTypes > m_outputTracks ATLAS_THREAD_SAFE {}
std::array< std::atomic< int >, SiCombinatorialTrackFinderData_xk::kNSeedTypes > m_extraTracks ATLAS_THREAD_SAFE {}
std::array< std::atomic< int >, SiCombinatorialTrackFinderData_xk::kNSeedTypes > m_bremTracks ATLAS_THREAD_SAFE {}
std::array< std::atomic< int >, SiCombinatorialTrackFinderData_xk::kNSeedTypes > m_seedsWithTrack ATLAS_THREAD_SAFE {}
std::array< std::atomic< double >, SiCombinatorialTrackFinderData_xk::kNSeedTypes > m_deSize ATLAS_THREAD_SAFE {}
std::vector< std::vector< double > > m_usedSeedsEta ATLAS_THREAD_SAFE
std::vector< std::vector< double > > m_seedsWithTracksEta ATLAS_THREAD_SAFE
std::vector< statAllTypesm_indexToEnum {kTwoClusters,kWrongInit,kWrongRoad,kNoTrack,kNotNewTrack,kBremAttempt}
std::unique_ptr< Trk::TrackParametersgetAtaPlane (MagField::AtlasFieldCache &fieldCache, SiTrackMakerEventData_xk &data, bool sss, const std::vector< const Trk::SpacePoint * > &SP, const EventContext &ctx) const
bool globalPositions (const Trk::SpacePoint &s0, const Trk::SpacePoint &s1, const Trk::SpacePoint &s2, double *p0, double *p1, double *p2) const
bool globalPosition (const Trk::SpacePoint &sp, const double *dir, double *p) const
 This is a refinement of the global position for strip space-points.
InDet::TrackQualityCuts setTrackQualityCuts (bool simpleTrack) const
bool newSeed (SiTrackMakerEventData_xk &data, const std::vector< const Trk::SpacePoint * > &Sp) const
bool isCaloCompatible (SiTrackMakerEventData_xk &data) const
bool isHadCaloCompatible (SiTrackMakerEventData_xk &data) const
double pTmin (double eta) const
MsgStream & dumpStatistics (MsgStream &out) const
MsgStream & dumpconditions (MsgStream &out) const
template<typename T, size_t N, size_t M>
void resetCounter (std::array< std::array< T, M >, N > &a) const
 helper for working with the stat arrays
template<typename T, size_t N>
void resetCounter (std::array< T, N > &a) const
static void globalDirections (const double *p0, const double *p1, const double *p2, double *d0, double *d1, double *d2)
 Here, we derive the local track direction in the space-points as the tangents to the estimated trajectory (assuming a circle in x-y and straight line in r-z)
static void detectorElementsSelection (SiTrackMakerEventData_xk &data, std::vector< const InDetDD::SiDetectorElement * > &DE)
static int kindSeed (const std::vector< const Trk::SpacePoint * > &Sp)
static int rapidity (const std::vector< const Trk::SpacePoint * > &Sp)
static bool isNewTrack (SiTrackMakerEventData_xk &data, Trk::Track *Tr)
static void clusterTrackMap (SiTrackMakerEventData_xk &data, Trk::Track *Tr)
static MsgStream & dumpevent (SiTrackMakerEventData_xk &data, MsgStream &out)

Detailed Description

InDet::SiTrackMaker_xk is algorithm which produce Trk::Track started from 3 space points information of SCT and Pixels in the road of InDetDD::SiDetectorElement* sorted in propagation order.

In AthenaMT, event dependent cache inside SiTrackMaker_xk is not preferred. SiTrackMakerEventData_xk class holds event dependent data for SiTrackMaker_xk. Its object is instantiated in SiSPSeededTrackFinder::execute.

Author
Igor..nosp@m.Gavr.nosp@m.ilenk.nosp@m.o@ce.nosp@m.rn.ch

Definition at line 63 of file SiTrackMaker_xk.h.

Member Enumeration Documentation

◆ kNStatEtaTypes

Enumerator
kUsedSeedsEta 
kSeedsWithTracksEta 

Definition at line 223 of file SiTrackMaker_xk.h.

◆ statAllTypes

Enumerator
kTotalInputSeeds 
kTotalUsedSeeds 
kTotalNoTrackPar 
kTotalBremSeeds 
kTwoClusters 
kWrongInit 
kWrongRoad 
kNoTrack 
kNotNewTrack 
kBremAttempt 
kOutputTracks 
kExtraTracks 
kBremTracks 
kDESize 
kSeedsWithTracks 

Definition at line 205 of file SiTrackMaker_xk.h.

205 {
213 kNoTrack,
219 kDESize,
221 };

Constructor & Destructor Documentation

◆ SiTrackMaker_xk() [1/3]

InDet::SiTrackMaker_xk::SiTrackMaker_xk ( const std::string & t,
const std::string & n,
const IInterface * p )

Definition at line 27 of file SiTrackMaker_xk.cxx.

29 : base_class(t, n, p)
30{
31
32}

◆ ~SiTrackMaker_xk()

virtual InDet::SiTrackMaker_xk::~SiTrackMaker_xk ( )
virtualdefault

◆ SiTrackMaker_xk() [2/3]

InDet::SiTrackMaker_xk::SiTrackMaker_xk ( )
privatedelete

◆ SiTrackMaker_xk() [3/3]

InDet::SiTrackMaker_xk::SiTrackMaker_xk ( const SiTrackMaker_xk & )
privatedelete

Member Function Documentation

◆ clusterTrackMap()

void InDet::SiTrackMaker_xk::clusterTrackMap ( SiTrackMakerEventData_xk & data,
Trk::Track * Tr )
staticprivate

Definition at line 1208 of file SiTrackMaker_xk.cxx.

1209{
1211 m = Tr->measurementsOnTrack()->begin(),
1212 me = Tr->measurementsOnTrack()->end ();
1213
1214 for (; m!=me; ++m) {
1215 const Trk::PrepRawData* prd = static_cast<const Trk::RIO_OnTrack*>(*m)->prepRawData();
1216 if (prd) data.clusterTrack().insert(std::make_pair(prd, Tr));
1217 }
1218}
char data[hepevt_bytes_allocation_ATLAS]
Definition HepEvt.cxx:11
DataModel_detail::const_iterator< DataVector > const_iterator
Standard const_iterator.
Definition DataVector.h:838
const DataVector< const MeasurementBase > * measurementsOnTrack() const
return a pointer to a vector of MeasurementBase (NOT including any that come from outliers).

◆ detectorElementsSelection()

void InDet::SiTrackMaker_xk::detectorElementsSelection ( SiTrackMakerEventData_xk & data,
std::vector< const InDetDD::SiDetectorElement * > & DE )
staticprivate

Definition at line 1085 of file SiTrackMaker_xk.cxx.

1087{
1088 std::vector<const InDetDD::SiDetectorElement*>::iterator d = DE.begin();
1089 while (d!=DE.end()) {
1090 if ((*d)->isPixel()) {
1091 if (!data.pix()) {
1092 d = DE.erase(d);
1093 continue;
1094 }
1095 } else if (!data.sct()) {
1096 d = DE.erase(d);
1097 continue;
1098 }
1099 ++d;
1100 }
1101}

◆ dump()

MsgStream & InDet::SiTrackMaker_xk::dump ( SiTrackMakerEventData_xk & data,
MsgStream & out ) const
override

Definition at line 399 of file SiTrackMaker_xk.cxx.

400{
401 out<<std::endl;
402 if (data.nprint()) return dumpevent(data, out);
403 return dumpconditions(out);
404}
static MsgStream & dumpevent(SiTrackMakerEventData_xk &data, MsgStream &out)
MsgStream & dumpconditions(MsgStream &out) const

◆ dumpconditions()

MsgStream & InDet::SiTrackMaker_xk::dumpconditions ( MsgStream & out) const
private

Definition at line 410 of file SiTrackMaker_xk.cxx.

411{
412
413 int n = 62-m_tracksfinder.type().size();
414 std::string s4; for (int i=0; i<n; ++i) s4.append(" "); s4.append("|");
415
416 std::string fieldmode[9] ={"NoField" ,"ConstantField","SolenoidalField",
417 "ToroidalField" ,"Grid3DField" ,"RealisticField" ,
418 "UndefinedField","AthenaField" , "?????" };
419
421
422 // Get AtlasFieldCache
423 SG::ReadCondHandle<AtlasFieldCacheCondObj> readHandle{m_fieldCondObjInputKey};
424 const AtlasFieldCacheCondObj* fieldCondObj{*readHandle};
425 if (fieldCondObj) {
426 MagField::AtlasFieldCache fieldCache;
427 fieldCondObj->getInitializedCache(fieldCache);
428 if (!fieldCache.solenoidOn()) fieldModeEnum = Trk::NoField;
429 }
430
431 Trk::MagneticFieldProperties fieldprop(fieldModeEnum);
432 int mode = fieldprop.magneticFieldMode();
433 if (mode<0 || mode>8 ) mode = 8;
434
435 n = 62-fieldmode[mode].size();
436 std::string s3; for (int i=0; i<n; ++i) s3.append(" "); s3.append("|");
437
438 n = 62-m_roadmaker.type().size();
439 std::string s6; for (int i=0; i<n; ++i) s6.append(" "); s6.append("|");
440
441 out<<"|----------------------------------------------------------------------"
442 <<"-------------------|"
443 <<std::endl;
444 out<<"| Tool for road builder | "<<m_roadmaker .type() <<s6<<std::endl;
445 out<<"} Tool for track finding | "<<m_tracksfinder.type()<<s4<<std::endl;
446 out<<"| Use association tool ? | "
447 <<std::setw(12)<<m_useassoTool
448 << std::endl;
449 out<<"| Magnetic field mode | "<<fieldmode[mode] <<s3<<std::endl;
450 out<<"| Use seeds filter ? | "<<std::setw(12)<<m_seedsfilter
451 <<" |"<<std::endl;
452 out<<"| Min pT of track (MeV) | "<<std::setw(12)<<std::setprecision(5)<<m_pTmin
453 <<" |"<<std::endl;
454 out<<"| Max Xi2 for cluster | "<<std::setw(12)<<std::setprecision(5)<<m_xi2max
455 <<" |"<<std::endl;
456 out<<"| Max Xi2 for outlayer | "<<std::setw(12)<<std::setprecision(5)<<m_xi2maxNoAdd
457 <<" |"<<std::endl;
458 out<<"| Max Xi2 for link | "<<std::setw(12)<<std::setprecision(5)<<m_xi2maxlink
459 <<" |"<<std::endl;
460 out<<"| Min number of clusters | "<<std::setw(12)<<m_nclusmin
461 <<" |"<<std::endl;
462 out<<"| Min number of wclusters | "<<std::setw(12)<<m_nwclusmin
463 <<" |"<<std::endl;
464 out<<"| Max number holes | "<<std::setw(12)<<m_nholesmax
465 <<" |"<<std::endl;
466 out<<"| Max holes gap | "<<std::setw(12)<<m_dholesmax
467 <<" |"<<std::endl;
468 out<<"|----------------------------------------------------------------------"
469 <<"-------------------|"
470 <<std::endl;
471 return out;
472}
void getInitializedCache(MagField::AtlasFieldCache &cache) const
get B field cache for evaluation as a function of 2-d or 3-d position.
BooleanProperty m_useassoTool
BooleanProperty m_seedsfilter
DoubleProperty m_xi2maxNoAdd
ToolHandle< InDet::ISiCombinatorialTrackFinder > m_tracksfinder
IntegerProperty m_dholesmax
ToolHandle< InDet::ISiDetElementsRoadMaker > m_roadmaker
DoubleProperty m_xi2maxlink
IntegerProperty m_nholesmax
Trk::MagneticFieldMode m_fieldModeEnum
IntegerProperty m_nwclusmin
SG::ReadCondHandleKey< AtlasFieldCacheCondObj > m_fieldCondObjInputKey
IntegerProperty m_nclusmin
bool solenoidOn() const
status of the magnets
MagneticFieldMode
MagneticFieldMode describing the field setup within a volume.
@ NoField
Field is set to 0., 0., 0.,.

◆ dumpevent()

MsgStream & InDet::SiTrackMaker_xk::dumpevent ( SiTrackMakerEventData_xk & data,
MsgStream & out )
staticprivate

Definition at line 479 of file SiTrackMaker_xk.cxx.

480{
481 out<<"|---------------------------------------------------------------------|"
482 <<std::endl;
483 out<<"| Number input seeds | "<<std::setw(12)<<data.inputseeds()
484 <<" |"<<std::endl;
485 out<<"| Number good seeds | "<<std::setw(12)<<data.goodseeds()
486 <<" |"<<std::endl;
487 out<<"| Number output tracks | "<<std::setw(12)<<data.findtracks()
488 <<" |"<<std::endl;
489 out<<"|---------------------------------------------------------------------|"
490 <<std::endl;
491 return out;
492}

◆ dumpStatistics()

MsgStream & InDet::SiTrackMaker_xk::dumpStatistics ( MsgStream & out) const
private

Definition at line 199 of file SiTrackMaker_xk.cxx.

200{
201
202 out<<"|----------------------|--------------|--------------|--------------|--------------|--------------|"
203 <<std::endl;
204 out<<"| Kind of seed | PPP | PPS | PSS | SSS | ALL |"
205 <<std::endl;
206 out<<"|----------------------|--------------|--------------|--------------|--------------|--------------|"
207 <<std::endl;
208 out<<"| Input seeds | "
209 <<std::setw(12)<<m_totalInputSeeds[0]<<" | "
210 <<std::setw(12)<<m_totalInputSeeds[1]<<" | "
211 <<std::setw(12)<<m_totalInputSeeds[2]<<" | "
212 <<std::setw(12)<<m_totalInputSeeds[3]<<" | "
213 <<std::setw(12)<<(m_totalInputSeeds[0]+m_totalInputSeeds[1]+m_totalInputSeeds[2]+m_totalInputSeeds[3])
214 <<" |"<<std::endl;
215
216 out<<"| No track parameters | "
217 <<std::setw(12)<<m_totalNoTrackPar[0]<<" | "
218 <<std::setw(12)<<m_totalNoTrackPar[1]<<" | "
219 <<std::setw(12)<<m_totalNoTrackPar[2]<<" | "
220 <<std::setw(12)<<m_totalNoTrackPar[3]<<" | "
221 <<std::setw(12)<<(m_totalNoTrackPar[0]+m_totalNoTrackPar[1]+ m_totalNoTrackPar[2]+m_totalNoTrackPar[3])
222 <<" |"<<std::endl;
223
224 out<<"| Used seeds | "
225 <<std::setw(12)<<m_totalUsedSeeds[0]<<" | "
226 <<std::setw(12)<<m_totalUsedSeeds[1]<<" | "
227 <<std::setw(12)<<m_totalUsedSeeds[2]<<" | "
228 <<std::setw(12)<<m_totalUsedSeeds[3]<<" | "
229 <<std::setw(12)<<(m_totalUsedSeeds[0]+m_totalUsedSeeds[1]+ m_totalUsedSeeds[2]+m_totalUsedSeeds[3])
230 <<" |"<<std::endl;
231
232 out<<"| Used seeds brem | "
233 <<std::setw(12)<<m_totalBremSeeds[0]<<" | "
234 <<std::setw(12)<<m_totalBremSeeds[1]<<" | "
235 <<std::setw(12)<<m_totalBremSeeds[2]<<" | "
236 <<std::setw(12)<<m_totalBremSeeds[3]<<" | "
237 <<std::setw(12)<<(m_totalBremSeeds[0]+m_totalBremSeeds[1]+m_totalBremSeeds[2]+m_totalBremSeeds[3])
238 <<" |"<<std::endl;
239
240 double tdetsize = 0.;
241 int goodseed = 0;
243 tdetsize+= m_deSize[i];
244 goodseed+= m_totalUsedSeeds[i];
245 if(m_totalUsedSeeds[i] > 0.) m_deSize[i]= m_deSize[i]/m_totalUsedSeeds[i];
246 }
247 if(goodseed > 0) tdetsize/=double(goodseed);
248
249 out<<"| Det elements in road | "
250 <<std::setw(12)<<std::setprecision(4)<<m_deSize[0]<<" | "
251 <<std::setw(12)<<std::setprecision(5)<<m_deSize[1]<<" | "
252 <<std::setw(12)<<std::setprecision(5)<<m_deSize[2]<<" | "
253 <<std::setw(12)<<std::setprecision(5)<<m_deSize[3]<<" | "
254 <<std::setw(12)<<std::setprecision(5)<<tdetsize
255 <<" |"<<std::endl;
256
257
258 out<<"| Two clusters on DE | "
259 <<std::setw(12)<<m_twoClusters[0]<<" | "
260 <<std::setw(12)<<m_twoClusters[1]<<" | "
261 <<std::setw(12)<<m_twoClusters[2]<<" | "
262 <<std::setw(12)<<m_twoClusters[3]<<" | "
263 <<std::setw(12)<<(m_twoClusters[0]+m_twoClusters[1]+m_twoClusters[2]+m_twoClusters[3])
264 <<" |"<<std::endl;
265
266 out<<"| Wrong DE road | "
267 <<std::setw(12)<<m_wrongRoad[0]<<" | "
268 <<std::setw(12)<<m_wrongRoad[1]<<" | "
269 <<std::setw(12)<<m_wrongRoad[2]<<" | "
270 <<std::setw(12)<<m_wrongRoad[3]<<" | "
271 <<std::setw(12)<<(m_wrongRoad[0]+m_wrongRoad[1]+m_wrongRoad[2]+m_wrongRoad[3])
272 <<" |"<<std::endl;
273
274 out<<"| Wrong initialization | "
275 <<std::setw(12)<<m_wrongInit[0]<<" | "
276 <<std::setw(12)<<m_wrongInit[1]<<" | "
277 <<std::setw(12)<<m_wrongInit[2]<<" | "
278 <<std::setw(12)<<m_wrongInit[3]<<" | "
279 <<std::setw(12)<<(m_wrongInit[0]+m_wrongInit[1]+m_wrongInit[2]+m_wrongInit[3])
280 <<" |"<<std::endl;
281
282 out<<"| Can not find track | "
283 <<std::setw(12)<<m_noTrack[0]<<" | "
284 <<std::setw(12)<<m_noTrack[1]<<" | "
285 <<std::setw(12)<<m_noTrack[2]<<" | "
286 <<std::setw(12)<<m_noTrack[3]<<" | "
287 <<std::setw(12)<<(m_noTrack[0]+m_noTrack[1]+m_noTrack[2]+m_noTrack[3])
288 <<" |"<<std::endl;
289
290 out<<"| It is not new track | "
291 <<std::setw(12)<<m_notNewTrack[0]<<" | "
292 <<std::setw(12)<<m_notNewTrack[1]<<" | "
293 <<std::setw(12)<<m_notNewTrack[2]<<" | "
294 <<std::setw(12)<<m_notNewTrack[3]<<" | "
295 <<std::setw(12)<<(m_notNewTrack[0]+m_notNewTrack[1]+m_notNewTrack[2]+m_notNewTrack[3])
296 <<" |"<<std::endl;
297
298 out<<"| Attempts brem model | "
299 <<std::setw(12)<<m_bremAttempt[0]<<" | "
300 <<std::setw(12)<<m_bremAttempt[1]<<" | "
301 <<std::setw(12)<<m_bremAttempt[2]<<" | "
302 <<std::setw(12)<<m_bremAttempt[3]<<" | "
303 <<std::setw(12)<<(m_bremAttempt[0]+m_bremAttempt[1]+m_bremAttempt[2]+m_bremAttempt[3])
304 <<" |"<<std::endl;
305
306 out<<"| Output tracks | "
307 <<std::setw(12)<<m_outputTracks[0]<<" | "
308 <<std::setw(12)<<m_outputTracks[1]<<" | "
309 <<std::setw(12)<<m_outputTracks[2]<<" | "
310 <<std::setw(12)<<m_outputTracks[3]<<" | "
311 <<std::setw(12)<<(m_outputTracks[0]+m_outputTracks[1]+m_outputTracks[2]+m_outputTracks[3])
312 <<" |"<<std::endl;
313
314 out<<"| Output extra tracks | "
315 <<std::setw(12)<<m_extraTracks[0]<<" | "
316 <<std::setw(12)<<m_extraTracks[1]<<" | "
317 <<std::setw(12)<<m_extraTracks[2]<<" | "
318 <<std::setw(12)<<m_extraTracks[3]<<" | "
319 <<std::setw(12)<<(m_extraTracks[0]+m_extraTracks[1]+m_extraTracks[2]+m_extraTracks[3])
320 <<" |"<<std::endl;
321
322 out<<"| Output tracks brem | "
323 <<std::setw(12)<<m_bremTracks[0]<<" | "
324 <<std::setw(12)<<m_bremTracks[1]<<" | "
325 <<std::setw(12)<<m_bremTracks[2]<<" | "
326 <<std::setw(12)<<m_bremTracks[3]<<" | "
327 <<std::setw(12)<<(m_bremTracks[0]+m_bremTracks[1]+m_bremTracks[2]+m_bremTracks[3])
328 <<" |"<<std::endl;
329
330 out<<"|----------------------|--------------|--------------|--------------|--------------|--------------|--------------|"
331 <<std::endl;
332 out<<"| Seeds with track | "
333 <<std::setw(12)<<m_seedsWithTrack[0]<<" | "
334 <<std::setw(12)<<m_seedsWithTrack[1]<<" | "
335 <<std::setw(12)<<m_seedsWithTrack[2]<<" | "
336 <<std::setw(12)<<m_seedsWithTrack[3]<<" | "
337 <<std::setw(12)<<(m_seedsWithTrack[0]+m_seedsWithTrack[1]+m_seedsWithTrack[2]+m_seedsWithTrack[3])
338 <<" | Number seeds |"<<std::endl;
339 out<<"|----------------------|--------------|--------------|--------------|--------------|--------------|--------------|"
340 <<std::endl;
341
342 std::vector<std::pair<int,std::string>> rapidityTablePrint;
343 //Defining the range values to be printed
344 rapidityTablePrint.emplace_back(0,std::string("| Track/Used 0.0-0.5 | "));
345 rapidityTablePrint.emplace_back(1,std::string("| 0.5-1.0 | "));
346 rapidityTablePrint.emplace_back(2,std::string("| 1.0-1.5 | "));
347 rapidityTablePrint.emplace_back(3,std::string("| 1.5-2.0 | "));
348 rapidityTablePrint.emplace_back(4,std::string("| 2.0-2.5 | "));
349 rapidityTablePrint.emplace_back(5,std::string("| 2.5-3.0 | "));
350 rapidityTablePrint.emplace_back(6,std::string("| 3.0-3.5 | "));
351 rapidityTablePrint.emplace_back(7,std::string("| 3.5-4.0 | "));
352
353
355
356 std::array<double,SiCombinatorialTrackFinderData_xk::kNSeedTypes+1> pu{};
357 pu.fill(0);
358
359 double totalUsedSeedsEta = 0.;
360
361 for(int j = 0; j != SiCombinatorialTrackFinderData_xk::kNSeedTypes; ++j){
362
363 if(m_usedSeedsEta[j][i]!=0.) pu[j]=m_seedsWithTracksEta[j][i]/m_usedSeedsEta[j][i];
364 totalUsedSeedsEta += m_usedSeedsEta[j][i];
365
366 }
367
368 if(totalUsedSeedsEta!=0.) {
369
370 for(int j = 0; j != SiCombinatorialTrackFinderData_xk::kNSeedTypes; ++j){
371
372 pu[SiCombinatorialTrackFinderData_xk::kNSeedTypes] += m_seedsWithTracksEta[j][i];
373
374 }
376 }
377
378 out<<rapidityTablePrint.at(i).second;
379 out << std::setw(12) << std::setprecision(4) << pu[0] << " | "
380 << std::setw(12) << std::setprecision(4) << pu[1] << " | "
381 << std::setw(12) << std::setprecision(4) << pu[2] << " | "
382 << std::setw(12) << std::setprecision(4) << pu[3] << " | "
383 << std::setw(12) << std::setprecision(4) << pu[4] << " | "
384 << std::setw(12)
385 << static_cast<int>(m_seedsWithTracksEta[0][i]) +
386 static_cast<int>(m_seedsWithTracksEta[1][i]) +
387 static_cast<int>(m_seedsWithTracksEta[2][i]) +
388 static_cast<int>(m_seedsWithTracksEta[3][i])
389 << " |" << std::endl;
390 }
391
392 out<<"|----------------------|--------------|--------------|--------------|--------------|--------------|--------------|"
393 <<std::endl;
394
395 return out;
396
397}

◆ endEvent()

void InDet::SiTrackMaker_xk::endEvent ( SiTrackMakerEventData_xk & data) const
overridevirtual

End event for track finder tool

correction to exclude memory fragmentation

Definition at line 602 of file SiTrackMaker_xk.cxx.

603{
605 m_tracksfinder->endEvent(data.combinatorialData());
606
608 data.clusterTrack().clear();
609
610 // end event for seed to track tool
611 if (m_seedsegmentsWrite) m_seedtrack->endEvent(data.conversionData());
612
613 // fill statistics
614 {
615 std::lock_guard<std::mutex> lock(m_counterMutex);
616
617 for(int K = 0; K != SiCombinatorialTrackFinderData_xk::kNSeedTypes; ++K) {
618 for(int r = 0; r != 8; ++r) {
619 m_usedSeedsEta[K][r] += data.summaryStatUsedInTrack()[kUsedSeedsEta][K][r];
620 m_seedsWithTracksEta[K][r] += data.summaryStatUsedInTrack()[kSeedsWithTracksEta][K][r];
621 }
622 }
623 }
624 for (int K = 0; K != SiCombinatorialTrackFinderData_xk::kNSeedTypes; ++K) {
625 m_totalInputSeeds[K] += data.summaryStatAll()[kTotalInputSeeds][K];
626 m_totalUsedSeeds[K] = m_totalUsedSeeds[K] + data.summaryStatAll()[kTotalUsedSeeds][K];
627 m_totalNoTrackPar[K] += data.summaryStatAll()[kTotalNoTrackPar][K];
628 m_totalBremSeeds[K] += data.summaryStatAll()[kTotalBremSeeds][K];
629 m_twoClusters[K] += data.summaryStatAll()[kTwoClusters][K];
630 m_wrongRoad[K] += data.summaryStatAll()[kWrongRoad][K];
631 m_wrongInit[K] += data.summaryStatAll()[kWrongInit][K];
632 m_noTrack[K] += data.summaryStatAll()[kNoTrack][K];
633 m_notNewTrack[K] += data.summaryStatAll()[kNotNewTrack][K];
634 m_bremAttempt[K] += data.summaryStatAll()[kBremAttempt][K];
635 m_outputTracks[K] += data.summaryStatAll()[kOutputTracks][K];
636 m_extraTracks[K] += data.summaryStatAll()[kExtraTracks][K];
637 m_bremTracks[K] += data.summaryStatAll()[kBremTracks][K];
638 m_seedsWithTrack[K] += data.summaryStatAll()[kSeedsWithTracks][K];
639 m_deSize[K] = m_deSize[K] + data.summaryStatAll()[kDESize][K];
640 }
641 // Print event information
642 //
643 if (msgLevel()<=0) {
644 data.nprint() = 1;
645 dump(data, msg(MSG::DEBUG));
646 }
647}
ToolHandle< InDet::ISeedToTrackConversionTool > m_seedtrack
BooleanProperty m_seedsegmentsWrite
MsgStream & dump(SiTrackMakerEventData_xk &data, MsgStream &out) const override
int r
Definition globals.cxx:22
MsgStream & msg
Definition testRead.cxx:32

◆ finalize()

StatusCode InDet::SiTrackMaker_xk::finalize ( )
overridevirtual

Definition at line 186 of file SiTrackMaker_xk.cxx.

187{
188 MsgStream &out = msg(MSG::INFO);
189 out << "::finalize() -- statistics:" <<std::endl;
190 dumpStatistics(out);
191 out<<endmsg;
192 return AlgTool::finalize();
193}
#define endmsg
MsgStream & dumpStatistics(MsgStream &out) const

◆ getAtaPlane()

std::unique_ptr< Trk::TrackParameters > InDet::SiTrackMaker_xk::getAtaPlane ( MagField::AtlasFieldCache & fieldCache,
SiTrackMakerEventData_xk & data,
bool sss,
const std::vector< const Trk::SpacePoint * > & SP,
const EventContext & ctx ) const
private

we need at least three space points on the seed.

for tracklets we select first, middle, and last spacepoint of the seed to improve pT estimate

for tracklets we select middle, 3rd-quarter, and last spacepoint of the seed

for tracklets we select first, penultimate, and last spacepoint of the seed

for tracklets we select first, second, and third spacepoint of the seed

get the first cluster on the first hit

and use the surface from this cluster as our reference plane

write the global positions into arrays. This includes an improved position estimate for strip spacepoints. If this improvement fails, the method can return false -> then we abort

translate second and third SP w.r.t first one

distance of second SP to first in transverse plane Also happens to be u-coordinate of second SP in conformal mapping

denominator for conformal mapping

coordinate system for conformal mapping - this is local x

u-coordinate of third SP in conformal mapping

v-coordinate of third SP in conformal mapping

A,B are slope and intercept of the straight line in the u,v plane connecting the three points.

Keep in mind that v1 == 0

From inserting A into linear equation. Note that Igor sneaks in a factor two here

Curvature estimate. (2R)²=(1+A²)/b² => 1/2R = b/sqrt(1+A²) = B / sqrt(1+A²).

estimate of the track dz/dr (1/tanTheta), corrected for curvature effects

local x of the surface in the global frame

local y of the surface in the global frame

centre of the surface in the global frame

location of the first SP w.r.t centre of the surface

local x and y - project onto local axes

silently switch off field if solenoid is off

if we are not running with "no field":

get the field at our first SP

field is in kiloTesla. So here we check for more than 0.1 Tesla

phi estimate

theta estimate

inverse transverse momentum estimate

if we have low field, use a straight-line estimate

note: Now no curvature correction

no pt estimate, assume min pt

treat absence of solenoid like the low-field case

apply the pt on the initial parameter estimate, with some margin

update qoverp

qoverp from qoverpt and theta

ref point = first SP

never done in main ATLAS tracking. Would otherwise check if the seed is compatible with a hadronic ROI

now we can return the initial track parameters we built, parameterised using the ref surface. Pass a nullptr for the covariance

Definition at line 879 of file SiTrackMaker_xk.cxx.

885{
887 if (theSeed.size() < 3) return nullptr;
888
889 std::vector<const Trk::SpacePoint*> SP;
890
892 if (m_trackletPoints == 1) {
893 unsigned int middleIdx = theSeed.size() == 3 ? 1 : theSeed.size()/2;
894 SP = {theSeed[0], theSeed[middleIdx], theSeed.back()};
895 }
896 // for tracklets we select last 3 spacepoints of the seed
897 else if (m_trackletPoints == 2) {
898 SP = {theSeed[theSeed.size() - 3], theSeed[theSeed.size() - 2], theSeed.back()};
899 }
901 else if (m_trackletPoints == 3) {
902 unsigned int middleIdx = theSeed.size() == 3 ? 0 : theSeed.size()/2;
903 unsigned int quarterIdx = theSeed.size() == 3 ? 1 : 3*theSeed.size()/4;
904 SP = {theSeed[middleIdx], theSeed[quarterIdx], theSeed.back()};
905 }
907 else if (m_trackletPoints == 4) {
908 SP = {theSeed[0], theSeed[theSeed.size() - 2], theSeed.back()};
909 }
911 else if (m_trackletPoints == 5) {
912 SP = {theSeed[0], theSeed[1], theSeed[2]};
913 }
914
915
917 const Trk::PrepRawData* cl = SP[0]->clusterList().first;
918 if (!cl) return nullptr;
920 const Trk::PlaneSurface* pla =
921 static_cast<const Trk::PlaneSurface*>(&cl->detectorElement()->surface());
922 if (!pla) return nullptr;
923
926 double p0[3],p1[3],p2[3];
927 if (!globalPositions(*(SP[0]),*(SP[1]),*(SP[2]),p0,p1,p2)) return nullptr;
928
930 double x0 = p0[0] ;
931 double y0 = p0[1] ;
932 double z0 = p0[2] ;
933 double x1 = p1[0]-x0;
934 double y1 = p1[1]-y0;
935 double x2 = p2[0]-x0;
936 double y2 = p2[1]-y0;
937 double z2 = p2[2]-z0;
938
941 double u1 = 1./sqrt(x1*x1+y1*y1) ;
943 double rn = x2*x2+y2*y2 ;
944 double r2 = 1./rn ;
946 double a = x1*u1 ;
947 double b = y1*u1 ;
949 double u2 = (a*x2+b*y2)*r2 ;
951 double v2 = (a*y2-b*x2)*r2 ;
954 double A = v2/(u2-u1) ;
955 double B = 2.*(v2-A*u2) ;
956 double C = B/sqrt(1.+A*A) ;
957 double T = z2*sqrt(r2)/(1.+.04*C*C*rn);
958 if(m_ITKGeometry){
959 T = std::abs(C) > 1.e-6 ? (z2*C)/asin(C*sqrt(rn)) : z2/sqrt(rn);
960 }
961
962 const Amg::Transform3D& Tp = pla->transform();
963
965 double Ax[3] = {Tp(0,0),Tp(1,0),Tp(2,0)};
967 double Ay[3] = {Tp(0,1),Tp(1,1),Tp(2,1)};
969 double D [3] = {Tp(0,3),Tp(1,3),Tp(2,3)};
971 double d[3] = {x0-D[0],y0-D[1],z0-D[2]};
973 data.par()[0] = d[0]*Ax[0]+d[1]*Ax[1]+d[2]*Ax[2];
974 data.par()[1] = d[0]*Ay[0]+d[1]*Ay[1]+d[2]*Ay[2];
975
978 if (!fieldCache.solenoidOn()) fieldModeEnum = Trk::NoField;
979
980 Trk::MagneticFieldProperties fieldprop(fieldModeEnum);
982 if (fieldprop.magneticFieldMode() > 0) {
983
985 double H[3],gP[3] ={x0,y0,z0};
986 fieldCache.getFieldZR(gP, H);
987
990 if (fabs(H[2])>.0001) {
992 data.par()[2] = atan2(b+a*A,a-b*A);
994 data.par()[3] = atan2(1.,T) ;
996 data.par()[5] = -C/(300.*H[2]) ;
997 } else {
1000 T = z2*sqrt(r2) ;
1001 data.par()[2] = atan2(y2,x2);
1002 data.par()[3] = atan2(1.,T) ;
1003 data.par()[5] = m_ITKGeometry ? 0.9/m_pTmin : 1./m_pTmin ;
1004 }
1005 }
1007 else {
1008 T = z2*sqrt(r2) ;
1009 data.par()[2] = atan2(y2,x2);
1010 data.par()[3] = atan2(1.,T) ;
1011 data.par()[5] = m_ITKGeometry ? 0.9/m_pTmin : 1./m_pTmin ;
1012 }
1013
1014 double pTm = pTmin(SP[0]->eta()); // all spacepoints should have approx. same eta
1015
1017 if (m_ITKGeometry){
1018 if(std::abs(data.par()[5])*pTm > 1) return nullptr;
1019 }
1020 else if(std::abs(data.par()[5])*m_pTmin > 1.1) return nullptr;
1021
1023 data.par()[4] = data.par()[5]/sqrt(1.+T*T);
1024 data.par()[6] = x0 ;
1025 data.par()[7] = y0 ;
1026 data.par()[8] = z0 ;
1027
1030 if (sss && !isHadCaloCompatible(data)) return nullptr;
1031
1034 std::unique_ptr<Trk::TrackParameters> T0 = pla->createUniqueTrackParameters(data.par()[0],
1035 data.par()[1],
1036 data.par()[2],
1037 data.par()[3],
1038 data.par()[4],
1039 std::nullopt);
1040
1041 if(m_ITKGeometry && m_tracksfinder->pTseed(data.combinatorialData(),*T0,SP,ctx) < pTm) return nullptr;
1042 else return T0;
1043
1044}
Scalar eta() const
pseudorapidity method
static Double_t a
static Double_t Tp(Double_t *t, Double_t *par)
#define H(x, y, z)
Definition MD5.cxx:114
bool globalPositions(const Trk::SpacePoint &s0, const Trk::SpacePoint &s1, const Trk::SpacePoint &s2, double *p0, double *p1, double *p2) const
BooleanProperty m_ITKGeometry
double pTmin(double eta) const
IntegerProperty m_trackletPoints
bool isHadCaloCompatible(SiTrackMakerEventData_xk &data) const
void getFieldZR(const double *ATH_RESTRICT xyz, double *ATH_RESTRICT bxyz, double *ATH_RESTRICT deriv=nullptr)
get B field valaue on the z-r plane at given position works only inside the solenoid.
virtual Surface::ChargedTrackParametersUniquePtr createUniqueTrackParameters(double l1, double l2, double phi, double theta, double qop, std::optional< AmgSymMatrix(5)> cov=std::nullopt) const override final
Use the Surface as a ParametersBase constructor, from local parameters - charged.
const Amg::Transform3D & transform() const
Returns HepGeom::Transform3D by reference.
struct color C
Eigen::Affine3d Transform3D
unsigned long long T
cl
print [x.__class__ for x in toList(dqregion.getSubRegions()) ]

◆ getTracks() [1/2]

std::list< Trk::Track * > InDet::SiTrackMaker_xk::getTracks ( const EventContext & ctx,
SiTrackMakerEventData_xk & data,
const std::vector< const Trk::SpacePoint * > & Sp ) const
overridevirtual

incremenet seed counter

0 or 3 typically (PPP or SSS), other numbers indicate PPS/PSS (1/2)

eta of the seed, rounded down to leading digit via int-conversion

more counter incrementation

prepare output list

if we run the SI track maker without using the Si, this becomes a trivial task...

check seed quality.

this checks if all of the clusters on the seed are already on one single existing track. If not, we consider the seed to be "good"

read the B-field cache

Get initial parameters estimation from our seed

if we failed to get the initial parameters, we bail out. Can happen in certain pathological cases (e.g. malformed strip hits), or if we would be running with calo-ROI strip seeds (we aren't)

otherwise, increment the 'good seeds' counter

Now, obtain a search road of detector elements. This is done by extrapolating our estimated starting parameters through the detector and collecting all detector elements reasonably close to the projected trajectory. This will populate the 'DE" list.

if we don't use all of pix and SCT, filter our list, erasing any that don't fit our requirements

if we did not find sufficient detector elements to fulfill the minimum cluster requirement, bail out. We will not be able to build a track satisfying the cuts.

update statistics tables - we have sufficient detector elements to have a chance of finding a track!

prepare a list of global positions

update another counter

Find possible list of tracks using space points space points information

Note: The branch below is the one taken in ATLAS default inside-out tracking for run-3

update stat tables

update the cluster-track-map to allow to filter any upcoming seeds with hits that are already taken

require sufficient free clusters on track

Definition at line 653 of file SiTrackMaker_xk.cxx.

655{
657 ++data.inputseeds();
658
660 int K = kindSeed(Sp);
662 int r = rapidity(Sp);
663
665 ++data.summaryStatAll()[kTotalInputSeeds][K];
666
668 std::list<Trk::Track*> tracks;
670 if (!data.pix() && !data.sct()) return tracks;
671
673 bool isGoodSeed{true};
676 if (m_seedsfilter) isGoodSeed=newSeed(data, Sp);
677 if (!m_seedsegmentsWrite && !isGoodSeed) {
678 return tracks;
679 }
680
681 // Get AtlasFieldCache
682 MagField::AtlasFieldCache fieldCache;
683
685 SG::ReadCondHandle<AtlasFieldCacheCondObj> readHandle{m_fieldCondObjInputKey, ctx};
686 const AtlasFieldCacheCondObj* fieldCondObj{*readHandle};
687 if (fieldCondObj == nullptr) {
688 ATH_MSG_ERROR("InDet::SiTrackMaker_xk::getTracks: Failed to retrieve AtlasFieldCacheCondObj with key " << m_fieldCondObjInputKey.key());
689 return tracks;
690 }
691 fieldCondObj->getInitializedCache (fieldCache);
692
694 std::unique_ptr<Trk::TrackParameters> Tp = nullptr;
695 Tp = getAtaPlane(fieldCache, data, false, Sp, ctx);
696 if (m_seedsegmentsWrite && Tp) {
697 m_seedtrack->executeSiSPSeedSegments(data.conversionData(), Tp.get(), isGoodSeed, Sp);
698 }
699
700 if (m_seedsegmentsWrite && !isGoodSeed) {
701 return tracks;
702 }
706 if (!Tp) {
707 ++data.summaryStatAll()[kTotalNoTrackPar][K];
708 return tracks;
709 }
710
712 ++data.goodseeds();
713
718
719 std::vector<const InDetDD::SiDetectorElement*> DE;
720
722 if (!m_cosmicTrack){
723 m_roadmaker->detElementsRoad(ctx, fieldCache, *Tp, Trk::alongMomentum, DE, data.roadMakerData());
724 }
725 else{
726 m_roadmaker->detElementsRoad(ctx, fieldCache, *Tp, Trk::oppositeMomentum, DE, data.roadMakerData());
727 }
728 } else {
729 int road_length = m_trigInDetRoadPredictorTool->getRoad(Sp, DE, ctx);
730 if(road_length == 0) {
731 return tracks;
732 }
733 }
734
736 if (!data.pix() || !data.sct()) detectorElementsSelection(data, DE);
737
740 if ( static_cast<int>(DE.size()) < m_nclusmin) {
741 return tracks;
742 }
743
745 data.summaryStatAll()[kDESize][K] += double(DE.size());
746 ++data.summaryStatAll()[kTotalUsedSeeds][K];
747
749 std::vector<Amg::Vector3D> Gp;
750
752 ++data.summaryStatUsedInTrack()[kUsedSeedsEta][K][r];
753
756 if (!m_useBremModel) {
758 Trk::Track* newTrack = m_trigInDetTrackFollowingTool->getTrack(Sp, DE, ctx);
759 if(newTrack != nullptr) tracks.push_back(newTrack);
760 }
761 else {
762 tracks = m_tracksfinder->getTracks(data.combinatorialData(), *Tp, Sp, Gp, DE, data.clusterTrack(),ctx);
763 }
764 } else if (!m_useCaloSeeds) {
765 ++data.summaryStatAll()[kTotalBremSeeds][K];
766 tracks = m_tracksfinder->getTracksWithBrem(data.combinatorialData(), *Tp, Sp, Gp, DE, data.clusterTrack(), false,ctx);
767 }
769 else if (isCaloCompatible(data)) {
770
771 ++data.summaryStatAll()[kTotalBremSeeds][K];
772 tracks = m_tracksfinder->getTracksWithBrem(data.combinatorialData(), *Tp, Sp, Gp, DE, data.clusterTrack(), true,ctx);
773 } else {
774 tracks = m_tracksfinder->getTracks (data.combinatorialData(), *Tp, Sp, Gp, DE, data.clusterTrack(),ctx);
775 }
776
778 std::array<bool,SiCombinatorialTrackFinderData_xk::kNCombStats> inf{0,0,0,0,0,0};
779 m_tracksfinder->fillStatistic(data.combinatorialData(),inf);
780 for (size_t p =0; p<inf.size(); ++p){
781 if(inf[p]) ++data.summaryStatAll()[m_indexToEnum[p]][K];
782 }
783
786 if (m_seedsfilter) {
787 std::list<Trk::Track*>::iterator t = tracks.begin();
788 while (t!=tracks.end()) {
790 if (!isNewTrack(data, *t)) {
791 delete (*t);
792 tracks.erase(t++);
793 } else {
794 if((*t)->info().trackProperties(Trk::TrackInfo::BremFit)) ++data.summaryStatAll()[kBremTracks][K];
795 clusterTrackMap(data, *t++);
796 }
797 }
798 }
799 data.findtracks() += tracks.size();
800
801 if(!tracks.empty()) {
802 ++data.summaryStatAll()[kSeedsWithTracks][K];
803 ++data.summaryStatUsedInTrack()[kSeedsWithTracksEta][K][r];
804 data.summaryStatAll()[kOutputTracks][K] += tracks.size();
805 data.summaryStatAll()[kExtraTracks][K] += (tracks.size()-1);
806 }
807
808 return tracks;
809
810}
#define ATH_MSG_ERROR(x)
bool newSeed(SiTrackMakerEventData_xk &data, const std::vector< const Trk::SpacePoint * > &Sp) const
std::unique_ptr< Trk::TrackParameters > getAtaPlane(MagField::AtlasFieldCache &fieldCache, SiTrackMakerEventData_xk &data, bool sss, const std::vector< const Trk::SpacePoint * > &SP, const EventContext &ctx) const
BooleanProperty m_useTrigInDetRoadPredictorTool
std::vector< statAllTypes > m_indexToEnum
BooleanProperty m_useBremModel
static void clusterTrackMap(SiTrackMakerEventData_xk &data, Trk::Track *Tr)
BooleanProperty m_cosmicTrack
BooleanProperty m_useTrigTrackFollowingTool
static int rapidity(const std::vector< const Trk::SpacePoint * > &Sp)
bool isCaloCompatible(SiTrackMakerEventData_xk &data) const
static int kindSeed(const std::vector< const Trk::SpacePoint * > &Sp)
BooleanProperty m_useCaloSeeds
ToolHandle< ITrigInDetRoadPredictorTool > m_trigInDetRoadPredictorTool
ToolHandle< ITrigInDetTrackFollowingTool > m_trigInDetTrackFollowingTool
static bool isNewTrack(SiTrackMakerEventData_xk &data, Trk::Track *Tr)
static void detectorElementsSelection(SiTrackMakerEventData_xk &data, std::vector< const InDetDD::SiDetectorElement * > &DE)
@ BremFit
A brem fit was performed on this track.
TStreamerInfo * inf
@ oppositeMomentum
@ alongMomentum

◆ getTracks() [2/2]

std::list< Trk::Track * > InDet::SiTrackMaker_xk::getTracks ( const EventContext & ctx,
SiTrackMakerEventData_xk & data,
const Trk::TrackParameters & Tp,
const std::vector< Amg::Vector3D > & Gp ) const
overridevirtual

Definition at line 816 of file SiTrackMaker_xk.cxx.

818{
819 ++data.inputseeds();
820 std::list<Trk::Track*> tracks;
821 if (!data.pix() && !data.sct()) return tracks;
822
823 ++data.goodseeds();
824
825 // Get AtlasFieldCache
826 MagField::AtlasFieldCache fieldCache;
827
828 SG::ReadCondHandle<AtlasFieldCacheCondObj> readHandle{m_fieldCondObjInputKey, ctx};
829 const AtlasFieldCacheCondObj* fieldCondObj{*readHandle};
830 if (fieldCondObj == nullptr) {
831 ATH_MSG_ERROR("InDet::SiTrackMaker_xk::getTracks: Failed to retrieve AtlasFieldCacheCondObj with key " << m_fieldCondObjInputKey.key());
832 return tracks;
833 }
834 fieldCondObj->getInitializedCache (fieldCache);
835
836 // Get detector elements road
837 //
838 std::vector<const InDetDD::SiDetectorElement*> DE;
839 if (!m_cosmicTrack) m_roadmaker->detElementsRoad(ctx, fieldCache, Tp,Trk::alongMomentum, DE, data.roadMakerData());
840 else m_roadmaker->detElementsRoad(ctx, fieldCache, Tp,Trk::oppositeMomentum,DE, data.roadMakerData());
841
842 if (!data.pix() || !data.sct()) detectorElementsSelection(data, DE);
843
844 if (static_cast<int>(DE.size()) < m_nclusmin) return tracks;
845
846 // Find possible list of tracks with trigger track parameters or global positions
847 //
848 std::vector<const Trk::SpacePoint*> Sp;
849
850 if (!m_useBremModel) {
851 tracks = m_tracksfinder->getTracks (data.combinatorialData(), Tp, Sp, Gp, DE, data.clusterTrack(),ctx);
852 } else if (!m_useCaloSeeds) {
853 tracks = m_tracksfinder->getTracksWithBrem(data.combinatorialData(), Tp, Sp, Gp, DE, data.clusterTrack(), false,ctx);
854 } else if (isCaloCompatible(data)) {
855 tracks = m_tracksfinder->getTracksWithBrem(data.combinatorialData(), Tp, Sp, Gp, DE, data.clusterTrack(), true,ctx);
856 } else {
857 tracks = m_tracksfinder->getTracks (data.combinatorialData(), Tp, Sp, Gp, DE, data.clusterTrack(),ctx);
858 }
859
860 if (m_seedsfilter) {
861 std::list<Trk::Track*>::iterator t = tracks.begin();
862 while (t!=tracks.end()) {
863 if (!isNewTrack(data, *t)) {
864 delete (*t);
865 tracks.erase(t++);
866 } else {
867 clusterTrackMap(data, *t++);
868 }
869 }
870 }
871 data.findtracks() += tracks.size();
872 return tracks;
873}

◆ globalDirections()

void InDet::SiTrackMaker_xk::globalDirections ( const double * p0,
const double * p1,
const double * p2,
double * d0,
double * d1,
double * d2 )
staticprivate

Here, we derive the local track direction in the space-points as the tangents to the estimated trajectory (assuming a circle in x-y and straight line in r-z)

transform transverse coordinates relative to the first SP

Now apply the same transform used in the seed maker: x,y --> u:=x/(x²+y²); v:=y/(x²+y²); in a frame centered around first SP

set x axis as direction unit vector in transverse plane, from first to second SP

reasoning: u = sqrt(d01)/d01, and v01 = 0

decompose the SP3-SP1 separation into components parallel...

and orthogonal to our new local x axis

squared distance third to first SP the transverse plane

u,v coordinates of third SP in the new frame

Now obtain A,B parameters from circle parameterisation in the new frame: V = (-x0/y0) x U + 1/(2y0) =: A x U + B. A ~= delta V / delta U, B = V - A x U, add a factor 2 for utility

v01 is 0, as the direction from SP 1 to 2 is the u-axis

Now we can resolve the related geometry. Note that A0, in the original frame, is related to the angle Phi1 between the tangent the central SP and the tendon between the first and second SP, alpha, as tan(alpha) = A0.

dz/dr along the seed = 1 / tan(theta)

1/sqrt(1+1/tan²(theta)) -> sin(theta)

cos theta (= sin(theta)/tan(theta))

multiply with cos(alpha) via A0

(a,b) parameterises the direction corresponding to the tendon between the first two SP. Now, go from there to the local tangents by removing or adding the angle 'alpha' from before, derived from A0. direction at first SP - rotated by - 1 x alpha w.r.t the tendon defining a and b The structure below comes from the formulae for sin / cos of a sum of two angles with (a,b) = r(cos phi / - sin phi)

px0: a sin theta cos alpha - b sin theta sin alpha

py0: b sin theta cos alpha + a sin theta sin alpha

pz0: cos theta

direction at second SP - rotated by + 1 x alpha w.r.t the tendon

px1: a sin theta cos alpha + b sin theta sin alpha

py1: b sin theta cos alpha - a sin theta sin alpha

pz1: cos theta

direction at third SP

px2: a sin theta cos alpha * C2 - b sin theta cos alpha S2

py2: b sin theta cos alpha * C2 + a sin theta cos alpha S2

Definition at line 1423 of file SiTrackMaker_xk.cxx.

1425{
1427 double x01 = p1[0]-p0[0] ;
1428 double y01 = p1[1]-p0[1] ;
1429 double x02 = p2[0]-p0[0] ;
1430 double y02 = p2[1]-p0[1] ;
1431
1435
1437 double d01 = x01*x01+y01*y01 ;
1438 double x1 = sqrt(d01) ;
1439 double u01 = 1./x1 ;
1440 double a = x01*u01 ;
1441 double b = y01*u01 ;
1443 double x2 = a*x02+b*y02 ;
1445 double y2 = a*y02-b*x02 ;
1446
1448 double d02 = x2*x2+y2*y2 ;
1450 double u02 = x2/d02 ;
1451 double v02 = y2/d02 ;
1456 double A0 = v02 /(u02-u01) ;
1457 double B0 = 2.*(v02-A0*u02) ;
1458
1464
1465 double C2 = (1.-B0*y2) ;
1466 double S2 = (A0+B0*x2) ;
1467
1468 double T = (p2[2]-p0[2])/sqrt(d02);
1469 double sinTheta = 1./sqrt(1.+T*T) ;
1470 double cosTheta = sinTheta*T ;
1471 double sinThetaCosAlpha = sinTheta / sqrt(1.+A0*A0) ;
1472 double Sa = sinThetaCosAlpha*a ;
1473 double Sb = sinThetaCosAlpha*b ;
1474
1479 d0[0] = Sa -Sb*A0;
1480 d0[1]= Sb +Sa*A0;
1481 d0[2]=cosTheta;
1482
1484 d1[0] = Sa +Sb*A0;
1485 d1[1]= Sb -Sa*A0;
1486 d1[2]=cosTheta;
1487
1489 d2[0] = Sa*C2-Sb*S2;
1490 d2[1]= Sb*C2+Sa*S2;
1491 d2[2]=cosTheta;
1492}
static const int B0
Definition AtlasPID.h:122
struct TBPatternUnitContext S2

◆ globalPosition()

bool InDet::SiTrackMaker_xk::globalPosition ( const Trk::SpacePoint & sp,
const double * dir,
double * p ) const
private

This is a refinement of the global position for strip space-points.

It uses the direction estimate to fix the hit position along the strip axis (non-sensitive direction).

pick up the two components of the space point

get the two ends of the strip in the global frame for each of the two clusters

get the "strip axis" in the global frame for each of the two clusters these define two planes in which the strips are sensitive

get the connection vector between the bottom ends of the two strips

divide max distance by the first strip length

Get the cross products of the direction vector and the strip axes. This gives us a normal representation of the planes containing both.

The two planes are still only defined up to a shift along the normal. Now, we fix this degree of freedom by requiring that the u-plane (perpendicular to momentum and second strip direction) should actually contain the second strip itself. To do so, we virtually 'hold' the u-plane at its intersection plane with the first strip, starting with it touching the very bottom of the first strip, and then move the intersection point along the strip. We can thus parameterise the shift by the distance we need to walk this way

Find the maximum distance we can go until we have walked past the full length of the first strip. Equivalent to the projection of the full strip length on the u-plane normal.

if no such component, bail out - no solution exists

The distance we need to walk to contain the second strip in the u-plane is is obtained by projecting the separation vector between the strip starting points onto the normal. (remember, we virtually start our shift with the u-plane intersecting the first strip at its starting point). We normalise the shifts to the maximum allowed shift.

this is playing the same game, but for shifting the v-plane w.r.t the version that intersects with the bottom of the second strip until we contain the first strip in it. du has the same length as dv (symmetry of the problem)

check that the result is valid. We want the final fixed planes to still intersect the respective other strip within some tolerance. Keep in mind that s=0 --> intersect start of strip, s == 1 --> intersect end of strip. We apply a tolerance beyond each end

now, update the position estimate as the intersection of the updated u-plane with the first strip.

Definition at line 1308 of file SiTrackMaker_xk.cxx.

1310{
1312 const Trk::PrepRawData* c0 = sp.clusterList().first;
1313 const Trk::PrepRawData* c1 = sp.clusterList().second;
1314
1315 const InDetDD::SiDetectorElement* de0 = static_cast<const InDet::SiCluster*>(c0)->detectorElement();
1316 const InDetDD::SiDetectorElement* de1 = static_cast<const InDet::SiCluster*>(c1)->detectorElement();
1317
1319 Amg::Vector2D localPos = c0->localPosition();
1320 std::pair<Amg::Vector3D,Amg::Vector3D> e0
1321 (de0->endsOfStrip(InDetDD::SiLocalPosition(localPos.y(),localPos.x(),0.)));
1322
1323 localPos = c1->localPosition();
1324 std::pair<Amg::Vector3D,Amg::Vector3D> e1
1325 (de1->endsOfStrip(InDetDD::SiLocalPosition(localPos.y(),localPos.x(),0.)));
1326
1329 double a0[3] = {e0.second.x()-e0.first.x(), e0.second.y()-e0.first.y(), e0.second.z()-e0.first.z()};
1330 double a1[3] = {e1.second.x()-e1.first.x(), e1.second.y()-e1.first.y(), e1.second.z()-e1.first.z()};
1332 double dr[3] = {e1.first .x()-e0.first.x(), e1.first .y()-e0.first.y(), e1.first .z()-e0.first.z()};
1333
1335 double d0 = m_distmax/sqrt(a0[0]*a0[0]+a0[1]*a0[1]+a0[2]*a0[2]);
1336
1339 double u[3] = {a1[1]*dir[2]-a1[2]*dir[1],a1[2]*dir[0]-a1[0]*dir[2],a1[0]*dir[1]-a1[1]*dir[0]};
1340 double v[3] = {a0[1]*dir[2]-a0[2]*dir[1],a0[2]*dir[0]-a0[0]*dir[2],a0[0]*dir[1]-a0[1]*dir[0]};
1341
1348
1351 double du = a0[0]*u[0]+a0[1]*u[1]+a0[2]*u[2];
1352
1354 if (du==0. ) return false;
1355
1362 double s0 = (dr[0]*u[0]+dr[1]*u[1]+dr[2]*u[2])/du;
1367 double s1 = (dr[0]*v[0]+dr[1]*v[1]+dr[2]*v[2])/du;
1368
1373 if (s0 < -d0 || s0 > 1.+d0 || s1 < -d0 || s1 > 1.+d0) return false;
1374
1377 p[0] = e0.first.x()+s0*a0[0];
1378 p[1] = e0.first.y()+s0*a0[1];
1379 p[2] = e0.first.z()+s0*a0[2];
1380
1381 return true;
1382}
static Double_t sp
static Double_t s0
std::pair< Amg::Vector3D, Amg::Vector3D > endsOfStrip(const Amg::Vector2D &position) const
Special method for SCT to retrieve the two ends of a "strip" Returned coordinates are in global frame...
double a0
Definition globals.cxx:27
Eigen::Matrix< double, 2, 1 > Vector2D
@ u
Enums for curvilinear frames.
Definition ParamDefs.h:77
double e0(const xAOD::CaloCluster &cluster)
return the uncorrected cluster energy in pre-sampler
double e1(const xAOD::CaloCluster &cluster)
return the uncorrected cluster energy in 1st sampling

◆ globalPositions()

bool InDet::SiTrackMaker_xk::globalPositions ( const Trk::SpacePoint & s0,
const Trk::SpacePoint & s1,
const Trk::SpacePoint & s2,
double * p0,
double * p1,
double * p2 ) const
private

first, fill the arrays with the global positions of the 3 points

for PPP seeds, we are done

for SSS, we need some extra work

try to refine the position estimate for strip SP using the cluster information in combination with the direction estimate

Definition at line 1266 of file SiTrackMaker_xk.cxx.

1269{
1270
1272 p0[0] = s0.globalPosition().x();
1273 p0[1] = s0.globalPosition().y();
1274 p0[2] = s0.globalPosition().z();
1275
1276 p1[0] = s1.globalPosition().x();
1277 p1[1] = s1.globalPosition().y();
1278 p1[2] = s1.globalPosition().z();
1279
1280 p2[0] = s2.globalPosition().x();
1281 p2[1] = s2.globalPosition().y();
1282 p2[2] = s2.globalPosition().z();
1283
1285 if (!s0.clusterList().second && !s1.clusterList().second && !s2.clusterList().second) return true;
1286
1288 double dir0[3],dir1[3],dir2[3];
1289
1290 globalDirections(p0,p1,p2,dir0,dir1,dir2);
1291
1294 if (s0.clusterList().second && !globalPosition(s0,dir0,p0)) return false;
1295 if (s1.clusterList().second && !globalPosition(s1,dir1,p1)) return false;
1296 if (s2.clusterList().second && !globalPosition(s2,dir2,p2)) return false;
1297
1298 return true;
1299}
static void globalDirections(const double *p0, const double *p1, const double *p2, double *d0, double *d1, double *d2)
Here, we derive the local track direction in the space-points as the tangents to the estimated trajec...
bool globalPosition(const Trk::SpacePoint &sp, const double *dir, double *p) const
This is a refinement of the global position for strip space-points.

◆ initialize()

StatusCode InDet::SiTrackMaker_xk::initialize ( )
overridevirtual

Get beam geometry

read the config string for the field mode

this one is the default

Get detector elements road maker tool

Get combinatorial track finder tool

Get trigger track following tool

Get trigger road predictor tool

Get seed to track conversion tool This is used if we want to write out the seeds for performance studies

flag for HI running

this is on by default in offline tracking

useSSSfilter is on by default in offline, m_useHClusSeed is off by default.

pt cut can never be below 20 MeV

initialize counters

Definition at line 38 of file SiTrackMaker_xk.cxx.

39{
42 if (not m_beamSpotKey.empty()) {
43 ATH_CHECK( m_beamSpotKey.initialize() );
44 }
45
47 if (m_fieldmode == "NoField") m_fieldModeEnum = Trk::NoField;
48 else if (m_fieldmode == "MapSolenoid") m_fieldModeEnum = Trk::FastField;
50
53 if ( m_roadmaker.retrieve().isFailure() ) {
54 ATH_MSG_FATAL( "Failed to retrieve tool " << m_roadmaker );
55 return StatusCode::FAILURE;
56 }
57 ATH_MSG_DEBUG( "Retrieved tool " << m_roadmaker );
58
61 if ( m_tracksfinder.retrieve().isFailure() ) {
62 ATH_MSG_FATAL( "Failed to retrieve tool " << m_tracksfinder );
63 return StatusCode::FAILURE;
64 }
65 ATH_MSG_DEBUG( "Retrieved tool " << m_tracksfinder );
66
70 if ( m_trigInDetTrackFollowingTool.retrieve().isFailure() ) {
71 ATH_MSG_FATAL( "Failed to retrieve tool " << m_trigInDetTrackFollowingTool );
72 return StatusCode::FAILURE;
73 }
74 ATH_MSG_DEBUG( "Retrieved tool " << m_trigInDetTrackFollowingTool );
75 } else {
77 }
78
82 if ( m_trigInDetRoadPredictorTool.retrieve().isFailure() ) {
83 ATH_MSG_FATAL( "Failed to retrieve tool " << m_trigInDetRoadPredictorTool );
84 return StatusCode::FAILURE;
85 }
86 ATH_MSG_DEBUG( "Retrieved tool " << m_trigInDetRoadPredictorTool );
87 } else {
89 }
90
95 if (m_seedtrack.retrieve().isFailure()) {
96 ATH_MSG_FATAL( "Failed to retrieve tool " << m_seedtrack );
97 return StatusCode::FAILURE;
98 }
99 ATH_MSG_DEBUG( "Retrieved tool " << m_seedtrack );
100 } else {
101 m_seedtrack.disable();
102 }
103
105 m_heavyion = false;
106
107
108 // TrackpatternRecoInfo preparation
109 //
110 if (m_patternName == "SiSpacePointsSeedMaker_Cosmic" ) {
112 } else if (m_patternName == "SiSpacePointsSeedMaker_HeavyIon" ) {
114 m_heavyion = true;
115 } else if (m_patternName == "SiSpacePointsSeedMaker_LowMomentum") {
117 } else if (m_patternName == "SiSpacePointsSeedMaker_BeamGas" ) {
119 } else if (m_patternName == "SiSpacePointsSeedMaker_ForwardTracks" ) {
121 } else if (m_patternName == "SiSpacePointsSeedMaker_LargeD0" ) {
123 } else if (m_patternName == "SiSpacePointsSeedMaker_ITkConversionTracks") {
125 } else {
126 m_trackinfo.setPatternRecognitionInfo(Trk::TrackInfo::SiSPSeededFinder );
127 }
128
133
135 if (m_pTmin < 20.) m_pTmin = 20.;
136
137 if (m_etabins.size()>0) {
138 if (m_ptbins.size() > (m_etabins.size()-1)) {
139 ATH_MSG_ERROR( "No. of cut values bigger than eta bins");
140 return StatusCode::FAILURE;
141 }
142
143 if (m_ptbins.size() < (m_etabins.size()-1)) {
144 ATH_MSG_DEBUG( "No. of cut values smaller than eta bins. Extending size..." );
145 m_ptbins.value().resize(m_etabins.size()-1, m_ptbins.value().back());
146 }
147
148 for (auto& pt : m_ptbins) {
149 if (pt < 20.) pt = 20.;
150 }
151 }
152
153
155 resetCounter(m_totalInputSeeds);
156 resetCounter(m_totalNoTrackPar);
157 resetCounter(m_totalUsedSeeds);
158 resetCounter(m_outputTracks);
159 resetCounter(m_totalBremSeeds);
160 resetCounter(m_bremTracks);
161 resetCounter(m_twoClusters);
162 resetCounter(m_wrongRoad);
163 resetCounter(m_wrongInit);
164 resetCounter(m_noTrack);
165 resetCounter(m_notNewTrack);
166 resetCounter(m_bremAttempt);
167 resetCounter(m_seedsWithTrack);
168 resetCounter(m_deSize);
170 m_seedsWithTracksEta.resize(SiCombinatorialTrackFinderData_xk::kNSeedTypes);
171 std::fill(m_usedSeedsEta.begin(), m_usedSeedsEta.end(),
173 std::fill(m_seedsWithTracksEta.begin(), m_seedsWithTracksEta.end(),
175
177 ATH_CHECK( m_fieldCondObjInputKey.initialize());
179 return StatusCode::SUCCESS;
180}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_FATAL(x)
#define ATH_MSG_DEBUG(x)
DoubleArrayProperty m_etabins
SG::ReadHandleKey< ROIPhiRZContainer > m_caloCluster
void resetCounter(std::array< std::array< T, M >, N > &a) const
helper for working with the stat arrays
StringProperty m_patternName
SG::ReadHandleKey< ROIPhiRZContainer > m_caloHad
BooleanProperty m_useSSSfilter
BooleanProperty m_useHClusSeed
SG::ReadCondHandleKey< InDet::BeamSpotData > m_beamSpotKey
DoubleArrayProperty m_ptbins
@ SiSpacePointsSeedMaker_Cosmic
Entries allowing to distinguish different seed makers.
@ SiSpacePointsSeedMaker_LargeD0
Large d0 for displaced vertex searches.
@ SiSpacePointsSeedMaker_ForwardTracks
Entries allowing to distinguish different seed makers.
@ SiSpacePointsSeedMaker_ITkConversionTracks
ITkConversion Track flag.
@ SiSPSeededFinder
Tracks from SiSPSeedFinder.
@ FastField
call the fast field access method of the FieldSvc
@ FullField
Field is set to be realistic, but within a given Volume.

◆ isCaloCompatible()

bool InDet::SiTrackMaker_xk::isCaloCompatible ( SiTrackMakerEventData_xk & data) const
private

Definition at line 1388 of file SiTrackMaker_xk.cxx.

1389{
1390 if (!data.caloClusterROIEM()) return false;
1391
1392
1393 double F = data.par()[2] ;
1394 double E = -log(tan(.5*data.par()[3])) ;
1395 double R = sqrt(data.par()[6]*data.par()[6]+data.par()[7]*data.par()[7]);
1396 double Z = data.par()[8] ;
1397 return data.caloClusterROIEM()->hasMatchingROI(F, E, R, Z, m_phiWidth, m_etaWidth);
1398}
#define F(x, y, z)
Definition MD5.cxx:112
double R(const INavigable4Momentum *p1, const double v_eta, const double v_phi)

◆ isHadCaloCompatible()

bool InDet::SiTrackMaker_xk::isHadCaloCompatible ( SiTrackMakerEventData_xk & data) const
private

Definition at line 1404 of file SiTrackMaker_xk.cxx.

1405{
1406 if (!data.caloClusterROIHad()) return false;
1407
1408 double F = data.par()[2] ;
1409 double E = -log(tan(.5*data.par()[3])) ;
1410 double R = sqrt(data.par()[6]*data.par()[6]+data.par()[7]*data.par()[7]);
1411 double Z = data.par()[8] ;
1412
1413 return data.caloClusterROIHad()->hasMatchingROI(F, E, R, Z, m_phiWidth, m_etaWidth);
1414}

◆ isNewTrack()

bool InDet::SiTrackMaker_xk::isNewTrack ( SiTrackMakerEventData_xk & data,
Trk::Track * Tr )
staticprivate

Definition at line 1224 of file SiTrackMaker_xk.cxx.

1225{
1226 const Trk::PrepRawData* prd [100];
1227 std::multimap<const Trk::PrepRawData*,const Trk::Track*>::const_iterator
1228 ti,t[100],te = data.clusterTrack().end();
1229
1230 int n = 0;
1231
1233 m = Tr->measurementsOnTrack()->begin(),
1234 me = Tr->measurementsOnTrack()->end ();
1235
1236 for (; m!=me; ++m) {
1237
1238 const Trk::PrepRawData* pr = static_cast<const Trk::RIO_OnTrack*>(*m)->prepRawData();
1239 if (pr) {
1240 prd[n] =pr;
1241 t [n] = data.clusterTrack().find(prd[n]);
1242 if (t[n]==te) return true;
1243 ++n;
1244 }
1245 }
1246
1247 if (!n) return true;
1248 int nclt = n;
1249
1250 for (int i=0; i!=n; ++i) {
1251 int nclmax = 0;
1252 for (ti=t[i]; ti!=te; ++ti) {
1253 if ( (*ti).first != prd[i] ) break;
1254 int ncl = (*ti).second->measurementsOnTrack()->size();
1255 if (ncl > nclmax) nclmax = ncl;
1256 }
1257 if (nclt > nclmax) return true;
1258 }
1259 return false;
1260}

◆ kindSeed()

int InDet::SiTrackMaker_xk::kindSeed ( const std::vector< const Trk::SpacePoint * > & Sp)
staticprivate

Definition at line 1181 of file SiTrackMaker_xk.cxx.

1182{
1183
1184 if(Sp.size()!=3) return 0;//correct handling of Pixel-only ITk tracklets
1185
1186 std::vector<const Trk::SpacePoint*>::const_iterator s=Sp.begin(),se=Sp.end();
1187
1188 int n = 0;
1189 for(; s!=se; ++s) {if((*s)->clusterList().second) ++n;}
1190 return n;
1191}

◆ newEvent()

void InDet::SiTrackMaker_xk::newEvent ( const EventContext & ctx,
SiTrackMakerEventData_xk & data,
bool PIX,
bool SCT ) const
overridevirtual

initialize beam position

propagate pixel / strip usage to the event data object

build a holder for the configured track quality cuts

Setup New event for track finder tool

Erase cluster to track association m_seedsfilter is true in the vast majority of all applications, so this is usually done

Erase statistic information

retrieve calo seeds for brem fit

retrieve hadronic seeds for SSS seed filter

if we want to write out the seeds, also call newEvent for the seed-to-track converter

Definition at line 498 of file SiTrackMaker_xk.cxx.

499{
500
502 data.xybeam()[0] = 0.;
503 data.xybeam()[1] = 0.;
504 if (not m_beamSpotKey.empty()) {
505 SG::ReadCondHandle<InDet::BeamSpotData> beamSpotHandle { m_beamSpotKey, ctx };
506 if (beamSpotHandle.isValid()) {
507 data.xybeam()[0] = beamSpotHandle->beamPos()[0];
508 data.xybeam()[1] = beamSpotHandle->beamPos()[1];
509 }
510 }
512 data.pix() = PIX and m_usePix;
513 data.sct() = SCT and m_useSct;
514
516 InDet::TrackQualityCuts trackquality = setTrackQualityCuts(false);
517
519 m_tracksfinder->newEvent(ctx, data.combinatorialData(), m_trackinfo, trackquality);
520
524 if (m_seedsfilter) data.clusterTrack().clear();
525
528 data.inputseeds() = 0;
529 data.goodseeds() = 0;
530 data.findtracks() = 0;
531 resetCounter(data.summaryStatAll());
532 resetCounter(data.summaryStatUsedInTrack());
533
536 SG::ReadHandle<ROIPhiRZContainer> calo_rois(m_caloCluster, ctx);
537 if (!calo_rois.isValid()) {
538 ATH_MSG_FATAL("Failed to get EM Calo cluster collection " << m_caloCluster );
539 }
540 data.setCaloClusterROIEM(*calo_rois);
541 }
542
545 SG::ReadHandle<ROIPhiRZContainer> calo_rois(m_caloHad, ctx);
546 if (!calo_rois.isValid()) {
547 ATH_MSG_FATAL("Failed to get Had Calo cluster collection " << m_caloHad );
548 }
549 data.setCaloClusterROIHad(*calo_rois);
550 }
552 if (m_seedsegmentsWrite) m_seedtrack->newEvent(data.conversionData(), m_trackinfo, m_patternName);
553}
@ SCT
Definition RegSelEnums.h:25
InDet::TrackQualityCuts setTrackQualityCuts(bool simpleTrack) const

◆ newSeed()

bool InDet::SiTrackMaker_xk::newSeed ( SiTrackMakerEventData_xk & data,
const std::vector< const Trk::SpacePoint * > & Sp ) const
private

counter for clusters on track

start by checking the first cluster - always needed

lookup if the cluster is already used by another track

add tracks consuming our seed space-points to the trackseed list.

increment cluster counter

now, prepare to check also the second cluster on any strip seed

if we don't have one, nothing to do

otherwise, same game as before. Note that a track consuming both clusters on a strip hit is counted twice into the map

incremenent counter again

check if at least on cluster is not already used by any track. This works since the multiset allows adding the same track multiple times If this is the case, we accept the seed.

in the case of HI reco, we accept any 3-cluster (PPP) seed.

Now we look for the track consuming the largest number of clusters

This is done by looping over all tracks using any of our clusters, and counting the appearance of each track in the multiset. If one single track contains all of the clusters (--> is included n times), we reject this seed.

loop over the list of tracks

if this is a new track, reset the counter

otherwise increment the counter. If this track has all clusters from the seed on it, reject the event

If we have no single track 'eating' all of our clusters, we accept the seed

Definition at line 1112 of file SiTrackMaker_xk.cxx.

1113{
1114 std::multiset<const Trk::Track*> trackseed;
1115 std::multimap<const Trk::PrepRawData*,const Trk::Track*>::const_iterator iter_clusterOnTrack,iter_clusterOnTrackEnd = data.clusterTrack().end();
1116
1118 size_t n = 0;
1119
1120 for (const Trk::SpacePoint* spacePoint : Sp) {
1121
1123 const Trk::PrepRawData* prd = spacePoint->clusterList().first;
1124
1126 for (iter_clusterOnTrack = data.clusterTrack().find(prd); iter_clusterOnTrack!=iter_clusterOnTrackEnd; ++iter_clusterOnTrack) {
1127 if ((*iter_clusterOnTrack).first!=prd) break;
1129 trackseed.insert((*iter_clusterOnTrack).second);
1130 }
1132 ++n;
1134 prd = spacePoint->clusterList().second;
1136 if (!prd) continue;
1138 for (iter_clusterOnTrack = data.clusterTrack().find(prd); iter_clusterOnTrack!=iter_clusterOnTrackEnd; ++iter_clusterOnTrack) {
1139 if ((*iter_clusterOnTrack).first!=prd) break;
1140 trackseed.insert((*iter_clusterOnTrack).second);
1141 }
1143 ++n;
1144 }
1148 if(trackseed.size() < n) return true;
1150 if( m_heavyion && n==3 ) return true;
1151
1153
1158 const Trk::Track* currentTrack {nullptr};
1159 size_t clustersOnCurrent = 1;
1161 for(const Trk::Track* track : trackseed) {
1163 if(track != currentTrack) {
1164 currentTrack = track;
1165 clustersOnCurrent = 1;
1166 continue ;
1167 }
1170 if(++clustersOnCurrent == n) return false;
1171 }
1173 return clustersOnCurrent!=n;
1174}

◆ newTrigEvent()

void InDet::SiTrackMaker_xk::newTrigEvent ( const EventContext & ctx,
SiTrackMakerEventData_xk & data,
bool PIX,
bool SCT ) const
overridevirtual

Definition at line 559 of file SiTrackMaker_xk.cxx.

560{
561 data.pix() = PIX && m_usePix;
562 data.sct() = SCT && m_useSct;
563 bool simpleTrack = true;
564
565 InDet::TrackQualityCuts trackquality = setTrackQualityCuts(simpleTrack);
566
567 // New event for track finder tool
568 //
569 m_tracksfinder->newEvent(ctx, data.combinatorialData(), m_trackinfo, trackquality);
570
571 // Erase cluster to track association
572 //
573 if (m_seedsfilter) data.clusterTrack().clear();
574
575 // Erase statistic information
576 //
577 data.inputseeds() = 0;
578 data.goodseeds() = 0;
579 data.findtracks() = 0;
582 data.summaryStatAll()[i][k] = 0.;
583 }
584 }
586 for (int k = 0; k != SiCombinatorialTrackFinderData_xk::kNSeedTypes; ++k) {
588 data.summaryStatUsedInTrack()[i][k][r] = 0.;
589 }
590 }
591 }
592
593 // Add possibility to write seed segment information
594 if (m_seedsegmentsWrite) m_seedtrack->newEvent(data.conversionData(), m_trackinfo, m_patternName);
595
596}

◆ operator=()

SiTrackMaker_xk & InDet::SiTrackMaker_xk::operator= ( const SiTrackMaker_xk & )
privatedelete

◆ pTmin()

double InDet::SiTrackMaker_xk::pTmin ( double eta) const
private

Definition at line 1496 of file SiTrackMaker_xk.cxx.

1497{
1498 if (m_ptbins.size() == 0) return m_pTmin;
1499 double aeta = std::abs(eta);
1500 for(int n = int(m_ptbins.size()-1); n>=0; --n) {
1501 if(aeta > m_etabins[n]) return m_ptbins[n];
1502 }
1503 return m_pTmin;
1504}

◆ rapidity()

int InDet::SiTrackMaker_xk::rapidity ( const std::vector< const Trk::SpacePoint * > & Sp)
staticprivate

Definition at line 1193 of file SiTrackMaker_xk.cxx.

1194{
1195 if(Sp.size() < 2) return 0;
1196
1197 Amg::Vector3D delta_sp = Sp[0]->globalPosition() - Sp[1]->globalPosition();
1198 float eta = 2*std::abs(delta_sp.eta());
1199 int n = int(eta);
1200 if(n > 7) n = 7;
1201 return n;
1202}
Eigen::Matrix< double, 3, 1 > Vector3D

◆ resetCounter() [1/2]

template<typename T, size_t N, size_t M>
void InDet::SiTrackMaker_xk::resetCounter ( std::array< std::array< T, M >, N > & a) const
inlineprivate

helper for working with the stat arrays

Definition at line 265 of file SiTrackMaker_xk.h.

265 {
266 for (auto & subarr : a) resetCounter(subarr);
267 }

◆ resetCounter() [2/2]

template<typename T, size_t N>
void InDet::SiTrackMaker_xk::resetCounter ( std::array< T, N > & a) const
inlineprivate

Definition at line 268 of file SiTrackMaker_xk.h.

268 {
269 std::fill(a.begin(),a.end(),0);
270 }

◆ setTrackQualityCuts()

InDet::TrackQualityCuts InDet::SiTrackMaker_xk::setTrackQualityCuts ( bool simpleTrack) const
private

Definition at line 1050 of file SiTrackMaker_xk.cxx.

1051{
1052 InDet::TrackQualityCuts trackquality;
1053 // Integer cuts
1054 //
1055 trackquality.setIntCut("MinNumberOfClusters" ,m_nclusmin );
1056 trackquality.setIntCut("MinNumberOfWClusters",m_nwclusmin );
1057 trackquality.setIntCut("MaxNumberOfHoles" ,m_nholesmax );
1058 trackquality.setIntCut("MaxHolesGap" ,m_dholesmax );
1059
1060 if (m_useassoTool) trackquality.setIntCut("UseAssociationTool",1);
1061 else trackquality.setIntCut("UseAssociationTool",0);
1062 if (m_cosmicTrack) trackquality.setIntCut("CosmicTrack" ,1);
1063 else trackquality.setIntCut("CosmicTrack" ,0);
1064 if (simpleTrack) trackquality.setIntCut("SimpleTrack" ,1);
1065 else trackquality.setIntCut("SimpleTrack" ,0);
1066 if (m_multitracks) trackquality.setIntCut("doMultiTracksProd" ,1);
1067 else trackquality.setIntCut("doMultiTracksProd" ,0);
1068
1069 // Double cuts
1070 //
1071 trackquality.setDoubleCut("pTmin" ,m_pTmin );
1072 trackquality.setDoubleCut("pTminBrem" ,m_pTminBrem );
1073 trackquality.setDoubleCut("MaxXi2forCluster" ,m_xi2max );
1074 trackquality.setDoubleCut("MaxXi2forOutlier" ,m_xi2maxNoAdd);
1075 trackquality.setDoubleCut("MaxXi2forSearch" ,m_xi2maxlink );
1076 trackquality.setDoubleCut("MaxXi2MultiTracks" ,m_xi2multitracks);
1077
1078 return trackquality;
1079}
BooleanProperty m_multitracks
DoubleProperty m_xi2multitracks
void setIntCut(const std::string &, int)
void setDoubleCut(const std::string &, double)

Member Data Documentation

◆ ATLAS_THREAD_SAFE [1/17]

std::array<std::atomic<int>,SiCombinatorialTrackFinderData_xk::kNSeedTypes> m_totalInputSeeds InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 186 of file SiTrackMaker_xk.h.

186{};

◆ ATLAS_THREAD_SAFE [2/17]

std::array<std::atomic<int>,SiCombinatorialTrackFinderData_xk::kNSeedTypes> m_totalUsedSeeds InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 187 of file SiTrackMaker_xk.h.

187{};

◆ ATLAS_THREAD_SAFE [3/17]

std::array<std::atomic<int>,SiCombinatorialTrackFinderData_xk::kNSeedTypes> m_totalNoTrackPar InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 188 of file SiTrackMaker_xk.h.

188{};

◆ ATLAS_THREAD_SAFE [4/17]

std::array<std::atomic<int>,SiCombinatorialTrackFinderData_xk::kNSeedTypes> m_totalBremSeeds InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 189 of file SiTrackMaker_xk.h.

189{};

◆ ATLAS_THREAD_SAFE [5/17]

std::array<std::atomic<int>,SiCombinatorialTrackFinderData_xk::kNSeedTypes> m_twoClusters InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 190 of file SiTrackMaker_xk.h.

190{};

◆ ATLAS_THREAD_SAFE [6/17]

std::array<std::atomic<int>,SiCombinatorialTrackFinderData_xk::kNSeedTypes> m_wrongRoad InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 191 of file SiTrackMaker_xk.h.

191{};

◆ ATLAS_THREAD_SAFE [7/17]

std::array<std::atomic<int>,SiCombinatorialTrackFinderData_xk::kNSeedTypes> m_wrongInit InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 192 of file SiTrackMaker_xk.h.

192{};

◆ ATLAS_THREAD_SAFE [8/17]

std::array<std::atomic<int>,SiCombinatorialTrackFinderData_xk::kNSeedTypes> m_noTrack InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 193 of file SiTrackMaker_xk.h.

193{};

◆ ATLAS_THREAD_SAFE [9/17]

std::array<std::atomic<int>,SiCombinatorialTrackFinderData_xk::kNSeedTypes> m_notNewTrack InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 194 of file SiTrackMaker_xk.h.

194{};

◆ ATLAS_THREAD_SAFE [10/17]

std::array<std::atomic<int>,SiCombinatorialTrackFinderData_xk::kNSeedTypes> m_bremAttempt InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 195 of file SiTrackMaker_xk.h.

195{};

◆ ATLAS_THREAD_SAFE [11/17]

std::array<std::atomic<int>,SiCombinatorialTrackFinderData_xk::kNSeedTypes> m_outputTracks InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 196 of file SiTrackMaker_xk.h.

196{};

◆ ATLAS_THREAD_SAFE [12/17]

std::array<std::atomic<int>,SiCombinatorialTrackFinderData_xk::kNSeedTypes> m_extraTracks InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 197 of file SiTrackMaker_xk.h.

197{};

◆ ATLAS_THREAD_SAFE [13/17]

std::array<std::atomic<int>,SiCombinatorialTrackFinderData_xk::kNSeedTypes> m_bremTracks InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 198 of file SiTrackMaker_xk.h.

198{};

◆ ATLAS_THREAD_SAFE [14/17]

std::array<std::atomic<int>,SiCombinatorialTrackFinderData_xk::kNSeedTypes> m_seedsWithTrack InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 199 of file SiTrackMaker_xk.h.

199{};

◆ ATLAS_THREAD_SAFE [15/17]

std::array<std::atomic<double>,SiCombinatorialTrackFinderData_xk::kNSeedTypes> m_deSize InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 200 of file SiTrackMaker_xk.h.

200{};

◆ ATLAS_THREAD_SAFE [16/17]

std::vector<std::vector<double> > m_usedSeedsEta InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE
mutableprivate

Definition at line 202 of file SiTrackMaker_xk.h.

◆ ATLAS_THREAD_SAFE [17/17]

std::vector<std::vector<double> > m_seedsWithTracksEta InDet::SiTrackMaker_xk::ATLAS_THREAD_SAFE
mutableprivate

Definition at line 203 of file SiTrackMaker_xk.h.

◆ m_beamSpotKey

SG::ReadCondHandleKey<InDet::BeamSpotData> InDet::SiTrackMaker_xk::m_beamSpotKey {this, "BeamSpotKey", "BeamSpotData", "SG key for beam spot"}
private

Definition at line 131 of file SiTrackMaker_xk.h.

131{this, "BeamSpotKey", "BeamSpotData", "SG key for beam spot"};

◆ m_caloCluster

SG::ReadHandleKey<ROIPhiRZContainer> InDet::SiTrackMaker_xk::m_caloCluster {this, "EMROIPhiRZContainer", ""}
private

Definition at line 133 of file SiTrackMaker_xk.h.

133{this, "EMROIPhiRZContainer", ""};

◆ m_caloHad

SG::ReadHandleKey<ROIPhiRZContainer> InDet::SiTrackMaker_xk::m_caloHad {this, "HadROIPhiRZContainer", ""}
private

Definition at line 134 of file SiTrackMaker_xk.h.

134{this, "HadROIPhiRZContainer", ""};

◆ m_cosmicTrack

BooleanProperty InDet::SiTrackMaker_xk::m_cosmicTrack {this, "CosmicTrack", false, "Is it cosmic track"}
private

Definition at line 145 of file SiTrackMaker_xk.h.

145{this, "CosmicTrack", false, "Is it cosmic track"};

◆ m_counterMutex

std::mutex InDet::SiTrackMaker_xk::m_counterMutex
mutableprivate

Definition at line 185 of file SiTrackMaker_xk.h.

◆ m_dholesmax

IntegerProperty InDet::SiTrackMaker_xk::m_dholesmax {this, "nHolesGapMax", 2, "Max holes gap"}
private

Definition at line 164 of file SiTrackMaker_xk.h.

164{this, "nHolesGapMax", 2, "Max holes gap"};

◆ m_distmax

DoubleProperty InDet::SiTrackMaker_xk::m_distmax {this, "MaxDistanceForSCTsp", 5.}
private

Definition at line 161 of file SiTrackMaker_xk.h.

161{this, "MaxDistanceForSCTsp", 5.};

◆ m_etabins

DoubleArrayProperty InDet::SiTrackMaker_xk::m_etabins {this, "etaBins", {}, "eta bins"}
private

Definition at line 170 of file SiTrackMaker_xk.h.

170{this, "etaBins", {}, "eta bins"};

◆ m_etaWidth

DoubleProperty InDet::SiTrackMaker_xk::m_etaWidth {this, "etaWidth", 0.3}
private

Definition at line 169 of file SiTrackMaker_xk.h.

169{this, "etaWidth", 0.3};

◆ m_fieldCondObjInputKey

SG::ReadCondHandleKey<AtlasFieldCacheCondObj> InDet::SiTrackMaker_xk::m_fieldCondObjInputKey {this, "AtlasFieldCacheCondObj", "fieldCondObj", "Name of the Magnetic Field conditions object key"}
private

Definition at line 132 of file SiTrackMaker_xk.h.

132{this, "AtlasFieldCacheCondObj", "fieldCondObj", "Name of the Magnetic Field conditions object key"};

◆ m_fieldmode

StringProperty InDet::SiTrackMaker_xk::m_fieldmode {this, "MagneticFieldMode", "MapSolenoid", "Mode of magnetic field"}
private

Definition at line 140 of file SiTrackMaker_xk.h.

140{this, "MagneticFieldMode", "MapSolenoid", "Mode of magnetic field"};

◆ m_fieldModeEnum

Trk::MagneticFieldMode InDet::SiTrackMaker_xk::m_fieldModeEnum {Trk::FullField}
private

Definition at line 178 of file SiTrackMaker_xk.h.

◆ m_heavyion

bool InDet::SiTrackMaker_xk::m_heavyion {false}
private

Definition at line 177 of file SiTrackMaker_xk.h.

177{false}; // Is it heavy ion events

◆ m_indexToEnum

std::vector<statAllTypes> InDet::SiTrackMaker_xk::m_indexToEnum {kTwoClusters,kWrongInit,kWrongRoad,kNoTrack,kNotNewTrack,kBremAttempt}
private

◆ m_ITKGeometry

BooleanProperty InDet::SiTrackMaker_xk::m_ITKGeometry {this, "ITKGeometry", false, "ITK geometry"}
private

Definition at line 151 of file SiTrackMaker_xk.h.

151{this, "ITKGeometry", false, "ITK geometry"};

◆ m_LRTmode

BooleanProperty InDet::SiTrackMaker_xk::m_LRTmode {this, "LRTMode", false}
private

Definition at line 155 of file SiTrackMaker_xk.h.

155{this, "LRTMode", false};

◆ m_multitracks

BooleanProperty InDet::SiTrackMaker_xk::m_multitracks {this, "doMultiTracksProd", false}
private

Definition at line 146 of file SiTrackMaker_xk.h.

146{this, "doMultiTracksProd", false};

◆ m_nclusmin

IntegerProperty InDet::SiTrackMaker_xk::m_nclusmin {this, "nClustersMin", 6, "Min number clusters"}
private

Definition at line 165 of file SiTrackMaker_xk.h.

165{this, "nClustersMin", 6, "Min number clusters"};

◆ m_nholesmax

IntegerProperty InDet::SiTrackMaker_xk::m_nholesmax {this, "nHolesMax", 2, "Max number holes"}
private

Definition at line 163 of file SiTrackMaker_xk.h.

163{this, "nHolesMax", 2, "Max number holes"};

◆ m_nwclusmin

IntegerProperty InDet::SiTrackMaker_xk::m_nwclusmin {this, "nWeightedClustersMin", 6, "Min umber weighted clusters(pix=2 sct=1)"}
private

Definition at line 166 of file SiTrackMaker_xk.h.

166{this, "nWeightedClustersMin", 6, "Min umber weighted clusters(pix=2 sct=1)"};

◆ m_patternName

StringProperty InDet::SiTrackMaker_xk::m_patternName {this, "TrackPatternRecoInfo", "SiSPSeededFinder", "Name of the pattern recognition"}
private

Definition at line 141 of file SiTrackMaker_xk.h.

141{this, "TrackPatternRecoInfo", "SiSPSeededFinder", "Name of the pattern recognition"};

◆ m_phiWidth

DoubleProperty InDet::SiTrackMaker_xk::m_phiWidth {this, "phiWidth", 0.3}
private

Definition at line 168 of file SiTrackMaker_xk.h.

168{this, "phiWidth", 0.3};

◆ m_ptbins

DoubleArrayProperty InDet::SiTrackMaker_xk::m_ptbins {this, "pTBins", {}, "pT bins"}
private

Definition at line 171 of file SiTrackMaker_xk.h.

171{this, "pTBins", {}, "pT bins"};

◆ m_pTmin

DoubleProperty InDet::SiTrackMaker_xk::m_pTmin {this, "pTmin", 500., "min pT"}
private

Definition at line 159 of file SiTrackMaker_xk.h.

159{this, "pTmin", 500., "min pT"};

◆ m_pTminBrem

DoubleProperty InDet::SiTrackMaker_xk::m_pTminBrem {this, "pTminBrem", 1000., "min pT for Brem mode"}
private

Definition at line 160 of file SiTrackMaker_xk.h.

160{this, "pTminBrem", 1000., "min pT for Brem mode"};

◆ m_roadmaker

ToolHandle<InDet::ISiDetElementsRoadMaker> InDet::SiTrackMaker_xk::m_roadmaker {this, "RoadTool", "InDet::SiDetElementsRoadMaker_xk"}
private

Definition at line 122 of file SiTrackMaker_xk.h.

122{this, "RoadTool", "InDet::SiDetElementsRoadMaker_xk"};

◆ m_seedsegmentsWrite

BooleanProperty InDet::SiTrackMaker_xk::m_seedsegmentsWrite {this, "SeedSegmentsWrite", false, "Call seed to track conversion"}
private

Definition at line 152 of file SiTrackMaker_xk.h.

152{this, "SeedSegmentsWrite", false, "Call seed to track conversion"};

◆ m_seedsfilter

BooleanProperty InDet::SiTrackMaker_xk::m_seedsfilter {this, "UseSeedFilter", true, "Use seed filter"}
private

Definition at line 139 of file SiTrackMaker_xk.h.

139{this, "UseSeedFilter", true, "Use seed filter"};

◆ m_seedtrack

ToolHandle<InDet::ISeedToTrackConversionTool> InDet::SiTrackMaker_xk::m_seedtrack {this, "SeedToTrackConversion", "InDet::SeedToTrackConversionTool"}
private

Definition at line 126 of file SiTrackMaker_xk.h.

126{this, "SeedToTrackConversion", "InDet::SeedToTrackConversionTool"};

◆ m_trackinfo

Trk::TrackInfo InDet::SiTrackMaker_xk::m_trackinfo
private

Definition at line 176 of file SiTrackMaker_xk.h.

◆ m_trackletPoints

IntegerProperty InDet::SiTrackMaker_xk::m_trackletPoints {this, "trackletPoints", 1, "Select which tracklet points to use"}
private

Definition at line 167 of file SiTrackMaker_xk.h.

167{this, "trackletPoints", 1, "Select which tracklet points to use"};

◆ m_tracksfinder

ToolHandle<InDet::ISiCombinatorialTrackFinder> InDet::SiTrackMaker_xk::m_tracksfinder {this, "CombinatorialTrackFinder", "InDet::SiCombinatorialTrackFinder_xk"}
private

Definition at line 123 of file SiTrackMaker_xk.h.

123{this, "CombinatorialTrackFinder", "InDet::SiCombinatorialTrackFinder_xk"};

◆ m_trigInDetRoadPredictorTool

ToolHandle<ITrigInDetRoadPredictorTool> InDet::SiTrackMaker_xk::m_trigInDetRoadPredictorTool {this, "TrigInDetRoadPredictorTool", "TrigInDetRoadPredictorTool_FTF"}
private

Definition at line 125 of file SiTrackMaker_xk.h.

125{this, "TrigInDetRoadPredictorTool", "TrigInDetRoadPredictorTool_FTF"};

◆ m_trigInDetTrackFollowingTool

ToolHandle<ITrigInDetTrackFollowingTool> InDet::SiTrackMaker_xk::m_trigInDetTrackFollowingTool {this, "TrigTrackFollowingTool", "TrigInDetTrackFollowingTool"}
private

Definition at line 124 of file SiTrackMaker_xk.h.

124{this, "TrigTrackFollowingTool", "TrigInDetTrackFollowingTool"};

◆ m_useassoTool

BooleanProperty InDet::SiTrackMaker_xk::m_useassoTool {this, "UseAssociationTool", false, "Use prd-track association tool"}
private

Definition at line 144 of file SiTrackMaker_xk.h.

144{this, "UseAssociationTool", false, "Use prd-track association tool"};

◆ m_useBremModel

BooleanProperty InDet::SiTrackMaker_xk::m_useBremModel {this, "useBremModel", false}
private

Definition at line 147 of file SiTrackMaker_xk.h.

147{this, "useBremModel", false};

◆ m_useCaloSeeds

BooleanProperty InDet::SiTrackMaker_xk::m_useCaloSeeds {this, "doCaloSeededBrem", false}
private

Definition at line 148 of file SiTrackMaker_xk.h.

148{this, "doCaloSeededBrem", false};

◆ m_useHClusSeed

BooleanProperty InDet::SiTrackMaker_xk::m_useHClusSeed {this, "doHadCaloSeedSSS", false, "Hadronic Calorimeter Seeds"}
private

Definition at line 150 of file SiTrackMaker_xk.h.

150{this, "doHadCaloSeedSSS", false, "Hadronic Calorimeter Seeds"};

◆ m_usePix

BooleanProperty InDet::SiTrackMaker_xk::m_usePix {this, "usePixel", true, "flags to set whether to use pixel/sct cluster, irrespective of what is in event"}
private

Definition at line 142 of file SiTrackMaker_xk.h.

142{this, "usePixel", true, "flags to set whether to use pixel/sct cluster, irrespective of what is in event"};

◆ m_useSct

BooleanProperty InDet::SiTrackMaker_xk::m_useSct {this, "useSCT", true}
private

Definition at line 143 of file SiTrackMaker_xk.h.

143{this, "useSCT", true};

◆ m_useSSSfilter

BooleanProperty InDet::SiTrackMaker_xk::m_useSSSfilter {this, "useSSSseedsFilter", true}
private

Definition at line 149 of file SiTrackMaker_xk.h.

149{this, "useSSSseedsFilter", true};

◆ m_useTrigInDetRoadPredictorTool

BooleanProperty InDet::SiTrackMaker_xk::m_useTrigInDetRoadPredictorTool {this, "useTrigInDetRoadPredictorTool", false, "Option to use TrigInDetRoadPredictorTool instead of ISiDetElementsRoadMaker"}
private

Definition at line 154 of file SiTrackMaker_xk.h.

154{this, "useTrigInDetRoadPredictorTool", false, "Option to use TrigInDetRoadPredictorTool instead of ISiDetElementsRoadMaker"};

◆ m_useTrigTrackFollowingTool

BooleanProperty InDet::SiTrackMaker_xk::m_useTrigTrackFollowingTool {this, "useTrigTrackFollowingTool", false, "Option to use TrigInDetTrackFollowingTool instead of SiCombinatorialTrackFinder_xk"}
private

Definition at line 153 of file SiTrackMaker_xk.h.

153{this, "useTrigTrackFollowingTool", false, "Option to use TrigInDetTrackFollowingTool instead of SiCombinatorialTrackFinder_xk"};

◆ m_xi2max

DoubleProperty InDet::SiTrackMaker_xk::m_xi2max {this, "Xi2max", 15., "max Xi2 for updators"}
private

Definition at line 156 of file SiTrackMaker_xk.h.

156{this, "Xi2max", 15., "max Xi2 for updators"};

◆ m_xi2maxlink

DoubleProperty InDet::SiTrackMaker_xk::m_xi2maxlink {this, "Xi2maxlink", 200., "max Xi2 for clusters"}
private

Definition at line 158 of file SiTrackMaker_xk.h.

158{this, "Xi2maxlink", 200., "max Xi2 for clusters"};

◆ m_xi2maxNoAdd

DoubleProperty InDet::SiTrackMaker_xk::m_xi2maxNoAdd {this, "Xi2maxNoAdd", 35., "max Xi2 for clusters"}
private

Definition at line 157 of file SiTrackMaker_xk.h.

157{this, "Xi2maxNoAdd", 35., "max Xi2 for clusters"};

◆ m_xi2multitracks

DoubleProperty InDet::SiTrackMaker_xk::m_xi2multitracks {this, "Xi2maxMultiTracks", 3., "max Xi2 for multi tracks"}
private

Definition at line 162 of file SiTrackMaker_xk.h.

162{this, "Xi2maxMultiTracks", 3., "max Xi2 for multi tracks"};

The documentation for this class was generated from the following files: