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