1 #ifndef ALIITSSAPTRACKER_H
2 #define ALIITSSAPTRACKER_H
7 #include <TStopwatch.h>
11 #include "AliExternalTrackParam.h"
12 #include "AliESDVertex.h"
14 //------- compilation options, comment out all for best performance ------
15 #define _TIMING_ // print timing info, use only in offline mode
16 //#define _CONTROLH_ // fill control histos, use only in offline mode
17 //#define _DEBUG_ // print debug info, use only in offline mode
18 //------------------------------------------------------------------------
24 class AliITSSAPTracker : public TObject
27 enum {kALrSPD1,kALrSPD2, kALrSDD1,kALrSDD2, kALrSSD1,kALrSSD2,kNLrActive};
28 enum {kLrBeamPime, kLrSPD1,kLrSPD2, kLrShield1, kLrSDD1,kLrSDD2, kLrShield2, kLrSSD1,kLrSSD2,
29 kMaxLrITS,kNLrPassive=kMaxLrITS-kNLrActive};
30 enum {kInvalidBit=BIT(14)};
40 typedef struct SPDtracklet SPDtracklet_t;
43 AliExternalTrackParam paramOut;
44 AliExternalTrackParam paramInw;
52 typedef struct ITStrack ITStrack_t;
55 virtual ~AliITSSAPTracker();
59 void Clear(Option_t *opt="");
60 void ClearTracklets();
61 void ClearTracks() {fTracks.clear();}
63 void SetSPDVertex(const AliESDVertex* v) {fSPDVertex = v;}
64 const AliESDVertex* GetSPDVertex() const {return fSPDVertex;}
65 void AddCluster(AliITSRecPoint* cl);
66 void SetBz(float v) {fBz = v;}
68 // methods for trackleting ---------------->>>
69 Bool_t FindTracklets();
70 Int_t AssociateClusterOfL2(int icl2);
71 Bool_t IsBlacklisted(int id1,int id2);
72 void Blacklist(int id1,int id2);
74 void SetPhiShift(float v=0.0045) {fPhiShift = v;}
75 void SetSigThetaTracklet(float v=0.025) {fSigThetaTracklet = v;}
76 void SetSigPhiTracklet(float v=0.08) {fSigPhiTracklet = v;}
77 void SetChi2CutTracklet(float v=1.5) {fChi2CutTracklet = v;}
79 Double_t GetClSystYErr2(Int_t lr) const {return fgkClSystYErr2[lr];}
80 Double_t GetClSystZErr2(Int_t lr) const {return fgkClSystZErr2[lr];}
82 int GetNTracklets() const {return (int)fTracklets.size();}
83 int GetNTracks() const {return fNTracks;}
84 void PrintTracklets() const;
85 void PrintTracklet(Int_t itr) const;
86 const AliITSSAPTracker::SPDtracklet_t& GetTracklet(int i) const {return fTracklets[i];}
87 // methods for trackleting ----------------<<<
89 // methods for track reconstruction ------->>>
90 Float_t GetMinPt() const {return fMinPt;}
91 void SetMinPt(Float_t v=0.3) {fMinPt = v<0.2 ? 0.2 : v;}
92 void CalcAuxTracking();
93 Bool_t CreateTrack(AliITSSAPTracker::ITStrack_t& track, AliITSSAPTracker::SPDtracklet_t& trlet);
94 void Tracklets2Tracks();
95 AliITSSAPLayer* GetLayer(int i) const {return (AliITSSAPLayer*)fLayers[i];}
96 Int_t GetActiveLayerID(int i) const {return fgkLr2Active[i];}
97 Float_t GetChi2TotCut(int ncl) const;
98 Bool_t CrossPassiveLayer(AliExternalTrackParam& track, Int_t lrID);
99 Bool_t FollowToLayer(AliITSSAPTracker::ITStrack_t& track, Int_t lrID);
100 Double_t GetXatLabRLin(AliExternalTrackParam& track, double r);
101 void CookLabel(AliITSSAPTracker::ITStrack_t& track);
102 void CookLabel(AliITSSAPTracker::SPDtracklet_t& tracklet);
103 void PrintTrack(const AliITSSAPTracker::ITStrack_t& track) const;
104 Bool_t IsObligatoryLayer(int lr) const {return !fSkipLayer[lr];}
105 Bool_t IsAcceptableTrack(const AliITSSAPTracker::ITStrack_t& track) const;
106 void PrintTracks() const;
107 Int_t GetTrackletMCTruth(AliITSSAPTracker::SPDtracklet_t& trlet) const;
109 Bool_t RefitInward(int itr);
110 void SetMaxMissedLayers(int n=0) { fMaxMissedLayers = n;}
111 Int_t GetMaxMissedLayers() const { return fMaxMissedLayers;}
112 AliITSSAPTracker::ITStrack_t& GetTrack(int i) const {return (ITStrack_t &)fTracks[i];}
113 // methods for track reconstruction -------<<<
115 // methods for vertex reconstruction ------>>>
116 Bool_t FitTrackVertex();
117 AliESDVertex& GetTrackVertex() const {return (AliESDVertex&)fTrackVertex;}
118 // methods for vertex reconstruction ------<<<
122 AliITSSAPLayer* fLayers[kNLrActive];
123 Bool_t fSkipLayer[kNLrActive]; //! skip layer
124 Int_t fNClusters[kNLrActive]; //! number of clusters per event
126 // for SPD trackleting ----------------- >>>
127 std::vector<bool> fSPD2Discard; //! status of SPD2 clusters during trackleting
128 std::vector<SPDtracklet_t> fTracklets; //! found tracklets
129 std::vector<int> fSPD1Tracklet; //! id+1 of traclet using this SPD1 cluster
130 TBits* fBlacklist; //! blacklisted combinations
131 Float_t fPhiShift; //! Phi shift reference value (at 0.5 T)
132 Float_t fSigThetaTracklet; //! sigTheta for tracklets
133 Float_t fSigPhiTracklet; //! sigPhi for tracklets
134 Float_t fChi2CutTracklet; //! cut on tracklet total chi2
135 Float_t fPhiShiftSc; //! dPhi offset to account for bending
136 Float_t fDThetaTrackletSc; //! max dTheta for tracklets with scaling from chi2 cut
137 Float_t fDPhiTrackletSc; //! max dPhi for tracklets with scaling from chi2 cut
138 Float_t fBz; //! Bz field in ITS
140 // auxilary precomputed stuff
141 Float_t fDPhiTol; //! tolerance on phi, accounting for bending
142 Float_t fDThSig2Inv; //! inverse of sigma dTheta
143 Float_t fDPhSig2Inv; //! inverse of sigma dPhi
144 // for SPD trackleting ----------------- <<<
146 // for track reconstruction ------------ >>>
147 Float_t fMinPt; //! user pt cutoff
148 Float_t fCurvMax; //! max curvature to reconstruct
149 Float_t fZSPD2CutMin; //! min Z of tracklet SPD2 to consider tracking
150 Float_t fZSPD2CutMax; //! max Z of tracklet SPD2 to consider tracking
151 Float_t fMaxChi2Tr2Cl; //! cut on cluster-to-track chi2
152 Float_t fAddErr2YspdVtx; //! additional error to Y of the SPD vertex in track fit
153 Float_t fAddErr2ZspdVtx; //! additional error to Z of the SPD vertex in track fit
154 Float_t fChi2TotCut[kNLrActive]; //! cut on total chi2 depending on track length
156 Float_t fNSigma2[kNLrActive]; //! N^2 sigma margin for cluster search
157 Float_t fYToler2[kNLrActive]; //! Y additional margin^2 for cluster search
158 Float_t fZToler2[kNLrActive]; //! Z additional margin^2 for cluster search
160 Float_t fMSDist[kNLrActive]; //! shift due to the MS for 1 GeV particle
161 Float_t fMSPhi[kNLrActive]; //! dphi due to the MS for 1 GeV particle
162 Float_t fTolPhiCrude[kNLrActive]; //! tolerance in dphi for particle of unknown momentum
163 Float_t fTolZCrude[kNLrActive]; //! tolerance in Z for particle of unknown momentum
164 Float_t fMissChi2Penalty; //! penalize missed clusters
165 Int_t fMaxMissedLayers; //! allow to miss at most this number of layers
166 Int_t fNTracks; //! n found tracks
167 std::vector<ITStrack_t> fTracks; //! found tracks container
168 AliESDVertex fTrackVertex; //! fitted track vertex
169 Bool_t fFitVertex; //! fit vertex with tracks
170 // for track reconstruction ------------ <<<
173 const AliESDVertex* fSPDVertex; //! external vertex
175 static const Float_t fgkRhoLITS[kMaxLrITS]; // <rho*L> for each material layer
176 static const Float_t fgkX2X0ITS[kMaxLrITS]; // <x/x0> for each material layer
177 static const Float_t fgkRLayITS[kMaxLrITS]; // radii of material layers
178 static const Float_t fgkRSpanITS[kMaxLrITS]; // half R span of the material layer
179 static const Float_t fgkZSpanITS[kMaxLrITS]; // half Z span of the material layer
180 static const Int_t fgkPassivLrITS[kNLrPassive]; // list of passive layer enums
181 static const Int_t fgkActiveLrITS[kNLrActive]; // list of active layer enums
182 static const Double_t fgkClSystYErr2[kNLrActive]; // syst error^2 for Y direction
183 static const Double_t fgkClSystZErr2[kNLrActive]; // syst error^2 for Y direction
185 static const Int_t fgkLr2Active[kMaxLrITS]; // conversion from LrID to ActiveLr ID
186 static const Int_t fgkLrDefBins[kNLrActive][2]; // default binning for cluster navigator
187 static const Int_t fgkDummyLabel; // dummy MC label
188 static const Float_t fgkDefMass; // default mass for tracking
192 enum {kSWTotal,kSWTracklets,kSWTracks,kSWVertex,kNSW};
194 const TStopwatch& GetStopwatch(int i) const {return fSW[i];}
195 const char* GetStopwatchName(int i) const {return fgkSWNames[i];}
197 static const char* fgkSWNames[kNSW];
198 TStopwatch fSW[kNSW];
204 TH2F *fHTrackletMC,*fHTrackletAll,*fHTrackletFake,*fHTrackMC,*fHTrackAll,*fHTrackFake,*fHVtxDiffXY
205 ,*fHVtxDiffXMlt,*fHVtxDiffYMlt,*fHVtxDiffZMlt
206 ,*fHVtxPullXMlt,*fHVtxPullYMlt,*fHVtxPullZMlt
208 ,*fHVtxMCSPDDiffXMlt,*fHVtxMCSPDDiffYMlt,*fHVtxMCSPDDiffZMlt
210 ,*fHVtxMCSPDPullXMlt,*fHVtxMCSPDPullYMlt,*fHVtxMCSPDPullZMlt
211 ,*fHChi2NDFvsPT,*fHChi2vsNC;
212 TH1F *fHVtxMltRef,*fHVtxOKMlt,*fHVtxDiffZ,*fHVtxMCSPDDiffZ;
215 void FillTrackingControlHistos(int lrID,int lbl,const AliExternalTrackParam* bestTr,
216 const double cpar[2],const double ccov[3],const AliITSRecPoint* bestCl);
218 Double_t* DefLogAx(double xMn,double xMx, int nbin);
220 void SaveHistos(const char* outFName="XXXITSTrackerControlH.root");
221 TObjArray* GetHistos() const {return (TObjArray*)&fArrHisto;}
222 enum {kHResidY,kHPullY,kHResidZ,kHPullZ,kHChi2Cl};
226 AliITSSAPTracker(const AliITSSAPTracker&);
227 AliITSSAPTracker& operator=(const AliITSSAPTracker&);
228 ClassDef(AliITSSAPTracker,0) // ITS SA primaries tracker/vertexer
232 //______________________________________________
233 inline Bool_t AliITSSAPTracker::IsBlacklisted(int id1,int id2)
235 // check if this combination is blacklisted
236 return fBlacklist->TestBitNumber(UInt_t(id1*fNClusters[0])+id2);
239 //______________________________________________
240 inline void AliITSSAPTracker::Blacklist(int id1,int id2)
242 // blacklist this combination
243 return fBlacklist->SetBitNumber(UInt_t(id1*fNClusters[0])+id2);
246 //______________________________________________
247 inline Float_t AliITSSAPTracker::GetChi2TotCut(int ncl) const
249 // return chi2 cut for given number of clusters. Min ncl=3
250 return fChi2TotCut[ncl-2];