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

InDet::SiCombinatorialTrackFinder_xk is algorithm which produce track-finding in the road of InDetDD::SiDetectorElement* of SCT and Pixels sorted in propagation order. More...

#include <SiCombinatorialTrackFinder_xk.h>

Inheritance diagram for InDet::SiCombinatorialTrackFinder_xk:
Collaboration diagram for InDet::SiCombinatorialTrackFinder_xk:

Public Member Functions

Standard tool methods
 SiCombinatorialTrackFinder_xk (const std::string &, const std::string &, const IInterface *)
virtual ~SiCombinatorialTrackFinder_xk ()=default
virtual StatusCode initialize () override
virtual StatusCode finalize () override
Main methods for local track finding
virtual const std::list< Trk::Track * > & getTracks (SiCombinatorialTrackFinderData_xk &data, const Trk::TrackParameters &, const std::vector< const Trk::SpacePoint * > &, const std::vector< Amg::Vector3D > &, std::vector< const InDetDD::SiDetectorElement * > &, const TrackQualityCuts &, const EventContext &ctx) const override
virtual const std::list< Trk::Track * > & getTracks (SiCombinatorialTrackFinderData_xk &data, const Trk::TrackParameters &, const std::vector< const Trk::SpacePoint * > &, const std::vector< Amg::Vector3D > &, std::vector< const InDetDD::SiDetectorElement * > &, std::multimap< const Trk::PrepRawData *, const Trk::Track * > &, const EventContext &ctx) const override
virtual const std::list< Trk::Track * > & getTracksWithBrem (SiCombinatorialTrackFinderData_xk &data, const Trk::TrackParameters &, const std::vector< const Trk::SpacePoint * > &, const std::vector< Amg::Vector3D > &, std::vector< const InDetDD::SiDetectorElement * > &, std::multimap< const Trk::PrepRawData *, const Trk::Track * > &, bool, const EventContext &ctx) const override
virtual double pTseed (SiCombinatorialTrackFinderData_xk &data, const Trk::TrackParameters &, const std::vector< const Trk::SpacePoint * > &, const EventContext &) const override
virtual void newEvent (const EventContext &ctx, SiCombinatorialTrackFinderData_xk &data) const override
virtual void newEvent (const EventContext &ctx, SiCombinatorialTrackFinderData_xk &data, Trk::TrackInfo, const TrackQualityCuts &) const override
virtual void endEvent (SiCombinatorialTrackFinderData_xk &data) const override
Print internal tool parameters and status
MsgStream & dump (SiCombinatorialTrackFinderData_xk &data, MsgStream &out) const override

Private Attributes

Tool handles
ToolHandle< IInDetConditionsToolm_pixelCondSummaryTool
ToolHandle< IInDetConditionsToolm_sctCondSummaryTool
ToolHandle< Trk::IBoundaryCheckToolm_boundaryCheckTool
PublicToolHandle< Trk::IPatternParametersPropagatorm_proptool
PublicToolHandle< Trk::IPatternParametersUpdatorm_updatortool
PublicToolHandle< Trk::IRIO_OnTrackCreatorm_riocreator
Data handles
SG::ReadHandleKey< InDet::PixelClusterContainer > m_pixcontainerkey
SG::ReadHandleKey< InDet::SCT_ClusterContainer > m_sctcontainerkey
SG::ReadCondHandleKey< InDet::SiDetElementBoundaryLinks_xk > m_boundaryPixelKey
SG::ReadCondHandleKey< InDet::SiDetElementBoundaryLinks_xk > m_boundarySCTKey
SG::ReadCondHandleKey< AtlasFieldCacheCondObjm_fieldCondObjInputKey {this, "AtlasFieldCacheCondObj", "fieldCondObj", "Name of the Magnetic Field conditions object key"}
SG::ReadHandleKey< InDet::SiDetectorElementStatusm_pixelDetElStatus {this, "PixelDetElStatus", "", "Key of SiDetectorElementStatus for Pixel"}
 Optional read handle to get status data to test whether a pixel detector element is good.
SG::ReadHandleKey< InDet::SiDetectorElementStatusm_sctDetElStatus {this, "SCTDetElStatus", "", "Key of SiDetectorElementStatus for SCT"}
 Optional read handle to get status data to test whether a SCT detector element is good.
Properties
BooleanProperty m_usePIX {this, "usePixel", true}
BooleanProperty m_useSCT {this, "useSCT", true}
BooleanProperty m_ITkGeometry {this, "ITkGeometry", false}
BooleanProperty m_doFastTracking {this, "doFastTracking", false}
StringProperty m_fieldmode {this, "MagneticFieldMode", "MapSolenoid", "Mode of magnetic field"}
DoubleProperty m_qualityCut {this, "TrackQualityCut", 9.3, "Simple track quality cut"}
FloatProperty m_minPtCut {this, "MinFinalPtCut", 100, "Cut on the pt of the final track. Must be >0 to avoid NANs when computing eta."}
float m_minPt2Cut =0
BooleanProperty m_writeHolesFromPattern {this, "writeHolesFromPattern", false,"Flag to activate writing hole info from the pattern recognition"}

Friends

class SiTrajectory_xk

Data members, which are updated in only initialize

enum  EStat_t {
  TwoCluster , WrongRoad , WrongInit , CantFindTrk ,
  NotNewTrk , BremAttempt , NumberOfStats , Success
}
 Array entries for data.statistic counter. More...
int m_outputlevel {0}
Trk::MagneticFieldProperties m_fieldprop
 Magnetic field properties.
EStat_t findTrack (SiCombinatorialTrackFinderData_xk &data, const Trk::TrackParameters &, const std::vector< const Trk::SpacePoint * > &, const std::vector< Amg::Vector3D > &, std::vector< const InDetDD::SiDetectorElement * > &, std::multimap< const Trk::PrepRawData *, const Trk::Track * > &, const EventContext &) const
Trk::TrackconvertToTrack (SiCombinatorialTrackFinderData_xk &data, const EventContext &ctx) const
Trk::TrackconvertToNextTrack (SiCombinatorialTrackFinderData_xk &data) const
void magneticFieldInit ()
void detectorElementLinks (std::vector< const InDetDD::SiDetectorElement * > &, std::vector< const InDet::SiDetElementBoundaryLink_xk * > &, const EventContext &ctx) const
MsgStream & dumpconditions (MsgStream &out) const
void initializeCombinatorialData (const EventContext &ctx, SiCombinatorialTrackFinderData_xk &data) const
virtual void fillStatistic (SiCombinatorialTrackFinderData_xk &data, std::array< bool, NumberOfStats > &information) const override
static void getTrackQualityCuts (SiCombinatorialTrackFinderData_xk &data, const TrackQualityCuts &)
static bool spacePointsToClusters (const std::vector< const Trk::SpacePoint * > &, std::vector< const InDet::SiCluster * > &, std::optional< std::reference_wrapper< std::vector< const InDetDD::SiDetectorElement * > > >=std::nullopt)
static MsgStream & dumpevent (SiCombinatorialTrackFinderData_xk &data, MsgStream &out)

Detailed Description

InDet::SiCombinatorialTrackFinder_xk is algorithm which produce track-finding in the road of InDetDD::SiDetectorElement* of SCT and Pixels sorted in propagation order.

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

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

Definition at line 60 of file SiCombinatorialTrackFinder_xk.h.

Member Enumeration Documentation

◆ EStat_t

Constructor & Destructor Documentation

◆ SiCombinatorialTrackFinder_xk()

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

Definition at line 34 of file SiCombinatorialTrackFinder_xk.cxx.

36 : base_class(t, n, p)
37{
38}

◆ ~SiCombinatorialTrackFinder_xk()

virtual InDet::SiCombinatorialTrackFinder_xk::~SiCombinatorialTrackFinder_xk ( )
virtualdefault

Member Function Documentation

◆ convertToNextTrack()

Trk::Track * InDet::SiCombinatorialTrackFinder_xk::convertToNextTrack ( SiCombinatorialTrackFinderData_xk & data) const
private

Definition at line 885 of file SiCombinatorialTrackFinder_xk.cxx.

