]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/SPECTRA/AliProtonAnalysis.h
Adding more QA histograms: eta vs phi vs Nclusters
[u/mrichter/AliRoot.git] / PWG2 / SPECTRA / AliProtonAnalysis.h
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 //
13 //    Origin: Panos Christakoglou | Panos.Christakoglou@cern.ch
14 //-------------------------------------------------------------------------
15
16 #include "TObject.h"
17 #include "TH1I.h"
18
19 class TF1;
20 class TH2D;
21 class TH1F;
22 class TList;
23
24 #include "AliPID.h"
25 #include "AliCFContainer.h"
26 class AliCFDataGrid;
27 class AliAODEvent;
28 class AliAODtrack;
29 class AliESDEvent;
30 class AliESDtrack;
31 class AliExternalTrackParam;
32 class AliStack;
33 class AliESDVertex;
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();
41
42   void SetEtaMode() {fAnalysisEtaMode = kTRUE;}
43
44   void UseTPCOnly() {fUseTPCOnly = kTRUE;}
45   void UseHybridTPC() {fUseHybridTPC = kTRUE;}
46   
47   void InitAnalysisHistograms(Int_t nbinsY, Float_t fLowY, Float_t fHighY,
48                               Int_t nbinsPt, Float_t fLowPt, Float_t fHighPt);
49   Bool_t ReadFromFile(const char* filename);
50   void Analyze(AliESDEvent *fESD, 
51                const AliESDVertex *vertex);
52   void Analyze(AliAODEvent *fAOD);
53   void Analyze(AliStack *stack, Bool_t iInclusive);
54   
55   AliCFContainer *GetProtonContainer() const {return fProtonContainer;}
56   AliCFContainer *GetAntiProtonContainer() const {return fAntiProtonContainer;}
57
58   TH2D *GetProtonYPtHistogram() const {return fHistYPtProtons;}
59   TH2D *GetAntiProtonYPtHistogram() const {return fHistYPtAntiProtons;}
60   TH1D *GetProtonYHistogram();
61   TH1D *GetAntiProtonYHistogram();
62   TH1D *GetProtonPtHistogram();
63   TH1D *GetAntiProtonPtHistogram();
64   TH1D *GetProtonCorrectedYHistogram();
65   TH1D *GetAntiProtonCorrectedYHistogram();
66   TH1D *GetProtonCorrectedPtHistogram();
67   TH1D *GetAntiProtonCorrectedPtHistogram();
68   
69   TH1D *GetYRatioHistogram();
70   TH1D *GetYRatioCorrectedHistogram(TH2D *gCorrectionMapProtons,
71                                     TH2D *gCorrectionMapAntiProtons);
72   TH1D *GetPtRatioHistogram();
73   TH1D *GetPtRatioCorrectedHistogram(TH2D *gCorrectionMapProtons,
74                                      TH2D *gCorrectionMapAntiProtons);
75   TH1D *GetYAsymmetryHistogram();
76   TH1D *GetPtAsymmetryHistogram();
77
78   TH1I *GetEventHistogram() const {return fHistEvents;}
79
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
84   //Cut functions
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;}
91   void    SetMinITSClusters(Int_t minITSClusters) {
92     fMinITSClusters = minITSClusters;
93     fMinITSClustersFlag = kTRUE;
94   }
95   void    SetMaxChi2PerITSCluster(Double_t maxChi2PerITSCluster) {
96     fMaxChi2PerITSCluster = maxChi2PerITSCluster;
97     fMaxChi2PerITSClusterFlag = kTRUE;
98   }
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   }
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;}
117   void    SetMaxSigmaToVertex(Double_t maxSigmaToVertex) {
118     fMaxSigmaToVertex = maxSigmaToVertex;
119     fMaxSigmaToVertexFlag = kTRUE;
120   }
121   void    SetMaxSigmaToVertexTPC(Double_t maxSigmaToVertex) {
122     fMaxSigmaToVertexTPC = maxSigmaToVertex;
123     fMaxSigmaToVertexTPCFlag = kTRUE;
124   }
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   }
145   void    SetITSRefit() {fITSRefitFlag = kTRUE;}
146   void    SetTPCRefit() {fTPCRefitFlag = kTRUE;}
147   void    SetESDpid() {fESDpidFlag = kTRUE;}
148   void    SetTPCpid() {fTPCpidFlag = kTRUE;}
149
150   //Prior probabilities
151   void SetPriorProbabilities(Double_t * const partFrac) {
152     for(Int_t i = 0; i < AliPID::kSPECIESN; i++) fPartFrac[i] = partFrac[i];} 
153   void SetPriorProbabilityFunctions(TF1 *const felectron, 
154                                     TF1 *const fmuon, 
155                                     TF1 *const fpion, 
156                                     TF1 *const fkaon, 
157                                     TF1 *const fproton) {
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
167   //interface to the correction framework
168   void Correct(Int_t step);
169   Bool_t ReadCorrectionContainer(const char* filename);
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;} 
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);}
182
183  private:
184   AliProtonAnalysis(const AliProtonAnalysis&); // Not implemented
185   AliProtonAnalysis& operator=(const AliProtonAnalysis&); // Not implemented
186
187   Bool_t   IsAccepted(AliESDEvent *esd,
188                       const AliESDVertex *vertex, 
189                       AliESDtrack *track);
190   Bool_t   IsInPhaseSpace(AliESDtrack *track);
191
192   Float_t  GetSigmaToVertex(AliESDtrack* esdTrack) const; 
193   Double_t Rapidity(Double_t Px, Double_t Py, Double_t Pz) const;
194   
195   Bool_t fAnalysisEtaMode; //run the analysis in eta or y
196
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
207   Double_t fMaxSigmaToVertexTPC; //max sigma to vertex cut
208   Double_t fMaxDCAXY, fMaxDCAXYTPC; //max DCA xy
209   Double_t fMaxDCAZ, fMaxDCAZTPC; //max DCA z
210   Double_t fMaxConstrainChi2; //max constrain chi2 - vertex
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
215   Bool_t fMaxSigmaToVertexTPCFlag; //shows if this cut is used or not
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
219   Bool_t fITSRefitFlag, fTPCRefitFlag; //shows if this cut is used or not
220   Bool_t fESDpidFlag, fTPCpidFlag; //shows if this cut is used or not
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
224   
225   //pid
226   Bool_t fFunctionProbabilityFlag; //flag: kTRUE if functions used
227   Double_t fPartFrac[10]; //prior probabilities
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
234   //Detectors
235   Bool_t fUseTPCOnly; //kTRUE if TPC only information is used
236   Bool_t fUseHybridTPC; //kTRUE if TPC info is used for momentum - PID and ITS for vertex & points
237
238   //Analysis containers
239   AliCFContainer *fProtonContainer; //container for protons
240   AliCFContainer *fAntiProtonContainer; //container for antiprotons
241   TH1I *fHistEvents; //event counter
242   TH2D *fHistYPtProtons; //Y-Pt of Protons
243   TH2D *fHistYPtAntiProtons; // Y-Pt of Antiprotons
244
245   //Corrections
246   TList *fEffGridListProtons; //list for the efficiency grid - protons 
247   TList *fCorrectionListProtons2D; //list for the 2d corrections 
248   TList *fEfficiencyListProtons1D; //list for the 1d efficiencies
249   TList *fCorrectionListProtons1D; //list for the 1d corrections 
250   TList *fEffGridListAntiProtons; //list for the efficiency grid - antiprotons 
251   TList *fCorrectionListAntiProtons2D; //list for the 2d corrections 
252   TList *fEfficiencyListAntiProtons1D; //list for the 1d efficiencies
253   TList *fCorrectionListAntiProtons1D; //list for the 1d corrections 
254   AliCFDataGrid *fCorrectProtons; //corrected data grid for protons
255   AliCFDataGrid *fCorrectAntiProtons; //corrected data grid for antiprotons
256
257   ClassDef(AliProtonAnalysis,0);
258 };
259
260 #endif