]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/esdTrackCuts/AliESDtrackCuts.h
1c758d625526deb150be328d42692e5827dd02ec
[u/mrichter/AliRoot.git] / PWG0 / esdTrackCuts / AliESDtrackCuts.h
1 //
2 //  Class for handling of ESD track cuts.
3 //
4 //  The class manages a number of track quality cuts, a
5 //  track-to-vertex cut and a number of kinematic cuts. Two methods
6 //  can be used to figure out if an ESD track survives the cuts:
7 //  AcceptTrack which takes a single AliESDtrack as argument and
8 //  returns kTRUE/kFALSE or GetAcceptedTracks which takes an AliESD
9 //  object and returns an TObjArray (of AliESDtracks) with the tracks
10 //  in the ESD that survived the cuts.
11 //
12 //
13 //  TODO: 
14 //  - add functionality to save and load cuts
15 //  - add different ways to make track to vertex cut
16 //  - add histograms for kinematic cut variables?
17 //  - upper and lower cuts for all (non-boolean) cuts
18 //  - update print method
19 //  - is there a smarter way to manage the cuts?
20 //  - put comments to each variable
21 //  - implement destructor !!!
22 //
23
24 #ifndef ALIESDTRACKCUTS_H
25 #define ALIESDTRACKCUTS_H
26
27 #include <TObject.h>
28 #include <TH2.h>
29
30 class AliESD;
31 class AliESDtrack;
32 class AliLog;
33
34 class AliESDtrackCuts : public TObject 
35 {
36
37 public:
38   AliESDtrackCuts();
39   virtual ~AliESDtrackCuts();
40
41   Bool_t AcceptTrack(AliESDtrack* esdTrack);
42   TObjArray* GetAcceptedTracks(AliESD* esd);
43   Int_t CountAcceptedTracks(AliESD* esd);
44
45   virtual Long64_t Merge(TCollection* list);
46   virtual void Copy(TObject &c) const;
47
48   //######################################################
49   // track quality cut setters  
50   void SetMinNClustersTPC(Int_t min=-1)          {fCutMinNClusterTPC=min;}
51   void SetMinNClustersITS(Int_t min=-1)          {fCutMinNClusterITS=min;}
52   void SetMaxChi2PerClusterTPC(Float_t max=1e99) {fCutMaxChi2PerClusterTPC=max;}
53   void SetMaxChi2PerClusterITS(Float_t max=1e99) {fCutMaxChi2PerClusterITS=max;}
54   void SetRequireTPCRefit(Bool_t b=kFALSE)       {fCutRequireTPCRefit=b;}
55   void SetRequireITSRefit(Bool_t b=kFALSE)       {fCutRequireITSRefit=b;}
56   void SetAcceptKingDaughters(Bool_t b=kFALSE)   {fCutAcceptKinkDaughters=b;}
57   void SetMaxCovDiagonalElements(Float_t c1=1e99, Float_t c2=1e99, Float_t c3=1e99, Float_t c4=1e99, Float_t c5=1e99) 
58     {fCutMaxC11=c1; fCutMaxC22=c2; fCutMaxC33=c3; fCutMaxC44=c4; fCutMaxC55=c5;}
59   
60   // track to vertex cut setters
61   void SetMinNsigmaToVertex(Float_t sigma=1e99)       {fCutNsigmaToVertex = sigma;}
62   void SetRequireSigmaToVertex(Bool_t b=kTRUE )       {fCutSigmaToVertexRequired = b;}
63   
64   // track kinmatic cut setters  
65   void SetPRange(Float_t r1=0, Float_t r2=1e99)       {fPMin=r1;   fPMax=r2;}
66   void SetPtRange(Float_t r1=0, Float_t r2=1e99)      {fPtMin=r1;  fPtMax=r2;}
67   void SetPxRange(Float_t r1=-1e99, Float_t r2=1e99)  {fPxMin=r1;  fPxMax=r2;}
68   void SetPyRange(Float_t r1=-1e99, Float_t r2=1e99)  {fPyMin=r1;  fPyMax=r2;}
69   void SetPzRange(Float_t r1=-1e99, Float_t r2=1e99)  {fPzMin=r1;  fPzMax=r2;}
70   void SetEtaRange(Float_t r1=-1e99, Float_t r2=1e99) {fEtaMin=r1; fEtaMax=r2;}
71   void SetRapRange(Float_t r1=-1e99, Float_t r2=1e99) {fRapMin=r1; fRapMax=r2;}
72
73   Float_t GetMinNsigmaToVertex() { return fCutNsigmaToVertex; } 
74
75   //######################################################
76   void SetHistogramsOn(Bool_t b=kFALSE) {fHistogramsOn = b;}
77   void DefineHistograms(Int_t color=1);
78   void SaveHistograms(Char_t* dir="track_selection");
79
80   Float_t GetSigmaToVertex(AliESDtrack* esdTrack);
81   
82   virtual void Print(const Option_t* = "") const;
83
84   // void SaveQualityCuts(Char_t* file)
85   // void LoadQualityCuts(Char_t* file)
86
87 protected:
88   void Init(); // sets everything to 0
89
90   enum { kNCuts = 21 };
91
92   //######################################################
93   // esd track quality cuts
94   static const Char_t* fgkCutNames[kNCuts]; //! names of cuts (for internal use)
95
96   Int_t   fCutMinNClusterTPC;         // min number of tpc clusters
97   Int_t   fCutMinNClusterITS;         // min number of its clusters  
98
99   Float_t fCutMaxChi2PerClusterTPC;   // max tpc fit chi2 per tpc cluster
100   Float_t fCutMaxChi2PerClusterITS;   // max its fit chi2 per its cluster
101
102   Float_t fCutMaxC11;                 // max cov. matrix diag. elements (res. y^2)
103   Float_t fCutMaxC22;                 // max cov. matrix diag. elements (res. z^2)
104   Float_t fCutMaxC33;                 // max cov. matrix diag. elements (res. sin(phi)^2)
105   Float_t fCutMaxC44;                 // max cov. matrix diag. elements (res. tan(theta_dip)^2)
106   Float_t fCutMaxC55;                 // max cov. matrix diag. elements (res. 1/pt^2)
107
108   Bool_t  fCutAcceptKinkDaughters;    // accepting kink daughters?
109   Bool_t  fCutRequireTPCRefit;        // require TPC refit
110   Bool_t  fCutRequireITSRefit;        // require ITS refit
111   
112   // track to vertex cut
113   Float_t fCutNsigmaToVertex;         // max number of estimated sigma from track-to-vertex
114   Bool_t  fCutSigmaToVertexRequired;  // cut track if sigma from track-to-vertex could not be calculated
115
116   // esd kinematics cuts
117   Float_t fPMin,   fPMax;             // definition of the range of the P
118   Float_t fPtMin,  fPtMax;            // definition of the range of the Pt
119   Float_t fPxMin,  fPxMax;            // definition of the range of the Px
120   Float_t fPyMin,  fPyMax;            // definition of the range of the Py
121   Float_t fPzMin,  fPzMax;            // definition of the range of the Pz
122   Float_t fEtaMin, fEtaMax;           // definition of the range of the eta
123   Float_t fRapMin, fRapMax;           // definition of the range of the y
124
125   //######################################################
126   // diagnostics histograms
127   Bool_t fHistogramsOn;               // histograms on/off
128
129   TH1F* fhNClustersITS[2];            //->
130   TH1F* fhNClustersTPC[2];            //->
131
132   TH1F* fhChi2PerClusterITS[2];       //->
133   TH1F* fhChi2PerClusterTPC[2];       //->
134
135   TH1F* fhC11[2];                     //->
136   TH1F* fhC22[2];                     //->
137   TH1F* fhC33[2];                     //->
138   TH1F* fhC44[2];                     //->
139   TH1F* fhC55[2];                     //->
140
141   TH1F* fhDXY[2];                     //->
142   TH1F* fhDZ[2];                      //->
143   TH2F* fhDXYvsDZ[2];                 //->
144
145   TH1F* fhDXYNormalized[2];           //->
146   TH1F* fhDZNormalized[2];            //->
147   TH2F* fhDXYvsDZNormalized[2];       //->
148
149   TH1F*  fhCutStatistics;             //-> statistics of what cuts the tracks did not survive
150   TH2F*  fhCutCorrelation;            //-> 2d statistics plot
151
152  private:
153
154   AliESDtrackCuts(const AliESDtrackCuts& pd);  // Copy Constructor
155   AliESDtrackCuts &operator=(const AliESDtrackCuts &c);
156   
157   
158   ClassDef(AliESDtrackCuts, 1)
159 };
160
161
162 #endif