]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/DPhi/AliAnalysisTaskDiHadron.h
moving AliAnalysisTaskDiHadron (Jason) to PWGCF
[u/mrichter/AliRoot.git] / PWGCF / Correlations / DPhi / AliAnalysisTaskDiHadron.h
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 //2- and 3-particle trigger particle correlation analysis
17 //Author: Jason Glyndwr Ulery, ulery@uni-frankfurt.de
18 //version 3.4, last revised 2010/08/15
19
20 #ifndef ALIANALYSISTASKDIHADRON_H
21 #define ALIANALYSISTASKDIHADRON_H
22
23
24 class TF1;
25 class TH1F;
26 class TH2F;
27 class TH3F;
28 class AliESDEvent;
29
30 #include "AliAnalysisTask.h"
31
32
33 class AliAnalysisTaskDiHadron : public AliAnalysisTask{
34   //Default constructor
35  public:
36
37
38   AliAnalysisTaskDiHadron(const char *name="AliAnalysisTaskDiHadron");
39
40   virtual ~AliAnalysisTaskDiHadron() {}
41   //virtual ~AliAnalysisTaskDiHadron();
42
43  
44
45   //AliAnalysisTask Functions
46  
47   virtual void ConnectInputData(Option_t *option);
48   virtual void CreateOutputObjects();
49   virtual void Exec(Option_t *option);
50   virtual void Terminate(Option_t *);
51  
52
53   void SetEfficiencies(Float_t EffFitPt, const TF1 *FitLow, const TF1 *FitHigh, Int_t NFitLowParam, Int_t NFitHighParam, Float_t *FitLowParam, Float_t *FitHighParam);
54   void SetBins(Int_t nBinPhi, Int_t nBinEta, Int_t nBinPhiEtaPhi, Int_t nBinPhiEtaEta, Int_t nBinPhi3, Int_t nBinEta3, Float_t dPhiMin, Float_t dPhiMax, Int_t NTPtBins, Int_t NMixBins, Int_t NCentBins,Int_t NAPtBins, Int_t NAPt3Bins, Int_t NVertexBins, Int_t NXEBin,Float_t *PtTrigArray, Float_t *PtAssocArray, Float_t *PtAssoc3Array1, Float_t *PtAssoc3Array2, Int_t *CentArrayMin, Int_t *CentArrayMax, Float_t *XEArray);
55   void SetOptions(Int_t fEfficiencyCorr, Int_t fDEBUG,Int_t fMCHistos);
56   void SetCuts(Int_t MinClutersTPC, Float_t MinClusterRatio, Float_t MaxTPCchi2, Int_t MinClustersITS, Float_t EtaCut, Float_t TrigEtaCut, Float_t NearPhiCut, Float_t XECut, Float_t MaxDCA, Float_t MaxDCAXY, Float_t MaxDCAZ, Int_t DCA2D, Int_t TPCRefit, Int_t ITSRefit, Int_t SPDCut, Float_t MinPtAssoc, Float_t MaxPtAssoc, Float_t VzCut, Int_t NIDs, const char *TrigIDArray);
57   Int_t CheckVertex(const AliESDEvent *rESD);
58   Int_t CheckTrigger(const AliESDEvent *rESD);
59   Int_t TrackCuts(const AliESDEvent *rESD, Float_t *rPt, Float_t *rEta, Float_t *rPhi, Short_t *rCharge, Float_t *rEff, Int_t **rPtAssoc3, Int_t *rNPtAssoc3, Int_t *rGoodTracks);
60   Int_t TrackCutsMC(AliMCEvent *rMC, Float_t *rPt, Float_t *rEta, Float_t *rPhi, Short_t *rCharge, Float_t *rEff, Int_t **rPtAssoc3, Int_t *rNPtAssoc3, Int_t *rGoodTracks);
61
62   
63  private:
64
65 //Maximum numbers for creating the arrays
66   enum{kNumberOfPtBins=20, 
67        kNumberOfCentBins=10, 
68        kNumberOfAPtBins=50,
69        kNumberOfApt3Bins=50,
70        kNumberOfTriggerIDs=10,
71        kNumberOfVertexBins=20,
72        kNumberOfXeBins=20,
73        kNumberOfEventsToMix=100};
74   // virtual void SetParameters();
75   AliESDEvent *fESD; //ESD data file
76   AliMCEvent *fMC;//Monty Carlo event file
77   TList *fOutput;//List used to store all the histograms for output
78   
79    //Cuts
80   
81   Int_t fMinClustersTPC;//Minimum Clusters in the TPC
82   Float_t fMinClusterRatio;//Ratio of MinClusters to Number of Possible Clusters in TPC
83   Float_t fMaxTPCchi2;//Max chi^2/ndf for TPC fit
84   Int_t fMinClustersITS;//Minimum Clusters in ITS
85   Float_t fEtaCut;//Associated particle eta range |eta|<fEtaCut
86   Float_t fTrigEtaCut;//Trigger particle eta range |eta|<fTrigEtaCut (should be <=fEtaCut)
87   Float_t fNearPhiCut;//Point in DeltaPhi seperating near and away sides
88   Float_t fXECut;//Point in DeltaPhi seperating near and away for XE (not finished)
89   Float_t fMaxDCA;//Max total DCA (if fDCA2D==0) (in cm)
90   Float_t fMaxDCAXY;//Max radial DCA (if fDCA2D==1)
91   Float_t fMaxDCAZ;//Max longitudional DCA (if fDCA2D==1)
92   Int_t fDCA2D;//1 cut on total DCA, 2 cut seperately on radial and logitudional, 3 pt dependent cut
93   Int_t fTPCRefit;//1 require TPCRefit 0 not
94   Int_t fITSRefit;//1 require ITSRefit 0 not
95   Int_t fSPDCut;//1 require a point in the SPD 0 not
96   Float_t fMinPtAssoc;//Minimum pT to look at 
97   Float_t fMaxPtAssoc;//Maximum pT to look at  
98   Float_t fVzCut;//z-vertex cut |z-vertex|<fVzCut in cm
99     Int_t fEfficiencyCorr;//Toggle correcting of efficiencies when filling histograms
100     Int_t fDEBUG;//if 1 extra printfs for debugging
101
102     //Binning
103     Int_t fnBinPhi;//Number of bins for DeltaPhi histograms
104     Int_t fnBinEta;//Number of bins for DeltaEta histograms
105     Int_t fnBinPhiEtaPhi;//Number of bins for DeltaPhi-DeltaEta histograms in Phi
106     Int_t fnBinPhiEtaEta;//Number of bins for DeltaPhi-DeltaEta histograms in Eta
107     Int_t fnBinPhi3;//Number of bins in 3-particle DeltaPhi-DeltaPhi (along each axis)
108     Int_t fnBinEta3;//Number of bins in 3-particle DeltaEta-DeltaEta (along each axis)
109     Float_t fPi;//Pi=3.14....
110     Float_t fdPhiMin;//lowest edge of the first bin in DeltaPhi plots
111     Float_t fdPhiMax;//maximum in DeltaPhi plots (should be 2*Pi+fdPhiMin)
112
113     //Parameters for settings
114     Int_t fNTPtBins;//Number of trigger pt bins
115     Int_t fNMix;//Number of events to mix
116     Int_t fNCentBins;//Number of centrality bins (they are allowed to overlap)
117     Int_t fNAPtBins;//Number of associated particle bins for 2-particle correlations (overlap not allowed)
118     Int_t fNAPt3Bins;//Number of associated particle bins for 3-particle correlations (overlap allowed)
119     Int_t fNVertexBins;//Number of bins in z-vertex for mixing
120     Int_t fNXEBins;//Number of bins for XE distribution (not finished)
121     Int_t fNIDs;//Number of trigger IDs (such as CINT1B)
122     Float_t fEffFitPt;//Pt cut used to divide efficiciency fits into low and high pt
123     Int_t fNFitLowParam;//Number of fit parameters for low pt efficeincy fit
124     Int_t fNFitHighParam;//Number of fit parameters for high pt efficeincy fit
125     Int_t fMCHistos;//Toggle whether to make Monty Carlo histograms (if made and ran on data will simply be empty not cause a crash)
126
127     TF1 *fFitLow;//Low pt efficiency fit function
128     TF1 *fFitHigh;//High pt efficiency fit function
129     Float_t *fFitLowParam; //fit parameters for low pt efficiency fit function
130     Float_t *fFitHighParam;//fit parameters for high pt efficiency fit function
131
132     Float_t *fPtTrigArray;//Array for trigger pt bins
133     Float_t *fPtAssocArray;//Array for associated pt bins (2-particle correlations)
134     //2 arrays to allow overlaping bins
135     Float_t *fPtAssoc3Array1;//lower bin edge for 3-particle correlations associated pt bins 
136     Float_t *fPtAssoc3Array2; //upper bin eged for 3-particle correlations associated pt bins
137     //2 arrays to allow oeverlaping bins 
138     Int_t *fCentArrayMin;//lower bin edge for centrality bins (based on SPD tracklet multiplicity)
139     Int_t *fCentArrayMax;//upper bin edge for centraily bins
140     Float_t *fXEArray;//XE bin array
141     char *fTrigIDArray;//array of trigger ids (such as CINT1B)
142     
143     Float_t fVertexArray[(kNumberOfVertexBins+1)];//array for dividing z vertex into bins for event mixing
144   
145   //Histograms
146     TH1F *fHistPt[kNumberOfCentBins][2];//Pt Distribution of tracks (no efficency correction)
147     TH1F *fHistPtEff[kNumberOfCentBins][2];//Pt Distribution of tracks (efficnecy corrected if set)
148     TH1F *fHistPtTrig[kNumberOfPtBins][kNumberOfCentBins][2];//Pt Distributions of tracks from triggered events
149     TH1F *fHistMult[2];//Multipuliciy Distributions
150     TH1F *fHistMultTrig[kNumberOfPtBins][2];//Multiplicity Distributions of triggered events
151
152     TH2F *fHistPhi[kNumberOfCentBins][2];//Phi Distributions of tracks
153     TH2F *fHistPhiTrig[kNumberOfPtBins][kNumberOfCentBins][2];//Phi Distributions fo tracks in triggered events
154     TH2F *fHistDeltaPhi[kNumberOfPtBins][kNumberOfCentBins][3][2];//DeltaPhi 2-particle distributions
155     TH2F *fHistDeltaPhiMix[kNumberOfPtBins][kNumberOfCentBins][3][2];//DeltaPhi 2-particle distributions from mixed events
156
157     TH2F *fHistPhiPt[kNumberOfCentBins][2];//Pt weighed phi distributions
158     TH2F *fHistPhiTrigPt[kNumberOfPtBins][kNumberOfCentBins][2];//Pt weighted phi distributions from triggered events
159     TH2F *fHistDeltaPhiPt[kNumberOfPtBins][kNumberOfCentBins][3][2];//DeltaPhi 2-partricle pt weighted distributions
160     TH2F *fHistDeltaPhiMixPt[kNumberOfPtBins][kNumberOfCentBins][3][2];//DeltaPhi 2-particle pt weigthed distributions from mixed events
161
162     TH2F *fHistEta[kNumberOfCentBins][2];//Eta distributions
163     TH2F *fHistEtaTrig[kNumberOfPtBins][kNumberOfCentBins][2];//Eta distributions from mixed events
164
165     TH2F *fHistEtaPt[kNumberOfCentBins][2];//Pt-weighted Eta distributions
166     TH2F *fHistEtaTrigPt[kNumberOfPtBins][kNumberOfCentBins][2];//Pt-weighted Eta diestibutions from mixed events
167
168     TH2F *fHistDeltaEtaN[kNumberOfPtBins][kNumberOfCentBins][3][2];//Near-side 2-particle DeltaEta distributions
169     TH2F *fHistDeltaEtaNMix[kNumberOfPtBins][kNumberOfCentBins][3][2];//Near-side 2-particle DeltaEta distributions from mxied events
170
171     TH2F *fHistDeltaEtaNPt[kNumberOfPtBins][kNumberOfCentBins][3][2];//Near-side pt weighted 2-particle DeltaEta distributions
172     TH2F *fHistDeltaEtaNMixPt[kNumberOfPtBins][kNumberOfCentBins][3][2];//Near-side pt-weighted 2-particle DeltaEta distributions from mixed events
173
174     TH2F *fHistDeltaEtaA[kNumberOfPtBins][kNumberOfCentBins][3][2];//Away-side DeltaEta distributions
175     TH2F *fHistDeltaEtaAMix[kNumberOfPtBins][kNumberOfCentBins][3][2];//Away-side DeltaEta distributions from mixed events
176
177     TH2F *fHistDeltaEtaAPt[kNumberOfPtBins][kNumberOfCentBins][3][2];//Away-side pt-weighted DeltaEta distributions from mixed events
178     TH2F *fHistDeltaEtaAMixPt[kNumberOfPtBins][kNumberOfCentBins][3][2];//Away-side pt-weigthed DeltaEta distributions from mixed events
179
180     TH1F *fHistNEvents[kNumberOfCentBins][2];//Number of events count
181     TH1F *fHistNTrigger[kNumberOfCentBins][2];//Number of triggers count
182     TH1F *fHistNTriggerPt[kNumberOfCentBins][2];//Pt-weighted number of triggers count
183     TH1F *fHistNMix[kNumberOfCentBins][2];//Number of mixed events count
184
185     TH3F *fHistPhiEta[kNumberOfCentBins][2];//Phi-Eta distributions
186     TH3F *fHistPhiEtaTrig[kNumberOfPtBins][kNumberOfCentBins][2];//Phi-Eta distributions in triggered events
187     TH3F *fHistDeltaPhiEta[kNumberOfPtBins][kNumberOfCentBins][2];//DeltaPhi-DeltaEta 2-particle distributions
188     TH3F *fHistDeltaPhiEtaMix[kNumberOfPtBins][kNumberOfCentBins][2];//DeltaPhi-DeltaEta 2-particle distributions from mixed events
189
190     TH3F *fHistPhiEtaPt[kNumberOfCentBins][2];//Pt-weighted Phi-Eta distributions
191     TH3F *fHistPhiEtaTrigPt[kNumberOfPtBins][kNumberOfCentBins][2];//Pt-weighted Phi-Eta distributions from triggered events
192     TH3F *fHistDeltaPhiEtaPt[kNumberOfPtBins][kNumberOfCentBins][2];//Pt-weighted DeltaPhi-DeltaEta 2-particle distributions
193     TH3F *fHistDeltaPhiEtaMixPt[kNumberOfPtBins][kNumberOfCentBins][2];//Pt-weigthed DeltaPhi-DeltaEta 2-particle distributions from mixed events
194
195     TH1F *fHistXEN[kNumberOfPtBins][kNumberOfCentBins][2];//XE Near-side (not finished)
196     TH1F *fHistXENMix[kNumberOfPtBins][kNumberOfCentBins][2];//XE near-side from mixed events (not finished)
197     TH1F *fHistXEA[kNumberOfPtBins][kNumberOfCentBins][2];//XE away-side (not finished)
198     TH1F *fHistXEAMix[kNumberOfPtBins][kNumberOfCentBins][2];//XE away-side from mxied events (not finished)
199
200   //Three-Particle Histograms
201     TH2F *fHistDeltaPhiPhi[kNumberOfPtBins][kNumberOfApt3Bins][kNumberOfCentBins][4][2];//DeltaPhi-DeltaPhi 3-particle distributions
202     TH2F *fHistDeltaPhiPhiMix[kNumberOfPtBins][kNumberOfApt3Bins][kNumberOfCentBins][4][2];//DeltaPhi-DeltaPhi 3-particle distributions from mixed events
203     TH2F *fHistDeltaPhiPhiSS[kNumberOfPtBins][kNumberOfApt3Bins][kNumberOfCentBins][4][2];//DeltaPhi-DeltaPhi 3-particle soft-soft background term (from mixing)
204     TH2F *fHistDeltaEtaEta[kNumberOfPtBins][kNumberOfApt3Bins][kNumberOfCentBins][4][2];//DeltaEta-DeltaEta 3-particle distributions
205     TH2F *fHistDeltaEtaEtaMix[kNumberOfPtBins][kNumberOfApt3Bins][kNumberOfCentBins][4][2];//DeltaEta-DeltaEta 3-particle distributions from Mixed Events
206     TH2F *fHistDeltaEtaEtaSS[kNumberOfPtBins][kNumberOfApt3Bins][kNumberOfCentBins][4][2];//DeltaEta-DeltaEta 3-particle distributions soft-soft background term (from mixing)
207
208
209  
210
211   //Arrays for event mixing (only need parameters of tracks passing cuts kept)
212     Float_t *fMPt[kNumberOfEventsToMix][kNumberOfCentBins][kNumberOfVertexBins][2];//Pt array (mixed events)
213     Float_t *fMPhi[kNumberOfEventsToMix][kNumberOfCentBins][kNumberOfVertexBins][2]; //Phi array (mixed events)
214     Int_t fMixTrack[kNumberOfEventsToMix][kNumberOfCentBins][kNumberOfVertexBins][2];//Number of tracks in the event stored (mixed events)
215     Float_t *fMEta[kNumberOfEventsToMix][kNumberOfCentBins][kNumberOfVertexBins][2];//Eta Array (mixed events)
216     Short_t *fMCharge[kNumberOfEventsToMix][kNumberOfCentBins][kNumberOfVertexBins][2];//Charge array (mixed events)
217     Float_t *fMEff[kNumberOfEventsToMix][kNumberOfCentBins][kNumberOfVertexBins][2];//Efficiency correction array (mixed events)
218     Int_t fMixPointer[kNumberOfCentBins][kNumberOfVertexBins][2];//which event should be replaced with newest mixed event (mixed events)
219     Int_t fMixEnd[kNumberOfCentBins][kNumberOfVertexBins][2]; //current depth of mixed event pool (mixed events)
220   //additional ones to speed up the 3-particle correlations
221   Short_t *fMPtAssoc3[kNumberOfEventsToMix][kNumberOfCentBins][kNumberOfVertexBins][2][10];//which 3-particle pt bin(s) the track is in, using a bit extra memory fixing the pt bins to 10, but makes the arrays much easier (mixed events)
222   Short_t *fMNPtAssoc3[kNumberOfEventsToMix][kNumberOfCentBins][kNumberOfVertexBins][2];//number of 3-particle pt bins the track is in (mixed events)
223   
224
225   //Arrays for Main Events
226   Float_t *ftPhi;//Phi array (current event)
227   Float_t *ftEta;//Eta array (current event)
228   Float_t *ftPt;//Pt Array (current event)
229   Short_t *ftCharge;//Charge Array (current event)
230   Float_t *ftEff;//Efficneicy Array (current event)
231   Int_t **ftPtAssoc3;//which 3-particle bins the track is in (current event)
232   Int_t *ftNPtAssoc3;//number of 3-particle bins the track is in (current event)
233
234   AliAnalysisTaskDiHadron(const AliAnalysisTaskDiHadron&);//not implimented
235   AliAnalysisTaskDiHadron& operator=(const AliAnalysisTaskDiHadron&);//not implimnted
236
237   ClassDef(AliAnalysisTaskDiHadron,1);
238   
239 };
240 #endif
241   
242