]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliESDtrackCuts.h
Make sure that before calling CreateOutputObjects for every task, the current directo...
[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 <TString.h>
25
26 #include "AliAnalysisCuts.h"
27
28 class AliESDEvent;
29 class AliESDtrack;
30 class AliLog;
31 class TTree;
32 class TH1;
33 class TH1F;
34 class TH2F;
35 class TF1;
36 class TCollection;
37 class TFormula;
38
39 class AliESDtrackCuts : public AliAnalysisCuts
40 {
41 public:
42   enum ITSClusterRequirement { kOff = 0, kNone, kAny, kFirst, kOnlyFirst, kSecond, kOnlySecond, kBoth };
43   enum Detector { kSPD = 0, kSDD, kSSD };
44   
45   AliESDtrackCuts(const Char_t* name = "AliESDtrackCuts", const Char_t* title = "");
46   virtual ~AliESDtrackCuts();
47
48   virtual Bool_t IsSelected(TObject* obj)
49        {return AcceptTrack((AliESDtrack*)obj);}
50   virtual Bool_t IsSelected(TList* /*list*/) {return kTRUE;}
51
52   Bool_t AcceptTrack(AliESDtrack* esdTrack);
53   TObjArray* GetAcceptedTracks(AliESDEvent* esd, Bool_t bTPC = kFALSE);
54   Int_t CountAcceptedTracks(AliESDEvent* esd);
55   
56   static Int_t GetReferenceMultiplicity(AliESDEvent* esd, Bool_t tpcOnly);
57
58   static AliESDtrack* GetTPCOnlyTrack(AliESDEvent* esd, Int_t iTrack);
59   
60   // Standard cut definitions
61   static AliESDtrackCuts* GetStandardTPCOnlyTrackCuts();
62   static AliESDtrackCuts* GetStandardITSTPCTrackCuts2009(Bool_t selPrimaries=kTRUE);
63
64
65   virtual Long64_t Merge(TCollection* list);
66   virtual void Copy(TObject &c) const;
67   AliESDtrackCuts(const AliESDtrackCuts& pd);  // Copy Constructor
68   AliESDtrackCuts &operator=(const AliESDtrackCuts &c);
69
70   //######################################################
71   // track quality cut setters  
72   void SetMinNClustersTPC(Int_t min=-1)          {fCutMinNClusterTPC=min;}
73   void SetMinNClustersITS(Int_t min=-1)          {fCutMinNClusterITS=min;}
74   void SetClusterRequirementITS(Detector det, ITSClusterRequirement req = kOff) { fCutClusterRequirementITS[det] = req; }
75   void SetMaxChi2PerClusterTPC(Float_t max=1e10) {fCutMaxChi2PerClusterTPC=max;}
76   void SetMaxChi2PerClusterITS(Float_t max=1e10) {fCutMaxChi2PerClusterITS=max;}
77   void SetRequireTPCRefit(Bool_t b=kFALSE)       {fCutRequireTPCRefit=b;}
78   void SetRequireTPCStandAlone(Bool_t b=kFALSE)  {fCutRequireTPCStandAlone=b;}
79   void SetRequireITSRefit(Bool_t b=kFALSE)       {fCutRequireITSRefit=b;}
80   void SetRequireITSStandAlone(Bool_t b,Bool_t rejectITSpureSA=kFALSE) {fCutRequireITSStandAlone = b; fCutRejectITSpureSA=rejectITSpureSA;}
81   void SetAcceptKinkDaughters(Bool_t b=kTRUE)    {fCutAcceptKinkDaughters=b;}
82   void SetAcceptSharedTPCClusters(Bool_t b=kTRUE){fCutAcceptSharedTPCClusters=b;}
83   void SetMaxFractionSharedTPCClusters(Float_t max=1e10) {fCutMaxFractionSharedTPCClusters=max;}
84   void SetMaxCovDiagonalElements(Float_t c1=1e10, Float_t c2=1e10, Float_t c3=1e10, Float_t c4=1e10, Float_t c5=1e10) 
85     {fCutMaxC11=c1; fCutMaxC22=c2; fCutMaxC33=c3; fCutMaxC44=c4; fCutMaxC55=c5;}
86   void SetMaxRel1PtUncertainty(Float_t max=1e10)      {fCutMaxRel1PtUncertainty=max;}
87
88   // track to vertex cut setters
89   void SetMaxNsigmaToVertex(Float_t sigma=1e10)       {fCutNsigmaToVertex = sigma; SetRequireSigmaToVertex(kTRUE);}
90   void SetRequireSigmaToVertex(Bool_t b=kTRUE)        {fCutSigmaToVertexRequired = b;}
91   void SetMaxDCAToVertexXY(Float_t dist=1e10)         {fCutMaxDCAToVertexXY = dist;}
92   void SetMaxDCAToVertexZ(Float_t dist=1e10)          {fCutMaxDCAToVertexZ = dist;}
93   void SetMinDCAToVertexXY(Float_t dist=0.)           {fCutMinDCAToVertexXY = dist;}
94   void SetMinDCAToVertexZ(Float_t dist=0.)            {fCutMinDCAToVertexZ = dist;}
95   void SetMaxDCAToVertexXYPtDep(const char *dist="");
96   void SetMaxDCAToVertexZPtDep(const char *dist=""); 
97   void SetMinDCAToVertexXYPtDep(const char *dist="");
98   void SetMinDCAToVertexZPtDep(const char *dist=""); 
99   void SetDCAToVertex2D(Bool_t b=kFALSE)              {fCutDCAToVertex2D = b;}
100
101
102   // getters
103
104   Int_t   GetMinNClusterTPC()        const   { return fCutMinNClusterTPC;}
105   Int_t   GetMinNClustersITS()       const   { return fCutMinNClusterITS;}
106   ITSClusterRequirement GetClusterRequirementITS(Detector det) const { return fCutClusterRequirementITS[det]; }
107   Float_t GetMaxChi2PerClusterTPC()  const   { return fCutMaxChi2PerClusterTPC;}
108   Float_t GetMaxChi2PerClusterITS()  const   { return fCutMaxChi2PerClusterITS;}
109   Bool_t  GetRequireTPCRefit()       const   { return fCutRequireTPCRefit;}
110   Bool_t  GetRequireTPCStandAlone()  const   { return fCutRequireTPCStandAlone;}
111   Bool_t  GetRequireITSRefit()       const   { return fCutRequireITSRefit;}
112   Bool_t  GetRequireITSStandAlone()  const   { return fCutRequireITSStandAlone; }
113   Bool_t  GetAcceptKinkDaughters()   const   { return fCutAcceptKinkDaughters;}
114   Bool_t  GetAcceptSharedTPCClusters()        const   {return fCutAcceptSharedTPCClusters;}
115   Float_t GetMaxFractionSharedTPCClusters()   const   {return fCutMaxFractionSharedTPCClusters;}
116   void    GetMaxCovDiagonalElements(Float_t& c1, Float_t& c2, Float_t& c3, Float_t& c4, Float_t& c5) 
117       {c1 = fCutMaxC11; c2 = fCutMaxC22; c3 = fCutMaxC33; c4 = fCutMaxC44; c5 = fCutMaxC55;}
118   Float_t GetMaxRel1PtUncertainty()  const   { return fCutMaxRel1PtUncertainty;}
119   Float_t GetMaxNsigmaToVertex()     const   { return fCutNsigmaToVertex;}
120   Float_t GetMaxDCAToVertexXY()      const   { return fCutMaxDCAToVertexXY;}
121   Float_t GetMaxDCAToVertexZ()       const   { return fCutMaxDCAToVertexZ;}
122   Float_t GetMinDCAToVertexXY()      const   { return fCutMinDCAToVertexXY;}
123   Float_t GetMinDCAToVertexZ()       const   { return fCutMinDCAToVertexZ;}
124   const char* GetMaxDCAToVertexXYPtDep() const   { return fCutMaxDCAToVertexXYPtDep;}
125   const char* GetMaxDCAToVertexZPtDep()  const   { return fCutMaxDCAToVertexZPtDep;}
126   const char* GetMinDCAToVertexXYPtDep() const   { return fCutMinDCAToVertexXYPtDep;}
127   const char* GetMinDCAToVertexZPtDep()  const   { return fCutMinDCAToVertexZPtDep;}
128   Bool_t  GetDCAToVertex2D()         const   { return fCutDCAToVertex2D;}
129   Bool_t  GetRequireSigmaToVertex( ) const   { return fCutSigmaToVertexRequired;}
130
131   void GetPRange(Float_t& r1, Float_t& r2)  const {r1=fPMin;   r2=fPMax;}
132   void GetPtRange(Float_t& r1, Float_t& r2) const {r1=fPtMin;  r2=fPtMax;}
133   void GetPxRange(Float_t& r1, Float_t& r2) const {r1=fPxMin;  r2=fPxMax;}
134   void GetPyRange(Float_t& r1, Float_t& r2) const {r1=fPyMin;  r2=fPyMax;}
135   void GetPzRange(Float_t& r1, Float_t& r2) const {r1=fPzMin;  r2=fPzMax;}
136   void GetEtaRange(Float_t& r1, Float_t& r2) const {r1=fEtaMin; r2=fEtaMax;}
137   void GetRapRange(Float_t& r1, Float_t& r2) const {r1=fRapMin; r2=fRapMax;}
138
139   // track kinmatic cut setters
140   void SetPRange(Float_t r1=0, Float_t r2=1e10)       {fPMin=r1;   fPMax=r2;}
141   void SetPtRange(Float_t r1=0, Float_t r2=1e10)      {fPtMin=r1;  fPtMax=r2;}
142   void SetPxRange(Float_t r1=-1e10, Float_t r2=1e10)  {fPxMin=r1;  fPxMax=r2;}
143   void SetPyRange(Float_t r1=-1e10, Float_t r2=1e10)  {fPyMin=r1;  fPyMax=r2;}
144   void SetPzRange(Float_t r1=-1e10, Float_t r2=1e10)  {fPzMin=r1;  fPzMax=r2;}
145   void SetEtaRange(Float_t r1=-1e10, Float_t r2=1e10) {fEtaMin=r1; fEtaMax=r2;}
146   void SetRapRange(Float_t r1=-1e10, Float_t r2=1e10) {fRapMin=r1; fRapMax=r2;}
147
148   //######################################################
149   void SetHistogramsOn(Bool_t b=kFALSE) {fHistogramsOn = b;}
150   void DefineHistograms(Int_t color=1);
151   virtual Bool_t LoadHistograms(const Char_t* dir = 0);
152   void SaveHistograms(const Char_t* dir = 0);
153   void DrawHistograms();
154
155   static Float_t GetSigmaToVertex(AliESDtrack* esdTrack);
156   
157   static void EnableNeededBranches(TTree* tree);
158
159   // void SaveQualityCuts(Char_t* file)
160   // void LoadQualityCuts(Char_t* file)
161
162         TH1F* GetDZNormalized(Int_t i) const { return fhDZNormalized[i]; }
163
164 protected:
165   void Init(); // sets everything to 0
166   Bool_t CheckITSClusterRequirement(ITSClusterRequirement req, Bool_t clusterL1, Bool_t clusterL2);
167   Bool_t CheckPtDepDCA(TString dist,Bool_t print=kFALSE) const;
168   void SetPtDepDCACuts(Double_t pt);
169
170   enum { kNCuts = 35 }; 
171
172   //######################################################
173   // esd track quality cuts
174   static const Char_t* fgkCutNames[kNCuts]; //! names of cuts (for internal use)
175
176   Int_t   fCutMinNClusterTPC;         // min number of tpc clusters
177   Int_t   fCutMinNClusterITS;         // min number of its clusters
178   
179   ITSClusterRequirement fCutClusterRequirementITS[3];  // detailed ITS cluster requirements for (SPD, SDD, SSD)
180
181   Float_t fCutMaxChi2PerClusterTPC;   // max tpc fit chi2 per tpc cluster
182   Float_t fCutMaxChi2PerClusterITS;   // max its fit chi2 per its cluster
183
184   Float_t fCutMaxC11;                 // max cov. matrix diag. elements (res. y^2)
185   Float_t fCutMaxC22;                 // max cov. matrix diag. elements (res. z^2)
186   Float_t fCutMaxC33;                 // max cov. matrix diag. elements (res. sin(phi)^2)
187   Float_t fCutMaxC44;                 // max cov. matrix diag. elements (res. tan(theta_dip)^2)
188   Float_t fCutMaxC55;                 // max cov. matrix diag. elements (res. 1/pt^2)
189
190   Float_t fCutMaxRel1PtUncertainty;   // max relative uncertainty of 1/pt
191
192   Bool_t  fCutAcceptKinkDaughters;    // accepting kink daughters?
193   Bool_t  fCutAcceptSharedTPCClusters;// accepting shared clusters in TPC?
194   Float_t fCutMaxFractionSharedTPCClusters; //Maximum fraction of shared clusters in TPC
195   Bool_t  fCutRequireTPCRefit;        // require TPC refit
196   Bool_t  fCutRequireTPCStandAlone;   // require TPC standalone tracks
197   Bool_t  fCutRequireITSRefit;        // require ITS refit
198   Bool_t  fCutRequireITSStandAlone;   // require ITS standalone tracks
199   Bool_t  fCutRejectITSpureSA;        // reject  ITS standalone tracks found using all ITS clusters
200
201   // track to vertex cut
202   Float_t fCutNsigmaToVertex;         // max number of estimated sigma from track-to-vertex
203   Bool_t  fCutSigmaToVertexRequired;  // cut track if sigma from track-to-vertex could not be calculated
204   Float_t fCutMaxDCAToVertexXY;       // track-to-vertex cut in max absolute distance in xy-plane
205   Float_t fCutMaxDCAToVertexZ;        // track-to-vertex cut in max absolute distance in z-plane
206   Float_t fCutMinDCAToVertexXY;       // track-to-vertex cut on min absolute distance in xy-plane
207   Float_t fCutMinDCAToVertexZ;        // track-to-vertex cut on min absolute distance in z-plane
208   // 
209   TString fCutMaxDCAToVertexXYPtDep;  // pt-dep track-to-vertex cut in max absolute distance in xy-plane
210   TString fCutMaxDCAToVertexZPtDep;   // pt-dep track-to-vertex cut in max absolute distance in z-plane
211   TString fCutMinDCAToVertexXYPtDep;  // pt-dep track-to-vertex cut on min absolute distance in xy-plane
212   TString fCutMinDCAToVertexZPtDep;   // pt-dep track-to-vertex cut on min absolute distance in z-plane
213
214   // only internal use, set via strings above
215   TFormula *f1CutMaxDCAToVertexXYPtDep;  // pt-dep track-to-vertex cut in max absolute distance in xy-plane
216   TFormula *f1CutMaxDCAToVertexZPtDep;   // pt-dep track-to-vertex cut in max absolute distance in z-plane
217   TFormula *f1CutMinDCAToVertexXYPtDep;  // pt-dep track-to-vertex cut on min absolute distance in xy-plane
218   TFormula *f1CutMinDCAToVertexZPtDep;   // pt-dep track-to-vertex cut on min absolute distance in z-plane
219
220   Bool_t  fCutDCAToVertex2D;          // if true a 2D DCA cut is made. Tracks are accepted if sqrt((DCAXY / fCutMaxDCAToVertexXY)^2 + (DCAZ / fCutMaxDCAToVertexZ)^2) < 1 AND sqrt((DCAXY / fCutMinDCAToVertexXY)^2 + (DCAZ / fCutMinDCAToVertexZ)^2) > 1
221
222   // esd kinematics cuts
223   Float_t fPMin,   fPMax;             // definition of the range of the P
224   Float_t fPtMin,  fPtMax;            // definition of the range of the Pt
225   Float_t fPxMin,  fPxMax;            // definition of the range of the Px
226   Float_t fPyMin,  fPyMax;            // definition of the range of the Py
227   Float_t fPzMin,  fPzMax;            // definition of the range of the Pz
228   Float_t fEtaMin, fEtaMax;           // definition of the range of the eta
229   Float_t fRapMin, fRapMax;           // definition of the range of the y
230
231   //######################################################
232   // diagnostics histograms
233   Bool_t fHistogramsOn;               // histograms on/off
234
235   TH1F* fhNClustersITS[2];            //->
236   TH1F* fhNClustersTPC[2];            //->
237
238   TH1F* fhChi2PerClusterITS[2];       //->
239   TH1F* fhChi2PerClusterTPC[2];       //->
240
241   TH1F* fhC11[2];                     //->
242   TH1F* fhC22[2];                     //->
243   TH1F* fhC33[2];                     //->
244   TH1F* fhC44[2];                     //->
245   TH1F* fhC55[2];                     //->
246
247   TH1F* fhRel1PtUncertainty[2];       //-> rel. uncertainty of 1/pt
248
249   TH1F* fhDXY[2];                     //->
250   TH1F* fhDZ[2];                      //->
251   TH1F* fhDXYDZ[2];                   //-> absolute distance sqrt(dxy**2 + dz**2) to vertex; if 2D cut is set, normalized to given values
252   TH2F* fhDXYvsDZ[2];                 //->
253
254   TH1F* fhDXYNormalized[2];           //->
255   TH1F* fhDZNormalized[2];            //->
256   TH2F* fhDXYvsDZNormalized[2];       //->
257   TH1F* fhNSigmaToVertex[2];          //->
258
259   TH1F* fhPt[2];                      //-> pt of esd tracks
260   TH1F* fhEta[2];                     //-> eta of esd tracks
261
262   TF1*  ffDTheoretical;               //-> theoretical distance to vertex normalized (2d gauss)
263
264   TH1F*  fhCutStatistics;             //-> statistics of what cuts the tracks did not survive
265   TH2F*  fhCutCorrelation;            //-> 2d statistics plot
266
267   ClassDef(AliESDtrackCuts, 11)
268 };
269
270
271 #endif