]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGCF/EBYE/LRC/AliLRCProcess.cxx
Compatibility with the Root trunk
[u/mrichter/AliRoot.git] / PWGCF / EBYE / LRC / AliLRCProcess.cxx
CommitLineData
41f63c69 1/**************************************************************************\r
2 * Author: Andrey Ivanov. *\r
3 * Contributors are mentioned in the code where appropriate. *\r
4 * *\r
5 * Permission to use, copy, modify and distribute this software and its *\r
6 * documentation strictly for non-commercial purposes is hereby granted *\r
7 * without fee, provided that the above copyright notice appears in all *\r
8 * copies and that both the copyright notice and this permission notice *\r
9 * appear in the supporting documentation. The authors make no claims *\r
10 * about the suitability of this software for any purpose. It is *\r
11 * provided "as is" without express or implied warranty. *\r
12 **************************************************************************/\r
13\r
d96e5666 14// This class is creatig TH2D histogramms for Nch - Nch , Nch - Pt , Pt - Pt \r
15// dirtributions for given ETA windows and some supplementary data for Long Range Correlation (LRC) analysis . \r
16// Class is designid to work with AliAnalysisTaskLRC\r
17\r
18// Author : Andrey Ivanov , St.Peterburg State University\r
19// Email: Andrey.Ivanov@cern.ch\r
20\r
d96e5666 21#include "Riostream.h"\r
22#include "AliLRCProcess.h"\r
41f63c69 23#include "TH1F.h"\r
24#include "TH1D.h"\r
25#include "TH2D.h"\r
26#include "TProfile.h"\r
27#include "TList.h"\r
28#include "TString.h"\r
29#include "TMath.h"\r
c64cb1f6 30\r
31using std::cout;\r
32using std::endl;\r
33\r
d96e5666 34ClassImp(AliLRCProcess)\r
35\r
41f63c69 36AliLRCProcess::AliLRCProcess():fIsEventOpend(kFALSE), fIsOnline(kFALSE), fDisplayInitOnDemandWarning(kTRUE), fEventCount(0),fStartForwardETA(0), fEndForwardETA(0), fStartForwardPhi(0),fEndForwardPhi(0),fStartBakwardETA(0), fEndBakwardETA(0),fStartBackwardPhi(0),fEndBackwardPhi(0),fHiPt(0),fLoPt(0),fHiMult(0),fLoMult(0),fMultBins(0),fPtBins(0),fSumPtFw(0), fSumPtBw(0), fSumPtBw2(0), fNchFw(0), fNchBw(0),fOutList(0), fShortDef(0),fOutputSlot(0), fHistPt(0),fHistEta(0),fHistNN(0),fHistPtN(0),fHistPtPt(0),fProfNberr(0),fProfNberrPtPt(0),fProfdPtB(0),fProfTestLRC(0),fHistPtForward(0),fHistEtaForward(0),fHistNchForward(0),fHistPhiForward(0),fHistPtBakward(0),fHistEtaBakward(0),fHistNchBakward(0),fHistPhiBakward(0){};\r
37\r
38AliLRCProcess::AliLRCProcess(Double_t _StartForwardETA,Double_t _EndForwardETA,Double_t _StartBakwardETA,Double_t _EndBakwardETA): fIsEventOpend(kFALSE), fIsOnline(kFALSE), fDisplayInitOnDemandWarning(kTRUE), fEventCount(0),fStartForwardETA(0), fEndForwardETA(0), fStartForwardPhi(0),fEndForwardPhi(0),fStartBakwardETA(0), fEndBakwardETA(0),fStartBackwardPhi(0),fEndBackwardPhi(0),fHiPt(0),fLoPt(0),fHiMult(0),fLoMult(0),fMultBins(0),fPtBins(0),fSumPtFw(0), fSumPtBw(0), fSumPtBw2(0),fNchFw(0), fNchBw(0),fOutList(0), fShortDef(0),fOutputSlot(0), fHistPt(0),fHistEta(0),fHistNN(0),fHistPtN(0),fHistPtPt(0),fProfNberr(0),fProfNberrPtPt(0),fProfdPtB(0),fProfTestLRC(0),fHistPtForward(0),fHistEtaForward(0),fHistNchForward(0),fHistPhiForward(0),fHistPtBakward(0),fHistEtaBakward(0),fHistNchBakward(0),fHistPhiBakward(0)\r
d96e5666 39{\r
40// Constructor with window setup makes ready-to-run processor\r
41fEventCount=0;\r
42\r
43SetETAWindows( _StartForwardETA, _EndForwardETA,_StartBakwardETA,_EndBakwardETA);\r
41f63c69 44SetHistPtRange(0.0,4.0,200);\r
45SetHistMultRange(0,100); \r
46SetForwardWindowPhi(0,2*TMath::Pi());\r
47SetBackwardWindowPhi(0,2*TMath::Pi());\r
c8ffda94 48}\r
49\r
d96e5666 50Bool_t AliLRCProcess::InitDataMembers()\r
51{\r
52// This method is actualy creating output histogramms\r
53// Thist method is to be called in CreateOutputObjects method of AliAnalysisTask\r
41f63c69 54//cout<<" # Init for "<<fShortDef<<" this="<<this<<"\n";\r
d96e5666 55if(fIsOnline)\r
56{\r
57Printf("Can't init data members more then one time! \n");\r
58return kFALSE;\r
59}\r
60fEventCount=0;\r
61fOutList = new TList();\r
62fOutList->SetName(fShortDef);\r
63\r
41f63c69 64 Double_t loMult,hiMult;\r
d96e5666 65 \r
41f63c69 66 loMult=fLoMult-0.5;\r
67 hiMult=fHiMult+0.5;\r
d96e5666 68\r
69 // Window statistics histograms\r
70 \r
71 // Forward\r
72 \r
73 fHistPtForward = new TH1D("fHistPtForward", "P_{T} distribution in Forward window", 100, 0.0, 5);\r
74 fHistPtForward->GetXaxis()->SetTitle("P_{T} (GeV/c)");\r
75 fHistPtForward->GetYaxis()->SetTitle("dN/dP_{T} (c/GeV)");\r
76 fHistPtForward->SetMarkerStyle(kFullCircle);\r
77 \r
78 \r
79 fHistEtaForward = new TH1D("fEtaForward", "#eta distribution in Forward window", 200, -2, 2);\r
80 fHistEtaForward->GetXaxis()->SetTitle("ETA");\r
81 fHistEtaForward->GetYaxis()->SetTitle("dN/ETA");\r
82 fHistEtaForward->SetMarkerStyle(kFullCircle);\r
83 \r
41f63c69 84 \r
85 fHistNchForward = new TH1D("fHistNchForward", "N_{ch} distribution in Forward window", fMultBins, loMult, hiMult);\r
d96e5666 86 fHistNchForward->GetXaxis()->SetTitle("N_{ch}");\r
87 fHistNchForward->GetYaxis()->SetTitle("dN/dN_{ch}");\r
88 fHistNchForward->SetMarkerStyle(kFullCircle);\r
41f63c69 89 \r
90 \r
91 fHistPhiForward = new TH1D("fPhiForward", "Phi distribution in Forward window", 36, 0, 2*TMath::Pi());\r
92 fHistPhiForward->GetXaxis()->SetTitle("Phi");\r
93 fHistPhiForward->GetYaxis()->SetTitle("dN/Phi");\r
94 fHistPhiForward->SetMarkerStyle(kFullCircle);\r
95 \r
96 \r
d96e5666 97 \r
98 // Bakward\r
99 \r
100 fHistPtBakward = new TH1D("fHistPtBakward", "P_{T} distribution in Bakward window", 100, 0.0, 5);\r
101 fHistPtBakward->GetXaxis()->SetTitle("P_{T} (GeV/c)");\r
102 fHistPtBakward->GetYaxis()->SetTitle("dN/dP_{T} (c/GeV)");\r
103 fHistPtBakward->SetMarkerStyle(kFullCircle);\r
104 \r
105 \r
106 fHistEtaBakward = new TH1D("fEtaBakward", "#eta distribution in Bakward window", 200, -2, 2);\r
107 fHistEtaBakward->GetXaxis()->SetTitle("ETA");\r
108 fHistEtaBakward->GetYaxis()->SetTitle("dN/ETA");\r
109 fHistEtaBakward->SetMarkerStyle(kFullCircle);\r
110 \r
111 \r
41f63c69 112 fHistNchBakward = new TH1D("fHistNchBakward", "N_{ch} distribution in Bakward window", fMultBins, loMult, hiMult);\r
d96e5666 113 fHistNchBakward->GetXaxis()->SetTitle("P_{T} (GeV/c)");\r
114 fHistNchBakward->GetYaxis()->SetTitle("dN/dP_{T} (c/GeV)");\r
115 fHistNchBakward->SetMarkerStyle(kFullCircle);\r
116\r
41f63c69 117 fHistPhiBakward = new TH1D("fPhiBakward", "#Phi distribution in Bakward window", 36, 0, 2*TMath::Pi());\r
118 fHistPhiBakward->GetXaxis()->SetTitle("Phi");\r
119 fHistPhiBakward->GetYaxis()->SetTitle("dN/Phi");\r
120 fHistPhiBakward->SetMarkerStyle(kFullCircle);\r
121 \r
d96e5666 122\r
123 //Overal statistics histograms\r
124 \r
125 fHistPt = new TH1F("fHistPt", "P_{T} distribution", 100, 0.0, 5.0);\r
126 fHistPt->GetXaxis()->SetTitle("P_{T} (GeV/c)");\r
127 fHistPt->GetYaxis()->SetTitle("dN/dP_{T} (c/GeV)");\r
128 fHistPt->SetMarkerStyle(kFullCircle);\r
129 \r
130 \r
131 fHistEta = new TH1F("fHistEta", "#eta distribution", 200, -2, 2);\r
132 fHistEta->GetXaxis()->SetTitle("ETA");\r
133 fHistEta->GetYaxis()->SetTitle("dN/ETA");\r
134 fHistEta->SetMarkerStyle(kFullCircle);\r
135 \r
136 \r
137 \r
138 // -------- LRC histograms\r
139 \r
41f63c69 140 fHistNN = new TH2D("NN","NN",fMultBins, loMult, hiMult,fMultBins, loMult, hiMult);\r
141 fHistPtN = new TH2D("PtN","PtN",fMultBins, loMult, hiMult,fPtBins,fLoPt,fHiPt);\r
d96e5666 142 fHistPtPt = new TH2D("PtPt","PtPt",fPtBins,fLoPt,fHiPt,fPtBins,fLoPt,fHiPt);\r
41f63c69 143 fProfNberr = new TProfile("nber","nber",fMultBins, loMult, hiMult);\r
144 fProfNberrPtPt = new TProfile("nberPtPt","nberPtPt",fPtBins,fLoPt,fHiPt);\r
145 fProfdPtB = new TProfile("dPtB","Overal multievent Pt_Backward (first bin) Pt_Backward^2 (sec. bin) ",16,0.5,16.5); \r
146 fProfTestLRC = new TProfile("TestLRC","Test LRC calculaion via TProfile",fMultBins, loMult, hiMult); \r
d96e5666 147\r
148\r
149 // ---------- Adding data members to output list\r
150 \r
151 // Adding overal statistics\r
152 \r
153 fOutList->Add(fHistPt);\r
154 fOutList->Add(fHistEta);\r
155 \r
156 //Adding LRC hists\r
157 \r
158 fOutList->Add(fHistNN);\r
159 fOutList->Add(fHistPtN);\r
160 fOutList->Add(fHistPtPt);\r
161 fOutList->Add(fProfNberr);\r
41f63c69 162 fOutList->Add(fProfNberrPtPt);\r
d96e5666 163 fOutList->Add(fProfdPtB);\r
164 fOutList->Add(fProfTestLRC);\r
165\r
166 \r
167 //Adding window statistics\r
168\r
169 fOutList->Add(fHistPtForward);\r
170 fOutList->Add(fHistEtaForward);\r
171 fOutList->Add(fHistNchForward);\r
172 fOutList->Add(fHistPtBakward);\r
173 fOutList->Add(fHistEtaBakward);\r
174 fOutList->Add(fHistNchBakward);\r
41f63c69 175 fOutList->Add(fHistPhiBakward);\r
176 fOutList->Add(fHistPhiForward);\r
d96e5666 177 \r
41f63c69 178 // Adding status to dPtB\r
179 \r
180 fProfdPtB->Fill(3 , fStartForwardETA);\r
181 fProfdPtB->Fill(4 , fEndForwardETA);\r
182 fProfdPtB->Fill(5 , fStartBakwardETA);\r
183 fProfdPtB->Fill(6 , fEndBakwardETA);\r
184 fProfdPtB->Fill(7 , fStartForwardPhi);\r
185 fProfdPtB->Fill(8 , fEndForwardPhi);\r
186 fProfdPtB->Fill(9 , fStartBackwardPhi);\r
187 fProfdPtB->Fill(10 , fEndBackwardPhi);\r
188 \r
189\r
d96e5666 190\r
191\r
192fIsOnline=kTRUE;\r
193return kTRUE;\r
194}\r
195AliLRCProcess::~AliLRCProcess()\r
196{\r
197//Destructor\r
198\r
199}\r
200\r
201// --------------------------------------- Setters ------------------\r
202 void AliLRCProcess::SetShortDef()\r
203 {\r
204 // Creating task and output container name \r
41f63c69 205 char str[80];\r
884bd1c4 206 snprintf(str,80, "TaskLRCw%3.1fto%3.1fvs%3.1fto%3.1f",fStartForwardETA,fEndForwardETA,fStartBakwardETA,fEndBakwardETA);\r
41f63c69 207\r
208 fShortDef= str;\r
d96e5666 209\r
210 }\r
211\r
212 void AliLRCProcess::SetForwardWindow(Double_t StartETA,Double_t EndETA)\r
213 {\r
214 //setter for the forward eta window\r
215 fStartForwardETA=StartETA;\r
216 fEndForwardETA=EndETA;\r
217 SetShortDef();\r
218 }\r
219 void AliLRCProcess::SetBackwardWindow(Double_t StartETA,Double_t EndETA)\r
220 {\r
221 //setter for the backward eta window\r
222 fStartBakwardETA=StartETA;\r
223 fEndBakwardETA=EndETA;\r
224 SetShortDef();\r
225 }\r
226 void AliLRCProcess::SetETAWindows(Double_t _StartForwardETA,Double_t _EndForwardETA,Double_t _StartBakwardETA,Double_t _EndBakwardETA)\r
227 {\r
228 //setter for the eta windows\r
229 fStartForwardETA=_StartForwardETA;\r
230 fEndForwardETA=_EndForwardETA;\r
231 fStartBakwardETA=_StartBakwardETA;\r
232 fEndBakwardETA=_EndBakwardETA;\r
233 SetShortDef();\r
234 }\r
235\r
41f63c69 236 void AliLRCProcess::SetHistPtRange(Double_t LoPt,Double_t HiPt,Int_t PtBins)\r
d96e5666 237 {\r
41f63c69 238// Sets Pt range and number of bins for Pt axis of histos\r
239 if(fIsOnline)\r
240 {\r
241 Printf("Can't change histos paramiters after InitDataMembers() was called! \n");\r
242 return ;\r
243 }\r
d96e5666 244 // Setter for Pt range and N bins in histos\r
245 fLoPt=LoPt;\r
246 fHiPt=HiPt;\r
247 fPtBins=PtBins; \r
248 }\r
41f63c69 249 void AliLRCProcess::SetHistMultRange(Int_t LoMult,Int_t HiMult,Int_t MultBins)\r
d96e5666 250 {\r
251 // Setter for multiplicity range and N bins in histos\r
41f63c69 252 if(fIsOnline)\r
253 {\r
254 Printf("Can't change histos paramiters after InitDataMembers() was called! \n");\r
255 return ;\r
256 }\r
d96e5666 257 fLoMult=LoMult;\r
258 fHiMult=HiMult;\r
259 if(!MultBins)\r
260 {\r
261 fMultBins=fHiMult-fLoMult+1;\r
262 }else\r
263 {\r
264 fMultBins=MultBins; \r
265 }\r
266 }\r
267\r
41f63c69 268 void AliLRCProcess::SetOutputSlotNumber(Int_t SlotNumber)\r
269 {\r
270 //Sets number of output slot for LRCProcessor\r
271 fOutputSlot=SlotNumber;\r
272 }\r
d96e5666 273\r
274//________________________________________________________________________\r
275\r
276\r
277\r
41f63c69 278TList* AliLRCProcess::CreateOutput() const\r
d96e5666 279{\r
280// Creates a link to output data TList\r
281return fOutList;\r
282}\r
283\r
41f63c69 284 TString AliLRCProcess::GetShortDef() const\r
d96e5666 285{\r
286return fShortDef;\r
287}\r
288\r
41f63c69 289Int_t AliLRCProcess::GetOutputSlotNumber() const\r
290{\r
291// Returns number of output slot for LRCProcessor\r
292return fOutputSlot;\r
293}\r
294\r
d96e5666 295void AliLRCProcess::StartEvent()\r
296{\r
297// Open new Event for track by track event import\r
41f63c69 298if(fIsEventOpend) // Check if trying to open event more than once !\r
d96e5666 299 {Printf("Event is already opened! Auto finishing ! \n");\r
c64cb1f6 300 cout<<fShortDef<<": event count = "<<fEventCount<<" "<<endl;\r
d96e5666 301 Printf("NchF = %i,NchB = %i \n",fNchFw,fNchBw);\r
302 \r
303 FinishEvent();\r
304 }\r
41f63c69 305if(!fIsOnline) // Autocreating histos if InitDataMembers was not called by hand\r
306{ \r
307Printf("InitDataMembers was not called by hand ! Autocreating histos...\n");\r
308InitDataMembers();\r
309}\r
310\r
d96e5666 311fNchFw=0;\r
312fSumPtFw=0;\r
313fNchBw=0;\r
314fSumPtBw=0;\r
41f63c69 315fSumPtBw2=0;\r
d96e5666 316\r
317fIsEventOpend=kTRUE;\r
318}\r
41f63c69 319void AliLRCProcess::AddTrackForward(Double_t Pt, Double_t Eta ,Double_t Phi)\r
320{\r
321// Imports track to the event directly to Forward window\r
322if(!fIsEventOpend)\r
323 {Printf("Event is not opened!\n");\r
324 return;}\r
325\r
326fNchFw++;\r
327fSumPtFw+=Pt;\r
328fHistPtForward->Fill(Pt);\r
329fHistEtaForward->Fill(Eta);\r
330fHistPhiForward->Fill(Phi);\r
331\r
332}\r
333void AliLRCProcess::AddTrackBackward(Double_t Pt, Double_t Eta ,Double_t Phi)\r
334{\r
335// Imports track to the event directly to Backward window\r
336if(!fIsEventOpend)\r
337 {Printf("Event is not opened!\n");\r
338 return;}\r
339 \r
340fNchBw++;\r
341fSumPtBw+=Pt;\r
342fSumPtBw2+=Pt*Pt;\r
343fProfdPtB->Fill(1,Pt); \r
344fProfdPtB->Fill(2,Pt*Pt);\r
345fHistPtBakward->Fill(Pt);\r
346fHistEtaBakward->Fill(Eta);\r
347fHistPhiBakward->Fill(Phi);\r
348}\r
349\r
350\r
351\r
352void AliLRCProcess::AddTrackPtEta(Double_t Pt, Double_t Eta ,Double_t Phi)\r
d96e5666 353{\r
354//Track by track event import : Imports track to the event \r
355\r
356if(!fIsEventOpend)\r
357 {Printf("Event is not opened!\n");\r
358 return;}\r
359\r
d96e5666 360 // Glabal trak data\r
41f63c69 361 fHistPt->Fill(Pt);\r
362 fHistEta->Fill(Eta);\r
363\r
364\r
d96e5666 365 //Forward window\r
41f63c69 366 if( (fStartForwardETA<Eta)&&(Eta<fEndForwardETA))\r
367 if( (fStartForwardPhi<Phi)&&(Phi<fEndForwardPhi))\r
368 AddTrackForward(Pt,Eta,Phi);\r
369\r
d96e5666 370 //Backward window\r
41f63c69 371 if((fStartBakwardETA<Eta)&&(Eta<fEndBakwardETA))\r
372 if((fStartBackwardPhi<Phi)&&(Phi<fEndBackwardPhi))\r
373 AddTrackBackward(Pt,Eta,Phi);\r
d96e5666 374\r
375}\r
376void AliLRCProcess::FinishEvent()\r
377{\r
378// Track by track event import : Close opened event and fill event summary histos\r
379\r
380if(!fIsEventOpend)\r
381 {Printf("Event is not opened!\n");\r
382 return;}\r
383\r
384 //Filling even-total data\r
385 \r
386 fHistNN->Fill(fNchFw,fNchBw);\r
387 \r
388 if(fNchBw!=0)\r
389 {\r
390 fSumPtBw=fSumPtBw/fNchBw;\r
391 fProfTestLRC->Fill(fNchFw,fSumPtBw);\r
392 fHistPtN->Fill(fNchFw,fSumPtBw);\r
393 fProfNberr->Fill(fNchFw,1.0/fNchBw);\r
394 \r
395 if(fNchFw!=0)\r
396 {\r
397 fSumPtFw=fSumPtFw/fNchFw;\r
398 fHistPtPt->Fill(fSumPtFw,fSumPtBw);\r
41f63c69 399 fProfNberrPtPt->Fill(fSumPtFw,1.0/fNchBw);\r
400 // dPtB for PtPt\r
401 fProfdPtB->Fill(15,fSumPtBw,fNchBw);\r
402 fProfdPtB->Fill(16,fSumPtBw2/fNchBw,fNchBw);\r
d96e5666 403 }\r
404 }\r
405 \r
406 \r
407 fHistNchForward->Fill(fNchFw);\r
408 fHistNchBakward->Fill(fNchBw);\r
409\r
410fEventCount++;\r
411fIsEventOpend=kFALSE;\r
412//cout<<fShortDef<<": event count = "<<fEventCount<<" ";\r
413// Printf("NchF = %i,NchB = %i",fNchFw,fNchBw);\r
414}\r
415\r