]>
Commit | Line | Data |
---|---|---|
734d2c12 | 1 | #ifndef ALIPROTONANALYSIS_H |
2 | #define ALIPROTONANALYSIS_H | |
3 | ||
4 | /* See cxx source for full Copyright notice */ | |
5 | ||
6 | ||
7 | /* $Id$ */ | |
8 | ||
9 | //------------------------------------------------------------------------- | |
10 | // Class AliProtonAnalysis | |
11 | // This is the class for the baryon (proton) analysis | |
12 | // | |
dd3fa486 | 13 | // Origin: Panos Christakoglou | Panos.Christakoglou@cern.ch |
734d2c12 | 14 | //------------------------------------------------------------------------- |
15 | ||
251e4034 | 16 | #include "TObject.h" |
3f6d0c08 | 17 | #include "TH1I.h" |
734d2c12 | 18 | |
aafecd8b | 19 | class TF1; |
251e4034 | 20 | class TH2D; |
21 | class TH1F; | |
9cd594db | 22 | class TList; |
3f6d0c08 | 23 | |
9cd594db | 24 | #include "AliPID.h" |
24421eb6 | 25 | #include "AliCFContainer.h" |
251e4034 | 26 | class AliCFDataGrid; |
ee4ca40d | 27 | class AliAODEvent; |
28 | class AliAODtrack; | |
734d2c12 | 29 | class AliESDEvent; |
30 | class AliESDtrack; | |
2b748670 | 31 | class AliExternalTrackParam; |
e4358d7f | 32 | class AliStack; |
6667f3a7 | 33 | class AliESDVertex; |
734d2c12 | 34 | |
35 | class AliProtonAnalysis : public TObject { | |
36 | public: | |
37 | AliProtonAnalysis(); | |
38 | AliProtonAnalysis(Int_t nbinsY, Float_t fLowY, Float_t fHighY, | |
39 | Int_t nbinsPt, Float_t fLowPt, Float_t fHighPt); | |
40 | virtual ~AliProtonAnalysis(); | |
2b748670 | 41 | |
dd3fa486 | 42 | void SetEtaMode() {fAnalysisEtaMode = kTRUE;} |
43 | ||
2b748670 | 44 | void UseTPCOnly() {fUseTPCOnly = kTRUE;} |
da169921 | 45 | void UseHybridTPC() {fUseHybridTPC = kTRUE;} |
734d2c12 | 46 | |
ef1a8dbd | 47 | void InitAnalysisHistograms(Int_t nbinsY, Float_t fLowY, Float_t fHighY, |
48 | Int_t nbinsPt, Float_t fLowPt, Float_t fHighPt); | |
2b748670 | 49 | Bool_t ReadFromFile(const char* filename); |
6667f3a7 | 50 | void Analyze(AliESDEvent *fESD, |
51 | const AliESDVertex *vertex); | |
ef1a8dbd | 52 | void Analyze(AliAODEvent *fAOD); |
3e23254a | 53 | void Analyze(AliStack *stack, Bool_t iInclusive); |
734d2c12 | 54 | |
9cd594db | 55 | AliCFContainer *GetProtonContainer() const {return fProtonContainer;} |
56 | AliCFContainer *GetAntiProtonContainer() const {return fAntiProtonContainer;} | |
251e4034 | 57 | |
9cd594db | 58 | TH2D *GetProtonYPtHistogram() const {return fHistYPtProtons;} |
59 | TH2D *GetAntiProtonYPtHistogram() const {return fHistYPtAntiProtons;} | |
734d2c12 | 60 | TH1D *GetProtonYHistogram(); |
61 | TH1D *GetAntiProtonYHistogram(); | |
62 | TH1D *GetProtonPtHistogram(); | |
63 | TH1D *GetAntiProtonPtHistogram(); | |
251e4034 | 64 | TH1D *GetProtonCorrectedYHistogram(); |
65 | TH1D *GetAntiProtonCorrectedYHistogram(); | |
66 | TH1D *GetProtonCorrectedPtHistogram(); | |
67 | TH1D *GetAntiProtonCorrectedPtHistogram(); | |
24421eb6 | 68 | |
734d2c12 | 69 | TH1D *GetYRatioHistogram(); |
3e23254a | 70 | TH1D *GetYRatioCorrectedHistogram(TH2D *gCorrectionMapProtons, |
71 | TH2D *gCorrectionMapAntiProtons); | |
734d2c12 | 72 | TH1D *GetPtRatioHistogram(); |
3e23254a | 73 | TH1D *GetPtRatioCorrectedHistogram(TH2D *gCorrectionMapProtons, |
74 | TH2D *gCorrectionMapAntiProtons); | |
734d2c12 | 75 | TH1D *GetYAsymmetryHistogram(); |
76 | TH1D *GetPtAsymmetryHistogram(); | |
3f6d0c08 | 77 | |
9cd594db | 78 | TH1I *GetEventHistogram() const {return fHistEvents;} |
41beb956 | 79 | |
3f6d0c08 | 80 | Int_t GetNumberOfAnalyzedEvents() {return (Int_t)fHistEvents->GetEntries();} |
81 | Bool_t PrintMean(TH1 *hist, Double_t edge); | |
82 | Bool_t PrintYields(TH1 *hist, Double_t edge); | |
83 | ||
734d2c12 | 84 | //Cut functions |
0008a5a6 | 85 | void SetPointOnITSLayer1() {fPointOnITSLayer1Flag = kTRUE;} |
86 | void SetPointOnITSLayer2() {fPointOnITSLayer2Flag = kTRUE;} | |
87 | void SetPointOnITSLayer3() {fPointOnITSLayer3Flag = kTRUE;} | |
88 | void SetPointOnITSLayer4() {fPointOnITSLayer4Flag = kTRUE;} | |
89 | void SetPointOnITSLayer5() {fPointOnITSLayer5Flag = kTRUE;} | |
90 | void SetPointOnITSLayer6() {fPointOnITSLayer6Flag = kTRUE;} | |
734d2c12 | 91 | void SetMinITSClusters(Int_t minITSClusters) { |
92 | fMinITSClusters = minITSClusters; | |
93 | fMinITSClustersFlag = kTRUE; | |
94 | } | |
734d2c12 | 95 | void SetMaxChi2PerITSCluster(Double_t maxChi2PerITSCluster) { |
96 | fMaxChi2PerITSCluster = maxChi2PerITSCluster; | |
97 | fMaxChi2PerITSClusterFlag = kTRUE; | |
98 | } | |
f24a70b0 | 99 | void SetMinTPCClusters(Int_t minTPCClusters) { |
100 | fMinTPCClusters = minTPCClusters; | |
101 | fMinTPCClustersFlag = kTRUE; | |
102 | } | |
103 | void SetMaxChi2PerTPCCluster(Double_t maxChi2PerTPCCluster) { | |
104 | fMaxChi2PerTPCCluster = maxChi2PerTPCCluster; | |
105 | fMaxChi2PerTPCClusterFlag = kTRUE; | |
106 | } | |
251e4034 | 107 | void SetMaxCov11(Double_t maxCov11) { |
108 | fMaxCov11 = maxCov11; fMaxCov11Flag = kTRUE;} | |
109 | void SetMaxCov22(Double_t maxCov22) { | |
110 | fMaxCov22 = maxCov22; fMaxCov22Flag = kTRUE;} | |
111 | void SetMaxCov33(Double_t maxCov33) { | |
112 | fMaxCov33 = maxCov33; fMaxCov33Flag = kTRUE;} | |
113 | void SetMaxCov44(Double_t maxCov44) { | |
114 | fMaxCov44 = maxCov44; fMaxCov44Flag = kTRUE;} | |
115 | void SetMaxCov55(Double_t maxCov55) { | |
116 | fMaxCov55 = maxCov55; fMaxCov55Flag = kTRUE;} | |
734d2c12 | 117 | void SetMaxSigmaToVertex(Double_t maxSigmaToVertex) { |
118 | fMaxSigmaToVertex = maxSigmaToVertex; | |
119 | fMaxSigmaToVertexFlag = kTRUE; | |
120 | } | |
f24a70b0 | 121 | void SetMaxSigmaToVertexTPC(Double_t maxSigmaToVertex) { |
122 | fMaxSigmaToVertexTPC = maxSigmaToVertex; | |
123 | fMaxSigmaToVertexTPCFlag = kTRUE; | |
124 | } | |
96f84c25 | 125 | void SetMaxDCAXY(Double_t maxDCAXY) { |
126 | fMaxDCAXY = maxDCAXY; | |
127 | fMaxDCAXYFlag = kTRUE; | |
128 | } | |
129 | void SetMaxDCAXYTPC(Double_t maxDCAXY) { | |
130 | fMaxDCAXYTPC = maxDCAXY; | |
131 | fMaxDCAXYTPCFlag = kTRUE; | |
132 | } | |
133 | void SetMaxDCAZ(Double_t maxDCAZ) { | |
134 | fMaxDCAZ = maxDCAZ; | |
135 | fMaxDCAZFlag = kTRUE; | |
136 | } | |
137 | void SetMaxDCAZTPC(Double_t maxDCAZ) { | |
138 | fMaxDCAZTPC = maxDCAZ; | |
139 | fMaxDCAZTPCFlag = kTRUE; | |
140 | } | |
141 | void SetMaxConstrainChi2(Double_t maxConstrainChi2) { | |
142 | fMaxConstrainChi2 = maxConstrainChi2; | |
143 | fMaxConstrainChi2Flag = kTRUE; | |
144 | } | |
734d2c12 | 145 | void SetITSRefit() {fITSRefitFlag = kTRUE;} |
146 | void SetTPCRefit() {fTPCRefitFlag = kTRUE;} | |
f24a70b0 | 147 | void SetESDpid() {fESDpidFlag = kTRUE;} |
148 | void SetTPCpid() {fTPCpidFlag = kTRUE;} | |
ef1a8dbd | 149 | |
734d2c12 | 150 | //Prior probabilities |
9cd594db | 151 | void SetPriorProbabilities(Double_t * const partFrac) { |
251e4034 | 152 | for(Int_t i = 0; i < AliPID::kSPECIESN; i++) fPartFrac[i] = partFrac[i];} |
9cd594db | 153 | void SetPriorProbabilityFunctions(TF1 *const felectron, |
154 | TF1 *const fmuon, | |
155 | TF1 *const fpion, | |
156 | TF1 *const fkaon, | |
157 | TF1 *const fproton) { | |
aafecd8b | 158 | fFunctionProbabilityFlag = kTRUE; |
159 | fElectronFunction = felectron; | |
160 | fMuonFunction = fmuon; | |
161 | fPionFunction = fpion; | |
162 | fKaonFunction = fkaon; | |
163 | fProtonFunction = fproton; | |
164 | } | |
165 | Double_t GetParticleFraction(Int_t i, Double_t p); | |
166 | ||
39f2a708 | 167 | //interface to the correction framework |
251e4034 | 168 | void Correct(Int_t step); |
39f2a708 | 169 | Bool_t ReadCorrectionContainer(const char* filename); |
9cd594db | 170 | TList *GetCorrectionListProtons2D() const {return fCorrectionListProtons2D;} |
171 | TList *GetEfficiencyListProtons1D() const {return fEfficiencyListProtons1D;} | |
172 | TList *GetCorrectionListProtons1D() const {return fCorrectionListProtons1D;} | |
173 | TList *GetCorrectionListAntiProtons2D() const {return fCorrectionListAntiProtons2D;} | |
174 | TList *GetEfficiencyListAntiProtons1D() const {return fEfficiencyListAntiProtons1D;} | |
175 | TList *GetCorrectionListAntiProtons1D() const {return fCorrectionListAntiProtons1D;} | |
24421eb6 | 176 | |
177 | //iStep=0->MC - iStep=1->Acceptance - iStep=2->Reconstruction - iStep=3->PID | |
178 | TH1D *GetUncorrectedProtonYHistogram(Int_t iStep) {return fProtonContainer->ShowProjection(0, iStep);} | |
179 | TH1D *GetUncorrectedProtonPtHistogram(Int_t iStep) {return fProtonContainer->ShowProjection(1, iStep);} | |
180 | TH1D *GetUncorrectedAntiProtonYHistogram(Int_t iStep) {return fAntiProtonContainer->ShowProjection(0, iStep);} | |
181 | TH1D *GetUncorrectedAntiProtonPtHistogram(Int_t iStep) {return fAntiProtonContainer->ShowProjection(1, iStep);} | |
39f2a708 | 182 | |
734d2c12 | 183 | private: |
8b8b0b7a | 184 | AliProtonAnalysis(const AliProtonAnalysis&); // Not implemented |
185 | AliProtonAnalysis& operator=(const AliProtonAnalysis&); // Not implemented | |
186 | ||
6667f3a7 | 187 | Bool_t IsAccepted(AliESDEvent *esd, |
188 | const AliESDVertex *vertex, | |
189 | AliESDtrack *track); | |
3e23254a | 190 | Bool_t IsInPhaseSpace(AliESDtrack *track); |
191 | ||
9cd594db | 192 | Float_t GetSigmaToVertex(AliESDtrack* esdTrack) const; |
193 | Double_t Rapidity(Double_t Px, Double_t Py, Double_t Pz) const; | |
734d2c12 | 194 | |
dd3fa486 | 195 | Bool_t fAnalysisEtaMode; //run the analysis in eta or y |
196 | ||
734d2c12 | 197 | Int_t fNBinsY; //number of bins in y |
198 | Float_t fMinY, fMaxY; //min & max value of y | |
199 | Int_t fNBinsPt; //number of bins in pT | |
200 | Float_t fMinPt, fMaxPt; //min & max value of pT | |
201 | ||
202 | //cuts | |
203 | Int_t fMinTPCClusters, fMinITSClusters; //min TPC & ITS clusters | |
204 | Double_t fMaxChi2PerTPCCluster, fMaxChi2PerITSCluster; //max chi2 per TPC & ITS cluster | |
205 | Double_t fMaxCov11, fMaxCov22, fMaxCov33, fMaxCov44, fMaxCov55; //max values of cov. matrix | |
206 | Double_t fMaxSigmaToVertex; //max sigma to vertex cut | |
f24a70b0 | 207 | Double_t fMaxSigmaToVertexTPC; //max sigma to vertex cut |
96f84c25 | 208 | Double_t fMaxDCAXY, fMaxDCAXYTPC; //max DCA xy |
209 | Double_t fMaxDCAZ, fMaxDCAZTPC; //max DCA z | |
210 | Double_t fMaxConstrainChi2; //max constrain chi2 - vertex | |
734d2c12 | 211 | Bool_t fMinTPCClustersFlag, fMinITSClustersFlag; //shows if this cut is used or not |
212 | Bool_t fMaxChi2PerTPCClusterFlag, fMaxChi2PerITSClusterFlag; //shows if this cut is used or not | |
213 | Bool_t fMaxCov11Flag, fMaxCov22Flag, fMaxCov33Flag, fMaxCov44Flag, fMaxCov55Flag; //shows if this cut is used or not | |
214 | Bool_t fMaxSigmaToVertexFlag; //shows if this cut is used or not | |
f24a70b0 | 215 | Bool_t fMaxSigmaToVertexTPCFlag; //shows if this cut is used or not |
96f84c25 | 216 | Bool_t fMaxDCAXYFlag, fMaxDCAXYTPCFlag; //shows if this cut is used or not |
217 | Bool_t fMaxDCAZFlag, fMaxDCAZTPCFlag; //shows if this cut is used or not | |
218 | Bool_t fMaxConstrainChi2Flag; //shows if this cut is used or not | |
734d2c12 | 219 | Bool_t fITSRefitFlag, fTPCRefitFlag; //shows if this cut is used or not |
f24a70b0 | 220 | Bool_t fESDpidFlag, fTPCpidFlag; //shows if this cut is used or not |
0008a5a6 | 221 | Bool_t fPointOnITSLayer1Flag, fPointOnITSLayer2Flag; //shows if this cut is used or not |
222 | Bool_t fPointOnITSLayer3Flag, fPointOnITSLayer4Flag; //shows if this cut is used or not | |
223 | Bool_t fPointOnITSLayer5Flag, fPointOnITSLayer6Flag; //shows if this cut is used or not | |
734d2c12 | 224 | |
225 | //pid | |
aafecd8b | 226 | Bool_t fFunctionProbabilityFlag; //flag: kTRUE if functions used |
ee4ca40d | 227 | Double_t fPartFrac[10]; //prior probabilities |
aafecd8b | 228 | TF1 *fElectronFunction; //momentum dependence of the prior probs |
229 | TF1 *fMuonFunction; //momentum dependence of the prior probs | |
230 | TF1 *fPionFunction; //momentum dependence of the prior probs | |
231 | TF1 *fKaonFunction; //momentum dependence of the prior probs | |
232 | TF1 *fProtonFunction; //momentum dependence of the prior probs | |
233 | ||
2b748670 | 234 | //Detectors |
235 | Bool_t fUseTPCOnly; //kTRUE if TPC only information is used | |
da169921 | 236 | Bool_t fUseHybridTPC; //kTRUE if TPC info is used for momentum - PID and ITS for vertex & points |
2b748670 | 237 | |
251e4034 | 238 | //Analysis containers |
239 | AliCFContainer *fProtonContainer; //container for protons | |
240 | AliCFContainer *fAntiProtonContainer; //container for antiprotons | |
3f6d0c08 | 241 | TH1I *fHistEvents; //event counter |
251e4034 | 242 | TH2D *fHistYPtProtons; //Y-Pt of Protons |
243 | TH2D *fHistYPtAntiProtons; // Y-Pt of Antiprotons | |
39f2a708 | 244 | |
245 | //Corrections | |
251e4034 | 246 | TList *fEffGridListProtons; //list for the efficiency grid - protons |
cdb3530f | 247 | TList *fCorrectionListProtons2D; //list for the 2d corrections |
248 | TList *fEfficiencyListProtons1D; //list for the 1d efficiencies | |
249 | TList *fCorrectionListProtons1D; //list for the 1d corrections | |
251e4034 | 250 | TList *fEffGridListAntiProtons; //list for the efficiency grid - antiprotons |
cdb3530f | 251 | TList *fCorrectionListAntiProtons2D; //list for the 2d corrections |
252 | TList *fEfficiencyListAntiProtons1D; //list for the 1d efficiencies | |
253 | TList *fCorrectionListAntiProtons1D; //list for the 1d corrections | |
251e4034 | 254 | AliCFDataGrid *fCorrectProtons; //corrected data grid for protons |
255 | AliCFDataGrid *fCorrectAntiProtons; //corrected data grid for antiprotons | |
256 | ||
734d2c12 | 257 | ClassDef(AliProtonAnalysis,0); |
258 | }; | |
259 | ||
260 | #endif |