]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ANALYSIS/AliESDtrackCuts.h
Moved from PWG0
[u/mrichter/AliRoot.git] / ANALYSIS / AliESDtrackCuts.h
CommitLineData
73318471 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//
22
23#ifndef ALIESDTRACKCUTS_H
24#define ALIESDTRACKCUTS_H
25
26#include <TF1.h>
27#include <TH2.h>
28#include "AliAnalysisCuts.h"
29
30class AliESD;
31class AliESDEvent;
32class AliESDtrack;
33class AliLog;
34class TTree;
35
36class AliESDtrackCuts : public AliAnalysisCuts
37{
38public:
39 AliESDtrackCuts(const Char_t* name = "AliESDtrackCuts", const Char_t* title = "");
40 virtual ~AliESDtrackCuts();
41 Bool_t IsSelected(TObject* obj)
42 {return AcceptTrack((AliESDtrack*)obj);}
43 Bool_t AcceptTrack(AliESDtrack* esdTrack);
44 TObjArray* GetAcceptedTracks(AliESD* esd);
45 Int_t CountAcceptedTracks(AliESD* esd);
46 TObjArray* GetAcceptedTracks(AliESDEvent* esd);
47 Int_t CountAcceptedTracks(AliESDEvent* esd);
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
70 // getters
71 Float_t GetMinNsigmaToVertex() { return fCutNsigmaToVertex;}
72 Bool_t GetRequireSigmaToVertex( ) { return fCutSigmaToVertexRequired;}
73
74 // track kinmatic cut setters
75 void SetPRange(Float_t r1=0, Float_t r2=1e10) {fPMin=r1; fPMax=r2;}
76 void SetPtRange(Float_t r1=0, Float_t r2=1e10) {fPtMin=r1; fPtMax=r2;}
77 void SetPxRange(Float_t r1=-1e10, Float_t r2=1e10) {fPxMin=r1; fPxMax=r2;}
78 void SetPyRange(Float_t r1=-1e10, Float_t r2=1e10) {fPyMin=r1; fPyMax=r2;}
79 void SetPzRange(Float_t r1=-1e10, Float_t r2=1e10) {fPzMin=r1; fPzMax=r2;}
80 void SetEtaRange(Float_t r1=-1e10, Float_t r2=1e10) {fEtaMin=r1; fEtaMax=r2;}
81 void SetRapRange(Float_t r1=-1e10, Float_t r2=1e10) {fRapMin=r1; fRapMax=r2;}
82
83 //######################################################
84 void SetHistogramsOn(Bool_t b=kFALSE) {fHistogramsOn = b;}
85 void DefineHistograms(Int_t color=1);
86 virtual Bool_t LoadHistograms(const Char_t* dir = 0);
87 void SaveHistograms(const Char_t* dir = 0);
88 void DrawHistograms();
89
90 Float_t GetSigmaToVertex(AliESDtrack* esdTrack);
91
92 static void EnableNeededBranches(TTree* tree);
93
94 // void SaveQualityCuts(Char_t* file)
95 // void LoadQualityCuts(Char_t* file)
96
97 TH1* GetDZNormalized(Int_t i) const { return fhDZNormalized[i]; }
98
99protected:
100 void Init(); // sets everything to 0
101
102 enum { kNCuts = 21 };
103
104 //######################################################
105 // esd track quality cuts
106 static const Char_t* fgkCutNames[kNCuts]; //! names of cuts (for internal use)
107
108 Int_t fCutMinNClusterTPC; // min number of tpc clusters
109 Int_t fCutMinNClusterITS; // min number of its clusters
110
111 Float_t fCutMaxChi2PerClusterTPC; // max tpc fit chi2 per tpc cluster
112 Float_t fCutMaxChi2PerClusterITS; // max its fit chi2 per its cluster
113
114 Float_t fCutMaxC11; // max cov. matrix diag. elements (res. y^2)
115 Float_t fCutMaxC22; // max cov. matrix diag. elements (res. z^2)
116 Float_t fCutMaxC33; // max cov. matrix diag. elements (res. sin(phi)^2)
117 Float_t fCutMaxC44; // max cov. matrix diag. elements (res. tan(theta_dip)^2)
118 Float_t fCutMaxC55; // max cov. matrix diag. elements (res. 1/pt^2)
119
120 Bool_t fCutAcceptKinkDaughters; // accepting kink daughters?
121 Bool_t fCutRequireTPCRefit; // require TPC refit
122 Bool_t fCutRequireITSRefit; // require ITS refit
123
124 // track to vertex cut
125 Float_t fCutNsigmaToVertex; // max number of estimated sigma from track-to-vertex
126 Bool_t fCutSigmaToVertexRequired; // cut track if sigma from track-to-vertex could not be calculated
127
128 // esd kinematics cuts
129 Float_t fPMin, fPMax; // definition of the range of the P
130 Float_t fPtMin, fPtMax; // definition of the range of the Pt
131 Float_t fPxMin, fPxMax; // definition of the range of the Px
132 Float_t fPyMin, fPyMax; // definition of the range of the Py
133 Float_t fPzMin, fPzMax; // definition of the range of the Pz
134 Float_t fEtaMin, fEtaMax; // definition of the range of the eta
135 Float_t fRapMin, fRapMax; // definition of the range of the y
136
137 //######################################################
138 // diagnostics histograms
139 Bool_t fHistogramsOn; // histograms on/off
140
141 TH1F* fhNClustersITS[2]; //->
142 TH1F* fhNClustersTPC[2]; //->
143
144 TH1F* fhChi2PerClusterITS[2]; //->
145 TH1F* fhChi2PerClusterTPC[2]; //->
146
147 TH1F* fhC11[2]; //->
148 TH1F* fhC22[2]; //->
149 TH1F* fhC33[2]; //->
150 TH1F* fhC44[2]; //->
151 TH1F* fhC55[2]; //->
152
153 TH1F* fhDXY[2]; //->
154 TH1F* fhDZ[2]; //->
155 TH2F* fhDXYvsDZ[2]; //->
156
157 TH1F* fhDXYNormalized[2]; //->
158 TH1F* fhDZNormalized[2]; //->
159 TH2F* fhDXYvsDZNormalized[2]; //->
160 TH1F* fhNSigmaToVertex[2]; //->
161
162 TH1F* fhPt[2]; //-> pt of esd tracks
163 TH1F* fhEta[2]; //-> eta of esd tracks
164
165 TF1* ffDTheoretical; //-> theoretical distance to vertex normalized (2d gauss)
166
167 TH1F* fhCutStatistics; //-> statistics of what cuts the tracks did not survive
168 TH2F* fhCutCorrelation; //-> 2d statistics plot
169
170 ClassDef(AliESDtrackCuts, 2)
171};
172
173
174#endif