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