886{
887 auto tsos = std::make_unique<Trk::TrackStates>(
888 data.trajectory().convertToNextTrackStateOnSurface());
889 if (tsos->empty()){
890 return nullptr;
891 }
892
893 // verify first track parameters
894 const Trk::TrackParameters *param = nullptr;
895 for (const Trk::TrackStateOnSurface *a_tsos : *tsos) {
896 const Trk::TrackParameters *param = a_tsos->trackParameters();
897 if (param) {
898 break;
899 }
900 }
901
902 if (param) {
903 auto momentum = param->momentum();
904 const auto pt2 = momentum.perp2();
905 // reject tracks with small pT
906 // The cut should be large enough otherwise eta computation of such tracks may yield NANs.
907 if (pt2 < m_minPt2Cut) {
908 ATH_MSG_WARNING( "Reject low pT track (pT = " << sqrt(pt2) << " < " << m_minPtCut.value() << ")");
909 return nullptr;
910 }
911 }
912 return new Trk::Track(data.trackinfo(),
913 std::move(tsos),
914 data.trajectory().convertToFitQuality());
915}
#define ATH_MSG_WARNING(x)
char data[hepevt_bytes_allocation_ATLAS]
Definition HepEvt.cxx:11
const Amg::Vector3D & momentum() const
Access method for the momentum.
ParametersBase< TrackParametersDim, Charged > TrackParameters

◆ convertToTrack()

Trk::Track * InDet::SiCombinatorialTrackFinder_xk::convertToTrack ( SiCombinatorialTrackFinderData_xk & data,
const EventContext & ctx ) const
private

Definition at line 839 of file SiCombinatorialTrackFinder_xk.cxx.

840{
841 const Trk::PatternTrackParameters *param = data.trajectory().firstParameters();
842 if (param) {
843 double pt = param->transverseMomentum();
844 // reject tracks with small pT
845 // The cut should be large enough otherwise eta computation of such tracks may yield NANs.
846 if (pt < m_minPtCut.value()) {
847 ATH_MSG_DEBUG( "Reject low pT track (pT = " << pt << " < " << m_minPtCut.value() << ")");
848 return nullptr;
849 }
850 }
851 if (!data.simpleTrack()) {
852 return new Trk::Track(
853 data.trackinfo(),
854 std::make_unique<Trk::TrackStates>(
855 data.trajectory().convertToTrackStateOnSurface(
856 data.cosmicTrack())),
857 data.trajectory().convertToFitQuality());
858 }
859
860 Trk::TrackInfo info = data.trackinfo();
861 info.setPatternRecognitionInfo(Trk::TrackInfo::SiSPSeededFinderSimple);
862 info.setParticleHypothesis(Trk::pion);
863 if (!data.flagToReturnFailedTrack()) {
864 return new Trk::Track(
865 info,
866 std::make_unique<Trk::TrackStates>(
867 data.trajectory().convertToSimpleTrackStateOnSurface(
868 data.cosmicTrack(), ctx)),
869 data.trajectory().convertToFitQuality());
870 } else {
871 return new Trk::Track(
872 info,
873 std::make_unique<Trk::TrackStates>(
874 data.trajectory()
875 .convertToSimpleTrackStateOnSurfaceForDisTrackTrigger(
876 data.cosmicTrack(), ctx)),
877 data.trajectory().convertToFitQuality());
878 }
879}
#define ATH_MSG_DEBUG(x)
double transverseMomentum() const
@ SiSPSeededFinderSimple
for tracks processed by the trigger version of the SiSPSeededFinder

◆ detectorElementLinks()

void InDet::SiCombinatorialTrackFinder_xk::detectorElementLinks ( std::vector< const InDetDD::SiDetectorElement * > & DE,
std::vector< const InDet::SiDetElementBoundaryLink_xk * > & DEL,
const EventContext & ctx ) const
private

Definition at line 989 of file SiCombinatorialTrackFinder_xk.cxx.

993{
994 const InDet::SiDetElementBoundaryLinks_xk* boundaryPixel{nullptr};
995 const InDet::SiDetElementBoundaryLinks_xk* boundarySCT{nullptr};
996 if (m_usePIX) {
997 SG::ReadCondHandle<InDet::SiDetElementBoundaryLinks_xk> boundaryPixelHandle(m_boundaryPixelKey,ctx);
998 boundaryPixel = *boundaryPixelHandle;
999 if (boundaryPixel==nullptr) {
1000 ATH_MSG_FATAL(m_boundaryPixelKey.fullKey() << " returns null pointer");
1001 }
1002 }
1003 if (m_useSCT) {
1004 SG::ReadCondHandle<InDet::SiDetElementBoundaryLinks_xk> boundarySCTHandle(m_boundarySCTKey,ctx);
1005 boundarySCT = *boundarySCTHandle;
1006 if (boundarySCT==nullptr) {
1007 ATH_MSG_FATAL(m_boundarySCTKey.fullKey() << " returns null pointer");
1008 }
1009 }
1010
1011 DEL.reserve(DE.size());
1012 for (const InDetDD::SiDetectorElement* d: DE) {
1013 IdentifierHash id = d->identifyHash();
1014 if (d->isPixel() && boundaryPixel && id < boundaryPixel->size()) DEL.push_back(&(*boundaryPixel)[id]);
1015 else if (d->isSCT() && boundarySCT && id < boundarySCT->size()) DEL.push_back(&(*boundarySCT)[id]);
1016 }
1017}
#define ATH_MSG_FATAL(x)
SG::ReadCondHandleKey< InDet::SiDetElementBoundaryLinks_xk > m_boundaryPixelKey
SG::ReadCondHandleKey< InDet::SiDetElementBoundaryLinks_xk > m_boundarySCTKey
std::vector< SiDetElementBoundaryLink_xk > SiDetElementBoundaryLinks_xk

◆ dump()

MsgStream & InDet::SiCombinatorialTrackFinder_xk::dump ( SiCombinatorialTrackFinderData_xk & data,
MsgStream & out ) const
override

Definition at line 119 of file SiCombinatorialTrackFinder_xk.cxx.

120{
121 if (not data.isInitialized()) initializeCombinatorialData(Gaudi::Hive::currentContext(), data);
122
123 out<<std::endl;
124 if (data.nprint()) return dumpevent(data, out);
125 return dumpconditions(out);
126}
static MsgStream & dumpevent(SiCombinatorialTrackFinderData_xk &data, MsgStream &out)
void initializeCombinatorialData(const EventContext &ctx, SiCombinatorialTrackFinderData_xk &data) const

◆ dumpconditions()

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

Definition at line 132 of file SiCombinatorialTrackFinder_xk.cxx.

133{
134 int n = 62-m_proptool.type().size();
135 std::string s1;
136 for (int i=0; i<n; ++i) s1.append(" ");
137 s1.append("|");
138
139 std::string fieldmode[9] ={"NoField" , "ConstantField", "SolenoidalField",
140 "ToroidalField" , "Grid3DField" , "RealisticField" ,
141 "UndefinedField", "AthenaField" , "?????" };
142
143 int mode = m_fieldprop.magneticFieldMode();
144 if (mode<0 || mode>8 ) mode = 8;
145
146 n = 62-fieldmode[mode].size();
147 std::string s3;
148 for (int i=0; i<n; ++i) s3.append(" ");
149 s3.append("|");
150
151 n = 62-m_updatortool.type().size();
152 std::string s4;
153 for (int i=0; i<n; ++i) s4.append(" ");
154 s4.append("|");
155
156 n = 62-m_riocreator.type().size();
157 std::string s5;
158 for (int i=0; i<n; ++i) s5.append(" ");
159 s5.append("|");
160
161 n = 62-m_pixcontainerkey.key().size();
162 std::string s7;
163 for (int i=0; i<n; ++i) s7.append(" ");
164 s7.append("|");
165
166 n = 62-m_sctcontainerkey.key().size();
167 std::string s8;
168 for (int i=0; i<n; ++i) s8.append(" ");
169 s8.append("|");
170
171 out<<"|----------------------------------------------------------------------"
172 <<"-------------------|"
173 <<std::endl;
174 if (m_usePIX) {
175 out<<"| Pixel clusters location | "<<m_pixcontainerkey.key() <<s7<<std::endl;
176 }
177 if (m_useSCT) {
178 out<<"| SCT clusters location | "<<m_sctcontainerkey.key() <<s8<<std::endl;
179 }
180 out<<"| Tool for propagation | "<<m_proptool .type()<<s1<<std::endl;
181 out<<"| Tool for updator | "<<m_updatortool.type()<<s4<<std::endl;
182 out<<"| Tool for rio on track | "<<m_riocreator .type()<<s5<<std::endl;
183 out<<"| Magnetic field mode | "<<fieldmode[mode] <<s3<<std::endl;
184 out<<"|----------------------------------------------------------------------"
185 <<"-------------------|"
186 <<std::endl;
187 return out;
188}
PublicToolHandle< Trk::IPatternParametersPropagator > m_proptool
SG::ReadHandleKey< InDet::SCT_ClusterContainer > m_sctcontainerkey
PublicToolHandle< Trk::IPatternParametersUpdator > m_updatortool
Trk::MagneticFieldProperties m_fieldprop
Magnetic field properties.
SG::ReadHandleKey< InDet::PixelClusterContainer > m_pixcontainerkey
PublicToolHandle< Trk::IRIO_OnTrackCreator > m_riocreator

