1 #ifndef ALIPROTONANALYSIS_H
2 #define ALIPROTONANALYSIS_H
4 /* See cxx source for full Copyright notice */
9 //-------------------------------------------------------------------------
10 // Class AliProtonAnalysis
11 // This is the class for the baryon (proton) analysis
13 // Origin: Panos Christakoglou, UOA-CERN, Panos.Christakoglou@cern.ch
14 //-------------------------------------------------------------------------
26 #include "AliCFContainer.h"
32 class AliExternalTrackParam;
35 class AliProtonAnalysis : public TObject {
38 AliProtonAnalysis(Int_t nbinsY, Float_t fLowY, Float_t fHighY,
39 Int_t nbinsPt, Float_t fLowPt, Float_t fHighPt);
40 virtual ~AliProtonAnalysis();
42 void UseTPCOnly() {fUseTPCOnly = kTRUE;}
44 void InitAnalysisHistograms(Int_t nbinsY, Float_t fLowY, Float_t fHighY,
45 Int_t nbinsPt, Float_t fLowPt, Float_t fHighPt);
46 Bool_t ReadFromFile(const char* filename);
47 void Analyze(AliESDEvent *fESD);
48 void Analyze(AliAODEvent *fAOD);
49 void Analyze(AliStack *stack);
51 AliCFContainer *GetProtonContainer() {return fProtonContainer;}
52 AliCFContainer *GetAntiProtonContainer() {return fAntiProtonContainer;}
54 TH2D *GetProtonYPtHistogram() {return fHistYPtProtons;}
55 TH2D *GetAntiProtonYPtHistogram() {return fHistYPtAntiProtons;}
56 TH1D *GetProtonYHistogram();
57 TH1D *GetAntiProtonYHistogram();
58 TH1D *GetProtonPtHistogram();
59 TH1D *GetAntiProtonPtHistogram();
60 TH1D *GetProtonCorrectedYHistogram();
61 TH1D *GetAntiProtonCorrectedYHistogram();
62 TH1D *GetProtonCorrectedPtHistogram();
63 TH1D *GetAntiProtonCorrectedPtHistogram();
65 TH1D *GetYRatioHistogram();
66 TH1D *GetPtRatioHistogram();
67 TH1D *GetYAsymmetryHistogram();
68 TH1D *GetPtAsymmetryHistogram();
70 TH1I *GetEventHistogram() {return fHistEvents;}
72 Int_t GetNumberOfAnalyzedEvents() {return (Int_t)fHistEvents->GetEntries();}
73 Bool_t PrintMean(TH1 *hist, Double_t edge);
74 Bool_t PrintYields(TH1 *hist, Double_t edge);
77 void SetMinITSClusters(Int_t minITSClusters) {
78 fMinITSClusters = minITSClusters;
79 fMinITSClustersFlag = kTRUE;
81 void SetMaxChi2PerITSCluster(Double_t maxChi2PerITSCluster) {
82 fMaxChi2PerITSCluster = maxChi2PerITSCluster;
83 fMaxChi2PerITSClusterFlag = kTRUE;
85 void SetMinTPCClusters(Int_t minTPCClusters) {
86 fMinTPCClusters = minTPCClusters;
87 fMinTPCClustersFlag = kTRUE;
89 void SetMaxChi2PerTPCCluster(Double_t maxChi2PerTPCCluster) {
90 fMaxChi2PerTPCCluster = maxChi2PerTPCCluster;
91 fMaxChi2PerTPCClusterFlag = kTRUE;
93 void SetMaxCov11(Double_t maxCov11) {
94 fMaxCov11 = maxCov11; fMaxCov11Flag = kTRUE;}
95 void SetMaxCov22(Double_t maxCov22) {
96 fMaxCov22 = maxCov22; fMaxCov22Flag = kTRUE;}
97 void SetMaxCov33(Double_t maxCov33) {
98 fMaxCov33 = maxCov33; fMaxCov33Flag = kTRUE;}
99 void SetMaxCov44(Double_t maxCov44) {
100 fMaxCov44 = maxCov44; fMaxCov44Flag = kTRUE;}
101 void SetMaxCov55(Double_t maxCov55) {
102 fMaxCov55 = maxCov55; fMaxCov55Flag = kTRUE;}
103 void SetMaxSigmaToVertex(Double_t maxSigmaToVertex) {
104 fMaxSigmaToVertex = maxSigmaToVertex;
105 fMaxSigmaToVertexFlag = kTRUE;
107 void SetMaxSigmaToVertexTPC(Double_t maxSigmaToVertex) {
108 fMaxSigmaToVertexTPC = maxSigmaToVertex;
109 fMaxSigmaToVertexTPCFlag = kTRUE;
111 void SetMaxDCAXY(Double_t maxDCAXY) {
112 fMaxDCAXY = maxDCAXY;
113 fMaxDCAXYFlag = kTRUE;
115 void SetMaxDCAXYTPC(Double_t maxDCAXY) {
116 fMaxDCAXYTPC = maxDCAXY;
117 fMaxDCAXYTPCFlag = kTRUE;
119 void SetMaxDCAZ(Double_t maxDCAZ) {
121 fMaxDCAZFlag = kTRUE;
123 void SetMaxDCAZTPC(Double_t maxDCAZ) {
124 fMaxDCAZTPC = maxDCAZ;
125 fMaxDCAZTPCFlag = kTRUE;
127 void SetMaxConstrainChi2(Double_t maxConstrainChi2) {
128 fMaxConstrainChi2 = maxConstrainChi2;
129 fMaxConstrainChi2Flag = kTRUE;
131 void SetITSRefit() {fITSRefitFlag = kTRUE;}
132 void SetTPCRefit() {fTPCRefitFlag = kTRUE;}
133 void SetESDpid() {fESDpidFlag = kTRUE;}
134 void SetTPCpid() {fTPCpidFlag = kTRUE;}
136 //Prior probabilities
137 void SetPriorProbabilities(Double_t *partFrac) {
138 for(Int_t i = 0; i < AliPID::kSPECIESN; i++) fPartFrac[i] = partFrac[i];}
139 void SetPriorProbabilityFunctions(TF1 *felectron, TF1 *fmuon, TF1 *fpion, TF1 *fkaon, TF1 *fproton) {
140 fFunctionProbabilityFlag = kTRUE;
141 fElectronFunction = felectron;
142 fMuonFunction = fmuon;
143 fPionFunction = fpion;
144 fKaonFunction = fkaon;
145 fProtonFunction = fproton;
147 Double_t GetParticleFraction(Int_t i, Double_t p);
149 //interface to the correction framework
150 void Correct(Int_t step);
151 Bool_t ReadCorrectionContainer(const char* filename);
152 TList *GetCorrectionListProtons2D() {return fCorrectionListProtons2D;}
153 TList *GetEfficiencyListProtons1D() {return fEfficiencyListProtons1D;}
154 TList *GetCorrectionListProtons1D() {return fCorrectionListProtons1D;}
155 TList *GetCorrectionListAntiProtons2D() {return fCorrectionListAntiProtons2D;}
156 TList *GetEfficiencyListAntiProtons1D() {return fEfficiencyListAntiProtons1D;}
157 TList *GetCorrectionListAntiProtons1D() {return fCorrectionListAntiProtons1D;}
159 //iStep=0->MC - iStep=1->Acceptance - iStep=2->Reconstruction - iStep=3->PID
160 TH1D *GetUncorrectedProtonYHistogram(Int_t iStep) {return fProtonContainer->ShowProjection(0, iStep);}
161 TH1D *GetUncorrectedProtonPtHistogram(Int_t iStep) {return fProtonContainer->ShowProjection(1, iStep);}
162 TH1D *GetUncorrectedAntiProtonYHistogram(Int_t iStep) {return fAntiProtonContainer->ShowProjection(0, iStep);}
163 TH1D *GetUncorrectedAntiProtonPtHistogram(Int_t iStep) {return fAntiProtonContainer->ShowProjection(1, iStep);}
166 AliProtonAnalysis(const AliProtonAnalysis&); // Not implemented
167 AliProtonAnalysis& operator=(const AliProtonAnalysis&); // Not implemented
169 Bool_t IsAccepted(AliESDtrack *track);
170 Float_t GetSigmaToVertex(AliESDtrack* esdTrack);
171 Double_t Rapidity(Double_t Px, Double_t Py, Double_t Pz);
173 Int_t fNBinsY; //number of bins in y
174 Float_t fMinY, fMaxY; //min & max value of y
175 Int_t fNBinsPt; //number of bins in pT
176 Float_t fMinPt, fMaxPt; //min & max value of pT
179 Int_t fMinTPCClusters, fMinITSClusters; //min TPC & ITS clusters
180 Double_t fMaxChi2PerTPCCluster, fMaxChi2PerITSCluster; //max chi2 per TPC & ITS cluster
181 Double_t fMaxCov11, fMaxCov22, fMaxCov33, fMaxCov44, fMaxCov55; //max values of cov. matrix
182 Double_t fMaxSigmaToVertex; //max sigma to vertex cut
183 Double_t fMaxSigmaToVertexTPC; //max sigma to vertex cut
184 Double_t fMaxDCAXY, fMaxDCAXYTPC; //max DCA xy
185 Double_t fMaxDCAZ, fMaxDCAZTPC; //max DCA z
186 Double_t fMaxConstrainChi2; //max constrain chi2 - vertex
187 Bool_t fMinTPCClustersFlag, fMinITSClustersFlag; //shows if this cut is used or not
188 Bool_t fMaxChi2PerTPCClusterFlag, fMaxChi2PerITSClusterFlag; //shows if this cut is used or not
189 Bool_t fMaxCov11Flag, fMaxCov22Flag, fMaxCov33Flag, fMaxCov44Flag, fMaxCov55Flag; //shows if this cut is used or not
190 Bool_t fMaxSigmaToVertexFlag; //shows if this cut is used or not
191 Bool_t fMaxSigmaToVertexTPCFlag; //shows if this cut is used or not
192 Bool_t fMaxDCAXYFlag, fMaxDCAXYTPCFlag; //shows if this cut is used or not
193 Bool_t fMaxDCAZFlag, fMaxDCAZTPCFlag; //shows if this cut is used or not
194 Bool_t fMaxConstrainChi2Flag; //shows if this cut is used or not
195 Bool_t fITSRefitFlag, fTPCRefitFlag; //shows if this cut is used or not
196 Bool_t fESDpidFlag, fTPCpidFlag; //shows if this cut is used or not
199 Bool_t fFunctionProbabilityFlag; //flag: kTRUE if functions used
200 Double_t fPartFrac[10]; //prior probabilities
201 TF1 *fElectronFunction; //momentum dependence of the prior probs
202 TF1 *fMuonFunction; //momentum dependence of the prior probs
203 TF1 *fPionFunction; //momentum dependence of the prior probs
204 TF1 *fKaonFunction; //momentum dependence of the prior probs
205 TF1 *fProtonFunction; //momentum dependence of the prior probs
208 Bool_t fUseTPCOnly; //kTRUE if TPC only information is used
210 //Analysis containers
211 AliCFContainer *fProtonContainer; //container for protons
212 AliCFContainer *fAntiProtonContainer; //container for antiprotons
213 TH1I *fHistEvents; //event counter
214 TH2D *fHistYPtProtons; //Y-Pt of Protons
215 TH2D *fHistYPtAntiProtons; // Y-Pt of Antiprotons
218 TList *fEffGridListProtons; //list for the efficiency grid - protons
219 TList *fCorrectionListProtons2D; //list for the 2d corrections
220 TList *fEfficiencyListProtons1D; //list for the 1d efficiencies
221 TList *fCorrectionListProtons1D; //list for the 1d corrections
222 TList *fEffGridListAntiProtons; //list for the efficiency grid - antiprotons
223 TList *fCorrectionListAntiProtons2D; //list for the 2d corrections
224 TList *fEfficiencyListAntiProtons1D; //list for the 1d efficiencies
225 TList *fCorrectionListAntiProtons1D; //list for the 1d corrections
226 AliCFDataGrid *fCorrectProtons; //corrected data grid for protons
227 AliCFDataGrid *fCorrectAntiProtons; //corrected data grid for antiprotons
229 ClassDef(AliProtonAnalysis,0);