]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliESDtrackCuts.h
Primary track quality selection for kink mothers.
[u/mrichter/AliRoot.git] / ANALYSIS / 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 AliESDEvent
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 histograms for kinematic cut variables?
16 //  - upper and lower cuts for all (non-boolean) cuts
17 //  - update print method
18 //  - put comments to each variable
19 //
20
21 #ifndef ALIESDTRACKCUTS_H
22 #define ALIESDTRACKCUTS_H
23
24 #include <TF1.h>
25 #include <TH2.h>
26 #include "AliAnalysisCuts.h"
27
28 class AliESDEvent;
29 class AliESDtrack;
30 class AliLog;
31 class TTree;
32
33 class AliESDtrackCuts : public AliAnalysisCuts
34 {
35 public:
36   AliESDtrackCuts(const Char_t* name = "AliESDtrackCuts", const Char_t* title = "");
37   virtual ~AliESDtrackCuts();
38
39   Bool_t IsSelected(TObject* obj)
40        {return AcceptTrack((AliESDtrack*)obj);}
41   Bool_t IsSelected(TList* /*list*/) {return kTRUE;}
42
43   Bool_t AcceptTrack(AliESDtrack* esdTrack);
44   TObjArray* GetAcceptedTracks(AliESDEvent* esd,Bool_t bTPC = kFALSE);
45   Int_t CountAcceptedTracks(AliESDEvent* esd);
46
47   static AliESDtrack* GetTPCOnlyTrack(AliESDEvent* esd, Int_t iTrack);
48
49   virtual Long64_t Merge(TCollection* list);
50   virtual void Copy(TObject &c) const;
51   AliESDtrackCuts(const AliESDtrackCuts& pd);  // Copy Constructor
52   AliESDtrackCuts &operator=(const AliESDtrackCuts &c);
53
54   //######################################################
55   // track quality cut setters  
56   void SetMinNClustersTPC(Int_t min=-1)          {fCutMinNClusterTPC=min;}
57   void SetMinNClustersITS(Int_t min=-1)          {fCutMinNClusterITS=min;}
58   void SetMaxChi2PerClusterTPC(Float_t max=1e10) {fCutMaxChi2PerClusterTPC=max;}
59   void SetMaxChi2PerClusterITS(Float_t max=1e10) {fCutMaxChi2PerClusterITS=max;}
60   void SetRequireTPCRefit(Bool_t b=kFALSE)       {fCutRequireTPCRefit=b;}
61   void SetRequireITSRefit(Bool_t b=kFALSE)       {fCutRequireITSRefit=b;}
62   void SetAcceptKingDaughters(Bool_t b=kFALSE)   {fCutAcceptKinkDaughters=b;}
63   void SetMaxCovDiagonalElements(Float_t c1=1e10, Float_t c2=1e10, Float_t c3=1e10, Float_t c4=1e10, Float_t c5=1e10) 
64     {fCutMaxC11=c1; fCutMaxC22=c2; fCutMaxC33=c3; fCutMaxC44=c4; fCutMaxC55=c5;}
65
66   // track to vertex cut setters
67   void SetMinNsigmaToVertex(Float_t sigma=1e10)       {fCutNsigmaToVertex = sigma;}
68   void SetRequireSigmaToVertex(Bool_t b=kTRUE )       {fCutSigmaToVertexRequired = b;}
69   void SetDCAToVertex(Float_t dist=1e10)              {fCutDCAToVertex = dist;}
70
71   // getters
72
73   Int_t   GetMinNClusterTPC()        const   { return fCutMinNClusterTPC;}
74   Int_t   GetMinNClustersITS()       const   { return fCutMinNClusterITS;}
75   Float_t GetMaxChi2PerClusterTPC()  const   { return fCutMaxChi2PerClusterTPC;}
76   Float_t GetMaxChi2PerClusterITS()  const   { return fCutMaxChi2PerClusterITS;}
77   Bool_t  GetRequireTPCRefit()       const   { return fCutRequireTPCRefit;}
78   Bool_t  GetRequireITSRefit()       const   { return fCutRequireITSRefit;}
79   Bool_t  GetAcceptKingDaughters()   const   { return fCutAcceptKinkDaughters;}
80   void    GetMaxCovDiagonalElements(Float_t& c1, Float_t& c2, Float_t& c3, Float_t& c4, Float_t& c5) 
81       {c1 = fCutMaxC11; c2 = fCutMaxC22; c3 = fCutMaxC33; c4 = fCutMaxC44; c5 = fCutMaxC55;}
82   Float_t GetMinNsigmaToVertex()     const   { return fCutNsigmaToVertex;}
83   Bool_t  GetRequireSigmaToVertex( ) const   { return fCutSigmaToVertexRequired;}
84   
85   // track kinmatic cut setters
86   void SetPRange(Float_t r1=0, Float_t r2=1e10)       {fPMin=r1;   fPMax=r2;}
87   void SetPtRange(Float_t r1=0, Float_t r2=1e10)      {fPtMin=r1;  fPtMax=r2;}
88   void SetPxRange(Float_t r1=-1e10, Float_t r2=1e10)  {fPxMin=r1;  fPxMax=r2;}
89   void SetPyRange(Float_t r1=-1e10, Float_t r2=1e10)  {fPyMin=r1;  fPyMax=r2;}
90   void SetPzRange(Float_t r1=-1e10, Float_t r2=1e10)  {fPzMin=r1;  fPzMax=r2;}
91   void SetEtaRange(Float_t r1=-1e10, Float_t r2=1e10) {fEtaMin=r1; fEtaMax=r2;}
92   void SetRapRange(Float_t r1=-1e10, Float_t r2=1e10) {fRapMin=r1; fRapMax=r2;}
93
94   //######################################################
95   void SetHistogramsOn(Bool_t b=kFALSE) {fHistogramsOn = b;}
96   void DefineHistograms(Int_t color=1);
97   virtual Bool_t LoadHistograms(const Char_t* dir = 0);
98   void SaveHistograms(const Char_t* dir = 0);
99   void DrawHistograms();
100
101   Float_t GetSigmaToVertex(AliESDtrack* esdTrack);
102   
103   static void EnableNeededBranches(TTree* tree);
104
105   // void SaveQualityCuts(Char_t* file)
106   // void LoadQualityCuts(Char_t* file)
107
108         TH1* GetDZNormalized(Int_t i) const { return fhDZNormalized[i]; }
109
110 protected:
111   void Init(); // sets everything to 0
112
113   enum { kNCuts = 22 };
114
115   //######################################################
116   // esd track quality cuts
117   static const Char_t* fgkCutNames[kNCuts]; //! names of cuts (for internal use)
118
119   Int_t   fCutMinNClusterTPC;         // min number of tpc clusters
120   Int_t   fCutMinNClusterITS;         // min number of its clusters
121
122   Float_t fCutMaxChi2PerClusterTPC;   // max tpc fit chi2 per tpc cluster
123   Float_t fCutMaxChi2PerClusterITS;   // max its fit chi2 per its cluster
124
125   Float_t fCutMaxC11;                 // max cov. matrix diag. elements (res. y^2)
126   Float_t fCutMaxC22;                 // max cov. matrix diag. elements (res. z^2)
127   Float_t fCutMaxC33;                 // max cov. matrix diag. elements (res. sin(phi)^2)
128   Float_t fCutMaxC44;                 // max cov. matrix diag. elements (res. tan(theta_dip)^2)
129   Float_t fCutMaxC55;                 // max cov. matrix diag. elements (res. 1/pt^2)
130
131   Bool_t  fCutAcceptKinkDaughters;    // accepting kink daughters?
132   Bool_t  fCutRequireTPCRefit;        // require TPC refit
133   Bool_t  fCutRequireITSRefit;        // require ITS refit
134
135   // track to vertex cut
136   Float_t fCutNsigmaToVertex;         // max number of estimated sigma from track-to-vertex
137   Bool_t  fCutSigmaToVertexRequired;  // cut track if sigma from track-to-vertex could not be calculated
138   Float_t fCutDCAToVertex;            // track-to-vertex cut in absolute distance
139
140   // esd kinematics cuts
141   Float_t fPMin,   fPMax;             // definition of the range of the P
142   Float_t fPtMin,  fPtMax;            // definition of the range of the Pt
143   Float_t fPxMin,  fPxMax;            // definition of the range of the Px
144   Float_t fPyMin,  fPyMax;            // definition of the range of the Py
145   Float_t fPzMin,  fPzMax;            // definition of the range of the Pz
146   Float_t fEtaMin, fEtaMax;           // definition of the range of the eta
147   Float_t fRapMin, fRapMax;           // definition of the range of the y
148
149   //######################################################
150   // diagnostics histograms
151   Bool_t fHistogramsOn;               // histograms on/off
152
153   TH1F* fhNClustersITS[2];            //->
154   TH1F* fhNClustersTPC[2];            //->
155
156   TH1F* fhChi2PerClusterITS[2];       //->
157   TH1F* fhChi2PerClusterTPC[2];       //->
158
159   TH1F* fhC11[2];                     //->
160   TH1F* fhC22[2];                     //->
161   TH1F* fhC33[2];                     //->
162   TH1F* fhC44[2];                     //->
163   TH1F* fhC55[2];                     //->
164
165   TH1F* fhDXY[2];                     //->
166   TH1F* fhDZ[2];                      //->
167   TH1F* fhDXYDZ[2];                   //-> absolute distance sqrt(dxy**2 + dz**2) to vertex
168   TH2F* fhDXYvsDZ[2];                 //->
169
170   TH1F* fhDXYNormalized[2];           //->
171   TH1F* fhDZNormalized[2];            //->
172   TH2F* fhDXYvsDZNormalized[2];       //->
173   TH1F* fhNSigmaToVertex[2];          //->
174
175   TH1F* fhPt[2];                      //-> pt of esd tracks
176   TH1F* fhEta[2];                     //-> eta of esd tracks
177
178   TF1*  ffDTheoretical;               //-> theoretical distance to vertex normalized (2d gauss)
179
180   TH1F*  fhCutStatistics;             //-> statistics of what cuts the tracks did not survive
181   TH2F*  fhCutCorrelation;            //-> 2d statistics plot
182
183   ClassDef(AliESDtrackCuts, 2)
184 };
185
186
187 #endif