◆ dumpevent()

MsgStream & InDet::SiCombinatorialTrackFinder_xk::dumpevent ( SiCombinatorialTrackFinderData_xk & data,
MsgStream & out )
staticprivate

Definition at line 194 of file SiCombinatorialTrackFinder_xk.cxx.

195{
196 out<<"|---------------------------------------------------------------------|"
197 <<std::endl;
198 out<<"| Min pT of track (MeV) | "<<std::setw(12)<<std::setprecision(5)<<data.pTmin()
199 <<" |"<<std::endl;
200 out<<"| Max Xi2 for cluster | "<<std::setw(12)<<std::setprecision(5)<<data.xi2max()
201 <<" |"<<std::endl;
202 out<<"| Max Xi2 for outlayer | "<<std::setw(12)<<std::setprecision(5)<<data.xi2maxNoAdd()
203 <<" |"<<std::endl;
204 out<<"| Max Xi2 for link | "<<std::setw(12)<<std::setprecision(5)<<data.xi2maxlink()
205 <<" |"<<std::endl;
206 out<<"| Min number of clusters | "<<std::setw(12)<<data.nclusmin()
207 <<" |"<<std::endl;
208 out<<"| Min number of wclusters | "<<std::setw(12)<<data.nwclusmin()
209 <<" |"<<std::endl;
210 out<<"| Max number holes | "<<std::setw(12)<<data.nholesmax()
211 <<" |"<<std::endl;
212 out<<"| Max holes gap | "<<std::setw(12)<<data.dholesmax()
213 <<" |"<<std::endl;
214 out<<"| Use PRD to track assoc.?| "<<std::setw(12)<<(data.PRDtoTrackMap() ? "yes" : "no ")
215 <<" |"<<std::endl;
216 out<<"|---------------------------------------------------------------------|"
217 <<std::endl;
218 out<<"| Number input seeds | "<<std::setw(12)<<data.inputseeds()
219 <<" |"<<std::endl;
220 out<<"| Number accepted seeds | "<<std::setw(12)<<data.goodseeds()
221 <<" |"<<std::endl;
222 out<<"| Number initial tracks | "<<std::setw(12)<<data.inittracks()
223 <<" |"<<std::endl;
224 out<<"| Number wrong DE roads | "<<std::setw(12)<<data.roadbug()
225 <<" |"<<std::endl;
226 out<<"| Number output tracks | "<<std::setw(12)<<data.findtracks()
227 <<" |"<<std::endl;
228 out<<"|---------------------------------------------------------------------|"
229 <<std::endl;
230 return out;
231}

◆ endEvent()

void InDet::SiCombinatorialTrackFinder_xk::endEvent ( SiCombinatorialTrackFinderData_xk & data) const
overridevirtual

Definition at line 303 of file SiCombinatorialTrackFinder_xk.cxx.

304{
305 if (not data.isInitialized()) initializeCombinatorialData(Gaudi::Hive::currentContext(), data);
306
307 // Print event information
308 //
309 if (m_outputlevel<=0) {
310 data.nprint() = 1;
311 dump(data, msg(MSG::DEBUG));
312 }
313}
MsgStream & dump(SiCombinatorialTrackFinderData_xk &data, MsgStream &out) const override
MsgStream & msg
Definition testRead.cxx:32

◆ fillStatistic()

void InDet::SiCombinatorialTrackFinder_xk::fillStatistic ( SiCombinatorialTrackFinderData_xk & data,
std::array< bool, NumberOfStats > & information ) const
overrideprivatevirtual

Definition at line 1126 of file SiCombinatorialTrackFinder_xk.cxx.

1127{
1128 for(int i=0; i!=NumberOfStats; ++i) information[i] = data.statistic()[i];
1129}

◆ finalize()

StatusCode InDet::SiCombinatorialTrackFinder_xk::finalize ( )
overridevirtual

Definition at line 110 of file SiCombinatorialTrackFinder_xk.cxx.

111{
112 return AlgTool::finalize();
113}

◆ findTrack()

InDet::SiCombinatorialTrackFinder_xk::EStat_t InDet::SiCombinatorialTrackFinder_xk::findTrack ( SiCombinatorialTrackFinderData_xk & data,
const Trk::TrackParameters & Tp,
const std::vector< const Trk::SpacePoint * > & Sp,
const std::vector< Amg::Vector3D > & Gp,
std::vector< const InDetDD::SiDetectorElement * > & DE,
std::multimap< const Trk::PrepRawData *, const Trk::Track * > & PT,
const EventContext & ctx ) const
private

init event data

populate a list of boundary links for the detector elements on our search road

Retrieve cached pointers to SG collections, or create the cache

Cluster list preparationn

in inside-out track finding, Sp.size() is typically 3 for TRT-seeded backtracking, it is 2 both applications go into this branch

returns false if two clusters are on the same detector element

use case if we have a set of global positions rather than space points to start from

use case if we have neither space-points nor global posittions, but track parameters to start from

Build initial trajectory

This will initialize the trajectory using the clusters we have and the parameter estimate

if the initialisation fails (indicating this is probably a bad attempt) and we are running with global positions instead of seeding

reset our cluster list

try again using the clusters from the track parameters only

if it worked now, set the quality flag to true

this can never happen?!

this can never happen either?!

if the last cluster on track is in the pixels, this is assumed to come from a pixel seed

max #iterations

Track finding

Strategy for pixel seeds

refine if needed

check if we found enough clusters

case of a strip seed or mixed PPS

first application of hit cut

backward smooting

apply hit cut again following smoothing step

refine if needed

quality cut

refine the hole cut

Definition at line 583 of file SiCombinatorialTrackFinder_xk.cxx.

