]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGCF/FEMTOSCOPY/AliFemtoUser/AliFemtoCorrFctnDEtaDPhi.cxx
increase double the energy range for the pT and E dependent histograms
[u/mrichter/AliRoot.git] / PWGCF / FEMTOSCOPY / AliFemtoUser / AliFemtoCorrFctnDEtaDPhi.cxx
CommitLineData
76ce4b5b 1////////////////////////////////////////////////////////////////////////////////
2// //
3// AliFemtoCorrFctnDEtaDPhi - A correlation function that analyzes //
4// two particle correlations with respect to the azimuthal angle (phi) //
5// and pseudorapidity (eta) difference //
6// //
7// Authors: Adam Kisiel Adam.Kisiel@cern.ch //
8// //
9////////////////////////////////////////////////////////////////////////////////
10
11#include "AliFemtoCorrFctnDEtaDPhi.h"
12#include "AliFemtoModelHiddenInfo.h"
13//#include "AliFemtoHisto.hh"
14#include <cstdio>
15#include <TMath.h>
16
17#ifdef __ROOT__
18ClassImp(AliFemtoCorrFctnDEtaDPhi)
19#endif
20
21#define PIH 1.57079632679489656
22#define PIT 6.28318530717958623
76ce4b5b 23
24//____________________________
25AliFemtoCorrFctnDEtaDPhi::AliFemtoCorrFctnDEtaDPhi(char* title, const int& aPhiBins=20, const int& aEtaBins=20):
26 AliFemtoCorrFctn(),
27 fDPhiDEtaNumerator(0),
28 fDPhiDEtaDenominator(0),
29 fDPhiNumerator(0),
30 fDPhiDenominator(0),
31 fDCosNumerator(0),
32 fDCosDenominator(0),
33 fDoPtAnalysis(0),
34 fDPhiPtNumerator(0),
35 fDPhiPtDenominator(0),
36 fDCosPtNumerator(0),
973a91f8 37 fDCosPtDenominator(0),
38 fPhi(0),
3b97e69c 39 fEta(0),
40 fYtYtNumerator(0),
d381642d 41 fYtYtDenominator(0),
42 fIfCorrection(kNone),
43 fPtCorrectionsNum(0),
44 fPtCorrectionsDen(0),
45 fEtaCorrectionsNum(0),
46 fEtaCorrectionsDen(0)
76ce4b5b 47{
4d019bab 48
49 fphiL = (-(int)(aPhiBins/4)+0.5)*2.*TMath::Pi()/aPhiBins;
50 fphiT = 2*TMath::Pi()+(-(int)(aPhiBins/4)+0.5)*2.*TMath::Pi()/aPhiBins;
51
76ce4b5b 52 // set up numerator
53 char tTitNumD[101] = "NumDPhiDEta";
54 strncat(tTitNumD,title, 100);
4d019bab 55 fDPhiDEtaNumerator = new TH2D(tTitNumD,title,aPhiBins,fphiL,fphiT,aEtaBins,-2.0,2.0);
76ce4b5b 56 // set up denominator
57 char tTitDenD[101] = "DenDPhiDEta";
58 strncat(tTitDenD,title, 100);
4d019bab 59 fDPhiDEtaDenominator = new TH2D(tTitDenD,title,aPhiBins,fphiL,fphiT,aEtaBins,-2.0,2.0);
76ce4b5b 60
61 // set up numerator
62 char tTitNumDPhi[101] = "NumDPhi";
63 strncat(tTitNumDPhi,title, 100);
64 fDPhiNumerator = new TH1D(tTitNumDPhi,title,aPhiBins*2,-0.5*TMath::Pi(),1.5*TMath::Pi());
65 // set up denominator
66 char tTitDenDPhi[101] = "DenDPhi";
67 strncat(tTitDenDPhi,title, 100);
68 fDPhiDenominator = new TH1D(tTitDenDPhi,title,aPhiBins*2,-0.5*TMath::Pi(),1.5*TMath::Pi());
69
70 // set up numerator
71 char tTitNumDCos[101] = "NumDCos";
72 strncat(tTitNumDCos,title, 100);
73 fDCosNumerator = new TH1D(tTitNumDCos,title,aPhiBins*2,-1.0,1.0);
74 // set up denominator
75 char tTitDenDCos[101] = "DenDCos";
76 strncat(tTitDenDCos,title, 100);
77 fDCosDenominator = new TH1D(tTitDenDCos,title,aPhiBins*2,-1.0,1.0);
78
973a91f8 79 char tTitPhi[101] = "Phi";
80 strncat(tTitPhi,title, 100);
81 fPhi = new TH1D(tTitPhi,title,90,-TMath::Pi(),TMath::Pi());
82
83 char tTitEta[101] = "Eta";
84 strncat(tTitEta,title, 100);
85 fEta = new TH1D(tTitEta,title,90,-1.2,1.2);
86
3b97e69c 87 // set up numerator
88 char tTitYtNum[101] = "NumYtYt";
89 strncat(tTitYtNum,title, 100);
90 fYtYtNumerator = new TH2D(tTitYtNum,title,aPhiBins,1,5,aEtaBins,1,5);
91 // set up denominator
92 char tTitYtYtDen[101] = "DenYtYt";
93 strncat(tTitYtYtDen,title, 100);
94 fYtYtDenominator = new TH2D(tTitYtYtDen,title,aPhiBins,1,5,aEtaBins,1,5);
95
96
d381642d 97 char tTitPtCorrectionsNum[101] = "NumpT1pT2EtaPhi";
98 strncat(tTitPtCorrectionsNum,title, 100);
99 char tTitPtCorrectionsDen[101] = "DenpT1pT2EtaPhi";
100 strncat(tTitPtCorrectionsDen,title, 100);
101
102 Int_t nbins[4] = {20,20,aPhiBins,aEtaBins};
103 Double_t xmin[4] = {0,0,-0.5*TMath::Pi(),-2.0};
104 Double_t xmax[4] = {4,4,1.5*TMath::Pi(),2.0};
105
106
107 fPtCorrectionsNum = new THnSparseF(tTitPtCorrectionsNum,title,4,nbins,xmin,xmax);
108 fPtCorrectionsDen = new THnSparseF(tTitPtCorrectionsDen,title,4,nbins,xmin,xmax);
109
110 char tTitEtaCorrectionsNum[101] = "NumEta1Eta2EtaPhi";
111 strncat(tTitEtaCorrectionsNum,title, 100);
112 char tTitEtaCorrectionsDen[101] = "DenEta1Eta2EtaPhi";
113 strncat(tTitEtaCorrectionsDen,title, 100);
114
115 Double_t xmineta[4] = {-1,1,-0.5*TMath::Pi(),-2.0};
116 Double_t xmaxeta[4] = {-1,1,1.5*TMath::Pi(),2.0};
117
118 fEtaCorrectionsNum = new THnSparseF(tTitEtaCorrectionsNum,title,4,nbins,xmineta,xmaxeta);
119 fEtaCorrectionsDen = new THnSparseF(tTitEtaCorrectionsDen,title,4,nbins,xmineta,xmaxeta);
120
121 // THnSparse(const char* name, const char* title, Int_t dim,
122 // const Int_t* nbins, const Double_t* xmin, const Double_t* xmax,
123 // Int_t chunksize);
973a91f8 124
76ce4b5b 125 // to enable error bar calculation...
126 fDPhiDEtaNumerator->Sumw2();
127 fDPhiDEtaDenominator->Sumw2();
128 fDPhiNumerator->Sumw2();
129 fDPhiDenominator->Sumw2();
130 fDCosNumerator->Sumw2();
131 fDCosDenominator->Sumw2();
973a91f8 132 fPhi->Sumw2();
133 fEta->Sumw2();
3b97e69c 134 fYtYtNumerator->Sumw2();
135 fYtYtDenominator->Sumw2();
d381642d 136 fPtCorrectionsNum->Sumw2();
137 fPtCorrectionsDen->Sumw2();
76ce4b5b 138}
139
140//____________________________
141AliFemtoCorrFctnDEtaDPhi::AliFemtoCorrFctnDEtaDPhi(const AliFemtoCorrFctnDEtaDPhi& aCorrFctn) :
142 AliFemtoCorrFctn(),
143 fDPhiDEtaNumerator(0),
144 fDPhiDEtaDenominator(0),
145 fDPhiNumerator(0),
146 fDPhiDenominator(0),
147 fDCosNumerator(0),
148 fDCosDenominator(0),
149 fDoPtAnalysis(0),
150 fDPhiPtNumerator(0),
151 fDPhiPtDenominator(0),
152 fDCosPtNumerator(0),
973a91f8 153 fDCosPtDenominator(0),
154 fPhi(0),
3b97e69c 155 fEta(0),
156 fYtYtNumerator(0),
d381642d 157 fYtYtDenominator(0),
158 fIfCorrection(kNone),
159 fPtCorrectionsNum(0),
160 fPtCorrectionsDen(0),
161 fEtaCorrectionsNum(0),
162 fEtaCorrectionsDen(0)
76ce4b5b 163{
164 // copy constructor
165 if (aCorrFctn.fDPhiDEtaNumerator)
166 fDPhiDEtaNumerator = new TH2D(*aCorrFctn.fDPhiDEtaNumerator);
167 else
168 fDPhiDEtaNumerator = 0;
169 if (aCorrFctn.fDPhiDEtaDenominator)
170 fDPhiDEtaDenominator = new TH2D(*aCorrFctn.fDPhiDEtaDenominator);
171 else
172 fDPhiDEtaDenominator = 0;
173
174 if (aCorrFctn.fDPhiNumerator)
175 fDPhiNumerator = new TH1D(*aCorrFctn.fDPhiNumerator);
176 else
177 fDPhiNumerator = 0;
178 if (aCorrFctn.fDPhiDenominator)
179 fDPhiDenominator = new TH1D(*aCorrFctn.fDPhiDenominator);
180 else
181 fDPhiDenominator = 0;
182
183 if (aCorrFctn.fDCosNumerator)
184 fDCosNumerator = new TH1D(*aCorrFctn.fDCosNumerator);
185 else
186 fDCosNumerator = 0;
187 if (aCorrFctn.fDCosDenominator)
188 fDCosDenominator = new TH1D(*aCorrFctn.fDCosDenominator);
189 else
190 fDCosDenominator = 0;
191
192 if (aCorrFctn.fDPhiPtNumerator)
193 fDPhiPtNumerator = new TH2D(*aCorrFctn.fDPhiPtNumerator);
194 else
195 fDPhiPtNumerator = 0;
196 if (aCorrFctn.fDPhiPtDenominator)
197 fDPhiPtDenominator = new TH2D(*aCorrFctn.fDPhiPtDenominator);
198 else
199 fDPhiPtDenominator = 0;
200
201 if (aCorrFctn.fDCosPtNumerator)
202 fDCosPtNumerator = new TH2D(*aCorrFctn.fDCosPtNumerator);
203 else
204 fDCosPtNumerator = 0;
205 if (aCorrFctn.fDCosPtDenominator)
206 fDCosPtDenominator = new TH2D(*aCorrFctn.fDCosPtDenominator);
207 else
208 fDCosPtDenominator = 0;
973a91f8 209 if (aCorrFctn.fPhi)
210 fPhi = new TH1D(*aCorrFctn.fPhi);
211 else
212 fPhi = 0;
213 if (aCorrFctn.fEta)
214 fEta = new TH1D(*aCorrFctn.fEta);
215 else
216 fEta = 0;
76ce4b5b 217
3b97e69c 218 if (aCorrFctn.fYtYtNumerator)
219 fYtYtNumerator = new TH2D(*aCorrFctn.fDPhiDEtaDenominator);
220 else
221 fYtYtNumerator = 0;
222
223 if (aCorrFctn.fYtYtDenominator)
224 fYtYtDenominator = new TH2D(*aCorrFctn.fDPhiDEtaDenominator);
225 else
226 fYtYtDenominator = 0;
d381642d 227
4d019bab 228 fphiL = aCorrFctn.fphiL;
229 fphiT = aCorrFctn.fphiT;
230
d381642d 231// if (aCorrFctn.fPtCorrectionsNum)
232// fPtCorrectionsNum = new THnSparseF(*aCorrFctn.fPtCorrectionsNum);
233// else
234// fPtCorrectionsNum = 0;
235
236// if (aCorrFctn.fPtCorrectionsDen)
237// fPtCorrectionsDen = new THnSparseF(*aCorrFctn.fPtCorrectionsDen);
238// else
239// fPtCorrectionsDen = 0;
240
241
242
76ce4b5b 243}
244//____________________________
245AliFemtoCorrFctnDEtaDPhi::~AliFemtoCorrFctnDEtaDPhi(){
246 // destructor
247 delete fDPhiDEtaNumerator;
248 delete fDPhiDEtaDenominator;
249 delete fDPhiNumerator;
250 delete fDPhiDenominator;
251 delete fDCosNumerator;
252 delete fDCosDenominator;
253 if (fDoPtAnalysis) {
254 delete fDPhiPtNumerator;
255 delete fDPhiPtDenominator;
256 delete fDCosPtNumerator;
257 delete fDCosPtDenominator;
258 }
973a91f8 259 delete fPhi;
260 delete fEta;
3b97e69c 261
262 delete fYtYtNumerator;
263 delete fYtYtDenominator;
d381642d 264
265 delete fPtCorrectionsNum;
266 delete fPtCorrectionsDen;
267 delete fEtaCorrectionsNum;
268 delete fEtaCorrectionsDen;
76ce4b5b 269}
270//_________________________
271AliFemtoCorrFctnDEtaDPhi& AliFemtoCorrFctnDEtaDPhi::operator=(const AliFemtoCorrFctnDEtaDPhi& aCorrFctn)
272{
273 // assignment operator
274 if (this == &aCorrFctn)
275 return *this;
276
277 if (aCorrFctn.fDPhiDEtaNumerator)
278 fDPhiDEtaNumerator = new TH2D(*aCorrFctn.fDPhiDEtaNumerator);
279 else
280 fDPhiDEtaNumerator = 0;
281 if (aCorrFctn.fDPhiDEtaDenominator)
282 fDPhiDEtaDenominator = new TH2D(*aCorrFctn.fDPhiDEtaDenominator);
283 else
284 fDPhiDEtaDenominator = 0;
285
286 if (aCorrFctn.fDPhiNumerator)
287 fDPhiNumerator = new TH1D(*aCorrFctn.fDPhiNumerator);
288 else
289 fDPhiNumerator = 0;
290 if (aCorrFctn.fDPhiDenominator)
291 fDPhiDenominator = new TH1D(*aCorrFctn.fDPhiDenominator);
292 else
293 fDPhiDenominator = 0;
294
295 if (aCorrFctn.fDCosNumerator)
296 fDCosNumerator = new TH1D(*aCorrFctn.fDCosNumerator);
297 else
298 fDCosNumerator = 0;
299 if (aCorrFctn.fDCosDenominator)
300 fDCosDenominator = new TH1D(*aCorrFctn.fDCosDenominator);
301 else
302 fDCosDenominator = 0;
303
304 if (aCorrFctn.fDPhiPtNumerator)
305 fDPhiPtNumerator = new TH2D(*aCorrFctn.fDPhiPtNumerator);
306 else
307 fDPhiPtNumerator = 0;
308 if (aCorrFctn.fDPhiPtDenominator)
309 fDPhiPtDenominator = new TH2D(*aCorrFctn.fDPhiPtDenominator);
310 else
311 fDPhiPtDenominator = 0;
312
313 if (aCorrFctn.fDCosPtNumerator)
314 fDCosPtNumerator = new TH2D(*aCorrFctn.fDCosPtNumerator);
315 else
316 fDCosPtNumerator = 0;
317 if (aCorrFctn.fDCosPtDenominator)
318 fDCosPtDenominator = new TH2D(*aCorrFctn.fDCosPtDenominator);
319 else
320 fDCosPtDenominator = 0;
973a91f8 321 if (aCorrFctn.fPhi)
322 fPhi = new TH1D(*aCorrFctn.fPhi);
323 else
324 fPhi = 0;
325 if (aCorrFctn.fEta)
326 fEta = new TH1D(*aCorrFctn.fEta);
327 else
328 fEta = 0;
76ce4b5b 329
3b97e69c 330 if (aCorrFctn.fYtYtNumerator)
331 fYtYtNumerator = new TH2D(*aCorrFctn.fDPhiDEtaDenominator);
332 else
333 fYtYtNumerator = 0;
334
335 if (aCorrFctn.fYtYtDenominator)
336 fYtYtDenominator = new TH2D(*aCorrFctn.fDPhiDEtaDenominator);
337 else
338 fYtYtDenominator = 0;
339
d381642d 340 fIfCorrection = kNone;
341
4d019bab 342 fphiL = aCorrFctn.fphiL;
343 fphiT = aCorrFctn.fphiT;
344
d381642d 345// if (aCorrFctn.fPtCorrectionsNum)
346// fPtCorrectionsNum = new THnSparseF(*aCorrFctn.fPtCorrectionsNum);
347// else
348// fPtCorrectionsNum = 0;
349
350// if (aCorrFctn.fPtCorrectionsDen)
351// fPtCorrectionsDen = new THnSparseF(*aCorrFctn.fPtCorrectionsDen);
352// else
353// fPtCorrectionsDen = 0;
354
355
356
76ce4b5b 357 return *this;
358}
359//_________________________
360void AliFemtoCorrFctnDEtaDPhi::Finish(){
361 // here is where we should normalize, fit, etc...
362 // we should NOT Draw() the histos (as I had done it below),
363 // since we want to insulate ourselves from root at this level
364 // of the code. Do it instead at root command line with browser.
365 // mShareNumerator->Draw();
366 //mShareDenominator->Draw();
367 //mRatio->Draw();
368
369}
370
371//____________________________
372AliFemtoString AliFemtoCorrFctnDEtaDPhi::Report(){
373 // create report
374 string stemp = "TPC Ncls Correlation Function Report:\n";
375 char ctemp[100];
376 snprintf(ctemp , 100, "Number of entries in numerator:\t%E\n",fDPhiDEtaNumerator->GetEntries());
377 stemp += ctemp;
378 snprintf(ctemp , 100, "Number of entries in denominator:\t%E\n",fDPhiDEtaDenominator->GetEntries());
379 stemp += ctemp;
380 // stemp += mCoulombWeight->Report();
381 AliFemtoString returnThis = stemp;
382 return returnThis;
383}
384//____________________________
385void AliFemtoCorrFctnDEtaDPhi::AddRealPair( AliFemtoPair* pair){
386 // add real (effect) pair
387 if (fPairCut)
388 if (!fPairCut->Pass(pair)) return;
389
973a91f8 390 /*double phi1 = pair->Track1()->Track()->P().Phi();
76ce4b5b 391 double phi2 = pair->Track2()->Track()->P().Phi();
392 double eta1 = pair->Track1()->Track()->P().PseudoRapidity();
973a91f8 393 double eta2 = pair->Track2()->Track()->P().PseudoRapidity();*/
394
395 double phi1 = pair->Track1()->FourMomentum().Phi();
396 double phi2 = pair->Track2()->FourMomentum().Phi();
397 double eta1 = pair->Track1()->FourMomentum().PseudoRapidity();
398 double eta2 = pair->Track2()->FourMomentum().PseudoRapidity();
76ce4b5b 399
400 double dphi = phi1 - phi2;
4d019bab 401 while (dphi<fphiL) dphi+=PIT;
402 while (dphi>fphiT) dphi-=PIT;
76ce4b5b 403
404 double deta = eta1 - eta2;
405
d381642d 406 double px1 = pair->Track1()->Track()->P().x();
407 double py1 = pair->Track1()->Track()->P().y();
408 //double pz1 = pair->Track1()->Track()->P().z();
76ce4b5b 409
d381642d 410 double px2 = pair->Track2()->Track()->P().x();
411 double py2 = pair->Track2()->Track()->P().y();
412 //double pz2 = pair->Track2()->Track()->P().z();
76ce4b5b 413
d381642d 414 double pt1 = TMath::Hypot(px1, py1);
415 double pt2 = TMath::Hypot(px2, py2);
76ce4b5b 416// double ptmin = pt1>pt2 ? pt2 : pt1;
417
418// double cosphi = (px1*px2 + py1*py2 + pz1*pz2)/
419// sqrt((px1*px1 + py1*py1 + pz1*pz1)*(px2*px2 + py2*py2 + pz2*pz2));
420
421 fDPhiDEtaNumerator->Fill(dphi, deta);
422
423 fDPhiNumerator->Fill(dphi);
424// fDCosNumerator->Fill(cosphi);
425
426 if (fDoPtAnalysis) {
427// fDPhiPtNumerator->Fill(dphi, ptmin);
428// fDCosPtNumerator->Fill(cosphi, ptmin);
429 }
430
973a91f8 431 fPhi->Fill(phi1);
432 fEta->Fill(eta1);
433
3b97e69c 434 double PionMass = 0.13956995;
d381642d 435 double yt1 = TMath::Log(sqrt(1+(pt1/PionMass)*(pt1/PionMass))+(pt1/PionMass));
436 double yt2 = TMath::Log(sqrt(1+(pt2/PionMass)*(pt2/PionMass))+(pt2/PionMass));
3b97e69c 437 fYtYtNumerator->Fill(yt1,yt2);
438
d381642d 439 if(fIfCorrection)
440 {
441 if(fIfCorrection == kPt){
442 Double_t val[] = {pt1,pt2,dphi,deta};
443 fPtCorrectionsNum->Fill(val);
444 }
445 if(fIfCorrection == kEta){
446 Double_t val[] = {eta1,eta2,dphi,deta};
447 fEtaCorrectionsNum->Fill(val);
448 }
449
450 }
451
76ce4b5b 452}
453//____________________________
454void AliFemtoCorrFctnDEtaDPhi::AddMixedPair( AliFemtoPair* pair){
455 // add mixed (background) pair
456 if (fPairCut)
457 if (!fPairCut->Pass(pair)) return;
458
973a91f8 459 /*double phi1 = pair->Track1()->Track()->P().Phi();
76ce4b5b 460 double phi2 = pair->Track2()->Track()->P().Phi();
461 double eta1 = pair->Track1()->Track()->P().PseudoRapidity();
973a91f8 462 double eta2 = pair->Track2()->Track()->P().PseudoRapidity();*/
463
464 double phi1 = pair->Track1()->FourMomentum().Phi();
465 double phi2 = pair->Track2()->FourMomentum().Phi();
466 double eta1 = pair->Track1()->FourMomentum().PseudoRapidity();
467 double eta2 = pair->Track2()->FourMomentum().PseudoRapidity();
76ce4b5b 468
469 double dphi = phi1 - phi2;
4d019bab 470 while (dphi<fphiL) dphi+=PIT;
471 while (dphi>fphiT) dphi-=PIT;
76ce4b5b 472
473 double deta = eta1 - eta2;
474
d381642d 475 double px1 = pair->Track1()->Track()->P().x();
476 double py1 = pair->Track1()->Track()->P().y();
477 //double pz1 = pair->Track1()->Track()->P().z();
76ce4b5b 478
d381642d 479 double px2 = pair->Track2()->Track()->P().x();
480 double py2 = pair->Track2()->Track()->P().y();
481 //double pz2 = pair->Track2()->Track()->P().z();
76ce4b5b 482
d381642d 483 double pt1 = TMath::Hypot(px1, py1);
484 double pt2 = TMath::Hypot(px2, py2);
76ce4b5b 485// double ptmin = pt1>pt2 ? pt2 : pt1;
486
487// double cosphi = (px1*px2 + py1*py2 + pz1*pz2)/
488// sqrt((px1*px1 + py1*py1 + pz1*pz1)*(px2*px2 + py2*py2 + pz2*pz2));
489
490 fDPhiDEtaDenominator->Fill(dphi, deta);
491
492 fDPhiDenominator->Fill(dphi);
493// fDCosDenominator->Fill(cosphi);
494
3b97e69c 495 //if (fDoPtAnalysis) {
76ce4b5b 496 // fDPhiPtDenominator->Fill(dphi, ptmin);
497 // fDCosPtDenominator->Fill(cosphi, ptmin);
3b97e69c 498 //}
499
500 double PionMass = 0.13956995;
d381642d 501 double yt1 = TMath::Log(sqrt(1+(pt1/PionMass)*(pt1/PionMass))+(pt1/PionMass));
502 double yt2 = TMath::Log(sqrt(1+(pt2/PionMass)*(pt2/PionMass))+(pt2/PionMass));
503 fYtYtDenominator->Fill(yt1,yt2);
504
505 if(fIfCorrection)
506 {
507 if(fIfCorrection == kPt){
508 Double_t val[] = {pt1,pt2,dphi,deta};
509 fPtCorrectionsDen->Fill(val);
510 }
511 if(fIfCorrection == kEta){
512 Double_t val[] = {eta1,eta2,dphi,deta};
513 fEtaCorrectionsDen->Fill(val);
514 }
515 }
516
76ce4b5b 517}
518
519
520void AliFemtoCorrFctnDEtaDPhi::WriteHistos()
521{
522 // Write out result histograms
523 fDPhiDEtaNumerator->Write();
524 fDPhiDEtaDenominator->Write();
973a91f8 525 /*fDPhiNumerator->Write();
76ce4b5b 526 fDPhiDenominator->Write();
527 fDCosNumerator->Write();
528 fDCosDenominator->Write();
529 if (fDoPtAnalysis) {
530 fDPhiPtNumerator->Write();
531 fDPhiPtDenominator->Write();
532 fDCosPtNumerator->Write();
533 fDCosPtDenominator->Write();
973a91f8 534 }*/
535 fPhi->Write();
536 fEta->Write();
3b97e69c 537
d381642d 538 if(fIfCorrection){
539 if(fIfCorrection==kPt){
540 fPtCorrectionsNum->Write();
541 fPtCorrectionsDen->Write();}
542 if(fIfCorrection==kEta){
543 fEtaCorrectionsNum->Write();
544 fEtaCorrectionsDen->Write();}
545 }
76ce4b5b 546}
547
548TList* AliFemtoCorrFctnDEtaDPhi::GetOutputList()
549{
550 // Prepare the list of objects to be written to the output
551 TList *tOutputList = new TList();
552
553 tOutputList->Add(fDPhiDEtaNumerator);
554 tOutputList->Add(fDPhiDEtaDenominator);
973a91f8 555 /*tOutputList->Add(fDPhiNumerator);
76ce4b5b 556 tOutputList->Add(fDPhiDenominator);
557 tOutputList->Add(fDCosNumerator);
558 tOutputList->Add(fDCosDenominator);
559 if (fDoPtAnalysis) {
560 tOutputList->Add(fDPhiPtNumerator);
561 tOutputList->Add(fDPhiPtDenominator);
562 tOutputList->Add(fDCosPtNumerator);
563 tOutputList->Add(fDCosPtDenominator);
973a91f8 564 }*/
565 tOutputList->Add(fPhi);
566 tOutputList->Add(fEta);
3b97e69c 567 tOutputList->Add(fYtYtNumerator);
568 tOutputList->Add(fYtYtDenominator);
d381642d 569
570 if(fIfCorrection){
571 if(fIfCorrection==kPt){
572 tOutputList->Add(fPtCorrectionsNum);
573 tOutputList->Add(fPtCorrectionsDen);
574 }
575 if(fIfCorrection==kEta){
576 tOutputList->Add(fEtaCorrectionsNum);
577 tOutputList->Add(fEtaCorrectionsDen);
578 }
579 }
76ce4b5b 580 return tOutputList;
581
582}
583
584void AliFemtoCorrFctnDEtaDPhi::SetDoPtAnalysis(int do2d)
585{
586 fDoPtAnalysis = do2d;
587
588 int aPhiBins = fDPhiDEtaNumerator->GetNbinsX();
589 const char *title = fDPhiDEtaNumerator->GetTitle();
590
591 // set up numerator
592 char tTitNumDPhiPt[101] = "NumDPhiPt";
593 strncat(tTitNumDPhiPt,title, 100);
4d019bab 594 fDPhiPtNumerator = new TH2D(tTitNumDPhiPt,title,aPhiBins*2,-0.5*TMath::Pi(),3./2.*TMath::Pi(), 30, 0.0, 3.0);
76ce4b5b 595 // set up denominator
596 char tTitDenDPhiPt[101] = "DenDPhiPt";
597 strncat(tTitDenDPhiPt,title, 100);
4d019bab 598 fDPhiPtDenominator = new TH2D(tTitDenDPhiPt,title,aPhiBins*2,-0.5*TMath::Pi(),3./2.*TMath::Pi(), 30, 0.0, 3.0);
76ce4b5b 599
600 // set up numerator
601 char tTitNumDCosPt[101] = "NumDCosPt";
602 strncat(tTitNumDCosPt,title, 100);
603 fDCosPtNumerator = new TH2D(tTitNumDCosPt,title,aPhiBins*2,-1.0,1.0, 30, 0.0, 3.0);
604 // set up denominator
605 char tTitDenDCosPt[101] = "DenDCosPt";
606 strncat(tTitDenDCosPt,title, 100);
607 fDCosPtDenominator = new TH2D(tTitDenDCosPt,title,aPhiBins*2,-1.0,1.0, 30, 0.0, 3.0);
608
609 fDPhiPtNumerator->Sumw2();
610 fDPhiPtDenominator->Sumw2();
611 fDCosPtNumerator->Sumw2();
612 fDCosPtDenominator->Sumw2();
613
614}
d381642d 615
616void AliFemtoCorrFctnDEtaDPhi::SetDoCorrections(CorrectionType doCorr)
617{
618 fIfCorrection = doCorr;
619}