]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HLT/ITS/trackingSAP/AliITSSAPTracker.h
first version of HLT ITS SA Primaries tracker/vertexer
[u/mrichter/AliRoot.git] / HLT / ITS / trackingSAP / AliITSSAPTracker.h
CommitLineData
c9598e1c 1#ifndef ALIITSSAPTRACKER_H
2#define ALIITSSAPTRACKER_H
3
4#include <TH1.h>
5#include <TH2.h>
6#include <TBits.h>
7#include <TStopwatch.h>
8#include <TObjArray.h>
9#include <algorithm>
10#include <vector>
11#include "AliExternalTrackParam.h"
12#include "AliESDVertex.h"
13
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//------------------------------------------------------------------------
19
20class AliITSRecPoint;
21class AliESDVertex;
22class AliITSSAPLayer;
23
24class AliITSSAPTracker : public TObject
25{
26 public:
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)};
31 //
32 struct SPDtracklet {
33 int id1;
34 int id2;
35 int label;
36 float dphi;
37 float dtht;
38 float chi2;
39 };
40 typedef struct SPDtracklet SPDtracklet_t;
41 //
42 struct ITStrack {
43 AliExternalTrackParam paramOut;
44 AliExternalTrackParam paramInw;
45 float chi2;
46 short ncl;
47 short nmiss;
48 int clID[6];
49 int label;
50 int trackletID;
51 };
52 typedef struct ITStrack ITStrack_t;
53 //
54 AliITSSAPTracker();
55 virtual ~AliITSSAPTracker();
56 //
57 void ProcessEvent();
58 void Init();
59 void Clear(Option_t *opt="");
60 void ClearTracklets();
61 void ClearTracks() {fTracks.clear();}
62 //
63 void SetSPDVertex(const AliESDVertex* v) {fSPDVertex = v;}
64 void AddCluster(AliITSRecPoint* cl);
65 void SetBz(float v) {fBz = v;}
66 //
67 // methods for trackleting ---------------->>>
68 Bool_t FindTracklets();
69 Int_t AssociateClusterOfL2(int icl2);
70 Bool_t IsBlacklisted(int id1,int id2);
71 void Blacklist(int id1,int id2);
72 //
73 void SetPhiShift(float v=0.0045) {fPhiShift = v;}
74 void SetSigThetaTracklet(float v=0.025) {fSigThetaTracklet = v;}
75 void SetSigPhiTracklet(float v=0.08) {fSigPhiTracklet = v;}
76 void SetChi2CutTracklet(float v=1.5) {fChi2CutTracklet = v;}
77 //
78 Double_t GetClSystYErr2(Int_t lr) const {return fgkClSystYErr2[lr];}
79 Double_t GetClSystZErr2(Int_t lr) const {return fgkClSystZErr2[lr];}
80 //
81 int GetNTracklets() const {return (int)fTracklets.size();}
82 int GetNTracks() const {return fNTracks;}
83 void PrintTracklets() const;
84 void PrintTracklet(Int_t itr) const;
85 // methods for trackleting ----------------<<<
86 //
87 // methods for track reconstruction ------->>>
88 Float_t GetMinPt() const {return fMinPt;}
89 void SetMinPt(Float_t v=0.3) {fMinPt = v<0.2 ? 0.2 : v;}
90 void CalcAuxTracking();
91 Bool_t CreateTrack(AliITSSAPTracker::ITStrack_t& track, AliITSSAPTracker::SPDtracklet_t& trlet);
92 void Tracklets2Tracks();
93 AliITSSAPLayer* GetLayer(int i) const {return (AliITSSAPLayer*)fLayers[i];}
94 Int_t GetActiveLayerID(int i) const {return fgkLr2Active[i];}
95 Float_t GetChi2TotCut(int ncl) const;
96 Bool_t CrossPassiveLayer(AliExternalTrackParam& track, Int_t lrID);
97 Bool_t FollowToLayer(AliITSSAPTracker::ITStrack_t& track, Int_t lrID);
98 Double_t GetXatLabRLin(AliExternalTrackParam& track, double r);
99 void CookLabel(AliITSSAPTracker::ITStrack_t& track);
100 void CookLabel(AliITSSAPTracker::SPDtracklet_t& tracklet);
101 void PrintTrack(const AliITSSAPTracker::ITStrack_t& track) const;
102 Bool_t IsObligatoryLayer(int lr) const {return !fSkipLayer[lr];}
103 Bool_t IsAcceptableTrack(const AliITSSAPTracker::ITStrack_t& track) const;
104 void PrintTracks() const;
105 Int_t GetTrackletMCTruth(AliITSSAPTracker::SPDtracklet_t& trlet) const;
106 void RefitInward();
107 Bool_t RefitInward(int itr);
108 AliITSSAPTracker::ITStrack_t& GetTrack(int i) const {return (ITStrack_t &)fTracks[i];}
109 // methods for track reconstruction -------<<<
110 //
111 // methods for vertex reconstruction ------>>>
112 Bool_t FitTrackVertex();
113 AliESDVertex& GetTrackVertex() const {return (AliESDVertex&)fTrackVertex;}
114 // methods for vertex reconstruction ------<<<
115 //
116 protected:
117 //
118 AliITSSAPLayer* fLayers[kNLrActive];
119 Bool_t fSkipLayer[kNLrActive]; //! skip layer
120 Int_t fNClusters[kNLrActive]; //! number of clusters per event
121 //
122 // for SPD trackleting ----------------- >>>
123 std::vector<bool> fSPD2Discard; //! status of SPD2 clusters during trackleting
124 std::vector<SPDtracklet_t> fTracklets; //! found tracklets
125 std::vector<int> fSPD1Tracklet; //! id+1 of traclet using this SPD1 cluster
126 TBits* fBlacklist; //! blacklisted combinations
127 Float_t fPhiShift; //! Phi shift reference value (at 0.5 T)
128 Float_t fSigThetaTracklet; //! sigTheta for tracklets
129 Float_t fSigPhiTracklet; //! sigPhi for tracklets
130 Float_t fChi2CutTracklet; //! cut on tracklet total chi2
131 Float_t fPhiShiftSc; //! dPhi offset to account for bending
132 Float_t fDThetaTrackletSc; //! max dTheta for tracklets with scaling from chi2 cut
133 Float_t fDPhiTrackletSc; //! max dPhi for tracklets with scaling from chi2 cut
134 Float_t fBz; //! Bz field in ITS
135 //
136 // auxilary precomputed stuff
137 Float_t fDPhiTol; //! tolerance on phi, accounting for bending
138 Float_t fDThSig2Inv; //! inverse of sigma dTheta
139 Float_t fDPhSig2Inv; //! inverse of sigma dPhi
140 // for SPD trackleting ----------------- <<<
141 //
142 // for track reconstruction ------------ >>>
143 Float_t fMinPt; //! user pt cutoff
144 Float_t fCurvMax; //! max curvature to reconstruct
145 Float_t fZSPD2CutMin; //! min Z of tracklet SPD2 to consider tracking
146 Float_t fZSPD2CutMax; //! max Z of tracklet SPD2 to consider tracking
147 Float_t fMaxChi2Tr2Cl; //! cut on cluster-to-track chi2
148 Float_t fAddErr2YspdVtx; //! additional error to Y of the SPD vertex in track fit
149 Float_t fAddErr2ZspdVtx; //! additional error to Z of the SPD vertex in track fit
150 Float_t fChi2TotCut[kNLrActive]; //! cut on total chi2 depending on track length
151 //
152 Float_t fNSigma2[kNLrActive]; //! N^2 sigma margin for cluster search
153 Float_t fYToler2[kNLrActive]; //! Y additional margin^2 for cluster search
154 Float_t fZToler2[kNLrActive]; //! Z additional margin^2 for cluster search
155
156 Float_t fMSDist[kNLrActive]; //! shift due to the MS for 1 GeV particle
157 Float_t fMSPhi[kNLrActive]; //! dphi due to the MS for 1 GeV particle
158 Float_t fTolPhiCrude[kNLrActive]; //! tolerance in dphi for particle of unknown momentum
159 Float_t fTolZCrude[kNLrActive]; //! tolerance in Z for particle of unknown momentum
160 Float_t fMissChi2Penalty; //! penalize missed clusters
161 Int_t fMaxMissedLayers; //! allow to miss at most this number of layers
162 Int_t fNTracks; //! n found tracks
163 std::vector<ITStrack_t> fTracks; //! found tracks container
164 AliESDVertex fTrackVertex; //! fitted track vertex
165 Bool_t fFitVertex; //! fit vertex with tracks
166 // for track reconstruction ------------ <<<
167 //
168
169 const AliESDVertex* fSPDVertex; //! external vertex
170
171 static const Float_t fgkRhoLITS[kMaxLrITS]; // <rho*L> for each material layer
172 static const Float_t fgkX2X0ITS[kMaxLrITS]; // <x/x0> for each material layer
173 static const Float_t fgkRLayITS[kMaxLrITS]; // radii of material layers
174 static const Float_t fgkRSpanITS[kMaxLrITS]; // half R span of the material layer
175 static const Float_t fgkZSpanITS[kMaxLrITS]; // half Z span of the material layer
176 static const Int_t fgkPassivLrITS[kNLrPassive]; // list of passive layer enums
177 static const Int_t fgkActiveLrITS[kNLrActive]; // list of active layer enums
178 static const Double_t fgkClSystYErr2[kNLrActive]; // syst error^2 for Y direction
179 static const Double_t fgkClSystZErr2[kNLrActive]; // syst error^2 for Y direction
180
181 static const Int_t fgkLr2Active[kMaxLrITS]; // conversion from LrID to ActiveLr ID
182 static const Int_t fgkLrDefBins[kNLrActive][2]; // default binning for cluster navigator
183 static const Int_t fgkDummyLabel; // dummy MC label
184 static const Float_t fgkDefMass; // default mass for tracking
185 //
186#ifdef _TIMING_
187 public:
188 enum {kSWTotal,kSWTracklets,kSWTracks,kSWVertex,kNSW};
189 void PrintTiming();
190 const TStopwatch& GetStopwatch(int i) const {return fSW[i];}
191 const char* GetStopwatchName(int i) const {return fgkSWNames[i];}
192 protected:
193 static const char* fgkSWNames[kNSW];
194 TStopwatch fSW[kNSW];
195#endif
196 //
197#ifdef _CONTROLH_
198 protected:
199 TObjArray fArrHisto;
200 TH2F *fHTrackletMC,*fHTrackletAll,*fHTrackletFake,*fHTrackMC,*fHTrackAll,*fHTrackFake,*fHVtxDiffXY
201 ,*fHVtxDiffXMlt,*fHVtxDiffYMlt,*fHVtxDiffZMlt
202 ,*fHVtxPullXMlt,*fHVtxPullYMlt,*fHVtxPullZMlt
203 ,*fHVtxMCSPDDiffXY
204 ,*fHVtxMCSPDDiffXMlt,*fHVtxMCSPDDiffYMlt,*fHVtxMCSPDDiffZMlt
205 ,*fHVtxMCSPDPullXY
206 ,*fHVtxMCSPDPullXMlt,*fHVtxMCSPDPullYMlt,*fHVtxMCSPDPullZMlt
207 ,*fHChi2NDFvsPT,*fHChi2vsNC;
208 TH1F *fHVtxMltRef,*fHVtxOKMlt,*fHVtxDiffZ,*fHVtxMCSPDDiffZ;
209 //
210 void FillRecoStat();
211 void FillTrackingControlHistos(int lrID,int lbl,const AliExternalTrackParam* bestTr,
212 const double cpar[2],const double ccov[3],const AliITSRecPoint* bestCl);
213 void BookHistos();
214 Double_t* DefLogAx(double xMn,double xMx, int nbin);
215 public:
216 void SaveHistos(const char* outFName="XXXITSTrackerControlH.root");
217 TObjArray* GetHistos() const {return (TObjArray*)&fArrHisto;}
218 enum {kHResidY,kHPullY,kHResidZ,kHPullZ,kHChi2Cl};
219#endif
220
221 private:
222 AliITSSAPTracker(const AliITSSAPTracker&);
223 AliITSSAPTracker& operator=(const AliITSSAPTracker&);
224 ClassDef(AliITSSAPTracker,0) // ITS SA primaries tracker/vertexer
225};
226
227
228//______________________________________________
229inline Bool_t AliITSSAPTracker::IsBlacklisted(int id1,int id2)
230{
231 // check if this combination is blacklisted
232 return fBlacklist->TestBitNumber(UInt_t(id1*fNClusters[0])+id2);
233}
234
235//______________________________________________
236inline void AliITSSAPTracker::Blacklist(int id1,int id2)
237{
238 // blacklist this combination
239 return fBlacklist->SetBitNumber(UInt_t(id1*fNClusters[0])+id2);
240}
241
242//______________________________________________
243inline Float_t AliITSSAPTracker::GetChi2TotCut(int ncl) const
244{
245 // return chi2 cut for given number of clusters. Min ncl=3
246 return fChi2TotCut[ncl-2];
247}
248
249#endif