591{
593 if (not data.isInitialized()) initializeCombinatorialData(ctx, data);
594
596 std::vector<const InDet::SiDetElementBoundaryLink_xk*> DEL;
597 detectorElementLinks(DE, DEL,ctx);
598
600 const InDet::PixelClusterContainer* p_pixcontainer = data.pixContainer();
601 if (m_usePIX && !p_pixcontainer) {
602 SG::ReadHandle<InDet::PixelClusterContainer> pixcontainer(m_pixcontainerkey,ctx);
603 if(!pixcontainer.isValid()){
604 //the pixel container is not valid
605 //we are supposed to use pixel so stop here
606 REPORT_ERROR (StatusCode::FAILURE)<< "Unable to dereference pix container, handle not valid";
607 return WrongInit;
608 }
609
610 p_pixcontainer = pixcontainer.ptr();
611 data.setPixContainer(p_pixcontainer);
612 }
613 const InDet::SCT_ClusterContainer* p_sctcontainer = data.sctContainer();
614 if (m_useSCT && !p_sctcontainer) {
615 SG::ReadHandle<InDet::SCT_ClusterContainer> sctcontainer(m_sctcontainerkey,ctx);
616
617 if(!sctcontainer.isValid()){
618 //the sct container is not valid
619 //we are supposed to use sct so stop here
620 REPORT_ERROR (StatusCode::FAILURE)<< "Unable to dereference sct container, handle not valid";
621 return WrongInit;
622 }
623
624 p_sctcontainer = sctcontainer.ptr();
625 data.setSctContainer(p_sctcontainer);
626 }
627
629 std::vector<const InDet::SiCluster*> Cl;
630 bool isTwoPointSeed = false;
631
635 if (Sp.size() > 1) {
636
638 if (!spacePointsToClusters(Sp,Cl)) {
639 return TwoCluster;
640 }
641 if (Sp.size()==2) isTwoPointSeed = true;
642 }
644 else if (Gp.size() > 2) {
645 if (!data.trajectory().globalPositionsToClusters(p_pixcontainer, p_sctcontainer, Gp, DEL, PT, Cl)) return TwoCluster;
646 } else {
648 if (!data.trajectory().trackParametersToClusters(p_pixcontainer, p_sctcontainer, Tp, DEL, PT, Cl, ctx)) return TwoCluster;
649 }
650 ++data.goodseeds();
651
653 bool Qr;
655 bool Q = data.trajectory().initialize(m_usePIX, m_useSCT, p_pixcontainer, p_sctcontainer, Tp, Cl, DEL, Qr,ctx);
656
659 if (!Q && Sp.size() < 2 && Gp.size() > 3) {
661 Cl.clear();
663 if (!data.trajectory().trackParametersToClusters(p_pixcontainer, p_sctcontainer, Tp, DEL, PT, Cl,ctx)) return TwoCluster;
664
665 if (!data.trajectory().initialize(m_usePIX, m_useSCT, p_pixcontainer, p_sctcontainer, Tp, Cl, DEL, Qr,ctx)) return TwoCluster;
667
668 Q = Qr = true;
669 }
670
672 if (!Qr){
673 ++data.roadbug();
674 return WrongRoad;
675 }
677 if (!Q) return WrongInit;
678
679 ++data.inittracks();
681 bool pixseed = data.trajectory().isLastPixel();
683 int itmax = 30;
684 if (!data.useFastTracking() and data.simpleTrack()) itmax = 10;
685 if (data.heavyIon()) itmax = 50;
686
687 //
688 bool toReturnFailedTrack = data.flagToReturnFailedTrack();
689 if( toReturnFailedTrack ) data.setResultCode(SiCombinatorialTrackFinderData_xk::ResultCode::Unrecoverable);
690
692 if (pixseed) {
693 if (!data.trajectory().forwardExtension (false,itmax,ctx)) return CantFindTrk;
694 if (!data.trajectory().backwardSmoother (false,ctx) ) return CantFindTrk;
695 if (!data.trajectory().backwardExtension(itmax,ctx) ) return CantFindTrk;
696 if (data.isITkGeometry() &&
697 (data.trajectory().nclusters() < data.nclusmin() ||
698 data.trajectory().ndf() < data.nwclusmin()) ) return CantFindTrk;
699
701 if(!data.useFastTracking() && data.trajectory().difference() > 0){
702 if (!data.trajectory().forwardFilter(ctx)) {
703 if( toReturnFailedTrack ) {
705 }
706 else {
707 return CantFindTrk;
708 }
709 }
710
711 if (!data.trajectory().backwardSmoother (false,ctx) ) {
712 if( toReturnFailedTrack ) {
714 }
715 else {
716 return CantFindTrk;
717 }
718 }
719 }
720
721 int na = data.trajectory().nclustersNoAdd();
723 if (data.trajectory().nclusters()+na < data.nclusmin() ||
724 data.trajectory().ndf() < data.nwclusmin()) {
725 if( toReturnFailedTrack ) {
727 }
728 else {
729 return CantFindTrk;
730 }
731 }
732 }
733
735 else { // Strategy for mixed seeds
736 if (!data.trajectory().backwardSmoother(isTwoPointSeed,ctx) ) return CantFindTrk;
737 if (!data.trajectory().backwardExtension(itmax,ctx) ) return CantFindTrk;
738 if (!data.trajectory().forwardExtension(true,itmax,ctx)) return CantFindTrk;
739
741 int na = data.trajectory().nclustersNoAdd();
742 if (data.trajectory().nclusters()+na < data.nclusmin() ||
743 data.trajectory().ndf() < data.nwclusmin()) return CantFindTrk;
745 if (!data.trajectory().backwardSmoother(false,ctx) ) return CantFindTrk;
746
748 na = data.trajectory().nclustersNoAdd();
749 if (data.trajectory().nclusters()+na < data.nclusmin() ||
750 data.trajectory().ndf() < data.nwclusmin()) {
751 if( toReturnFailedTrack ) {
753 }
754 else {
755 return CantFindTrk;
756 }
757 }
758
760 if (data.trajectory().difference() > 0) {
761 if (!data.trajectory().forwardFilter(ctx)) {
762 if( toReturnFailedTrack ) {
764 }
765 else {
766 return CantFindTrk;
767 }
768 }
769 if (!data.trajectory().backwardSmoother (false, ctx)) {
770 if( toReturnFailedTrack ) {
772 }
773 else {
774 return CantFindTrk;
775 }
776 }
777 }
778 }
779
781 if (data.trajectory().qualityOptimization() < (m_qualityCut*data.nclusmin())) {
782 if( toReturnFailedTrack ) {
784 }
785 else {
786 return CantFindTrk;
787 }
788 }
789
790 if (data.trajectory().pTfirst () < data.pTmin () &&
791 data.trajectory().nclusters() < data.nclusmin()) {
792 if( toReturnFailedTrack ) {
794 }
795 else {
796 return CantFindTrk;
797 }
798 }
799
800 if (data.trajectory().nclusters() < data.nclusminb() ||
801 data.trajectory().ndf () < data.nwclusmin()) {
802 if( toReturnFailedTrack ) {
804 }
805 else {
806 return CantFindTrk;
807 }
808 }
809
812 data.trajectory().updateHoleSearchResult();
813 if (!data.trajectory().getHoleSearchResult().passPatternHoleCut) {
814 if( toReturnFailedTrack ) {
816 }
817 else {
818 return CantFindTrk;
819 }
820 }
821 }
822
823 if( toReturnFailedTrack ) {
825 return CantFindTrk;
826 }
827 else {
829 }
830 }
831
832 return Success;
833}
#define REPORT_ERROR(SC)
Report an error.
static Double_t Tp(Double_t *t, Double_t *par)
static bool spacePointsToClusters(const std::vector< const Trk::SpacePoint * > &, std::vector< const InDet::SiCluster * > &, std::optional< std::reference_wrapper< std::vector< const InDetDD::SiDetectorElement * > > >=std::nullopt)
void detectorElementLinks(std::vector< const InDetDD::SiDetectorElement * > &, std::vector< const InDet::SiDetElementBoundaryLink_xk * > &, const EventContext &ctx) const
Trk::PrepRawDataContainer< SCT_ClusterCollection > SCT_ClusterContainer
Trk::PrepRawDataContainer< PixelClusterCollection > PixelClusterContainer

◆ getTrackQualityCuts()

void InDet::SiCombinatorialTrackFinder_xk::getTrackQualityCuts ( SiCombinatorialTrackFinderData_xk & data,
const TrackQualityCuts & Cuts )
staticprivate

Definition at line 1023 of file SiCombinatorialTrackFinder_xk.cxx.

1025{
1026 // Integer cuts
1027 //
1028 int intCut = 0;
1029
1030 if (!Cuts.getIntCut("CosmicTrack" , intCut)) intCut = 0;
1031 data.setCosmicTrack(intCut);
1032
1033 if (!Cuts.getIntCut("MinNumberOfClusters" , intCut)) intCut = 7;
1034 data.setNclusmin(intCut);
1035
1036 data.setNclusminb(std::max(3, data.nclusmin()-1));
1037
1038 if (!Cuts.getIntCut("MinNumberOfWClusters", intCut)) intCut = 7;
1039 data.setNwclusmin(intCut);
1040
1041 if (!Cuts.getIntCut("MaxNumberOfHoles" , intCut)) intCut = 2;
1042 if(!data.cosmicTrack() && intCut>2) intCut = 2;
1043 data.setNholesmax(intCut);
1044
1045 if (!Cuts.getIntCut("MaxHolesGap" , intCut)) intCut = 2;
1046 if(!data.cosmicTrack() && intCut>2) intCut = 2;
1047 if(intCut > data.nholesmax()) intCut = data.nholesmax();
1048 data.setDholesmax(intCut);
1049
1050 data.tools().setHolesClusters(data.nholesmax(), data.dholesmax(),
1051 data.nclusmin());
1052
1053 if (!Cuts.getIntCut("UseAssociationTool" , intCut)) intCut = 0;
1054 data.tools().setAssociation(intCut);
1055
1056 if (!Cuts.getIntCut("SimpleTrack" , intCut)) intCut = 0;
1057 data.setSimpleTrack(bool(intCut));
1058
1059 // Double cuts
1060 //
1061 double doubleCut = 0.;
1062
1063 if (!Cuts.getDoubleCut("pTmin" , doubleCut)) doubleCut = 500.;
1064 data.setPTmin(doubleCut);
1065
1066 if (!Cuts.getDoubleCut("pTminBrem" , doubleCut)) doubleCut = 1000.;
1067 data.setPTminBrem(doubleCut);
1068
1069 if (!Cuts.getDoubleCut("MaxXi2forCluster" , doubleCut)) doubleCut = 9.;
1070 data.setXi2max(doubleCut);
1071
1072 if (!Cuts.getDoubleCut("MaxXi2forOutlier" , doubleCut)) doubleCut = 25.;
1073 if (!data.cosmicTrack() && doubleCut > 25.) doubleCut = 25.;
1074 if (doubleCut <= data.xi2max()) doubleCut = data.xi2max()+5.;
1075 data.setXi2maxNoAdd(doubleCut);
1076
1077 if (!Cuts.getDoubleCut("MaxXi2forSearch" , doubleCut)) doubleCut = 100.;
1078 data.setXi2maxlink(doubleCut);
1079
1080 data.tools().setXi2pTmin(data.xi2max(), data.xi2maxNoAdd(),
1081 data.xi2maxlink(), data.pTmin());
1082
1083 if (!Cuts.getIntCut ("doMultiTracksProd", intCut)) intCut = 0;
1084 if (!Cuts.getDoubleCut("MaxXi2MultiTracks", doubleCut)) doubleCut = 7.;
1085 if (!data.cosmicTrack() && doubleCut > 7.) doubleCut = 7.;
1086 data.tools().setMultiTracks(intCut, doubleCut);
1087
1088 data.trajectory().setParameters();
1089}

