]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/EBYE/PIDFluctuation/task/AliEbyEPidRatioPhy.cxx
Update from PR Analysis: drathee
[u/mrichter/AliRoot.git] / PWGCF / EBYE / PIDFluctuation / task / AliEbyEPidRatioPhy.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: ALICE Offline.                                                 *
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
17 //=========================================================================//
18 //             AliEbyE Analysis for Particle Ratio Fluctuation             //
19 //                   Deepika Rathee  | Satyajit Jena                       //
20 //                   drathee@cern.ch | sjena@cern.ch                       //
21 //                  Date: Wed Jul  9 18:38:30 CEST 2014                    //
22 //          New approch to find particle ratio to reduce memory            //
23 //                             (Test Only)                                 //
24 //=========================================================================//
25
26 #include "TMath.h"
27 #include "TAxis.h"
28 #include "TProfile.h" 
29 #include "TProfile2D.h" 
30 #include "TH2D.h" 
31 #include "TH3D.h" 
32
33 #include "AliStack.h"
34 #include "AliMCEvent.h"
35 #include "AliMCParticle.h"
36 #include "AliESDEvent.h"
37 #include "AliESDtrackCuts.h"
38
39 #include "AliCentrality.h"
40 #include "AliTracker.h"
41
42 #include "AliAODEvent.h"
43 #include "AliAODTrack.h"
44 #include "AliAODMCParticle.h"
45
46 #include "AliEbyEPidRatioPhy.h"
47
48 using namespace std;
49
50 ClassImp(AliEbyEPidRatioPhy)
51 //________________________________________________________________________
52 AliEbyEPidRatioPhy::AliEbyEPidRatioPhy() : 
53   AliEbyEPidRatioBase("Dist", "Dist"),
54   fOutList(NULL),
55   fOrder(8),
56   fNNp(6),
57   fNp(NULL),
58   fNpPt(NULL),
59   fNMCNp(7),
60   fMCNp(NULL),
61   fMCNpPt(NULL),
62   fRedFactp(NULL),
63   fPtBinHist(NULL) {
64   AliLog::SetClassDebugLevel("AliEbyEPidRatioPhy",10);
65 }
66
67 //________________________________________________________________________
68 AliEbyEPidRatioPhy::~AliEbyEPidRatioPhy() {
69   // Destructor
70
71  
72   for (Int_t ii = 0; ii < fNNp; ++ii) {
73     for (Int_t kk = 0; kk < 2; ++kk)
74       if (fNp[ii][kk]) delete[] fNp[ii][kk];
75     if (fNp[ii]) delete[] fNp[ii];
76   }
77   if (fNp) delete[] fNp;
78
79
80  for (Int_t ii = 0; ii < fNMCNp; ++ii) {
81     for (Int_t kk = 0; kk < 2; ++kk)
82       if (fMCNp[ii][kk]) delete[] fMCNp[ii][kk];
83     if (fMCNp[ii]) delete[] fMCNp[ii];
84   }
85   if (fMCNp) delete[] fMCNp;
86
87   for (Int_t ii = 0; ii <= fOrder; ++ii) 
88     if (fRedFactp[ii]) delete[] fRedFactp[ii];
89   if (fRedFactp) delete[] fRedFactp;
90
91   if(fPtBinHist) delete fPtBinHist;
92
93  for (Int_t ii = 0; ii < fNNp; ++ii) {
94     for (Int_t kk = 0; kk < 2; ++kk)
95       for (Int_t jj = 0; jj < 2; ++jj) {
96         if (fNpPt[ii][kk][jj]) delete[] fNpPt[ii][kk][jj];
97         if(fNpPt[ii][kk]) delete[] fNpPt[ii][kk];
98       }
99     if (fNpPt[ii]) delete[] fNpPt[ii];
100   }
101   if (fNpPt) delete[] fNpPt;
102
103
104
105  for (Int_t ii = 0; ii < fNNp; ++ii) {
106     for (Int_t kk = 0; kk < 2; ++kk)
107       for (Int_t jj = 0; jj < 2; ++jj) {
108         if (fMCNpPt[ii][kk][jj]) delete[] fMCNpPt[ii][kk][jj];
109         if(fMCNpPt[ii][kk]) delete[] fMCNpPt[ii][kk];
110       }
111     if (fMCNpPt[ii]) delete[] fMCNpPt[ii];
112   }
113   if (fMCNpPt) delete[] fMCNpPt;
114
115
116
117
118
119   return;
120 }
121
122 //________________________________________________________________________
123 void AliEbyEPidRatioPhy::Process() {
124   ProcessTracks();
125   if (fIsMC)
126     ProcessParticles();
127   return;
128 }
129
130 //________________________________________________________________________
131 void AliEbyEPidRatioPhy::Init() {
132   fNp = new Int_t**[fNNp];
133   for (Int_t ii = 0 ; ii < fNNp; ++ii) {
134     fNp[ii] = new Int_t*[4];
135     for (Int_t kk = 0 ; kk < 4; ++kk)
136       fNp[ii][kk] = new Int_t[2];
137   }
138
139  fNpPt = new Int_t***[fNNp];
140   for (Int_t ii = 0 ; ii < fNNp; ++ii) {
141     fNpPt[ii] = new Int_t**[4];
142     for (Int_t kk = 0 ; kk < 4; ++kk) {
143       fNpPt[ii][kk] = new Int_t*[2];
144       for (Int_t ll = 0 ; ll < 2; ++ll)
145         fNpPt[ii][kk][ll] = new Int_t[AliEbyEPidRatioHelper::fgkfHistNBinsPt];
146     }
147   }
148
149   fMCNp = new Int_t**[fNMCNp];
150   for (Int_t ii = 0 ; ii < fNMCNp; ++ii) {
151     fMCNp[ii] = new Int_t*[4];
152     for (Int_t kk = 0 ; kk < 4; ++kk)
153       fMCNp[ii][kk] = new Int_t[2];
154   }
155   
156  fMCNpPt = new Int_t***[fNMCNp];
157   for (Int_t ii = 0 ; ii < fNMCNp; ++ii) {
158     fMCNpPt[ii] = new Int_t**[4];
159     for (Int_t kk = 0 ; kk < 4; ++kk) {
160       fMCNpPt[ii][kk] = new Int_t*[2];
161       for (Int_t ll = 0 ; ll < 2; ++ll)
162         fMCNpPt[ii][kk][ll] = new Int_t[AliEbyEPidRatioHelper::fgkfHistNBinsPt];
163     }
164   }
165
166   fRedFactp = new Double_t*[fOrder+1];
167   for (Int_t ii = 0 ; ii <= fOrder; ++ii)
168     fRedFactp[ii] = new Double_t[2];
169 }
170
171 //________________________________________________________________________
172 void AliEbyEPidRatioPhy::Reset() {
173   for (Int_t ii = 0; ii < fNNp; ++ii) 
174     for (Int_t kk = 0; kk < 4; ++kk)
175       for (Int_t jj = 0; jj < 2; ++jj)
176         fNp[ii][kk][jj] = 0;
177   for (Int_t ii = 0; ii < fNMCNp; ++ii) 
178     for (Int_t kk = 0; kk < 4; ++kk)
179       for (Int_t jj = 0; jj < 2; ++jj)
180         fMCNp[ii][kk][jj] = 0;
181   for (Int_t ii = 0; ii <= fOrder; ++ii) 
182     for (Int_t jj = 0; jj < 2; ++jj)
183       fRedFactp[ii][jj] = 1.;
184
185   for (Int_t ii = 0; ii < fNNp; ++ii) 
186     for (Int_t kk = 0; kk < 4; ++kk)
187       for (Int_t jj = 0; jj < 2; ++jj)
188         for (Int_t ll = 0; ll < AliEbyEPidRatioHelper::fgkfHistNBinsPt; ++ll)
189           fNpPt[ii][kk][jj][ll] = 0;
190
191  for (Int_t ii = 0; ii < fNMCNp; ++ii) 
192     for (Int_t kk = 0; kk < 4; ++kk)
193       for (Int_t jj = 0; jj < 2; ++jj)
194         for (Int_t ll = 0; ll < AliEbyEPidRatioHelper::fgkfHistNBinsPt; ++ll)
195           fMCNpPt[ii][kk][jj][ll] = 0;
196
197 }
198
199 //________________________________________________________________________
200 void AliEbyEPidRatioPhy::CreateHistograms() {
201
202   Float_t etaRange[2];
203   fESDTrackCuts->GetEtaRange(etaRange[0],etaRange[1]);
204
205   Float_t ptRange[2];
206   fESDTrackCuts->GetPtRange(ptRange[0],ptRange[1]);
207   TString sTitle("");
208   fPtBinHist = new TH1F("hPtBinHist","Make the pT Bins",AliEbyEPidRatioHelper::fgkfHistNBinsPt, AliEbyEPidRatioHelper::fgkfHistRangePt[0], AliEbyEPidRatioHelper::fgkfHistRangePt[1]);
209   if (fIsRatio) AddHistSetRatio("Ratio",       Form("%s, #it{p}_{T} [%.1f,%.1f]", sTitle.Data(), ptRange[0], ptRange[1]));
210   
211   AddHistSetCent("Phy",       Form("%s, #it{p}_{T} [%.1f,%.1f]", sTitle.Data(), ptRange[0], ptRange[1]));
212   // AddHistSetCent("PhyTPC",    Form("%s, #it{p}_{T} [%.1f,%.1f]", sTitle.Data(), ptRange[0], fHelper->GetMinPtForTOFRequired()));
213   // AddHistSetCent("PhyTOF",    Form("%s, #it{p}_{T} [%.1f,%.1f]", sTitle.Data(), fHelper->GetMinPtForTOFRequired(), ptRange[1]));
214   if (fIsPtBin) AddHistSetCentPt("PhyBin",    Form("%s, #it{p}_{T} [%.1f,%.1f]", sTitle.Data(), ptRange[0], ptRange[1]));  
215
216
217 #if USE_PHI
218   AddHistSetCent("Phyphi",    Form("%s,#it{p}_{T} [%.1f,%.1f], #varphi [%.2f,%.2f]", sTitle.Data(), ptRange[0], ptRange[1], fHelper->GetPhiMin(), fHelper->GetPhiMax()));
219   //  AddHistSetCent("PhyTPCphi", Form("%s, #it{p}_{T} [%.1f,%.1f], #varphi [%.2f,%.2f]",sTitle.Data(), ptRange[0], fHelper->GetMinPtForTOFRequired(), fHelper->GetPhiMin(), fHelper->GetPhiMax()));
220   // AddHistSetCent("PhyTOFphi", Form("%s, #it{p}_{T} [%.1f,%.1f], #varphi [%.2f,%.2f]",sTitle.Data(), fHelper->GetMinPtForTOFRequired(), ptRange[1], fHelper->GetPhiMin(), fHelper->GetPhiMax()));
221   if (fIsPtBin) AddHistSetCentPt("PhyBinPhi", Form("%s, #it{p}_{T} [%.1f,%.1f], #varphi [%.2f,%.2f]", 
222                                      sTitle.Data(), ptRange[0], ptRange[1], fHelper->GetPhiMin(), fHelper->GetPhiMax()));
223 #endif
224
225   if (fIsMC) {
226     TString sMCTitle("");
227   
228     if (fIsRatio) AddHistSetRatio("MCRatio",       Form("%s, #it{p}_{T} [%.1f,%.1f]", sTitle.Data(), ptRange[0], ptRange[1]));
229
230     AddHistSetCent("MC",      Form("%s", sTitle.Data()));
231     AddHistSetCent("MCpt",    Form("%s, #it{p}_{T} [%.1f,%.1f]", sMCTitle.Data(), ptRange[0], ptRange[1]));
232     // AddHistSetCent("MCTPC",   Form("%s, #it{p}_{T} [%.1f,%.1f]", sMCTitle.Data(), ptRange[0], fHelper->GetMinPtForTOFRequired()));
233     //AddHistSetCent("MCTOF",   Form("%s, #it{p}_{T} [%.1f,%.1f]", sMCTitle.Data(), fHelper->GetMinPtForTOFRequired(), ptRange[1]));
234     if (fIsPtBin) AddHistSetCentPt("MCBin",    Form("%s, #it{p}_{T} [%.1f,%.1f]", sMCTitle.Data(), ptRange[0], ptRange[1])); 
235 #if USE_PHI
236     AddHistSetCent("MCPhi",   Form("%s, #it{p}_{T} [%.1f,%.1f], #varphi [%.2f,%.2f]",sMCTitle.Data(), ptRange[0], ptRange[1], fHelper->GetPhiMin(), fHelper->GetPhiMax()));
237     //AddHistSetCent("MCTPCphi",Form("%s, #it{p}_{T} [%.1f,%.1f], #varphi [%.2f,%.2f]",sMCTitle.Data(), ptRange[0], fHelper->GetMinPtForTOFRequired(), fHelper->GetPhiMin(), fHelper->GetPhiMax()));
238     //AddHistSetCent("MCTOFphi",Form("%s, #it{p}_{T} [%.1f,%.1f], #varphi [%.2f,%.2f]",sMCTitle.Data(), fHelper->GetMinPtForTOFRequired(), ptRange[1], fHelper->GetPhiMin(), fHelper->GetPhiMax()));
239  if (fIsPtBin) AddHistSetCentPt("MCBinPhi", Form("%s, #it{p}_{T} [%.1f,%.1f], #varphi [%.2f,%.2f]", 
240                                         sMCTitle.Data(), ptRange[0], ptRange[1], fHelper->GetPhiMin(), fHelper->GetPhiMax()));
241 #endif
242     
243   }
244   
245   return;
246 }
247
248 //________________________________________________________________________
249 Int_t AliEbyEPidRatioPhy::ProcessTracks() {
250   Float_t etaRange[2];
251   fESDTrackCuts->GetEtaRange(etaRange[0],etaRange[1]);
252   
253   Float_t ptRange[2];
254   fESDTrackCuts->GetPtRange(ptRange[0],ptRange[1]);
255   // -- Track Loop
256   for (Int_t idxTrack = 0; idxTrack < fNTracks; ++idxTrack) {
257     AliVTrack *track = (fESD) ? static_cast<AliVTrack*>(fESD->GetTrack(idxTrack)) : static_cast<AliVTrack*>(fAOD->GetTrack(idxTrack)); 
258     // -- Check if track is accepted for basic parameters
259     if (!fHelper->IsTrackAcceptedBasicCharged(track))
260       continue;
261     // -- Check if accepted - ESD
262     if (fESD && !fESDTrackCuts->AcceptTrack(dynamic_cast<AliESDtrack*>(track)))
263       continue;
264     // -- Check if accepted - AOD
265     if (fAOD){
266       AliAODTrack * trackAOD = dynamic_cast<AliAODTrack*>(track);
267       
268       if (!trackAOD) {
269         AliError("Pointer to dynamic_cast<AliAODTrack*>(track) = ZERO");
270         continue;
271       }
272       if (!trackAOD->TestFilterBit(fAODtrackCutBit))
273         continue;
274       if(!(track->Pt() > ptRange[0] && track->Pt() <= ptRange[1] && TMath::Abs(track->Eta()) <= etaRange[1]))
275         continue;
276     }
277
278     if (!fHelper->IsTrackAcceptedDCA(track))
279       continue;
280     
281     Int_t iPid = 0;
282     Double_t pid[3];
283     if      (fHelper->IsTrackAcceptedPID(track, pid, (AliPID::kPion)))   iPid = 1;
284     else if (fHelper->IsTrackAcceptedPID(track, pid, (AliPID::kKaon)))   iPid = 2;
285     else if (fHelper->IsTrackAcceptedPID(track, pid, (AliPID::kProton))) iPid = 3;
286     else iPid = 0;
287    
288     //cout << idxTrack <<" --- PHY ---- " << iPid  << endl;
289
290
291     Double_t yP;
292     if (iPid != 0 && !fHelper->IsTrackAcceptedRapidity(track, yP, iPid))
293       continue;
294     
295     Int_t idxPart = (track->Charge() < 0) ? 0 : 1;
296     // -- in pt Range
297     fNp[0][0][idxPart] += 1;
298     if(iPid != 0) fNp[0][iPid][idxPart] += 1;
299
300     Int_t idxPt = fPtBinHist->FindBin(track->Pt()) -1;
301
302     if (idxPt>=0 && idxPt < AliEbyEPidRatioHelper::fgkfHistNBinsPt) {
303       fNpPt[0][0][idxPart][idxPt] += 1;
304       if(iPid != 0) fNpPt[0][iPid][idxPart][idxPt] += 1;
305     }
306
307     // -- in TPC pt Range
308     /*    if (track->Pt() <= fHelper->GetMinPtForTOFRequired()) {
309       fNp[1][0][idxPart] += 1;
310       if(iPid != 0) fNp[1][iPid][idxPart] += 1;
311     }
312     // -- in TPC+TOF pt Range
313     if (track->Pt() > fHelper->GetMinPtForTOFRequired()){
314       fNp[2][0][idxPart] += 1;
315       if(iPid != 0) fNp[2][iPid][idxPart] += 1;
316       }*/
317
318 #if USE_PHI
319     if(!fHelper->IsTrackAcceptedPhi(track))
320       continue;
321     
322     // -- in pt Range
323     fNp[3][iPid][idxPart] += 1;
324     if(iPid != 0) fNp[3][0][idxPart] += 1;
325
326     if (idxPt>=0 && idxPt < AliEbyEPidRatioHelper::fgkfHistNBinsPt) {
327       fNpPt[1][0][idxPart][idxPt] += 1;
328       if(iPid != 0) fNpPt[1][iPid][idxPart][idxPt] += 1;
329     }
330
331     // -- in TPC pt Range
332     /*   if (track->Pt() <= fHelper->GetMinPtForTOFRequired()) {
333       fNp[4][0][idxPart] += 1;
334       if(iPid != 0)fNp[4][iPid][idxPart] += 1;
335     }
336     // -- in TPC+TOF pt Range
337     if (track->Pt() > fHelper->GetMinPtForTOFRequired()) {
338       fNp[5][0][idxPart] += 1;
339       if(iPid != 0) fNp[5][iPid][idxPart] += 1;*/
340     }
341 #endif
342   } // for (Int_t idxTrack = 0; idxTrack < fESD->GetNumberOfTracks(); ++idxTrack) {
343  
344 FillHistSetCent("Phy",        0, kFALSE);
345 if (fIsRatio) FillHistSetRatio("Ratio",        0, kFALSE);
346 //  FillHistSetCent("PhyTPC",     1, kFALSE);
347 //  FillHistSetCent("PhyTOF",     2, kFALSE);
348 if (fIsPtBin)  FillHistSetCentPt("PhyBin",    0, kFALSE);
349 #if USE_PHI
350   FillHistSetCent("Phyphi",     3, kFALSE);
351 //  FillHistSetCent("PhyTPCphi",  4, kFALSE);
352 // FillHistSetCent("PhyTOFphi",  5, kFALSE);
353 if (fIsPtBin) FillHistSetCentPt("PhyBinPhi", 1, kFALSE);
354  #endif
355
356   return 0;
357 }
358
359 //________________________________________________________________________
360 Int_t AliEbyEPidRatioPhy::ProcessParticles() {
361   // -- Process primary particles from the stack and fill histograms
362   
363   Float_t etaRange[2];
364   fESDTrackCuts->GetEtaRange(etaRange[0],etaRange[1]);
365
366   Float_t ptRange[2];
367   fESDTrackCuts->GetPtRange(ptRange[0],ptRange[1]);
368
369   for (Int_t idxMC = 0; idxMC < fStack->GetNprimary(); ++idxMC) {
370     AliVParticle* particle = (fESD) ? fMCEvent->GetTrack(idxMC) : NULL;
371
372     if (!particle) 
373       continue;
374     if (!fHelper->IsParticleAcceptedBasicCharged(particle, idxMC))
375       continue;
376     
377     Int_t iPid = 0;  
378     if      (TMath::Abs(particle->PdgCode()) ==  211) iPid = 1; // pion
379     else if (TMath::Abs(particle->PdgCode()) ==  321) iPid = 2; // kaon
380     else if (TMath::Abs(particle->PdgCode()) == 2212) iPid = 3; // proton
381     else    iPid = 0;
382     
383     Double_t yMC;
384     if ((iPid != 0) && !fHelper->IsParticleAcceptedRapidity(particle, yMC, iPid))
385       continue;
386
387     // -- Check eta window -- for charged particles
388     if ((iPid == 0) && TMath::Abs(particle->Eta()) > etaRange[1])
389       continue;
390     
391     Int_t idxPart = (particle->PdgCode() < 0) ? 0 : 1;
392
393     fMCNp[0][0][idxPart]    += 1.;        
394     if(iPid != 0) fMCNp[0][iPid][idxPart] += 1.;        
395     
396     // -- Check main pt window
397     if (!(particle->Pt() > ptRange[0] && particle->Pt() <= ptRange[1]))
398       continue;
399     
400     // -- in pt Range
401     fMCNp[1][0][idxPart]    += 1.;        
402     if(iPid != 0)fMCNp[1][iPid][idxPart] += 1.;        
403     
404
405     Int_t idxPt = fPtBinHist->FindBin(particle->Pt()) -1;
406     
407     if (idxPt>=0 && idxPt < AliEbyEPidRatioHelper::fgkfHistNBinsPt) {
408       fMCNpPt[0][0][idxPart][idxPt] += 1;
409       if(iPid != 0) fMCNpPt[0][iPid][idxPart][idxPt] += 1;
410     }
411
412
413     // -- in TPC pt Range
414     /*   if (particle->Pt() <= fHelper->GetMinPtForTOFRequired()) {
415       fMCNp[2][0][idxPart]    += 1;
416       if(iPid != 0)fMCNp[2][iPid][idxPart] += 1;
417     }
418     // -- in TPC+TOF pt Range
419     if (particle->Pt() > fHelper->GetMinPtForTOFRequired()) {
420       fMCNp[3][0][idxPart]    += 1;
421       if(iPid != 0) fMCNp[3][iPid][idxPart] += 1;
422       }*/
423    
424 #if USE_PHI
425     if(!fHelper->IsParticleAcceptedPhi(particle))
426       continue;
427     // idxPhi = 1;
428
429     if (idxPt>=0 && idxPt < AliEbyEPidRatioHelper::fgkfHistNBinsPt) {
430       fMCNpPt[1][0][idxPart][idxPt] += 1;
431       if(iPid != 0) fMCNpPt[1][iPid][idxPart][idxPt] += 1;
432     }
433
434     // -- in pt Range
435     fMCNp[4][0][idxPart]    += 1;
436     if(iPid != 0)fMCNp[4][iPid][idxPart] += 1;
437     
438     /*  // -- in TPC pt Range
439     if (particle->Pt() <= fHelper->GetMinPtForTOFRequired()) {
440       fMCNp[5][0][idxPart]    += 1;
441       if(iPid != 0)fMCNp[5][iPid][idxPart] += 1;
442     }
443     
444     // -- in TPC+TOF pt Range
445     if (particle->Pt() > fHelper->GetMinPtForTOFRequired()){
446       fMCNp[6][0][idxPart]    += 1;
447       if(iPid != 0)fMCNp[6][iPid][idxPart] += 1;
448       }*/
449 #endif
450   } // for (Int_t idxMC = 0; idxMC < nPart; ++idxMC) {
451   
452   FillHistSetCent("MC",        0, kTRUE);
453   if (fIsRatio) FillHistSetRatio("MCRatio",     0, kTRUE);
454   FillHistSetCent("MCpt",      1, kTRUE);
455
456   //  FillHistSetCent("MCTPC",     2, kTRUE);
457   // FillHistSetCent("MCTOF",     3, kTRUE);
458   if (fIsPtBin) FillHistSetCentPt("MCBin",    0, kTRUE);
459
460 #if USE_PHI
461   FillHistSetCent("MCphi",     4, kTRUE);
462   if (fIsPtBin)  FillHistSetCentPt("MCBinPhi",   1, kTRUE);
463   
464   // FillHistSetCent("MCTPCphi",  5, kTRUE);
465   // FillHistSetCent("MCTOFphi",  6, kTRUE);
466
467   
468 #endif
469
470   return 0;
471 }
472
473 //________________________________________________________________________
474 void  AliEbyEPidRatioPhy::AddHistSetCent(const Char_t *name, const Char_t *title)  {
475   TString sName(name);
476   TString sTitle(title);
477  
478   Float_t etaRange[2];
479
480   fESDTrackCuts->GetEtaRange(etaRange[0],etaRange[1]);
481
482   //TList *list[4];
483   fOutList->Add(new TList);
484   TList *list =  static_cast<TList*>(fOutList->Last());
485   list->SetName(Form("f%s", name));
486   list->SetOwner(kTRUE);
487   
488
489   Int_t nBinsCent         =  AliEbyEPidRatioHelper::fgkfHistNBinsCent;
490   Double_t centBinRange[] = {AliEbyEPidRatioHelper::fgkfHistRangeCent[0], AliEbyEPidRatioHelper::fgkfHistRangeCent[1]};
491
492  
493   for (Int_t iPid = 0; iPid < 4; ++iPid) {
494     // fOutList->Add(new TList);
495     // list[iPid] = static_cast<TList*>(fOutList->Last());
496     // list[iPid]->SetName(Form("f%s_%s", name,AliEbyEPidRatioHelper::fgkPidName[iPid]));
497     // list[iPid]->SetOwner(kTRUE);
498     
499     TString sNetTitle(Form("%s - %s", AliEbyEPidRatioHelper::fgkPidLatex[iPid][1], AliEbyEPidRatioHelper::fgkPidLatex[iPid][0]));
500
501     sTitle = (iPid != 0 ) ? Form("|y| < %.1f", fHelper->GetRapidityMax()) : Form(" |#eta|<%.1f", etaRange[1]);
502
503     for (Int_t idx = 1; idx <= fOrder; ++idx) {
504       list->Add(new TProfile(Form("fProf%s%sNet%dM", AliEbyEPidRatioHelper::fgkPidName[iPid],name, idx), 
505                              Form("(%s)^{%d} : %s;Centrality(100);(%s)^{%d}",sNetTitle.Data(), idx, sTitle.Data(), sNetTitle.Data(), idx),
506                              100,-0.5,99.5));
507     }
508     
509     for (Int_t ii = 0; ii <= fOrder; ++ii) {
510       for (Int_t kk = 0; kk <= fOrder; ++kk) {
511         list->Add(new TProfile(Form("fProf%s%sNetF%02d%02d", AliEbyEPidRatioHelper::fgkPidName[iPid], name, ii, kk),
512                                Form("f_{%02d%02d} : %s;Centrality(100);f_{%02d%02d}", ii, kk, sTitle.Data(), ii, kk),
513                                100,-0.5,99.5));
514       }
515     }
516   
517   }  
518
519   /* fOutList->Add(new TList);
520      TList *list_nu = static_cast<TList*>(fOutList->Last());
521      list_nu->SetName(Form("f%s_nu", name));
522      list_nu->SetOwner(kTRUE);
523   */
524
525   for (Int_t iPhy = 0; iPhy < 46; ++iPhy) { 
526     list->Add(new TProfile(Form("fProf%sNu%02d",name,iPhy),Form("Physics Variable for index %d | %s ; Centrality;",iPhy,name),100,-0.5,99.5));
527   }
528   
529    
530   for (Int_t iPid = 0; iPid < 4; ++iPid) {
531     TString sNetTitle(Form("%s - %s", AliEbyEPidRatioHelper::fgkPidLatex[iPid][1], AliEbyEPidRatioHelper::fgkPidLatex[iPid][0]));
532     sTitle = (iPid != 0 ) ? Form(" |y|<%.1f", fHelper->GetRapidityMax()) : Form(" |#eta| < %.1f", etaRange[1]);
533
534     for (Int_t idx = 1; idx <= fOrder; ++idx) {
535       list->Add(new TProfile(Form("fProfBin%s%sNet%dM", AliEbyEPidRatioHelper::fgkPidName[iPid],name, idx), 
536                              Form("(%s)^{%d} : %s;Centrality(11);(%s)^{%d}", sNetTitle.Data(), idx, sTitle.Data(), sNetTitle.Data(), idx),
537                              nBinsCent, centBinRange[0], centBinRange[1]));
538     }
539     
540     for (Int_t ii = 0; ii <= fOrder; ++ii) {
541       for (Int_t kk = 0; kk <= fOrder; ++kk) {
542         list->Add(new TProfile(Form("fProfBin%s%sNetF%02d%02d", AliEbyEPidRatioHelper::fgkPidName[iPid], name, ii, kk),
543                                Form("f_{%02d%02d} : %s;Centrality(11);f_{%02d%02d}", ii, kk, sTitle.Data(), ii, kk),
544                                nBinsCent, centBinRange[0], centBinRange[1]));
545       }
546     }
547   
548   }  
549     
550
551   for (Int_t iPhy = 0; iPhy < 46; ++iPhy) { 
552     list->Add(new TProfile(Form("fProfBin%sNu%02d",name,iPhy),Form("Physics Variable for index %d | %s ; Centrality;",iPhy,name),nBinsCent, centBinRange[0], centBinRange[1]));
553   }
554   
555   return;
556 }
557
558
559
560 //________________________________________________________________________
561 void  AliEbyEPidRatioPhy::AddHistSetRatio(const Char_t *name, const Char_t *title)  {
562   TString sName(name);
563   TString sTitle(title);
564  
565   Float_t etaRange[2];
566   fESDTrackCuts->GetEtaRange(etaRange[0],etaRange[1]);
567
568   //TList *list[4];
569   fOutList->Add(new TList);
570   TList *list =  static_cast<TList*>(fOutList->Last());
571   list->SetName(Form("f%s", name));
572   list->SetOwner(kTRUE);
573   
574   Int_t nRbin     = 10000;
575   Double_t mRat[] = {0,1};
576
577   Int_t nBinsCent         =  AliEbyEPidRatioHelper::fgkfHistNBinsCent;
578   Double_t centBinRange[] = {AliEbyEPidRatioHelper::fgkfHistRangeCent[0], AliEbyEPidRatioHelper::fgkfHistRangeCent[1]};
579
580   TString xyz = Form("|y| < %.1f",fHelper->GetRapidityMax()); 
581
582   list->Add(new TH2F(Form("fHistRatioKPi%s",name), 
583                      Form("(%s %s) : K/#pi;Centrality(11);K/#pi", xyz.Data(), sTitle.Data()),
584                      nBinsCent, centBinRange[0], centBinRange[1], nRbin,mRat[0],mRat[1]));
585   
586   list->Add(new TH2F(Form("fHistRatioKpPip%s",name), 
587                      Form("(%s %s) : K^{+}/#pi^{+};Centrality(11);K^{+}/#pi^{+}", xyz.Data(), sTitle.Data()),
588                      nBinsCent, centBinRange[0], centBinRange[1], nRbin,mRat[0],mRat[1]));
589   
590   list->Add(new TH2F(Form("fHistRatioKmPip%s",name), 
591                      Form("(%s %s) : K^{-}/#pi^{+};Centrality(11);K^{-}/#pi^{+}", xyz.Data(), sTitle.Data()),
592                      nBinsCent, centBinRange[0], centBinRange[1], nRbin,mRat[0],mRat[1]));
593   
594   list->Add(new TH2F(Form("fHistRatioKmPim%s",name), 
595                      Form("(%s %s) : K^{-}/#pi^{-};Centrality(11);K^{-}/#pi^{-}", xyz.Data(), sTitle.Data()),
596                      nBinsCent, centBinRange[0], centBinRange[1], nRbin,mRat[0],mRat[1]));
597
598
599
600  list->Add(new TH2F(Form("fHistRatioPK%s",name), 
601                      Form("(%s %s) : P/K;Centrality(11);P/K", xyz.Data(), sTitle.Data()),
602                      nBinsCent, centBinRange[0], centBinRange[1], nRbin,mRat[0],mRat[1]));
603   
604   list->Add(new TH2F(Form("fHistRatioPpKp%s",name), 
605                      Form("(%s %s) : P/K^{+};Centrality(11);P/K^{+}", xyz.Data(), sTitle.Data()),
606                      nBinsCent, centBinRange[0], centBinRange[1], nRbin,mRat[0],mRat[1]));
607   
608   list->Add(new TH2F(Form("fHistRatioPmKp%s",name), 
609                      Form("(%s %s) : #bar{P}/K^{+};Centrality(11);#bar{P}/K^{+}", xyz.Data(), sTitle.Data()),
610                      nBinsCent, centBinRange[0], centBinRange[1], nRbin,mRat[0],mRat[1]));
611   
612   list->Add(new TH2F(Form("fHistRatioPmKm%s",name), 
613                      Form("(%s %s) : #bar{P}/K^{-};Centrality(11);#bar{P}/K^{-}", xyz.Data(), sTitle.Data()),
614                      nBinsCent, centBinRange[0], centBinRange[1], nRbin,mRat[0],mRat[1]));
615
616
617
618  list->Add(new TH2F(Form("fHistRatioPPi%s",name), 
619                      Form("(%s %s) : P/#pi;Centrality(11);K/#pi", xyz.Data(), sTitle.Data()),
620                      nBinsCent, centBinRange[0], centBinRange[1], nRbin,mRat[0],mRat[1]));
621   
622   list->Add(new TH2F(Form("fHistRatioPpPip%s",name), 
623                      Form("(%s %s) : P/#pi^{+};Centrality(11);P/#pi^{+}", xyz.Data(), sTitle.Data()),
624                      nBinsCent, centBinRange[0], centBinRange[1], nRbin,mRat[0],mRat[1]));
625   
626   list->Add(new TH2F(Form("fHistRatioPmPip%s",name), 
627                      Form("(%s %s) : #bar{P}/#pi^{+};Centrality(11);#bar{P}/#pi^{+}", xyz.Data(), sTitle.Data()),
628                      nBinsCent, centBinRange[0], centBinRange[1], nRbin,mRat[0],mRat[1]));
629   
630   list->Add(new TH2F(Form("fHistRatioPmPim%s",name), 
631                      Form("(%s %s) : #bar{P}/#pi^{-};Centrality(11);#bar{P}/#pi^{-}", xyz.Data(), sTitle.Data()),
632                      nBinsCent, centBinRange[0], centBinRange[1], nRbin,mRat[0],mRat[1]));
633   return;
634 }
635
636 //________________________________________________________________________
637 void AliEbyEPidRatioPhy::FillHistSetCent(const Char_t *name, Int_t idx, Bool_t isMC)  {
638   /*
639     printf(" !! %2d !! ", idx);
640     for (Int_t iPid = 0; iPid < 4; ++iPid) {
641     printf("%7d %7d ", fNp[idx][iPid][0] , fNp[idx][iPid][1]);
642     }
643     Printf("");
644   */
645   
646   Int_t ***np = (isMC) ? fMCNp : fNp;
647   
648     /*    
649           printf(" ** %2d ** ", idx);
650           for (Int_t iPid = 0; iPid < 4; ++iPid) {
651           printf("%7d %7d ", np[idx][iPid][0] , np[idx][iPid][1]);
652           }
653           Printf("");
654     */
655   
656     // TList *list[4];
657   
658   Float_t centralityBin = fHelper->GetCentralityBin();
659   Float_t centralityPer = fHelper->GetCentralityPercentile();//fHelper->GetCentralityBin();
660   
661   TList *list = static_cast<TList*>(fOutList->FindObject(Form("f%s",name)));
662   
663   
664
665   for (Int_t iPid = 0; iPid < 4; ++iPid) {
666     Int_t deltaNp = np[idx][iPid][1]-np[idx][iPid][0];  
667     Double_t delta = 1.;
668     for (Int_t idxOrder = 1; idxOrder <= fOrder; ++idxOrder) {
669       delta *= deltaNp;
670       (static_cast<TProfile*>(list->FindObject(Form("fProfBin%s%sNet%dM", AliEbyEPidRatioHelper::fgkPidName[iPid], name, idxOrder))))->Fill(centralityBin, delta);
671       (static_cast<TProfile*>(list->FindObject(Form("fProf%s%sNet%dM", AliEbyEPidRatioHelper::fgkPidName[iPid], name, idxOrder))))->Fill(centralityPer, delta);
672     }
673     
674     for (Int_t idxOrder = 0; idxOrder <= fOrder; ++ idxOrder) {
675       fRedFactp[idxOrder][0]  = 1.;
676       fRedFactp[idxOrder][1]  = 1.;
677     }
678     
679     for (Int_t idxOrder = 1; idxOrder <= fOrder; ++ idxOrder) {
680       fRedFactp[idxOrder][0]  = fRedFactp[idxOrder-1][0]  * Double_t(np[idx][iPid][0]-(idxOrder-1));
681       fRedFactp[idxOrder][1]  = fRedFactp[idxOrder-1][1]  * Double_t(np[idx][iPid][1]-(idxOrder-1));
682     }
683     
684     for (Int_t ii = 0; ii <= fOrder; ++ii) {   // ii -> p    -> n1
685       for (Int_t kk = 0; kk <= fOrder; ++kk) { // kk -> pbar -> n2
686         Double_t fik = fRedFactp[ii][1] * fRedFactp[kk][0];   // n1 *n2 -> p * pbar
687         (static_cast<TProfile*>(list->FindObject(Form("fProfBin%s%sNetF%02d%02d", AliEbyEPidRatioHelper::fgkPidName[iPid], name, ii, kk))))->Fill(centralityBin, fik);
688         (static_cast<TProfile*>(list->FindObject(Form("fProf%s%sNetF%02d%02d", AliEbyEPidRatioHelper::fgkPidName[iPid], name, ii, kk))))->Fill(centralityPer, fik);
689       }
690     }
691   }
692  
693   Double_t a[6][4]; Double_t b[22];
694   for (Int_t iPid = 0; iPid < 4; ++iPid) {
695     a[0][iPid] = np[idx][iPid][1]+np[idx][iPid][0];       // 0  n+ + n-
696     a[1][iPid] = np[idx][iPid][1];                        // 1  n+
697     a[2][iPid] = np[idx][iPid][0];                        // 2  n-
698     a[3][iPid] = np[idx][iPid][1]*np[idx][iPid][0];       // 3  n+ . n-
699     a[4][iPid] = np[idx][iPid][1]*(np[idx][iPid][1]-1);   // 4  n+ (n+ - 1)
700     a[5][iPid] = np[idx][iPid][0]*(np[idx][iPid][0]-1);   // 5  n- (n- - 1)
701   }
702   
703   b[0]  = a[0][0]*a[0][2];       // 24 N   K
704   b[1]  = a[0][1]*a[0][2];       // 25 Pi  K
705   b[2]  = a[1][1]*a[1][2];       // 26 pi+ k+
706   b[3]  = a[1][1]*a[2][2];       // 27 pi+ k-
707   b[4]  = a[2][1]*a[1][2];       // 28 pi- k+  
708   b[5]  = a[2][1]*a[2][2];       // 29 pi- k-
709   
710   b[6]  = a[0][0]*a[0][3];       // 30 N   P
711   b[7]  = a[0][2]*a[0][3];       // 31 K   P
712   b[8]  = a[1][2]*a[1][3];       // 32 k+  p+
713   b[9]  = a[1][2]*a[2][3];       // 33 k+  p-
714   b[10] = a[2][2]*a[1][3];       // 34 k-  p+
715   b[11] = a[2][2]*a[2][3];       // 35 k-  p-
716   
717   b[12] = a[0][0]*a[0][1];       // 36 N  Pi
718   b[13] = a[0][3]*a[0][1];       // 37 P  Pi
719   b[14] = a[1][3]*a[1][1];       // 38 p+ pi+
720   b[15] = a[1][3]*a[2][1];       // 39 p+ pi-
721   b[16] = a[2][3]*a[1][1];       // 40 p- pi+
722   b[17] = a[2][3]*a[2][1];       // 41 p- pi-
723   
724   b[18] = a[0][0]*(a[0][0] - 1); // 42 N ( N - 1 )
725   b[19] = a[0][1]*(a[0][1] - 1); // 43 Pi( Pi- 1 )
726   b[20] = a[0][2]*(a[0][1] - 1); // 44 K ( K - 1 )
727   b[21] = a[0][3]*(a[0][3] - 1); // 45 P ( P - 1 )
728   // TList *list_nu = static_cast<TList*>(fOutList->FindObject(Form("f%s_nu",name)));
729   Int_t k = 0;
730   for (Int_t i = 0; i < 6; i++) {
731     for (Int_t j = 0; j < 4; j++) {
732       (static_cast<TProfile*>(list->FindObject(Form("fProfBin%sNu%02d", name,k))))->Fill(centralityBin,a[i][j]); 
733       (static_cast<TProfile*>(list->FindObject(Form("fProf%sNu%02d", name,k))))->Fill(centralityPer,a[i][j]); 
734       k++;
735     }
736   }
737
738   for (Int_t j = 0; j < 22; j++) {
739     (static_cast<TProfile*>(list->FindObject(Form("fProfBin%sNu%02d", name,j+24))))->Fill(centralityBin,b[j]); 
740     (static_cast<TProfile*>(list->FindObject(Form("fProf%sNu%02d", name,j+24))))->Fill(centralityPer,b[j]); 
741   }
742   
743   return;
744 }
745
746
747 //________________________________________________________________________
748 void AliEbyEPidRatioPhy::FillHistSetRatio(const Char_t *name, Int_t idx, Bool_t isMC)  {
749    
750   Int_t ***np = (isMC) ? fMCNp : fNp;
751   Float_t centralityBin = fHelper->GetCentralityBin();
752   
753   TList *list = static_cast<TList*>(fOutList->FindObject(Form("f%s",name)));
754     
755   Double_t KPi   = -1; if(np[idx][1][1]+np[idx][1][0] != 0 ) KPi = (np[idx][2][1]+np[idx][2][0])/(np[idx][1][1]+np[idx][1][0]);
756   Double_t KpPip = -1;
757   Double_t KmPip = -1;if (np[idx][1][1] != 0 ) { KpPip  = (np[idx][2][1])/(np[idx][1][1]); KmPip =  (np[idx][2][0])/(np[idx][1][1]); }
758   Double_t KmPim = -1;if (np[idx][1][0] != 0) KmPim = (np[idx][2][0])/(np[idx][1][0]);
759   
760   (static_cast<TProfile*>(list->FindObject(Form("fHistRatioKPi%s",name))))->Fill(centralityBin, KPi);
761   (static_cast<TProfile*>(list->FindObject(Form("fHistRatioKpPip%s",name))))->Fill(centralityBin, KpPip);
762   (static_cast<TProfile*>(list->FindObject(Form("fHistRatioKmPip%s",name))))->Fill(centralityBin, KmPip);
763   (static_cast<TProfile*>(list->FindObject(Form("fHistRatioKmPim%s",name))))->Fill(centralityBin, KmPim);
764
765   Double_t PK   = -1; if(np[idx][2][1]+np[idx][2][0] != 0 ) PK = (np[idx][3][1]+np[idx][3][0])/(np[idx][2][1]+np[idx][2][0]);
766   Double_t PpKp = -1;
767   Double_t PmKp = -1;if (np[idx][2][1] != 0 ) { PpKp  = (np[idx][3][1])/(np[idx][2][1]); PmKp =  (np[idx][3][0])/(np[idx][2][1]); }
768   Double_t PmKm = -1;if (np[idx][2][0] != 0) PmKm = (np[idx][3][0])/(np[idx][2][0]);
769
770   (static_cast<TProfile*>(list->FindObject(Form("fHistRatioPK%s",name))))->Fill(centralityBin, PK);
771   (static_cast<TProfile*>(list->FindObject(Form("fHistRatioPpKp%s",name))))->Fill(centralityBin, PpKp);
772   (static_cast<TProfile*>(list->FindObject(Form("fHistRatioPmKp%s",name))))->Fill(centralityBin, PmKp);
773   (static_cast<TProfile*>(list->FindObject(Form("fHistRatioPmKm%s",name))))->Fill(centralityBin, PmKm);
774
775   Double_t PPi   = -1; if(np[idx][1][1]+np[idx][1][0] != 0 ) PPi = (np[idx][3][1]+np[idx][3][0])/(np[idx][1][1]+np[idx][1][0]);
776   Double_t PpPip = -1;
777   Double_t PmPip = -1;if (np[idx][1][1] != 0 ) { PpPip  = (np[idx][3][1])/(np[idx][1][1]); PmPip =  (np[idx][3][0])/(np[idx][1][1]); }
778   Double_t PmPim = -1;if (np[idx][1][0] != 0) PmPim = (np[idx][3][0])/(np[idx][1][0]);
779
780   (static_cast<TProfile*>(list->FindObject(Form("fHistRatioPPi%s",name))))->Fill(centralityBin, PPi);
781   (static_cast<TProfile*>(list->FindObject(Form("fHistRatioPpPip%s",name))))->Fill(centralityBin, PpPip);
782   (static_cast<TProfile*>(list->FindObject(Form("fHistRatioPmPip%s",name))))->Fill(centralityBin, PmPip);
783   (static_cast<TProfile*>(list->FindObject(Form("fHistRatioPmPim%s",name))))->Fill(centralityBin, PmPim);
784   
785   return;
786 }
787
788 /*
789 //________________________________________________________________________
790 void  AliEbyEPidRatioPhy::AddHistSetCentPt(const Char_t *name, const Char_t *title)  {
791   TString sName(name);
792   TString sTitle(title);
793   
794   Float_t etaRange[2];
795   fESDTrackCuts->GetEtaRange(etaRange[0],etaRange[1]);
796
797   //TList *list[4];
798   fOutList->Add(new TList);
799   TList *list =  static_cast<TList*>(fOutList->Last());
800   list->SetName(Form("f%s", name));
801   list->SetOwner(kTRUE);
802   
803
804   Int_t nBinsCent         =  AliEbyEPidRatioHelper::fgkfHistNBinsCent;
805   Double_t centBinRange[] = {AliEbyEPidRatioHelper::fgkfHistRangeCent[0], AliEbyEPidRatioHelper::fgkfHistRangeCent[1]};
806   Int_t nBinsPt           =  AliEbyEPidRatioHelper::fgkfHistNBinsPt;
807   Double_t ptBinRange[]   = {-0.5, AliEbyEPidRatioHelper::fgkfHistNBinsPt - 0.5};
808
809   for (Int_t idxPt  = 0; idxPt < AliEbyEPidRatioHelper::fgkfHistNBinsPt; ++idxPt) { 
810     for (Int_t iPid = 0; iPid < 4; ++iPid) {
811       TString sNetTitle(Form("%s - %s", AliEbyEPidRatioHelper::fgkPidLatex[iPid][1], AliEbyEPidRatioHelper::fgkPidLatex[iPid][0]));
812       sTitle = (iPid != 0 ) ? Form("|y| < %.1f", fHelper->GetRapidityMax()) : Form(" |#eta|<%.1f", etaRange[1]);
813       for (Int_t idx = 1; idx <= fOrder; ++idx) {
814         list->Add(new TProfile(Form("fProf%s%sNetPt%02d%02dM", AliEbyEPidRatioHelper::fgkPidName[iPid],name, idxPt, idx), 
815                                Form("(%s)^{%d} Pt: %s;Centrality(100);(%s)^{%d}",
816                                     sNetTitle.Data(), idx, sTitle.Data(), sNetTitle.Data(), idx),
817                                  100,-0.5,99.5));
818       }
819     
820       for (Int_t ii = 0; ii <= fOrder; ++ii) {
821         for (Int_t kk = 0; kk <= fOrder; ++kk) {
822           list->Add(new TProfile(Form("fProf%s%sNetPt%02dF%02d%02d", AliEbyEPidRatioHelper::fgkPidName[iPid], name, idxPt,ii, kk),
823                                    Form("f_{%02d%02d} Pt: %s;Centrality(100);f_{%02d%02d}", ii, kk, sTitle.Data(), ii, kk),
824                                    100,-0.5,99.5));
825         }
826       }
827       
828     }
829     
830     for (Int_t iPid = 0; iPid < 4; ++iPid) {
831       TString sNetTitle(Form("%s - %s", AliEbyEPidRatioHelper::fgkPidLatex[iPid][1], AliEbyEPidRatioHelper::fgkPidLatex[iPid][0]));
832       sTitle = (iPid != 0 ) ? Form("|y| < %.1f", fHelper->GetRapidityMax()) : Form(" |#eta|<%.1f", etaRange[1]);
833       for (Int_t idx = 1; idx <= fOrder; ++idx) {
834         list->Add(new TProfile(Form("fProfBin%s%sNetPt%02d%02dM", AliEbyEPidRatioHelper::fgkPidName[iPid],name, idxPt,idx), 
835                                  Form("(%s)^{%d} Pt: %s;Centrality(11);(%s)^{%d}", sNetTitle.Data(), idx, sTitle.Data(), sNetTitle.Data(), idx),
836                                  nBinsCent, centBinRange[0], centBinRange[1]));
837       }
838       
839       for (Int_t ii = 0; ii <= fOrder; ++ii) {
840         for (Int_t kk = 0; kk <= fOrder; ++kk) {
841           list->Add(new TProfile(Form("fProfBin%s%sNetPt%02dF%02d%02d", AliEbyEPidRatioHelper::fgkPidName[iPid], name, idxPt,ii, kk),
842                                    Form("f_{%02d%02d} Pt : %s;Centrality(11);f_{%02d%02d}", ii, kk, sTitle.Data(), ii, kk),
843                                    nBinsCent, centBinRange[0], centBinRange[1]));
844         }
845       }
846     }  
847   }
848   
849   return;
850 }
851
852 //________________________________________________________________________
853 void AliEbyEPidRatioPhy::FillHistSetCentPt(const Char_t *name, Int_t idx, Bool_t isMC)  {
854  
855   Int_t ****np = (isMC) ? fMCNpPt : fNpPt;
856   
857   Float_t centralityBin = fHelper->GetCentralityBin();
858   Float_t centralityPer = fHelper->GetCentralityPercentile();
859   
860   TList *list = static_cast<TList*>(fOutList->FindObject(Form("f%s",name)));
861   
862   
863   for (Int_t idxPt  = 0; idxPt < AliEbyEPidRatioHelper::fgkfHistNBinsPt; ++idxPt) {
864     for (Int_t iPid = 0; iPid < 4; ++iPid) {
865       Int_t deltaNp = np[idx][iPid][1][idxPt]-np[idx][iPid][0][idxPt];  
866       Double_t delta = 1.;
867       for (Int_t idxOrder = 1; idxOrder <= fOrder; ++idxOrder) {
868         delta *= deltaNp;
869         (static_cast<TProfile*>(list->FindObject(Form("fProfBin%s%sNetPt%02d%02dM", AliEbyEPidRatioHelper::fgkPidName[iPid], name, idxPt, idxOrder))))->Fill(centralityBin, delta);
870         (static_cast<TProfile*>(list->FindObject(Form("fProf%s%sNetPt%02d%02dM", AliEbyEPidRatioHelper::fgkPidName[iPid], name, idxPt,idxOrder))))->Fill(centralityPer, delta);
871       }
872     
873       for (Int_t idxOrder = 0; idxOrder <= fOrder; ++ idxOrder) {
874         fRedFactp[idxOrder][0]  = 1.;
875         fRedFactp[idxOrder][1]  = 1.;
876       }
877       
878       for (Int_t idxOrder = 1; idxOrder <= fOrder; ++ idxOrder) {
879         fRedFactp[idxOrder][0]  = fRedFactp[idxOrder-1][0] * Double_t(np[idx][iPid][0][idxPt]-(idxOrder-1));
880         fRedFactp[idxOrder][1]  = fRedFactp[idxOrder-1][1] * Double_t(np[idx][iPid][1][idxPt]-(idxOrder-1));
881       }
882       
883       for (Int_t ii = 0; ii <= fOrder; ++ii) {   // ii -> p    -> n1
884         for (Int_t kk = 0; kk <= fOrder; ++kk) { // kk -> pbar -> n2
885           Double_t fik = fRedFactp[ii][1] * fRedFactp[kk][0];   // n1 *n2 -> p * pbar
886           (static_cast<TProfile*>(list->FindObject(Form("fProfBin%s%sNetPt%02dF%02d%02d", AliEbyEPidRatioHelper::fgkPidName[iPid], name, idxPt,ii, kk))))->Fill(centralityBin, fik);
887           (static_cast<TProfile*>(list->FindObject(Form("fProf%s%sNetPt%02dF%02d%02d", AliEbyEPidRatioHelper::fgkPidName[iPid], name, idxPt,ii, kk))))->Fill(centralityPer, fik);
888         }
889       }
890     }
891   }// for (Int_t idxPt  = 0; idxPt < AliEbyEPidRatioHelper::fgkfHistNBinsPt; ++idxPt) {
892   
893   return;
894 }
895
896 */
897
898
899
900 //________________________________________________________________________
901 void  AliEbyEPidRatioPhy::AddHistSetCentPt(const Char_t *name, const Char_t *title)  {
902   TString sName(name);
903   TString sTitle(title);
904   
905   Float_t etaRange[2];
906   fESDTrackCuts->GetEtaRange(etaRange[0],etaRange[1]);
907
908   //TList *list[4];
909   fOutList->Add(new TList);
910   TList *list =  static_cast<TList*>(fOutList->Last());
911   list->SetName(Form("f%s", name));
912   list->SetOwner(kTRUE);
913   
914
915   Int_t nBinsCent         =  AliEbyEPidRatioHelper::fgkfHistNBinsCent;
916   Double_t centBinRange[] = {AliEbyEPidRatioHelper::fgkfHistRangeCent[0], AliEbyEPidRatioHelper::fgkfHistRangeCent[1]};
917   Int_t nBinsPt           =  AliEbyEPidRatioHelper::fgkfHistNBinsPt;
918   Double_t ptBinRange[]   = {-0.5, AliEbyEPidRatioHelper::fgkfHistNBinsPt - 0.5};
919  
920     for (Int_t iPid = 0; iPid < 4; ++iPid) {
921     TString sNetTitle(Form("%s - %s", AliEbyEPidRatioHelper::fgkPidLatex[iPid][1], AliEbyEPidRatioHelper::fgkPidLatex[iPid][0]));
922     sTitle = (iPid != 0 ) ? Form("|y| < %.1f", fHelper->GetRapidityMax()) : Form(" |#eta|<%.1f", etaRange[1]);
923     for (Int_t idx = 1; idx <= fOrder; ++idx) {
924       list->Add(new TProfile2D(Form("fProf%s%sNetPt%dM", AliEbyEPidRatioHelper::fgkPidName[iPid],name, idx), 
925                              Form("(%s)^{%d} Pt: %s;Centrality(100);(%s)^{%d}",sNetTitle.Data(), idx, sTitle.Data(), sNetTitle.Data(), idx),
926                                100,-0.5,99.5, nBinsPt, ptBinRange[0], ptBinRange[1]));
927     }
928     
929     for (Int_t ii = 0; ii <= fOrder; ++ii) {
930       for (Int_t kk = 0; kk <= fOrder; ++kk) {
931         list->Add(new TProfile2D(Form("fProf%s%sNetPtF%02d%02d", AliEbyEPidRatioHelper::fgkPidName[iPid], name, ii, kk),
932                                Form("f_{%02d%02d} Pt: %s;Centrality(100);f_{%02d%02d}", ii, kk, sTitle.Data(), ii, kk),
933                                100,-0.5,99.5, nBinsPt, ptBinRange[0], ptBinRange[1]));
934       }
935     }
936   
937     }
938   
939   for (Int_t iPid = 0; iPid < 4; ++iPid) {
940     TString sNetTitle(Form("%s - %s", AliEbyEPidRatioHelper::fgkPidLatex[iPid][1], AliEbyEPidRatioHelper::fgkPidLatex[iPid][0]));
941     sTitle = (iPid != 0 ) ? Form("|y| < %.1f", fHelper->GetRapidityMax()) : Form(" |#eta|<%.1f", etaRange[1]);
942     for (Int_t idx = 1; idx <= fOrder; ++idx) {
943       list->Add(new TProfile2D(Form("fProfBin%s%sNetPt%dM", AliEbyEPidRatioHelper::fgkPidName[iPid],name, idx), 
944                                Form("(%s)^{%d} Pt: %s;Centrality(11);(%s)^{%d}", sNetTitle.Data(), idx, sTitle.Data(), sNetTitle.Data(), idx),
945                                nBinsCent, centBinRange[0], centBinRange[1],nBinsPt, ptBinRange[0], ptBinRange[1]));
946     }
947     
948     for (Int_t ii = 0; ii <= fOrder; ++ii) {
949       for (Int_t kk = 0; kk <= fOrder; ++kk) {
950         list->Add(new TProfile2D(Form("fProfBin%s%sNetPtF%02d%02d", AliEbyEPidRatioHelper::fgkPidName[iPid], name, ii, kk),
951                                  Form("f_{%02d%02d} Pt : %s;Centrality(11);f_{%02d%02d}", ii, kk, sTitle.Data(), ii, kk),
952                                  nBinsCent, centBinRange[0], centBinRange[1],nBinsPt, ptBinRange[0], ptBinRange[1]));
953       }
954     }
955   }  
956    
957   return;
958 }
959
960 //________________________________________________________________________
961 void AliEbyEPidRatioPhy::FillHistSetCentPt(const Char_t *name, Int_t idx, Bool_t isMC)  {
962  
963   Int_t ****np = (isMC) ? fMCNpPt : fNpPt;
964   
965   Float_t centralityBin = fHelper->GetCentralityBin();
966   Float_t centralityPer = fHelper->GetCentralityPercentile();
967   
968   TList *list = static_cast<TList*>(fOutList->FindObject(Form("f%s",name)));
969   
970   
971   for (Int_t idxPt  = 0; idxPt < AliEbyEPidRatioHelper::fgkfHistNBinsPt; ++idxPt) {
972     for (Int_t iPid = 0; iPid < 4; ++iPid) {
973       Int_t deltaNp = np[idx][iPid][1][idxPt]-np[idx][iPid][0][idxPt];  
974       Double_t delta = 1.;
975       for (Int_t idxOrder = 1; idxOrder <= fOrder; ++idxOrder) {
976         delta *= deltaNp;
977         (static_cast<TProfile2D*>(list->FindObject(Form("fProfBin%s%sNetPt%dM", AliEbyEPidRatioHelper::fgkPidName[iPid], name, idxOrder))))->Fill(centralityBin,idxPt, delta);
978         (static_cast<TProfile2D*>(list->FindObject(Form("fProf%s%sNetPt%dM", AliEbyEPidRatioHelper::fgkPidName[iPid], name, idxOrder))))->Fill(centralityPer,idxPt, delta);
979       }
980     
981       for (Int_t idxOrder = 0; idxOrder <= fOrder; ++ idxOrder) {
982         fRedFactp[idxOrder][0]  = 1.;
983         fRedFactp[idxOrder][1]  = 1.;
984       }
985       
986       for (Int_t idxOrder = 1; idxOrder <= fOrder; ++ idxOrder) {
987         fRedFactp[idxOrder][0]  = fRedFactp[idxOrder-1][0] * Double_t(np[idx][iPid][0][idxPt]-(idxOrder-1));
988         fRedFactp[idxOrder][1]  = fRedFactp[idxOrder-1][1] * Double_t(np[idx][iPid][1][idxPt]-(idxOrder-1));
989       }
990       
991       for (Int_t ii = 0; ii <= fOrder; ++ii) {   // ii -> p    -> n1
992         for (Int_t kk = 0; kk <= fOrder; ++kk) { // kk -> pbar -> n2
993           Double_t fik = fRedFactp[ii][1] * fRedFactp[kk][0];   // n1 *n2 -> p * pbar
994           (static_cast<TProfile2D*>(list->FindObject(Form("fProfBin%s%sNetPtF%02d%02d", AliEbyEPidRatioHelper::fgkPidName[iPid], name, ii, kk))))->Fill(centralityBin,idxPt, fik);
995           (static_cast<TProfile2D*>(list->FindObject(Form("fProf%s%sNetPtF%02d%02d", AliEbyEPidRatioHelper::fgkPidName[iPid], name, ii, kk))))->Fill(centralityPer,idxPt, fik);
996         }
997       }
998     }
999   }// for (Int_t idxPt  = 0; idxPt < AliEbyEPidRatioHelper::fgkfHistNBinsPt; ++idxPt) {
1000   
1001   return;
1002 }
1003
1004