]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/SPECTRA/ChargedHadrons/dNdPt/AlidNdPtAnalysisPbPbAOD.h
modified dNdPt/AlidNdPtAnalysisPbPbAOD.{h,cxx}
[u/mrichter/AliRoot.git] / PWGLF / SPECTRA / ChargedHadrons / dNdPt / AlidNdPtAnalysisPbPbAOD.h
1 #ifndef AlidNdPtAnalysisPbPbAOD_H
2 #define AlidNdPtAnalysisPbPbAOD_H
3
4
5 //------------------------------------------------------------------------------
6 // AlidNdPtAnalysisPbPbAOD class used for dNdPt analysis in PbPb collision
7 // via AODs 
8 // 
9 // Author: P. Luettig, 15.05.2013
10 //------------------------------------------------------------------------------
11
12 class iostream;
13
14 class TObject;
15 class TFile;
16 class TCint;
17 class THnSparse;
18
19 #include "AliAnalysisTaskSE.h"
20
21
22 #include "TList.h"
23 #include "TFile.h"
24 #include "TH1.h"
25 #include "TH2.h"
26 #include "TH3.h"
27 #include "THnSparse.h"
28 #include "TClonesArray.h"
29
30 #include "TParticlePDG.h"
31 #include "TDatabasePDG.h"
32
33 #include "AliLog.h"
34 #include "AliCentrality.h"
35 #include "AliAODEvent.h"
36 #include "AliVEvent.h"
37
38 #include "AliInputEventHandler.h"
39 #include "AliAODInputHandler.h"
40 #include "AliAnalysisManager.h"
41 #include "AliMCEventHandler.h"
42 #include "AliAODMCHeader.h"
43 #include "AliAODMCParticle.h"
44 #include "AliGenHijingEventHeader.h"
45 #include "AliGenPythiaEventHeader.h"
46
47 #include "TSystem.h"
48 #include "TROOT.h"
49
50 class AlidNdPtAnalysisPbPbAOD : public AliAnalysisTaskSE {
51   public :
52     AlidNdPtAnalysisPbPbAOD(); 
53     AlidNdPtAnalysisPbPbAOD(const char *name);
54     ~AlidNdPtAnalysisPbPbAOD();
55     
56     virtual void UserCreateOutputObjects();
57     virtual void UserExec(Option_t *option);
58     virtual void Terminate(Option_t *);
59     
60     // Set binning for Histograms (if not set default binning is used)
61     void SetBinsMult(Int_t nbins, Double_t* edges)              { Printf("[I] Setting Mult Bins"); fMultNbins = nbins; fBinsMult = GetArrayClone(nbins,edges); }
62     void SetBinsPt(Int_t nbins, Double_t* edges)                { Printf("[I] Setting pT Bins"); fPtNbins = nbins; fBinsPt = GetArrayClone(nbins,edges); }
63     void SetBinsPtCorr(Int_t nbins, Double_t* edges)            { Printf("[I] Setting pTcorr Bins"); fPtCorrNbins = nbins; fBinsPtCorr = GetArrayClone(nbins,edges); }
64     void SetBinsEta(Int_t nbins, Double_t* edges)               { Printf("[I] Setting Eta Bins"); fEtaNbins = nbins; fBinsEta = GetArrayClone(nbins,edges); }
65     void SetBinsZv(Int_t nbins, Double_t* edges)                { Printf("[I] Setting Zv Bins"); fZvNbins = nbins; fBinsZv= GetArrayClone(nbins,edges); }
66     void SetBinsCentrality(Int_t nbins, Double_t* edges)        { Printf("[I] Setting Cent Bins"); fCentralityNbins = nbins; fBinsCentrality = GetArrayClone(nbins,edges); }
67     
68     
69     // set event cut variables
70     void SetCutMaxZVertex( Double_t d)                                  { dCutMaxZVertex = d; }
71     
72     // set track kinematic cut parameters
73     void SetCutPtRange(Double_t ptmin, Double_t ptmax)                  { dCutPtMin = ptmin; dCutPtMax = ptmax; }
74     void SetCutEtaRange(Double_t etamin, Double_t etamax)               { dCutEtaMin = etamin; dCutEtaMax = etamax; }
75     
76     // set track quality cut parameters
77     void SetCutRequireTPCRefit(Bool_t *b)                               { bCutRequireTPCRefit = b; } 
78     void SetCutMinNCrossedRowsTPC(Double_t d)                           { dCutMinNumberOfCrossedRows = d; }    
79     void SetCutMinRatioCrossedRowsOverFindableClustersTPC(Double_t d)   { dCutMinRatioCrossedRowsOverFindableClustersTPC = d; }
80     void SetCutMaxChi2PerClusterTPC(Double_t d)                         { dCutMaxChi2PerClusterTPC = d; }
81     void SetCutMaxFractionSharedTPCClusters(Double_t d)                 { dCutMaxFractionSharedTPCClusters = d; }
82     void SetCutMaxDCAToVertexZ(Double_t d)                              { dCutMaxDCAToVertexZ = d; }
83     void SetCutMaxDCAToVertexXY(Double_t d)                             { dCutMaxDCAToVertexXY = d; }
84     void SetCutRequireITSRefit(Bool_t *b)                               { bCutRequireITSRefit = b; } 
85     void SetCutMaxChi2PerClusterITS(Double_t d)                         { dCutMaxChi2PerClusterITS = d; }
86     void SetCutDCAToVertex2D(Bool_t *b)                                 { dCutDCAToVertex2D = b; } 
87     void SetCutRequireSigmaToVertex(Bool_t *b)                          { dCutRequireSigmaToVertex = b; } 
88     void SetCutMaxDCAToVertexXYPtDep(Double_t d0, Double_t d1, Double_t d2)
89     {
90       dCutMaxDCAToVertexXYPtDepPar0 = d0;
91       dCutMaxDCAToVertexXYPtDepPar1 = d1;
92       dCutMaxDCAToVertexXYPtDepPar2 = d2;
93     }
94     void SetCutAcceptKinkDaughters(Bool_t *b)                           { bCutAcceptKinkDaughters = b; } 
95     void SetCutMaxChi2TPCConstrainedGlobal(Double_t d)                  { dCutMaxChi2TPCConstrainedGlobal = d; }
96     
97     // getter for qualtiy track cuts
98     Double_t GetCutMinNCrossedRowsTPC()                                 { return dCutMinNumberOfCrossedRows; }
99     
100     THnSparseF *GetHistZvPtEtaCent() const { return hnZvPtEtaCent; }
101     TH1F *GetHistEventStatistics() const { return hEventStatistics; }
102     
103     const char * GetParticleName(Int_t pdg);
104     
105     AliGenHijingEventHeader* GetHijingEventHeader(AliAODMCHeader *header);
106     AliGenPythiaEventHeader* GetPythiaEventHeader(AliAODMCHeader *header);
107     
108 //     Int_t IsMCSecondary(AliAODMCParticle *part, TClonesArray *arrayMC);
109     Bool_t IsTrackAccepted(AliAODTrack *tr);
110     Bool_t IsMCTrackAccepted(AliAODMCParticle *part);
111     
112     Bool_t IsHijingParticle(const AliAODMCParticle *part, AliGenHijingEventHeader* hijingGenHeader);
113     Bool_t IsPythiaParticle(const AliAODMCParticle *part, AliGenPythiaEventHeader* pythiaGenHeader);
114     
115     static Double_t* GetArrayClone(Int_t n, Double_t* source);
116     
117   private :
118     
119     // Output List
120     TList       *fOutputList;
121     
122     // Histograms
123     TH1F        *hPt; // simple pT histogramm
124     TH1F        *hMCPt; // simple pT truth histogramm
125     THnSparseF  *hnZvPtEtaCent; //-> Zv:Pt:Eta:Cent
126     THnSparseF  *hnMCRecPrimZvPtEtaCent; //-> MC Zv:Pt:Eta:Cent
127     THnSparseF  *hnMCGenZvPtEtaCent; //-> MC Zv:Pt:Eta:Cent
128     THnSparseF  *hnMCRecSecZvPtEtaCent; //-> MC Zv:Pt:Eta:Cent, only secondaries
129     TH1F        *hEventStatistics; // contains statistics of number of events after each cut
130     TH1F        *hEventStatisticsCentrality; // contains number of events vs centrality, events need to have a track in kinematic range
131     TH1F        *hAllEventStatisticsCentrality; // contains number of events vs centrality, events need to be triggered
132     THnSparseF  *hnZvMultCent; // Zv:Mult:Cent
133     TH1F        *hTriggerStatistics; // contains number of events per trigger
134     TH1F        *hMCTrackPdgCode; // contains statistics of pdg codes of tracks
135     TH1F        *hMCTrackStatusCode; // contains statistics of status codes of tracks
136     TH1F        *hCharge; // charge distribution in data
137     TH1F        *hMCCharge; // charge distribution in MC
138     TH2F        *hMCPdgPt; // PDGvs PT for MC Particles
139     TH1F        *hMCHijingPrim; // number of particles, which are Hijing particles and primaries
140     TH1F        *hAccNclsTPC; //control histo: number of clusters in TPC for accepted tracks
141     TH1F        *hAccCrossedRowsTPC; //control histo: number of crossed rows in TPC for accepted tracks
142     TH2F        *hDCAPtAll; //control histo: DCA vs pT for all reconstructed tracks
143     TH2F        *hDCAPtAccepted; //control histo: DCA vs pT for all accepted reco tracks
144     TH2F        *hMCDCAPtSecondary; //control histo: DCA vs pT for all accepted reco track, which are secondaries (using MC info)
145     TH2F        *hMCDCAPtPrimary; //control histo: DCA vs pT for all accepted reco track, which are primaries (using MC info)
146     
147     
148     // global variables
149     Bool_t bIsMonteCarlo;
150      
151     // event cut variables
152     Double_t dCutMaxZVertex;
153     
154     // track kinematic cut variables
155     Double_t dCutPtMin;
156     Double_t dCutPtMax;
157     Double_t dCutEtaMin;
158     Double_t dCutEtaMax;
159     
160     // track quality cut variables
161     Bool_t      bCutRequireTPCRefit;
162     Double_t    dCutMinNumberOfCrossedRows;
163     Double_t    dCutMinRatioCrossedRowsOverFindableClustersTPC;
164     Double_t    dCutMaxChi2PerClusterTPC;
165     Double_t    dCutMaxFractionSharedTPCClusters;
166     Double_t    dCutMaxDCAToVertexZ;
167     Double_t    dCutMaxDCAToVertexXY;
168     Bool_t      bCutRequireITSRefit;
169     Double_t    dCutMaxChi2PerClusterITS;
170     Bool_t      dCutDCAToVertex2D;
171     Bool_t      dCutRequireSigmaToVertex;
172     Double_t    dCutMaxDCAToVertexXYPtDepPar0;
173     Double_t    dCutMaxDCAToVertexXYPtDepPar1;
174     Double_t    dCutMaxDCAToVertexXYPtDepPar2;
175     Bool_t      bCutAcceptKinkDaughters;
176     Double_t    dCutMaxChi2TPCConstrainedGlobal;
177     
178     //binning for THNsparse
179     Int_t fMultNbins;
180     Int_t fPtNbins;
181     Int_t fPtCorrNbins;
182     Int_t fEtaNbins;
183     Int_t fZvNbins;
184     Int_t fCentralityNbins;
185     Double_t* fBinsMult; //[fMultNbins]
186     Double_t* fBinsPt; //[fPtNbins]
187     Double_t* fBinsPtCorr; //[fPtCorrNbins]
188     Double_t* fBinsEta; //[fEtaNbins]
189     Double_t* fBinsZv; //[fZvNbins]
190     Double_t* fBinsCentrality; //[fCentralityNbins]
191     
192     AlidNdPtAnalysisPbPbAOD(const AlidNdPtAnalysisPbPbAOD&); // not implemented
193     AlidNdPtAnalysisPbPbAOD& operator=(const AlidNdPtAnalysisPbPbAOD&); // not implemented  
194     
195     ClassDef(AlidNdPtAnalysisPbPbAOD,2); // has to be at least 1, otherwise not streamable...
196 };
197
198 #endif