◆ getTracks() [1/2]

const std::list< Trk::Track * > & InDet::SiCombinatorialTrackFinder_xk::getTracks ( SiCombinatorialTrackFinderData_xk & data,
const Trk::TrackParameters & Tp,
const std::vector< const Trk::SpacePoint * > & Sp,
const std::vector< Amg::Vector3D > & Gp,
std::vector< const InDetDD::SiDetectorElement * > & DE,
const TrackQualityCuts & Cuts,
const EventContext & ctx ) const
overridevirtual

turn off brem fit & electron flags

remove existing tracks

try to find the tracks

if we didn't find anything, bail out

sort in step order

Definition at line 319 of file SiCombinatorialTrackFinder_xk.cxx.

327{
328
329 if (not data.isInitialized()) initializeCombinatorialData(ctx, data);
330
331 data.statistic().fill(false);
332
334 data.tools().setBremNoise(false, false);
336 data.tracks().erase(data.tracks().begin(), data.tracks().end());
337
338 ++data.inputseeds();
339 if (!m_usePIX && !m_useSCT) {
340 return data.tracks();
341 }
342
343 // Get track quality cuts information and write them into the event data... again?
345 std::multimap<const Trk::PrepRawData*, const Trk::Track*> PT;
346
348 EStat_t FT = findTrack(data, Tp, Sp, Gp, DE, PT, ctx);
349
351 if(FT!=Success) {
352 data.statistic()[FT] = true;
353 if( ! data.flagToReturnFailedTrack() ) return data.tracks();
354 }
355
357 data.trajectory().sortStep();
358
359
360 Trk::Track* t = convertToTrack(data,ctx);
361 if (t) {
362 ++data.findtracks();
363 if (m_writeHolesFromPattern) data.addPatternHoleSearchOutcome(t,data.trajectory().getHoleSearchResult());
364 data.tracks().push_back(t);
365 }
366
367 if (!data.tools().multiTrack() ||
368 data.simpleTrack() ||
369 Sp.size()<=2 ||
370 data.cosmicTrack() ||
371 data.trajectory().pTfirst() < data.tools().pTmin()) return data.tracks();
372
373 while ((t=convertToNextTrack(data))) {
374 ++data.findtracks();
375 if (m_writeHolesFromPattern) data.addPatternHoleSearchOutcome(t,data.trajectory().getHoleSearchResult());
376 data.tracks().push_back(t);
377 }
378 return data.tracks();
379}
Trk::Track * convertToNextTrack(SiCombinatorialTrackFinderData_xk &data) const
static void getTrackQualityCuts(SiCombinatorialTrackFinderData_xk &data, const TrackQualityCuts &)
Trk::Track * convertToTrack(SiCombinatorialTrackFinderData_xk &data, const EventContext &ctx) const
EStat_t findTrack(SiCombinatorialTrackFinderData_xk &data, const Trk::TrackParameters &, const std::vector< const Trk::SpacePoint * > &, const std::vector< Amg::Vector3D > &, std::vector< const InDetDD::SiDetectorElement * > &, std::multimap< const Trk::PrepRawData *, const Trk::Track * > &, const EventContext &) const
H5::PredType PT
Definition H5Traits.cxx:15

◆ getTracks() [2/2]

const std::list< Trk::Track * > & InDet::SiCombinatorialTrackFinder_xk::getTracks ( SiCombinatorialTrackFinderData_xk & data,
const Trk::TrackParameters & Tp,
const std::vector< const Trk::SpacePoint * > & Sp,
const std::vector< Amg::Vector3D > & Gp,
std::vector< const InDetDD::SiDetectorElement * > & DE,
std::multimap< const Trk::PrepRawData *, const Trk::Track * > & PT,
const EventContext & ctx ) const
overridevirtual

Definition at line 385 of file SiCombinatorialTrackFinder_xk.cxx.

393{
394
395 if (not data.isInitialized()) initializeCombinatorialData(ctx, data);
396
397 data.tools().setBremNoise(false, false);
398 data.tracks().erase(data.tracks().begin(), data.tracks().end());
399
400 data.statistic().fill(false);
401 ++data.inputseeds();
402 if (!m_usePIX && !m_useSCT) {
403 return data.tracks();
404 }
405
406 EStat_t FT = findTrack(data, Tp, Sp, Gp, DE, PT,ctx);
407 if(FT!=Success) {
408 data.statistic()[FT] = true;
409 if( ! data.flagToReturnFailedTrack() || std::abs(data.resultCode()) < SiCombinatorialTrackFinderData_xk::ResultCodeThreshold::RecoverableForDisTrk ) return data.tracks();
410 }
411 if (!data.trajectory().isNewTrack(PT))
412 {
413 data.statistic()[NotNewTrk] = true;
414 return data.tracks();
415 }
416 data.trajectory().sortStep();
417
418 if(data.useFastTracking()) {
419 if(!data.trajectory().filterWithPreciseClustersError(ctx)) {
420 data.statistic()[CantFindTrk] = true;
421 return data.tracks();
422 }
423 }
424
425 // Trk::Track production
426 //
427 Trk::Track* t = convertToTrack(data,ctx);
428 if (t==nullptr) return data.tracks(); // @TODO should one check if convertToNextTrack would yield anything ?
429 if (m_writeHolesFromPattern) data.addPatternHoleSearchOutcome(t,data.trajectory().getHoleSearchResult());
430
431 ++data.findtracks();
432 data.tracks().push_back(t);
433
434 if (!data.tools().multiTrack() ||
435 data.simpleTrack() ||
436 Sp.size()<=2 ||
437 data.cosmicTrack() ||
438 data.trajectory().pTfirst() < data.tools().pTmin()) return data.tracks();
439
440 while ((t=convertToNextTrack(data))) {
441 if (m_writeHolesFromPattern) data.addPatternHoleSearchOutcome(t,data.trajectory().getHoleSearchResult());
442 ++data.findtracks();
443 data.tracks().push_back(t);
444 }
445 return data.tracks();
446}

◆ getTracksWithBrem()

const std::list< Trk::Track * > & InDet::SiCombinatorialTrackFinder_xk::getTracksWithBrem ( SiCombinatorialTrackFinderData_xk & data,
const Trk::TrackParameters & Tp,
const std::vector< const Trk::SpacePoint * > & Sp,
const std::vector< Amg::Vector3D > & Gp,
std::vector< const InDetDD::SiDetectorElement * > & DE,
std::multimap< const Trk::PrepRawData *, const Trk::Track * > & PT,
bool isCaloCompatible,
const EventContext & ctx ) const
overridevirtual

Definition at line 453 of file SiCombinatorialTrackFinder_xk.cxx.

