]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/SPECTRA/AntiprotonToProton/AliProtonAnalysis.h
Updates to Trains. create a job-script to help
[u/mrichter/AliRoot.git] / PWGLF / SPECTRA / AntiprotonToProton / 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 class AliPID;
26 #include "AliCFContainer.h"
27 class AliCFDataGrid;
28 class AliAODEvent;
29 class AliAODtrack;
30 class AliESDEvent;
31 class AliESDtrack;
32 class AliExternalTrackParam;
33 class AliStack;
34 class AliESDVertex;
35 class AliProtonAnalysisBase;
36
37 class AliProtonAnalysis : public TObject {
38  public:
39   enum {
40     kStepSurvived        = 0,
41     kStepIsPrimary       = 1,
42     kStepIdentified      = 2,
43     kStepInPhaseSpace    = 3,
44     kNSteps = 4
45   };
46   AliProtonAnalysis();
47   AliProtonAnalysis(Int_t nbinsY, Float_t fLowY, Float_t fHighY,
48                     Int_t nbinsPt, Float_t fLowPt, Float_t fHighPt);
49   AliProtonAnalysis(Int_t nbinsY, Double_t *gY,
50                     Int_t nbinsPt, Double_t *gPt);
51   virtual ~AliProtonAnalysis();
52
53   void SetBaseAnalysis(AliProtonAnalysisBase * const baseAnalysis) {
54     fProtonAnalysisBase = baseAnalysis;}
55   AliProtonAnalysisBase *GetProtonAnalysisBaseObject() const {
56     return fProtonAnalysisBase;}
57
58   void InitAnalysisHistograms(Int_t nbinsY, Float_t fLowY, Float_t fHighY,
59                               Int_t nbinsPt, Float_t fLowPt, Float_t fHighPt);
60   void InitAnalysisHistograms(Int_t nbinsY, Double_t *gY,
61                               Int_t nbinsPt, Double_t *gPt);
62   Bool_t ReadFromFile(const char* filename);
63   void Analyze(AliESDEvent *fESD, 
64                const AliESDVertex *vertex);
65   void Analyze(AliAODEvent *fAOD);
66   void Analyze(AliStack *stack, Bool_t iInclusive);
67   
68   //QA for real data
69   void InitQA();
70   void FillQA(AliESDEvent *esd,
71               const AliESDVertex *vertex, 
72               AliESDtrack* track);
73   TList *GetQAList() {return fGlobalQAList;}
74
75   AliCFContainer *GetProtonContainer() const {return fProtonContainer;}
76   AliCFContainer *GetAntiProtonContainer() const {return fAntiProtonContainer;}
77
78   TH2D *GetProtonYPtHistogram() const {return fHistYPtProtons;}
79   TH2D *GetAntiProtonYPtHistogram() const {return fHistYPtAntiProtons;}
80   TH2D *GetProtonYPtCorrectedHistogram() const {
81     return fHistYPtProtonsCorrected;}
82   TH2D *GetAntiProtonCorrectedYPtHistogram() const {
83     return fHistYPtAntiProtonsCorrected;}
84   TH1D *GetProtonYHistogram();
85   TH1D *GetAntiProtonYHistogram();
86   TH1D *GetProtonPtHistogram();
87   TH1D *GetAntiProtonPtHistogram();
88   TH1D *GetProtonCorrectedYHistogram();
89   TH1D *GetAntiProtonCorrectedYHistogram();
90   TH1D *GetProtonCorrectedPtHistogram();
91   TH1D *GetAntiProtonCorrectedPtHistogram();
92
93   TH1F *GetEventStatistics() {return fHistEventStats;}
94
95   TList *GetYRatioHistogramsInPtBins();
96   TH1D *GetYRatioHistogram();
97   TH1D *GetYRatioCorrectedHistogram();
98   TH1D *GetPtRatioHistogram();
99   TH1D *GetPtRatioCorrectedHistogram();
100   TH1D *GetYAsymmetryHistogram();
101   TH1D *GetPtAsymmetryHistogram();
102   TH2D *GetProtonsAbsorptionMaps() {return fHistEfficiencyYPtProtons;}
103   TH2D *GetAntiProtonsAbsorptionMaps() {return fHistEfficiencyYPtAntiProtons;}
104
105   TH1I *GetEventHistogram() const {return fHistEvents;}
106
107   Int_t   GetNumberOfAnalyzedEvents() {return (Int_t)fHistEventStats->GetBinContent(5);}
108   Bool_t  PrintMean(TH1 *hist, Double_t edge);
109   Bool_t  PrintYields(TH1 *hist, Double_t edge); 
110
111   //interface to the correction framework
112   void Correct();
113   void Correct(Int_t step);
114   Bool_t ReadCorrectionContainer(const char* filename);
115
116   void SetCorrectionMapForFeedDown(const char* filename);
117   TH2D *GetCorrectionMapForFeedDownProtons() {
118     return fHistYPtCorrectionForFeedDownProtons;}
119   TH2D *GetCorrectionMapForFeedDownAntiProtons() {
120     return fHistYPtCorrectionForFeedDownAntiProtons;}
121
122   void SetCorrectionMapForCuts(const char* filename);
123   TH2D *GetCorrectionMapForCutsProtons() {
124     return fHistYPtCorrectionForCutsProtons;}
125   TH2D *GetCorrectionMapForCutsAntiProtons() {
126     return fHistYPtCorrectionForCutsAntiProtons;}
127
128   void SetCorrectionMapForTracking(const char* filename);
129   TH2D *GetCorrectionMapForTrackingProtons() {
130     return fHistYPtCorrectionForTrackingProtons;}
131   TH2D *GetCorrectionMapForTrackingAntiProtons() {
132     return fHistYPtCorrectionForTrackingAntiProtons;}
133
134
135   void SetCorrectionMapForSecondaries(const char* filename);
136   TH2D *GetCorrectionMapForSecondaries() {
137     return fHistYPtCorrectionForSecondaries;}
138
139   void SetCorrectionMapForCrossSection(const char* filename);
140   TH2D *GetProtonsCorrectionMapForCrossSection() {
141     return fHistCorrectionForCrossSectionYPtProtons;}
142   TH2D *GetAntiProtonsCorrectionMapForCrossSection() {
143     return fHistCorrectionForCrossSectionYPtAntiProtons;}
144
145   //iStep=0->MC - iStep=1->Acceptance - iStep=2->Reconstruction - iStep=3->PID
146   TH1D  *GetUncorrectedProtonYHistogram(Int_t iStep) {return fProtonContainer->ShowProjection(0, iStep);}
147   TH1D  *GetUncorrectedProtonPtHistogram(Int_t iStep) {return fProtonContainer->ShowProjection(1, iStep);}
148   TH1D  *GetUncorrectedAntiProtonYHistogram(Int_t iStep) {return fAntiProtonContainer->ShowProjection(0, iStep);}
149   TH1D  *GetUncorrectedAntiProtonPtHistogram(Int_t iStep) {return fAntiProtonContainer->ShowProjection(1, iStep);}
150
151  private:
152   AliProtonAnalysis(const AliProtonAnalysis&); // Not implemented
153   AliProtonAnalysis& operator=(const AliProtonAnalysis&); // Not implemented
154   
155   AliProtonAnalysisBase *fProtonAnalysisBase;//base analysis object
156
157   Int_t fNBinsY; //number of bins in y or eta
158   Double_t fMinY, fMaxY; //min & max value of y or eta
159   Int_t fNBinsPt;  //number of bins in pT
160   Double_t fMinPt, fMaxPt; //min & max value of pT
161   
162   //Analysis containers
163   AliCFContainer *fProtonContainer; //container for protons
164   AliCFContainer *fAntiProtonContainer; //container for antiprotons
165   TH1I *fHistEvents; //event counter
166   TH2D *fHistYPtProtons; //Y-Pt of Protons
167   TH2D *fHistYPtAntiProtons; // Y-Pt of Antiprotons
168   TH2D *fHistYPtProtonsCorrected; //Y-Pt of Protons (corrected)
169   TH2D *fHistYPtAntiProtonsCorrected; // Y-Pt of Antiprotons (corrected)
170   TH1F *fHistEventStats;//Event statistics
171   TList *fYRatioInPtBinsList;//TList of the eta dependent ratios for each pT bin
172
173   //Corrections
174   TList *fEffGridListProtons; //list for the efficiency grid - protons 
175   TList *fCorrectionListProtons2D; //list for the 2d corrections 
176   TList *fEfficiencyListProtons1D; //list for the 1d efficiencies
177   TList *fCorrectionListProtons1D; //list for the 1d corrections 
178   TList *fEffGridListAntiProtons; //list for the efficiency grid - antiprotons 
179   TList *fCorrectionListAntiProtons2D; //list for the 2d corrections 
180   TList *fEfficiencyListAntiProtons1D; //list for the 1d efficiencies
181   TList *fCorrectionListAntiProtons1D; //list for the 1d corrections 
182   AliCFDataGrid *fCorrectProtons; //corrected data grid for protons
183   AliCFDataGrid *fCorrectAntiProtons; //corrected data grid for antiprotons
184   TH2D *fHistEfficiencyYPtProtons;//efficiency 2D for the corrections - protons
185   TH2D *fHistEfficiencyYPtAntiProtons;//efficiency 2D for the corrections - antiprotons
186   TH2D *fHistCorrectionForCrossSectionYPtProtons;//correction for the proper cross-section - 2D protons
187   TH2D *fHistCorrectionForCrossSectionYPtAntiProtons;//correction for the proper cross-section - 2D antiprotons
188   Bool_t fHistCorrectionForCrossSectionFlag;//correct for cross-section
189   TH2D *fHistYPtCorrectionForCutsProtons;//correction factors for the cut efficiency (protons)
190   TH2D *fHistYPtCorrectionForCutsAntiProtons;//correction factors for the cut efficiency (antiprotons)
191   Bool_t fCorrectForCutsFlag;//correct for the cut efficiency
192   TH2D *fHistYPtCorrectionForTrackingProtons;//correction factors for the tracking efficiency (protons)
193   TH2D *fHistYPtCorrectionForTrackingAntiProtons;//correction factors for the tracking efficiency (antiprotons)
194   Bool_t fCorrectForTrackingFlag;//correct for the tracking efficiency
195   TH2D *fHistYPtCorrectionForFeedDownProtons;//correction factors for the feed-down contamination (protons)
196   TH2D *fHistYPtCorrectionForFeedDownAntiProtons;//correction factors for the feed-down contamination (antiprotons)
197   Bool_t fCorrectForFeedDownFlag;//correct for cut efficiency
198   TH2D *fHistYPtCorrectionForSecondaries;//correction factors for the background protons
199   Bool_t fCorrectForSecondariesFlag;//correct for secondaries
200
201   //QA lists
202   TList *fGlobalQAList; //global list
203   TList *fQA2DList; //QA 2D list
204   TList *fQAProtonsAcceptedList; //accepted protons
205   TList *fQAProtonsRejectedList; //rejected protons
206   TList *fQAAntiProtonsAcceptedList; //accepted antiprotons
207   TList *fQAAntiProtonsRejectedList; //rejected antiprotons
208
209   ClassDef(AliProtonAnalysis,1);
210 };
211
212 #endif