]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/AliITSUCATracker.h
142e68766993dc3f1fd11af40e73c4d2a19d0382
[u/mrichter/AliRoot.git] / ITS / UPGRADE / AliITSUCATracker.h
1 //-*- Mode: C++ -*-
2 // ************************************************************************
3 // This file is property of and copyright by the ALICE ITSU Project       *
4 // ALICE Experiment at CERN, All rights reserved.                         *
5 // See cxx source for full Copyright notice                               *
6 //                                                                        *
7 //*************************************************************************
8
9 #ifndef ALIITSUCATRACKER_H
10 #define ALIITSUCATRACKER_H
11
12 #define _TUNING_
13
14 #include <vector>
15
16 #include <TClonesArray.h>
17 #include "AliITSUCACell.h"
18 #include "AliITSUCATrackingStation.h"
19 typedef struct AliITSUCATrackingStation::ClsInfo ClsInfo_t;
20
21 #include "AliITSUTrackerGlo.h"
22
23 class AliITSUReconstructor;
24 class AliITSURecoDet;
25 class AliITSURecoLayer;
26 class AliITSUTrackCooked;
27
28 class TTree;
29 class AliCluster;
30 class AliESDEvent;
31
32 #ifdef _TUNING_
33 #include <TH1F.h>
34 #endif
35
36 using std::vector;
37
38 //__________________________________________________________________________________________________
39 class Doublets {
40 public:
41   Doublets(int xx = 0, int yy = 0, float tL = 0.f, float ph = 0.f)
42   : x((unsigned short)xx)
43   , y((unsigned short)yy)
44   , tanL(tL)
45   , phi(ph) {}
46   unsigned short x,y;
47   float tanL, phi;
48 };
49
50 //__________________________________________________________________________________________________
51 class AliITSUCATracker : public AliITSUTrackerGlo {
52 public:
53   AliITSUCATracker(AliITSUReconstructor* rec=0);
54   virtual ~AliITSUCATracker();
55
56   // These functions must be implemented
57   Int_t Clusters2Tracks(AliESDEvent *event);
58   Int_t PropagateBack(AliESDEvent *event);
59   Int_t RefitInward(AliESDEvent *event);
60   Int_t LoadClusters(TTree *ct);
61   void UnloadClusters();
62   AliCluster *GetCluster(Int_t index) const;
63
64   // Possibly, other public functions
65   Double_t GetMaterialBudget(const double* p0, const double* p1, double& x2x0, double& rhol) const;
66   Bool_t   GetSAonly() const { return fSAonly; }
67   void     SetChi2Cut(float cut) { fChi2Cut = cut; }
68   void     SetPhiCut(float cut) { fPhiCut = cut; }
69   void     SetSAonly(Bool_t sa=kTRUE) { fSAonly=sa; }
70   void     SetZCut(float cut) { fZCut = cut; }
71
72 #ifdef _TUNING_
73   bool                            fGood;
74   TH1F *                          fGoodCombChi2[5];
75   TH1F *                          fFakeCombChi2[5];
76   TH1F *                          fGoodCombN[4];
77   TH1F *                          fFakeCombN[4];
78   TH1F *                          fGDZ[6];
79   TH1F *                          fGDXY[6];
80   TH1F *                          fFDZ[6];
81   TH1F *                          fFDXY[6];
82   TH1F *                          fGDCAZ[5];
83   TH1F *                          fGDCAXY[5];
84   TH1F *                          fFDCAZ[5];
85   TH1F *                          fFDCAXY[5];
86   TH1F *                          fTan;
87   TH1F *                          fTanF;
88   TH1F *                          fPhi;
89   TH1F *                          fPhiF;
90   TH1F *                          fNEntries;
91   void ResetHistos();
92 #endif
93   //
94 protected:
95   bool   CellParams(int l, ClsInfo_t* c1, ClsInfo_t* c2, ClsInfo_t* c3, float &curv, float np[3]);
96   void   CellsTreeTraversal(vector<AliITSUCARoad> &roads, const int &iD, const int &doubl);
97   void   FindTracksCA(int iteration);
98   void   MakeCells(int iteration);
99   Bool_t RefitAt(Double_t xx, AliITSUTrackCooked *t, const AliITSUTrackCooked *c);
100   void   SetCuts(int it);
101   void   SetLabel(AliITSUTrackCooked &t, Float_t wrong);
102   
103 private:
104   AliITSUCATracker(const AliITSUCATracker&);
105   AliITSUCATracker &operator=(const AliITSUCATracker &tr);
106
107   // Data members
108
109   // classes for interfacing the geometry, materials etc
110   // Internal tracker arrays, layers, modules, etc
111   AliITSUCATrackingStation        fLayer[7];
112   vector<bool>                    fUsedClusters[7];
113   Float_t                         fChi2Cut;
114   Float_t                         fPhiCut;
115   Float_t                         fZCut;
116   vector<Doublets>                fDoublets[6];
117   vector<AliITSUCACell>           fCells[5];
118   TClonesArray                   *fCandidates[4];
119   Bool_t                          fSAonly;             // kTRUE if the standalone tracking only
120
121   
122   // Cuts
123   float fCPhi;
124   float fCDTanL;
125   float fCDPhi;
126   float fCZ;
127   float fCDCAz[5];
128   float fCDCAxy[5];
129   float fCDN[4];
130   float fCDP[4];
131   float fCDZ[6];
132
133   static const Double_t           fgkChi2Cut;      // chi2 cut during track merging
134   static const int                fgkNumberOfIterations;
135   static const float              fgkR[7];
136   //
137   ClassDef(AliITSUCATracker,2)   //ITSU stand-alone tracker
138 };
139
140 #endif // ALIITSUCATRACKER_H