462{
463
464 if (not data.isInitialized()) initializeCombinatorialData(ctx, data);
465
466 data.statistic().fill(false);
467
468 // Old information
469 //
470 int mult = 0;
471 if (data.tools().multiTrack()) mult = 1;
472 double Xi2m = data.tools().xi2multi();
473
474 data.tools().setBremNoise(false, true);
475 data.tracks().erase(data.tracks().begin(), data.tracks().end());
476
477 ++data.inputseeds();
478 if (!m_usePIX && !m_useSCT) {
479 return data.tracks();
480 }
481
482 EStat_t FT = findTrack(data,Tp,Sp,Gp,DE,PT,ctx);
483
484 bool Q = (FT==Success);
485 if (Q){
486 Q = data.trajectory().isNewTrack(PT);
487 }
488
489 int na = 0;
490 if (Q) {
491 data.trajectory().sortStep();
492
493 // Trk::Track production
494 //
495 Trk::TrackInfo oldinfo = data.trackinfo();
496 if (isCaloCompatible) data.trackinfo().setPatternRecognitionInfo(Trk::TrackInfo::TrackInCaloROI);
497
498 data.tools().setMultiTracks(0, Xi2m);
499 Trk::Track* t = convertToTrack(data,ctx);
500 data.trackinfo() = oldinfo;
501 data.tools().setMultiTracks(mult,Xi2m);
502
503 if (t==nullptr) return data.tracks(); // @TODO should one check whether the next findTrack call would yield something ?
504 if (m_writeHolesFromPattern) data.addPatternHoleSearchOutcome(t,data.trajectory().getHoleSearchResult());
505 ++data.findtracks();
506 data.tracks().push_back(t);
507 na = data.trajectory().nclusters();
508 if (na >=12 && !data.trajectory().nclustersNoAdd()) return data.tracks();
509
510 if (data.trajectory().pTfirst() < data.pTminBrem()) return data.tracks();
511 }
512 if ((*Sp.begin())->clusterList().second) return data.tracks();
513
514 // Repeat track finding using electron noise model
515 //
516 data.statistic()[BremAttempt] = true;
517
518 data.tools().setBremNoise(true,true);
519 FT = findTrack(data, Tp, Sp, Gp, DE, PT,ctx);
520
521 if (FT!=Success) {
522 data.statistic()[FT] = true;
523 return data.tracks();
524 }
525
526 if (!data.trajectory().isNewTrack(PT)) { data.statistic()[NotNewTrk] = true; return data.tracks(); }
527
528 int nb = data.trajectory().nclusters();
529 if (nb <= na) return data.tracks();
530
531 data.trajectory().sortStep();
532
533 if(data.useFastTracking()) {
534 if(!data.trajectory().filterWithPreciseClustersError(ctx)) {
535 data.statistic()[CantFindTrk] = true;
536 return data.tracks();
537 }
538 }
539
540 // Trk::Track production
541 //
542 Trk::TrackInfo oldinfo = data.trackinfo();
544 data.trackinfo().setTrackProperties(Trk::TrackInfo::BremFitSuccessful);
545 if (isCaloCompatible) data.trackinfo().setPatternRecognitionInfo(Trk::TrackInfo::TrackInCaloROI);
546
547 data.tools().setMultiTracks(0, Xi2m);
548 Trk::Track* t = convertToTrack(data, ctx);
549 data.trackinfo() = oldinfo;
550 data.tools().setMultiTracks(mult, Xi2m);
551
552 if (t==nullptr) return data.tracks();
553 if (m_writeHolesFromPattern) data.addPatternHoleSearchOutcome(t,data.trajectory().getHoleSearchResult());
554
555 ++data.findtracks();
556 data.tracks().push_back(t);
557 return data.tracks();
558}
void setTrackProperties(const TrackProperties &properties)
Methods setting the properties of track.
@ BremFit
A brem fit was performed on this track.
@ BremFitSuccessful
A brem fit was performed on this track and this fit was successful.

◆ initialize()

StatusCode InDet::SiCombinatorialTrackFinder_xk::initialize ( )
overridevirtual

Definition at line 44 of file SiCombinatorialTrackFinder_xk.cxx.

45{
46 // Get RungeKutta propagator tool
47 //
48 if ( m_proptool.retrieve().isFailure() ) {
49 ATH_MSG_FATAL("Failed to retrieve tool " << m_proptool);
50 return StatusCode::FAILURE;
51 }
52 ATH_MSG_DEBUG("Retrieved tool " << m_proptool);
53
54 // Get updator tool
55 //
56 if ( m_updatortool.retrieve().isFailure() ) {
57 ATH_MSG_FATAL("Failed to retrieve tool " << m_updatortool);
58 return StatusCode::FAILURE;
59 }
60 ATH_MSG_DEBUG("Retrieved tool " << m_updatortool);
61
62 // Get RIO_OnTrack creator
63 //
64 if ( m_riocreator.retrieve().isFailure() ) {
65 ATH_MSG_FATAL("Failed to retrieve tool " << m_riocreator);
66 return StatusCode::FAILURE;
67 }
68 ATH_MSG_DEBUG("Retrieved tool " << m_riocreator);
69
70 // disable pixel/SCT conditions summary tool: pixel/SCT are not used, the status event data is used and not being validated
71 ATH_CHECK( m_pixelCondSummaryTool.retrieve( DisableTool{!m_usePIX || (!m_pixelDetElStatus.empty() && !VALIDATE_STATUS_ARRAY_ACTIVATED)} ) );
72 ATH_CHECK( m_sctCondSummaryTool.retrieve( DisableTool{!m_useSCT || (!m_sctDetElStatus.empty() && !VALIDATE_STATUS_ARRAY_ACTIVATED)} ) );
73 //
74 // Get InDetBoundaryCheckTool
75 if ( m_boundaryCheckTool.retrieve().isFailure() ) {
76 ATH_MSG_FATAL("Failed to retrieve tool " << m_boundaryCheckTool);
77 return StatusCode::FAILURE;
78 }
79 ATH_MSG_DEBUG("Retrieved tool " << m_boundaryCheckTool);
80
81
82
83 // Setup callback for magnetic field
84 //
86
87 // Get output print level
88 //
89 m_outputlevel = msg().level()-MSG::DEBUG;
90
91 ATH_CHECK( m_pixcontainerkey.initialize (m_usePIX) );
93
94 ATH_CHECK( m_sctcontainerkey.initialize (m_useSCT) );
95 ATH_CHECK( m_boundarySCTKey.initialize (m_useSCT) );
96
97 // initialize conditions object key for field cache
98 //
99 ATH_CHECK( m_fieldCondObjInputKey.initialize() );
100 ATH_CHECK( m_pixelDetElStatus.initialize( !m_pixelDetElStatus.empty() && m_usePIX) );
101 ATH_CHECK( m_sctDetElStatus.initialize( !m_sctDetElStatus.empty() && m_useSCT) );
102 m_minPt2Cut = std::pow(m_minPtCut.value(),2);
103 return StatusCode::SUCCESS;
104}
#define ATH_CHECK
Evaluate an expression and check for errors.
ToolHandle< Trk::IBoundaryCheckTool > m_boundaryCheckTool
SG::ReadCondHandleKey< AtlasFieldCacheCondObj > m_fieldCondObjInputKey
ToolHandle< IInDetConditionsTool > m_pixelCondSummaryTool
SG::ReadHandleKey< InDet::SiDetectorElementStatus > m_sctDetElStatus
Optional read handle to get status data to test whether a SCT detector element is good.
SG::ReadHandleKey< InDet::SiDetectorElementStatus > m_pixelDetElStatus
Optional read handle to get status data to test whether a pixel detector element is good.
ToolHandle< IInDetConditionsTool > m_sctCondSummaryTool

◆ initializeCombinatorialData()

void InDet::SiCombinatorialTrackFinder_xk::initializeCombinatorialData ( const EventContext & ctx,
SiCombinatorialTrackFinderData_xk & data ) const
private

Add conditions object to SiCombinatorialTrackFinderData to be able to access the field cache for each new event Get conditions object for field cache

Must have set fieldCondObj BEFORE calling setTools because fieldCondObj is used there

Definition at line 1091 of file SiCombinatorialTrackFinder_xk.cxx.

