]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/FLOW/AliFlowCommon/AliFlowEventSimple.cxx
Adding possibility to use custom ExB correction map
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliFlowCommon / AliFlowEventSimple.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"
17#include "TObjArray.h"
26c4cbb9 18#include "TFile.h"
03a02aca 19#include "TList.h"
f1d945a1 20#include "TMath.h"
26c4cbb9 21#include "TH1F.h"
22#include "TH1D.h"
23#include "TProfile.h"
c076fda8 24#include "TBrowser.h"
f1d945a1 25#include "AliFlowVector.h"
26#include "AliFlowTrackSimple.h"
27#include "AliFlowEventSimple.h"
28
26c4cbb9 29/**************************************
30 * AliFlowEventSimple: A simple event *
31 * for flow analysis *
32 * *
33 * authors: Naomi van der Kolk *
34 * (kolk@nikhef.nl) *
35 * Ante Bilandzic *
36 * (anteb@nikhef.nl) *
37 * ***********************************/
38
f1d945a1 39ClassImp(AliFlowEventSimple)
40
41//-----------------------------------------------------------------------
42
46bec39c 43AliFlowEventSimple::AliFlowEventSimple():
44 fTrackCollection(NULL),
45 fNumberOfTracks(0),
c076fda8 46 fEventNSelTracksIntFlow(0),
47 fNumberOfTracksWrap(NULL),
a12990bb 48 fEventNSelTracksIntFlowWrap(NULL),
49 fMCReactionPlaneAngleWrap(NULL)
46bec39c 50{
51 cout << "AliFlowEventSimple: Default constructor to be used only by root for io" << endl;
52}
53
54//-----------------------------------------------------------------------
55
56AliFlowEventSimple::AliFlowEventSimple(Int_t aLenght):
57 fTrackCollection(NULL),
58 fNumberOfTracks(0),
c076fda8 59 fEventNSelTracksIntFlow(0),
60 fNumberOfTracksWrap(NULL),
a12990bb 61 fEventNSelTracksIntFlowWrap(NULL),
62 fMCReactionPlaneAngleWrap(NULL)
f1d945a1 63{
64 //constructor
e35ddff0 65 fTrackCollection = new TObjArray(aLenght) ;
f1d945a1 66}
67
68//-----------------------------------------------------------------------
69
e35ddff0 70AliFlowEventSimple::AliFlowEventSimple(const AliFlowEventSimple& anEvent):
bc6b015e 71 TObject(),
e35ddff0 72 fTrackCollection(anEvent.fTrackCollection),
73 fNumberOfTracks(anEvent.fNumberOfTracks),
c076fda8 74 fEventNSelTracksIntFlow(anEvent.fEventNSelTracksIntFlow),
a12990bb 75 fMCReactionPlaneAngle(anEvent.fMCReactionPlaneAngle),
c076fda8 76 fNumberOfTracksWrap(anEvent.fNumberOfTracksWrap),
a12990bb 77 fEventNSelTracksIntFlowWrap(anEvent.fEventNSelTracksIntFlowWrap),
78 fMCReactionPlaneAngleWrap(anEvent.fMCReactionPlaneAngleWrap)
f1d945a1 79{
80 //copy constructor
f1d945a1 81}
82
83//-----------------------------------------------------------------------
84
e35ddff0 85AliFlowEventSimple& AliFlowEventSimple::operator=(const AliFlowEventSimple& anEvent)
f1d945a1 86{
26c4cbb9 87 *fTrackCollection = *anEvent.fTrackCollection ;
e35ddff0 88 fNumberOfTracks = anEvent.fNumberOfTracks;
89 fEventNSelTracksIntFlow = anEvent.fEventNSelTracksIntFlow;
a12990bb 90 fMCReactionPlaneAngle = anEvent.fMCReactionPlaneAngle;
c076fda8 91 fNumberOfTracksWrap = anEvent.fNumberOfTracksWrap;
92 fEventNSelTracksIntFlowWrap = anEvent.fEventNSelTracksIntFlowWrap;
a12990bb 93 fMCReactionPlaneAngleWrap=anEvent.fMCReactionPlaneAngleWrap;
c076fda8 94
f1d945a1 95 return *this;
f1d945a1 96}
97
f1d945a1 98//-----------------------------------------------------------------------
99
100AliFlowEventSimple::~AliFlowEventSimple()
101{
102 //destructor
c076fda8 103 if (fTrackCollection) fTrackCollection->Delete(); delete fTrackCollection;
104 if (fNumberOfTracksWrap) delete fNumberOfTracksWrap;
105 if (fEventNSelTracksIntFlowWrap) delete fEventNSelTracksIntFlowWrap;
a12990bb 106 if (fMCReactionPlaneAngleWrap) delete fMCReactionPlaneAngleWrap;
f1d945a1 107}
108
109//-----------------------------------------------------------------------
110
111AliFlowTrackSimple* AliFlowEventSimple::GetTrack(Int_t i)
112{
113 //get track i from collection
e35ddff0 114 AliFlowTrackSimple* pTrack = (AliFlowTrackSimple*)TrackCollection()->At(i) ;
115 return pTrack;
f1d945a1 116}
117
118//-----------------------------------------------------------------------
03a02aca 119AliFlowVector AliFlowEventSimple::GetQ(Int_t n, TList *weightsList, Bool_t usePhiWeights, Bool_t usePtWeights, Bool_t useEtaWeights)
f1d945a1 120{
ae733b3b 121 // calculate Q-vector in harmonic n without weights (default harmonic n=2)
e35ddff0 122 Double_t dQX = 0.;
123 Double_t dQY = 0.;
124 AliFlowVector vQ;
125 vQ.Set(0.,0.);
9825d4a9 126
127 Int_t iOrder = n;
ae733b3b 128 Double_t iUsedTracks = 0;
26c4cbb9 129 Double_t dPhi=0.;
ae733b3b 130 Double_t dPt=0.;
131 Double_t dEta=0.;
26c4cbb9 132
133 AliFlowTrackSimple* pTrack = NULL;
134
135 Int_t nBinsPhi=0;
136 Double_t dBinWidthPt=0.;
ae733b3b 137 Double_t dPtMin=0.;
26c4cbb9 138 Double_t dBinWidthEta=0.;
ae733b3b 139 Double_t dEtaMin=0.;
26c4cbb9 140
ae733b3b 141 Double_t wPhi=1.; // weight Phi
142 Double_t wPt=1.; // weight Pt
143 Double_t wEta=1.; // weight Eta
26c4cbb9 144
03a02aca 145 TH1F *phiWeights = NULL;
ae733b3b 146 TH1D *ptWeights = NULL;
03a02aca 147 TH1D *etaWeights = NULL;
148
ae733b3b 149 Double_t dSumOfWeightsToPower2 = 0.; // sum_{i=1}^{n} pow((wPhi*wPt*wEta)_i, 2)
150 Double_t dSumOfWeightsToPower3 = 0.; // sum_{i=1}^{n} pow((wPhi*wPt*wEta)_i, 3)
151 Double_t dSumOfWeightsToPower4 = 0.; // sum_{i=1}^{n} pow((wPhi*wPt*wEta)_i, 4)
152 Double_t dSumOfWeightsToPower5 = 0.; // sum_{i=1}^{n} pow((wPhi*wPt*wEta)_i, 5)
153 Double_t dSumOfWeightsToPower6 = 0.; // sum_{i=1}^{n} pow((wPhi*wPt*wEta)_i, 6)
154 Double_t dSumOfWeightsToPower7 = 0.; // sum_{i=1}^{n} pow((wPhi*wPt*wEta)_i, 7)
155 Double_t dSumOfWeightsToPower8 = 0.; // sum_{i=1}^{n} pow((wPhi*wPt*wEta)_i, 8)
156
03a02aca 157 if(weightsList)
26c4cbb9 158 {
03a02aca 159 if(usePhiWeights)
26c4cbb9 160 {
03a02aca 161 phiWeights = dynamic_cast<TH1F *>(weightsList->FindObject("phi_weights"));
162 if(phiWeights) nBinsPhi = phiWeights->GetNbinsX();
163 }
164 if(usePtWeights)
26c4cbb9 165 {
03a02aca 166 ptWeights = dynamic_cast<TH1D *>(weightsList->FindObject("pt_weights"));
167 if(ptWeights)
168 {
169 dBinWidthPt = ptWeights->GetBinWidth(1); // assuming that all bins have the same width
ae733b3b 170 dPtMin = (ptWeights->GetXaxis())->GetXmin();
03a02aca 171 }
172 }
173 if(useEtaWeights)
26c4cbb9 174 {
03a02aca 175 etaWeights = dynamic_cast<TH1D *>(weightsList->FindObject("eta_weights"));
176 if(etaWeights)
177 {
178 dBinWidthEta = etaWeights->GetBinWidth(1); // assuming that all bins have the same width
ae733b3b 179 dEtaMin = (etaWeights->GetXaxis())->GetXmin();
03a02aca 180 }
181 }
182 } // end of if(weightsList)
26c4cbb9 183
03a02aca 184 // loop over tracks
26c4cbb9 185 for(Int_t i=0;i<fNumberOfTracks;i++)
186 {
187 pTrack = (AliFlowTrackSimple*)TrackCollection()->At(i);
188 if(pTrack)
189 {
190 if(pTrack->UseForIntegratedFlow())
f1d945a1 191 {
26c4cbb9 192 dPhi = pTrack->Phi();
193 dPt = pTrack->Pt();
194 dEta = pTrack->Eta();
ae733b3b 195
196 // determine Phi weight: (to be improved, I should here only access it + the treatment of gaps in the if statement)
197 if(phiWeights && nBinsPhi)
26c4cbb9 198 {
ae733b3b 199 wPhi = phiWeights->GetBinContent(1+(Int_t)(TMath::Floor(dPhi*nBinsPhi/TMath::TwoPi())));
26c4cbb9 200 }
ae733b3b 201 // determine v'(pt) weight:
202 if(ptWeights && dBinWidthPt)
26c4cbb9 203 {
ae733b3b 204 wPt=ptWeights->GetBinContent(1+(Int_t)(TMath::Floor((dPt-dPtMin)/dBinWidthPt)));
26c4cbb9 205 }
ae733b3b 206 // determine v'(eta) weight:
207 if(etaWeights && dBinWidthEta)
26c4cbb9 208 {
ae733b3b 209 wEta=etaWeights->GetBinContent(1+(Int_t)(TMath::Floor((dEta-dEtaMin)/dBinWidthEta)));
03a02aca 210 }
ae733b3b 211
212 // building up the weighted Q-vector:
26c4cbb9 213 dQX += wPhi*wPt*wEta*TMath::Cos(iOrder*dPhi);
214 dQY += wPhi*wPt*wEta*TMath::Sin(iOrder*dPhi);
ae733b3b 215
216 // weighted multiplicity:
217 iUsedTracks+=wPhi*wPt*wEta;
218
219 // weights raised to various powers are summed up:
220 dSumOfWeightsToPower2+=pow(wPhi*wPt*wEta, 2);
221 dSumOfWeightsToPower3+=pow(wPhi*wPt*wEta, 3);
222 dSumOfWeightsToPower4+=pow(wPhi*wPt*wEta, 4);
223 dSumOfWeightsToPower5+=pow(wPhi*wPt*wEta, 5);
224 dSumOfWeightsToPower6+=pow(wPhi*wPt*wEta, 6);
225 dSumOfWeightsToPower7+=pow(wPhi*wPt*wEta, 7);
226 dSumOfWeightsToPower8+=pow(wPhi*wPt*wEta, 8);
03a02aca 227
ae733b3b 228 } // end of if (pTrack->UseForIntegratedFlow())
229 } // end of if (pTrack)
26c4cbb9 230 else {cerr << "no particle!!!"<<endl;}
ae733b3b 231 } // loop over particles
26c4cbb9 232
e35ddff0 233 vQ.Set(dQX,dQY);
234 vQ.SetMult(iUsedTracks);
ae733b3b 235 vQ.SetSumOfWeightsToPower2(dSumOfWeightsToPower2);
236 vQ.SetSumOfWeightsToPower3(dSumOfWeightsToPower3);
237 vQ.SetSumOfWeightsToPower4(dSumOfWeightsToPower4);
238 vQ.SetSumOfWeightsToPower5(dSumOfWeightsToPower5);
239 vQ.SetSumOfWeightsToPower6(dSumOfWeightsToPower6);
240 vQ.SetSumOfWeightsToPower7(dSumOfWeightsToPower7);
241 vQ.SetSumOfWeightsToPower8(dSumOfWeightsToPower8);
26c4cbb9 242
e35ddff0 243 return vQ;
f1d945a1 244
5fef318d 245}
246
c076fda8 247//-----------------------------------------------------------------------
248void AliFlowEventSimple::Print(Option_t *option) const
249{
250 // -*-*-*-*-*Print some global quantities for this histogram collection class *-*-*-*-*-*-*-*
251 // ===============================================
252 // printf( "TH1.Print Name = %s, Entries= %d, Total sum= %g\n",GetName(),Int_t(fEntries),GetSumOfWeights());
a12990bb 253 printf( "Class.Print Name = %s, Total number of tracks= %d, Number of selected tracks= %d, MC EventPlaneAngle= %d",
254 GetName(),fNumberOfTracks, fEventNSelTracksIntFlow, fMCReactionPlaneAngle );
c076fda8 255
256 if (fTrackCollection) {
257 fTrackCollection->Print(option);
258 }
259 else {
260 printf( "Empty track collection \n");
261 }
262}
263
264//-----------------------------------------------------------------------
265 void AliFlowEventSimple::Browse(TBrowser *b)
266{
267 if (!b) return;
268 if (!fNumberOfTracksWrap) {
269 fNumberOfTracksWrap = new TParameter<int>("fNumberOfTracks", fNumberOfTracks);
270 b->Add(fNumberOfTracksWrap);
271 }
272 if (!fEventNSelTracksIntFlowWrap) {
273 fEventNSelTracksIntFlowWrap = new TParameter<int>("fEventNSelTracksIntFlow", fEventNSelTracksIntFlow);
274 b->Add(fEventNSelTracksIntFlowWrap);
275 }
a12990bb 276 if (!fMCReactionPlaneAngleWrap) {
277 fMCReactionPlaneAngleWrap = new TParameter<double>(" fMCReactionPlaneAngle", fMCReactionPlaneAngle);
278 b->Add( fMCReactionPlaneAngleWrap);
279 }
c076fda8 280 if (fTrackCollection) b->Add(fTrackCollection,"AliFlowTracksSimple");
281}
282
283
5fef318d 284
26c4cbb9 285
286
287
288
289
290
291
292