]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/FLOW/AliFlowCommonHist.cxx
Removed a cout statement
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliFlowCommonHist.cxx
CommitLineData
f1d945a1 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
f1d945a1 16#include "Riostream.h" //needed as include
17#include "AliFlowCommonConstants.h" //needed as include
18#include "AliFlowCommonHist.h"
19#include "AliFlowEventSimple.h"
20#include "AliFlowTrackSimple.h"
21
22#include "TString.h"
23#include "TProfile.h"
24#include "TMath.h" //needed as include
25#include "AliFlowVector.h"
26
27class TH1F;
28class TH1D;
29
30// AliFlowCommonHist:
31//
32// Description: Class to organise common histograms for Flow Analysis
33
34// authors: N. van der Kolk (kolk@nikhef.nl) and A. Bilandzic (anteb@nikhef.nl)
35
36
37ClassImp(AliFlowCommonHist)
38
39//-----------------------------------------------------------------------
40
41 AliFlowCommonHist::AliFlowCommonHist(TString input):
8232a5ec 42// fEvent(0),
f1d945a1 43 fTrack(0),
44 fHistMultOrig(0),
45 fHistMultInt(0),
46 fHistMultDiff(0),
47 fHistPtInt(0),
48 fHistPtDiff(0),
49 fHistPhiInt(0),
50 fHistPhiDiff(0),
51 fHistEtaInt(0),
52 fHistEtaDiff(0),
53 fHistProMeanPtperBin(0),
54 fHistQ(0)
55 {
56 //constructor creating histograms
57 Int_t fNbinsMult = AliFlowCommonConstants::GetNbinsMult();
58 Int_t fNbinsPt = AliFlowCommonConstants::GetNbinsPt();
59 Int_t fNbinsPhi = AliFlowCommonConstants::GetNbinsPhi();
60 Int_t fNbinsEta = AliFlowCommonConstants::GetNbinsEta();
61 Int_t fNbinsQ = AliFlowCommonConstants::GetNbinsQ();
62 TString name;
63
64 Double_t fMultMin = AliFlowCommonConstants::GetMultMin();
65 Double_t fMultMax = AliFlowCommonConstants::GetMultMax();
66 Double_t fPtMin = AliFlowCommonConstants::GetPtMin();
67 Double_t fPtMax = AliFlowCommonConstants::GetPtMax();
68 Double_t fPhiMin = AliFlowCommonConstants::GetPhiMin();
69 Double_t fPhiMax = AliFlowCommonConstants::GetPhiMax();
70 Double_t fEtaMin = AliFlowCommonConstants::GetEtaMin();
71 Double_t fEtaMax = AliFlowCommonConstants::GetEtaMax();
72 Double_t fQMin = AliFlowCommonConstants::GetQMin();
73 Double_t fQMax = AliFlowCommonConstants::GetQMax();
74
75 cout<<"The settings for the common histograms are as follows:"<<endl;
76 cout<<"Multiplicity: "<<fNbinsMult<<" bins between "<<fMultMin<<" and "<<fMultMax<<endl;
77 cout<<"Pt: "<<fNbinsPt<<" bins between "<<fPtMin<<" and "<<fPtMax<<endl;
78 cout<<"Phi: "<<fNbinsPhi<<" bins between "<<fPhiMin<<" and "<<fPhiMax<<endl;
79 cout<<"Eta: "<<fNbinsEta<<" bins between "<<fEtaMin<<" and "<<fEtaMax<<endl;
80 cout<<"Q: "<<fNbinsQ<<" bins between "<<fQMin<<" and "<<fQMax<<endl;
81
82 //Multiplicity
83 name = "Control_Flow_OrigMult_";
84 name +=input;
85 fHistMultOrig = new TH1F(name.Data(), name.Data(),fNbinsMult, fMultMin, fMultMax);
86 fHistMultOrig ->SetXTitle("Original Multiplicity");
87 fHistMultOrig ->SetYTitle("Counts");
88
89 name = "Control_Flow_MultInt_";
90 name +=input;
91 fHistMultInt = new TH1F(name.Data(), name.Data(),fNbinsMult, fMultMin, fMultMax);
92 fHistMultInt ->SetXTitle("Multiplicity for integrated flow");
93 fHistMultInt ->SetYTitle("Counts");
94
95 name = "Control_Flow_MultDiff_";
96 name +=input;
97 fHistMultDiff = new TH1F(name.Data(), name.Data(),fNbinsMult, fMultMin, fMultMax);
98 fHistMultDiff ->SetXTitle("Multiplicity for differential flow");
99 fHistMultDiff ->SetYTitle("Counts");
100
101 //Pt
102 name = "Control_Flow_PtInt_";
103 name +=input;
104 fHistPtInt = new TH1F(name.Data(), name.Data(),fNbinsPt, fPtMin, fPtMax);
105 fHistPtInt ->SetXTitle("Pt (GeV/c) for integrated flow");
106 fHistPtInt ->SetYTitle("Counts");
107
108 name = "Control_Flow_PtDiff_";
109 name +=input;
110 fHistPtDiff = new TH1F(name.Data(), name.Data(),fNbinsPt, fPtMin, fPtMax);
111 //binning has to be the same as for fHistProVPt! use to get Nprime!
112 fHistPtDiff ->SetXTitle("Pt (GeV/c) for differential flow");
113 fHistPtDiff ->SetYTitle("Counts");
114
115 //Phi
116 name = "Control_Flow_PhiInt_";
117 name +=input;
118 fHistPhiInt = new TH1F(name.Data(), name.Data(),fNbinsPhi, fPhiMin, fPhiMax);
119 fHistPhiInt ->SetXTitle("Phi for integrated flow");
120 fHistPhiInt ->SetYTitle("Counts");
121
122 name = "Control_Flow_PhiDiff_";
123 name +=input;
124 fHistPhiDiff = new TH1F(name.Data(), name.Data(),fNbinsPhi, fPhiMin, fPhiMax);
125 fHistPhiDiff ->SetXTitle("Phi for differential flow");
126 fHistPhiDiff ->SetYTitle("Counts");
127
128 //Eta
129 name = "Control_Flow_EtaInt_";
130 name +=input;
131 fHistEtaInt = new TH1F(name.Data(), name.Data(),fNbinsEta, fEtaMin, fEtaMax);
132 fHistEtaInt ->SetXTitle("Eta for integrated flow");
133 fHistEtaInt ->SetYTitle("Counts");
134
135 name = "Control_Flow_EtaDiff_";
136 name +=input;
137 fHistEtaDiff = new TH1F(name.Data(), name.Data(),fNbinsEta, fEtaMin, fEtaMax);
138 fHistEtaDiff ->SetXTitle("Eta for differential flow");
139 fHistEtaDiff ->SetYTitle("Counts");
140
141 //Mean Pt per pt bin
142 name = "Control_FlowPro_MeanPtperBin_";
143 name +=input;
144 fHistProMeanPtperBin = new TProfile(name.Data(), name.Data(),fNbinsPt,fPtMin,fPtMax);
145 fHistProMeanPtperBin ->SetXTitle("Pt");
146 fHistProMeanPtperBin ->SetYTitle("<Pt>");
147
148 //Q vector
149 name = "Control_Flow_Q_";
150 name +=input;
151 fHistQ = new TH1F(name.Data(), name.Data(),fNbinsQ, fQMin, fQMax);
152 fHistQ ->SetXTitle("Qvector/Mult");
153 fHistQ ->SetYTitle("Counts");
154}
155
156
157//-----------------------------------------------------------------------
158
159AliFlowCommonHist::~AliFlowCommonHist()
160{
161 //deletes histograms
162 delete fHistMultOrig;
163 delete fHistMultInt;
164 delete fHistMultDiff;
165 delete fHistPtInt;
166 delete fHistPtDiff;
167 delete fHistPhiInt;
168 delete fHistPhiDiff;
169 delete fHistEtaInt;
170 delete fHistEtaDiff;
171 delete fHistProMeanPtperBin;
172 delete fHistQ;
173}
174
175
176//-----------------------------------------------------------------------
177
8232a5ec 178Bool_t AliFlowCommonHist::FillControlHistograms(AliFlowEventSimple* Event)
f1d945a1 179{
180 //Fills the control histograms
8232a5ec 181 if (!Event){
f1d945a1 182 cout<<"##### FillControlHistograms: FlowEvent pointer null"<<endl;
183 return kFALSE;
184 }
185
186 Double_t fPt, fPhi, fEta;
187
188
189 //fill the histograms
8232a5ec 190 Int_t fNumberOfTracks = Event->NumberOfTracks();
f1d945a1 191 fHistMultOrig->Fill(fNumberOfTracks);
192
8232a5ec 193 AliFlowVector fQ = Event->GetQ();
f1d945a1 194 //weight by the Multiplicity
195 Double_t fQX = fQ.X()/fQ.GetMult();
196 Double_t fQY = fQ.Y()/fQ.GetMult();
197 fQ.Set(fQX,fQY);
198 fHistQ->Fill(fQ.Mod());
199
200 Int_t fMultInt = 0;
201 Int_t fMultDiff = 0;
202 for (Int_t i=0;i<fNumberOfTracks;i++) {
8232a5ec 203 fTrack = Event->GetTrack(i);
f1d945a1 204 if (fTrack ) {
205 if (fTrack->UseForIntegratedFlow()){
206 fPt = fTrack->Pt();
207 fHistPtInt->Fill(fPt);
208 fPhi = fTrack->Phi();
209 if (fPhi<0.) fPhi+=2*TMath::Pi();
210 fHistPhiInt->Fill(fPhi);
211 fEta = fTrack->Eta();
212 fHistEtaInt->Fill(fEta);
213 fMultInt++;
214 }
215 if (fTrack->UseForDifferentialFlow()){
216 fPt = fTrack->Pt();
217 fHistPtDiff->Fill(fPt);
218 fPhi = fTrack->Phi();
219 if (fPhi<0.) fPhi+=2*TMath::Pi();
220 fHistPhiDiff->Fill(fPhi);
221 fEta = fTrack->Eta();
222 fHistEtaDiff->Fill(fEta);
223 fHistProMeanPtperBin->Fill(fPt,fPt);
224 fMultDiff++;
225 }
226 } //track
227 } //loop over tracks
228
229 fHistMultInt->Fill(fMultInt);
230 fHistMultDiff->Fill(fMultDiff);
231
232 return kTRUE;
233}
234
235//-----------------------------------------------------------------------
236
237Double_t AliFlowCommonHist::GetEntriesInPtBin(Int_t fBin)
238{
239 //get entries in bin fBin from fHistPtDiff
240 Double_t fEntries = fHistPtDiff->GetBinContent(fBin);
241
242 return fEntries;
243
244}
245
246//-----------------------------------------------------------------------
247
248Double_t AliFlowCommonHist::GetMeanPt(Int_t fBin)
249{
250 //Get entry from bin fBin from fHistProMeanPtperBin
251 Double_t fMeanPt = fHistProMeanPtperBin->GetBinContent(fBin);
252
253 return fMeanPt;
254
0683b7d5 255}
256
257