1091 {
1092
1095 SG::ReadCondHandle<AtlasFieldCacheCondObj> readHandle{m_fieldCondObjInputKey, ctx};
1096 const AtlasFieldCacheCondObj* fieldCondObj{*readHandle};
1097 if (fieldCondObj == nullptr) {
1098 std::string msg = "InDet::SiCombinatorialTrackFinder_xk::initializeCombinatorialData: Failed to retrieve AtlasFieldCacheCondObj with key " + m_fieldCondObjInputKey.key();
1099 throw(std::runtime_error(msg));
1100 }
1101 data.setFieldCondObj(fieldCondObj);
1102
1104 data.setTools(&*m_proptool,
1105 &*m_updatortool,
1106 &*m_riocreator,
1109 &m_fieldprop,
1111 if (!m_pixelDetElStatus.empty()) {
1112 SG::ReadHandle<InDet::SiDetectorElementStatus> pixelDetElStatus( m_pixelDetElStatus, ctx);
1113 data.setPixelDetectorElementStatus( pixelDetElStatus.cptr() );
1114 }
1115 if (!m_sctDetElStatus.empty()) {
1116 SG::ReadHandle<InDet::SiDetectorElementStatus> sctDetElStatus( m_sctDetElStatus, ctx);
1117 data.setSCTDetectorElementStatus( sctDetElStatus.cptr() );
1118 }
1119
1120 // Set the ITk Geometry setup
1121 data.setITkGeometry(m_ITkGeometry);
1122 // Set the ITk Fast Tracking setup
1123 data.setFastTracking(m_doFastTracking);
1124}
#define VALIDATE_STATUS_ARRAY_ACTIVATED

◆ magneticFieldInit()

void InDet::SiCombinatorialTrackFinder_xk::magneticFieldInit ( )
private

Definition at line 921 of file SiCombinatorialTrackFinder_xk.cxx.

922{
923 // Build MagneticFieldProperties
924 //
925 if (m_fieldmode == "NoField" ) m_fieldprop = Trk::MagneticFieldProperties(Trk::NoField );
926 else if (m_fieldmode == "MapSolenoid") m_fieldprop = Trk::MagneticFieldProperties(Trk::FastField);
927 else m_fieldprop = Trk::MagneticFieldProperties(Trk::FullField);
928}
@ FastField
call the fast field access method of the FieldSvc
@ NoField
Field is set to 0., 0., 0.,.
@ FullField
Field is set to be realistic, but within a given Volume.

◆ newEvent() [1/2]

void InDet::SiCombinatorialTrackFinder_xk::newEvent ( const EventContext & ctx,
SiCombinatorialTrackFinderData_xk & data ) const
overridevirtual

Definition at line 237 of file SiCombinatorialTrackFinder_xk.cxx.

238{
239 if (not data.isInitialized()) initializeCombinatorialData(ctx, data);
240
241 // Erase statistic information
242 //
243 data.inputseeds() = 0;
244 data.goodseeds() = 0;
245 data.inittracks() = 0;
246 data.findtracks() = 0;
247 data.roadbug() = 0;
248
249 // Set track info
250 //
251 data.trackinfo().setPatternRecognitionInfo(Trk::TrackInfo::SiSPSeededFinder);
252 data.setCosmicTrack(0);
253
254 // Add conditions object to SiCombinatorialTrackFinderData to be able to access the field cache for each new event
255 // Get conditions object for field cache
256 SG::ReadCondHandle<AtlasFieldCacheCondObj> readHandle{m_fieldCondObjInputKey, ctx};
257 const AtlasFieldCacheCondObj* fieldCondObj{*readHandle};
258 if (fieldCondObj == nullptr) {
259 std::string msg = "InDet::SiCombinatorialTrackFinder_xk::newEvent: Failed to retrieve AtlasFieldCacheCondObj with key " + m_fieldCondObjInputKey.key();
260 throw(std::runtime_error(msg));
261 }
262 data.setFieldCondObj(fieldCondObj);
263}
@ SiSPSeededFinder
Tracks from SiSPSeedFinder.

◆ newEvent() [2/2]

void InDet::SiCombinatorialTrackFinder_xk::newEvent ( const EventContext & ctx,
SiCombinatorialTrackFinderData_xk & data,
Trk::TrackInfo info,
const TrackQualityCuts & Cuts ) const
overridevirtual

Get track quality cuts information from argument and write it into the event data

update pattern recognition flags in the event data based on track info arg

Definition at line 269 of file SiCombinatorialTrackFinder_xk.cxx.

271{
272
273 if (not data.isInitialized()) {
274 //Check if to use PRDAssociation before initializing all the tools
275 int useasso;
276 if(!Cuts.getIntCut ("UseAssociationTool" ,useasso )) useasso = 0;
277 data.tools().setAssociation(useasso);
278
280 }
281
282 newEvent(ctx, data);
283 data.trackinfo() = info;
284
288
289 data.setHeavyIon(false);
290 data.setCosmicTrack(0);
293 data.setCosmicTrack(1);
295 data.setHeavyIon(true);
296
297}
virtual void newEvent(const EventContext &ctx, SiCombinatorialTrackFinderData_xk &data) const override
@ SiSpacePointsSeedMaker_Cosmic
Entries allowing to distinguish different seed makers.

◆ pTseed()

double InDet::SiCombinatorialTrackFinder_xk::pTseed ( SiCombinatorialTrackFinderData_xk & data,
const Trk::TrackParameters & Tp,
const std::vector< const Trk::SpacePoint * > & Sp,
const EventContext & ctx ) const
overridevirtual

Definition at line 564 of file SiCombinatorialTrackFinder_xk.cxx.

569{
570 std::vector<const InDet::SiCluster*> Cl;
571 std::vector<const InDetDD::SiDetectorElement*> DE;
572 if(!spacePointsToClusters(Sp,Cl,DE)) return 0.;
573
574 std::vector<const InDet::SiDetElementBoundaryLink_xk*> DEL;
575 detectorElementLinks(DE,DEL,ctx);
576 return data.trajectory().pTseed(Tp,Cl,DEL,ctx);
577}

◆ spacePointsToClusters()

bool InDet::SiCombinatorialTrackFinder_xk::spacePointsToClusters ( const std::vector< const Trk::SpacePoint * > & Sp,
std::vector< const InDet::SiCluster * > & Sc,
std::optional< std::reference_wrapper< std::vector< const InDetDD::SiDetectorElement * > > > DE = std::nullopt )
staticprivate

loop over all SP

get the first cluster on an SP

add to list

for strips, also make sure to pick up the second one!

Detector elments test

here we reject cases where two subsequent clusters are on the same detector element

Definition at line 934 of file SiCombinatorialTrackFinder_xk.cxx.

936{
937 Sc.reserve(Sp.size());
939 for (const Trk::SpacePoint* s : Sp) {
941 const Trk::PrepRawData* p = s->clusterList().first;
942 if (p) {
944 const InDet::SiCluster* c = static_cast<const InDet::SiCluster*>(p);
945 if (c){
946 Sc.push_back(c);
947 }
948 }
950 p = s->clusterList().second;
951 if (p) {
952 const InDet::SiCluster* c = static_cast<const InDet::SiCluster*>(p);
953 if (c){
954 Sc.push_back(c);
955 }
956 }
957 }
958
960 std::vector<const InDet::SiCluster*>::iterator cluster = Sc.begin(), nextCluster, endClusters = Sc.end();
961
964 if (DE) {
965 DE->get().reserve(Sc.size());
966 }
967 for (; cluster != endClusters; ++cluster) {
968
969 const InDetDD::SiDetectorElement* de = (*cluster)->detectorElement();
970
971 nextCluster = cluster;
972 ++nextCluster;
973 for (; nextCluster != endClusters; ++nextCluster) {
974 if (de == (*nextCluster)->detectorElement()){
975 return false;
976 }
977 }
978 if (DE) {
979 DE->get().push_back(de);
980 }
981 }
982 return true;
983}

◆ SiTrajectory_xk

friend class SiTrajectory_xk
friend

Definition at line 64 of file SiCombinatorialTrackFinder_xk.h.

Member Data Documentation

◆ m_boundaryCheckTool

ToolHandle<Trk::IBoundaryCheckTool> InDet::SiCombinatorialTrackFinder_xk::m_boundaryCheckTool
private
Initial value:
{this, "BoundaryCheckTool",
"InDet::InDetBoundaryCheckTool", "Boundary checking tool for detector sensitivities"}

Definition at line 149 of file SiCombinatorialTrackFinder_xk.h.

149 {this, "BoundaryCheckTool",
150 "InDet::InDetBoundaryCheckTool", "Boundary checking tool for detector sensitivities"};

◆ m_boundaryPixelKey

SG::ReadCondHandleKey<InDet::SiDetElementBoundaryLinks_xk> InDet::SiCombinatorialTrackFinder_xk::m_boundaryPixelKey
private
Initial value:
{this, "PixelDetElementBoundaryLinks_xk",
"PixelDetElementBoundaryLinks_xk", "Key of InDet::SiDetElementBoundaryLinks_xk for Pixel"}

Definition at line 166 of file SiCombinatorialTrackFinder_xk.h.

166 {this, "PixelDetElementBoundaryLinks_xk",
167 "PixelDetElementBoundaryLinks_xk", "Key of InDet::SiDetElementBoundaryLinks_xk for Pixel"};

◆ m_boundarySCTKey

SG::ReadCondHandleKey<InDet::SiDetElementBoundaryLinks_xk> InDet::SiCombinatorialTrackFinder_xk::m_boundarySCTKey
private
Initial value:
{this, "SCT_DetElementBoundaryLinks_xk",
"SCT_DetElementBoundaryLinks_xk", "Key of InDet::SiDetElementBoundaryLinks_xk for SCT"}

Definition at line 168 of file SiCombinatorialTrackFinder_xk.h.

168 {this, "SCT_DetElementBoundaryLinks_xk",
169 "SCT_DetElementBoundaryLinks_xk", "Key of InDet::SiDetElementBoundaryLinks_xk for SCT"};

◆ m_doFastTracking

BooleanProperty InDet::SiCombinatorialTrackFinder_xk::m_doFastTracking {this, "doFastTracking", false}
private

Definition at line 191 of file SiCombinatorialTrackFinder_xk.h.

191{this, "doFastTracking", false};

◆ m_fieldCondObjInputKey

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

Definition at line 171 of file SiCombinatorialTrackFinder_xk.h.

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

◆ m_fieldmode

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

Definition at line 192 of file SiCombinatorialTrackFinder_xk.h.

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

◆ m_fieldprop

Trk::MagneticFieldProperties InDet::SiCombinatorialTrackFinder_xk::m_fieldprop
private

Magnetic field properties.

Definition at line 202 of file SiCombinatorialTrackFinder_xk.h.

◆ m_ITkGeometry

BooleanProperty InDet::SiCombinatorialTrackFinder_xk::m_ITkGeometry {this, "ITkGeometry", false}
private

Definition at line 190 of file SiCombinatorialTrackFinder_xk.h.

190{this, "ITkGeometry", false};

◆ m_minPt2Cut

float InDet::SiCombinatorialTrackFinder_xk::m_minPt2Cut =0
private

Definition at line 195 of file SiCombinatorialTrackFinder_xk.h.

◆ m_minPtCut

FloatProperty InDet::SiCombinatorialTrackFinder_xk::m_minPtCut {this, "MinFinalPtCut", 100, "Cut on the pt of the final track. Must be >0 to avoid NANs when computing eta."}
private

Definition at line 194 of file SiCombinatorialTrackFinder_xk.h.

194{this, "MinFinalPtCut", 100, "Cut on the pt of the final track. Must be >0 to avoid NANs when computing eta."};

◆ m_outputlevel

int InDet::SiCombinatorialTrackFinder_xk::m_outputlevel {0}
private

Definition at line 201 of file SiCombinatorialTrackFinder_xk.h.

201{0};

◆ m_pixcontainerkey

SG::ReadHandleKey<InDet::PixelClusterContainer> InDet::SiCombinatorialTrackFinder_xk::m_pixcontainerkey
private
Initial value:
{this, "PixelClusterContainer",
"PixelClusters"}

Definition at line 162 of file SiCombinatorialTrackFinder_xk.h.

162 {this, "PixelClusterContainer",
163 "PixelClusters"};

◆ m_pixelCondSummaryTool

ToolHandle<IInDetConditionsTool> InDet::SiCombinatorialTrackFinder_xk::m_pixelCondSummaryTool
private
Initial value:
{this, "PixelSummaryTool",
"PixelConditionsSummaryTool"}

Definition at line 145 of file SiCombinatorialTrackFinder_xk.h.

145 {this, "PixelSummaryTool",
146 "PixelConditionsSummaryTool"};

◆ m_pixelDetElStatus

SG::ReadHandleKey<InDet::SiDetectorElementStatus> InDet::SiCombinatorialTrackFinder_xk::m_pixelDetElStatus {this, "PixelDetElStatus", "", "Key of SiDetectorElementStatus for Pixel"}
private

Optional read handle to get status data to test whether a pixel detector element is good.

If set to e.g. PixelDetectorElementStatus the event data will be used instead of the pixel conditions summary tool.

Definition at line 177 of file SiCombinatorialTrackFinder_xk.h.

178{this, "PixelDetElStatus", "", "Key of SiDetectorElementStatus for Pixel"};

◆ m_proptool

PublicToolHandle<Trk::IPatternParametersPropagator> InDet::SiCombinatorialTrackFinder_xk::m_proptool
private
Initial value:
{this, "PropagatorTool",
"Trk::RungeKuttaPropagator/InDetPropagator"}

Definition at line 152 of file SiCombinatorialTrackFinder_xk.h.

152 {this, "PropagatorTool",
153 "Trk::RungeKuttaPropagator/InDetPropagator"};

◆ m_qualityCut

DoubleProperty InDet::SiCombinatorialTrackFinder_xk::m_qualityCut {this, "TrackQualityCut", 9.3, "Simple track quality cut"}
private

Definition at line 193 of file SiCombinatorialTrackFinder_xk.h.

193{this, "TrackQualityCut", 9.3, "Simple track quality cut"};

◆ m_riocreator

PublicToolHandle<Trk::IRIO_OnTrackCreator> InDet::SiCombinatorialTrackFinder_xk::m_riocreator
private
Initial value:
{this, "RIOonTrackTool",
"Trk::RIO_OnTrackCreator/RIO_OnTrackCreator"}

Definition at line 156 of file SiCombinatorialTrackFinder_xk.h.

156 {this, "RIOonTrackTool",
157 "Trk::RIO_OnTrackCreator/RIO_OnTrackCreator"};

◆ m_sctCondSummaryTool

ToolHandle<IInDetConditionsTool> InDet::SiCombinatorialTrackFinder_xk::m_sctCondSummaryTool
private
Initial value:
{this, "SctSummaryTool",
"InDetSCT_ConditionsSummaryTool/SCT_ConditionsSummaryTool", "Tool to retrieve SCT Conditions summary"}

Definition at line 147 of file SiCombinatorialTrackFinder_xk.h.

147 {this, "SctSummaryTool",
148 "InDetSCT_ConditionsSummaryTool/SCT_ConditionsSummaryTool", "Tool to retrieve SCT Conditions summary"};

◆ m_sctcontainerkey

SG::ReadHandleKey<InDet::SCT_ClusterContainer> InDet::SiCombinatorialTrackFinder_xk::m_sctcontainerkey
private
Initial value:
{this, "SCT_ClusterContainer",
"SCT_Clusters"}

Definition at line 164 of file SiCombinatorialTrackFinder_xk.h.

164 {this, "SCT_ClusterContainer",
165 "SCT_Clusters"};

◆ m_sctDetElStatus

SG::ReadHandleKey<InDet::SiDetectorElementStatus> InDet::SiCombinatorialTrackFinder_xk::m_sctDetElStatus {this, "SCTDetElStatus", "", "Key of SiDetectorElementStatus for SCT"}
private

Optional read handle to get status data to test whether a SCT detector element is good.

If set to e.g. SCTDetectorElementStatus the event data will be used instead of the SCT conditions summary tool.

Definition at line 183 of file SiCombinatorialTrackFinder_xk.h.

184{this, "SCTDetElStatus", "", "Key of SiDetectorElementStatus for SCT"};

◆ m_updatortool

PublicToolHandle<Trk::IPatternParametersUpdator> InDet::SiCombinatorialTrackFinder_xk::m_updatortool
private
Initial value:
{this, "UpdatorTool",
"Trk::KalmanUpdator_xk/InDetPatternUpdator"}

Definition at line 154 of file SiCombinatorialTrackFinder_xk.h.

154 {this, "UpdatorTool",
155 "Trk::KalmanUpdator_xk/InDetPatternUpdator"};

◆ m_usePIX

BooleanProperty InDet::SiCombinatorialTrackFinder_xk::m_usePIX {this, "usePixel", true}
private

Definition at line 188 of file SiCombinatorialTrackFinder_xk.h.

188{this, "usePixel", true};

◆ m_useSCT

BooleanProperty InDet::SiCombinatorialTrackFinder_xk::m_useSCT {this, "useSCT", true}
private

Definition at line 189 of file SiCombinatorialTrackFinder_xk.h.

189{this, "useSCT", true};

◆ m_writeHolesFromPattern

BooleanProperty InDet::SiCombinatorialTrackFinder_xk::m_writeHolesFromPattern {this, "writeHolesFromPattern", false,"Flag to activate writing hole info from the pattern recognition"}
private

Definition at line 196 of file SiCombinatorialTrackFinder_xk.h.

196{this, "writeHolesFromPattern", false,"Flag to activate writing hole info from the pattern recognition"};

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