]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/hfe/macros/DrawEfficiencyCentralitySource.C
Updates to run with deltas (L. Cunqueiro)
[u/mrichter/AliRoot.git] / PWGHF / hfe / macros / DrawEfficiencyCentralitySource.C
1 TH1I *GetNEvents(const Char_t *testfile,const Char_t *plus="");
2 void CorrectFromTheWidth(TH1D *h1);
3 TList *GetResults(const Char_t *testfile,const Char_t *plus="");
4 TObject* GetSpectrum(AliCFContainer *c, Int_t step);
5 TObject* GetEfficiency(AliCFContainer *c, Int_t step = 9);
6 TObject* GetEfficiency(AliCFContainer *c, Int_t step, Int_t step0);
7 AliCFContainer *GetContainerCentralitySource(AliCFContainer *container,Int_t bincentrality,Int_t source);
8 AliCFContainer *GetContainerSourceAsFunctionOfCentrality(AliCFContainer *container,Int_t source);
9
10 void DrawEfficiencyCentralitySource(Int_t bincentrality,Int_t source,const char *testfile = "/lustre/alice/train/V005.PbPb/2010-11-29_2209.4266/mergedPeriods/data/PbPb/LHC10h.pass1/HFEtask_PbPb.root");
11
12
13 void DrawEfficiencyCentralitySource(Int_t bincentrality,Int_t source,const char *testfile) {
14   
15   //
16   // source: 0 (charm), 1 (beauty), 2 (gamma), 3 (others)
17   // centrality:  from 0 to 100 with 0-5, 5-10, 10-15...
18   //
19   //
20
21   gStyle->SetPalette(1);
22   gStyle->SetOptStat(1111);
23   gStyle->SetPadBorderMode(0);
24   gStyle->SetCanvasColor(10);
25   gStyle->SetPadLeftMargin(0.13);
26   gStyle->SetPadRightMargin(0.13);
27
28   ///////////////////////////////////
29   // Take the stuff
30   //////////////////////////////////
31
32
33   TList *results = GetResults(testfile,"");
34   if(!results){
35     printf("No output objects: Calculation will terminate here\n");
36     return;
37   }
38
39   AliHFEcontainer *containerhfe = (AliHFEcontainer *) results->FindObject("trackContainer");
40   if(!containerhfe) {
41     printf("No hfe container \n");
42     return;
43   }
44
45   // 0: pt, 1: eta, 2: phi, 3: charge, 4: source, 5: centrality
46   AliCFContainer *sumcontaineresdd = containerhfe->MakeMergedCFContainer("sumesd","sumesd","MCTrackCont:recTrackContReco");
47   AliCFContainer *sumcontainermcc = containerhfe->MakeMergedCFContainer("summc","summc","MCTrackCont:recTrackContMC");
48
49   if(!sumcontaineresdd) {
50     printf("No container sum esd\n");
51     return;
52   }
53
54   if(!sumcontainermcc) {
55     printf("No container sum mc\n");
56     return;
57   }
58
59   // 0: pt, 1: eta, 2: phi for a given centrality and source
60   //AliCFContainer *sumcontaineresd = GetContainerSourceCentrality(sumcontaineresdd,bincentrality,source);
61   //AliCFContainer *sumcontainermc = GetContainerSourceCentrality(sumcontainermcc,bincentrality,source);
62
63   // 0: pt, 1: eta, 2: phi, 3: centrality for a given centrality and source
64   AliCFContainer *sumcontaineresd = GetContainerSourceAsFunctionOfCentrality(sumcontaineresdd,source);
65   AliCFContainer *sumcontainermc = GetContainerSourceAsFunctionOfCentrality(sumcontainermcc,source);
66
67  
68
69   Int_t nSteps = sumcontaineresd->GetNStep();
70   printf("In total %d steps\n",nSteps);
71
72   AliCFDataGrid *dataGrida = (AliCFDataGrid *) GetSpectrum(sumcontaineresdd,nSteps-1);
73  
74   TH1D *spectrumcentrality = (TH1D *) dataGrida->Project(5);
75   TH1D *spectrumpt = (TH1D *) dataGrida->Project(0);
76   TH2D *spectrumptc = (TH2D *) dataGrida->Project(5,0);
77
78   TAxis *xaxis = spectrumptc->GetXaxis();
79   Int_t bin0 = xaxis->FindBin(0.0);
80   Int_t bin5_s = xaxis->FindBin(4.99);
81   Int_t bin30 = xaxis->FindBin(30.0);
82   Int_t bin40_s = xaxis->FindBin(39.9);
83   Int_t bin70 = xaxis->FindBin(70.0);
84   Int_t bin80_s = xaxis->FindBin(79.9);
85
86   printf("Bin 0 %d\n",bin0);
87   printf("Bin 5 %d\n",bin5_s);
88   printf("Bin 30 %d\n",bin30);
89   printf("Bin 40 %d\n",bin40_s);
90   printf("Bin 70 %d\n",bin70);
91   printf("Bin 80 %d\n",bin80_s);
92   
93   TH1D *spectrumcentrality_0_5 = spectrumptc->ProjectionY("centrality_0_5",bin0,bin5_s);
94   TH1D *spectrumcentrality_30_40 = spectrumptc->ProjectionY("centrality_30_40",bin30,bin40_s);
95   TH1D *spectrumcentrality_70_80 = spectrumptc->ProjectionY("centrality_70_80",bin70,bin80_s);
96
97   Int_t numberOfEvents = (Int_t) containerhfe->GetNumberOfEvents();
98   
99   printf("Number of events for a %d after Event cut\n",numberOfEvents);
100
101   ////////////////////////////////
102   // Input after ITS&TPC refit
103   ///////////////////////////////
104   TCanvas * canvascpt = new TCanvas("RawSpectrumCentrality","RawSpectrumCentrality",1000,700);
105   canvascpt->Divide(2,1);
106   canvascpt->cd(1);
107   spectrumpt->SetTitle("");
108   spectrumpt->SetStats(0);
109   spectrumpt->SetLineColor(kBlue);
110   spectrumpt->SetMarkerColor(kBlue);
111   spectrumpt->SetMarkerStyle(25);
112   //
113   spectrumcentrality_0_5->SetTitle("");
114   spectrumcentrality_0_5->SetStats(0);
115   spectrumcentrality_0_5->SetLineColor(kRed);
116   spectrumcentrality_0_5->SetMarkerColor(kRed);
117   spectrumcentrality_0_5->SetMarkerStyle(26);
118   //
119   spectrumcentrality_30_40->SetTitle("");
120   spectrumcentrality_30_40->SetStats(0);
121   spectrumcentrality_30_40->SetLineColor(kMagenta);
122   spectrumcentrality_30_40->SetMarkerColor(kBlack);
123   spectrumcentrality_30_40->SetMarkerStyle(27);
124   //
125   spectrumcentrality_70_80->SetTitle("");
126   spectrumcentrality_70_80->SetStats(0);
127   spectrumcentrality_70_80->SetLineColor(kBlue);
128   spectrumcentrality_70_80->SetMarkerColor(kBlue);
129   spectrumcentrality_70_80->SetMarkerStyle(28);
130   //
131   spectrumpt->Draw();
132   spectrumcentrality_0_5->Draw("same");
133   spectrumcentrality_30_40->Draw("same");
134   spectrumcentrality_70_80->Draw("same");
135   TLegend *leg_different_centralities = new TLegend(0.4,0.6,0.89,0.89);
136   leg_different_centralities->AddEntry(spectrumpt,"Minimum-bias","p");
137   leg_different_centralities->AddEntry(spectrumcentrality_0_5,"0_5","p");
138   leg_different_centralities->AddEntry(spectrumcentrality_30_40,"30_40","p");
139   leg_different_centralities->AddEntry(spectrumcentrality_70_80,"70_80","p");
140   leg_different_centralities->Draw("same");
141
142   canvascpt->cd(2);
143   spectrumptc->Draw("colz");
144   
145   /////////////////////////////////////
146   // Take efficiencies
147   /////////////////////////////////////
148   
149   AliCFEffGrid  *efficiencystepkineITSTPC  = (AliCFEffGrid*)  GetEfficiency(sumcontaineresd,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepRecKineITSTPC,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepRecNoCut);
150   AliCFEffGrid  *efficiencystepPrim        = (AliCFEffGrid*)  GetEfficiency(sumcontaineresd,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepRecPrim,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepRecKineITSTPC);
151   AliCFEffGrid  *efficiencystepHFEcutsITS  = (AliCFEffGrid*)  GetEfficiency(sumcontaineresd,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsITS,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepRecPrim);
152   AliCFEffGrid  *efficiencystepHFEcutsTRD  = (AliCFEffGrid*)  GetEfficiency(sumcontaineresd,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsITS);
153   AliCFEffGrid  *efficiencystepPIDTOF      = (AliCFEffGrid*)  GetEfficiency(sumcontaineresd,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD+1,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD);
154   AliCFEffGrid  *efficiencystepPIDTPC      = (AliCFEffGrid*)  GetEfficiency(sumcontaineresd,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD+2,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD+1);
155
156   AliCFEffGrid  *efficiencystepPIDall      = (AliCFEffGrid*)  GetEfficiency(sumcontaineresd,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD+2,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD);
157
158   AliCFEffGrid  *efficiencystepall      = (AliCFEffGrid*)  GetEfficiency(sumcontaineresd,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD+2,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepRecNoCut);
159
160   //
161
162   AliCFEffGrid  *efficiencystepMCkineITSTPC  = (AliCFEffGrid*)  GetEfficiency(sumcontainermc,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepRecKineITSTPC,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepRecNoCut);
163   AliCFEffGrid  *efficiencystepMCPrim        = (AliCFEffGrid*)  GetEfficiency(sumcontainermc,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepRecPrim,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepRecKineITSTPC);
164   AliCFEffGrid  *efficiencystepMCHFEcutsITS  = (AliCFEffGrid*)  GetEfficiency(sumcontainermc,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsITS,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepRecPrim);
165   AliCFEffGrid  *efficiencystepMCHFEcutsTRD  = (AliCFEffGrid*)  GetEfficiency(sumcontainermc,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsITS);
166   AliCFEffGrid  *efficiencystepMCPIDTOF      = (AliCFEffGrid*)  GetEfficiency(sumcontainermc,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD+1,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD);
167   AliCFEffGrid  *efficiencystepMCPIDTPC      = (AliCFEffGrid*)  GetEfficiency(sumcontainermc,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD+2,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD+1);
168
169   AliCFEffGrid  *efficiencystepMCPIDall      = (AliCFEffGrid*)  GetEfficiency(sumcontainermc,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD+2,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD);
170
171   AliCFEffGrid  *efficiencystepMCall      = (AliCFEffGrid*)  GetEfficiency(sumcontainermc,AliHFEcuts::kNcutStepsMCTrack+AliHFEcuts::kStepHFEcutsTRD+2,0);
172
173   
174   ///////////////////////////////////////
175   // Plot 1D
176   ///////////////////////////////////////
177
178   /////////////////////
179   // Different cuts
180   /////////////////////
181   TCanvas * canvascomparison = new TCanvas("ITSTPCrefitStepESD","ITSTPCrefitStepESD",1000,700);
182   canvascomparison->Divide(2,2);
183  
184   canvascomparison->cd(1);
185   TH2D* h_effsteponlykineITSTPC1Donly = (TH2D *) efficiencystepkineITSTPC->Project(0,3);
186   h_effsteponlykineITSTPC1Donly->SetTitle("");
187   h_effsteponlykineITSTPC1Donly->SetStats(0);
188   h_effsteponlykineITSTPC1Donly->SetLineColor(kBlue);
189   h_effsteponlykineITSTPC1Donly->SetMarkerColor(kBlue);
190   h_effsteponlykineITSTPC1Donly->SetMarkerStyle(25);
191   h_effsteponlykineITSTPC1Donly->Draw("colz");
192    
193   canvascomparison->cd(2);
194   TH2D* h_effsteponlyPrim1Donly =( TH2D *) efficiencystepPrim->Project(0,3);
195   h_effsteponlyPrim1Donly->SetTitle("");
196   h_effsteponlyPrim1Donly->SetStats(0);
197   h_effsteponlyPrim1Donly->SetLineColor(kBlue);
198   h_effsteponlyPrim1Donly->SetMarkerColor(kBlue);
199   h_effsteponlyPrim1Donly->SetMarkerStyle(25);
200   h_effsteponlyPrim1Donly->Draw("colz");
201
202   
203   canvascomparison->cd(3);
204   TH2D* h_effsteponlyHFEcutsITS1Donly = (TH2D *) efficiencystepHFEcutsITS->Project(0,3);
205   h_effsteponlyHFEcutsITS1Donly->SetTitle("");
206   h_effsteponlyHFEcutsITS1Donly->SetStats(0);
207   h_effsteponlyHFEcutsITS1Donly->SetLineColor(kBlue);
208   h_effsteponlyHFEcutsITS1Donly->SetMarkerColor(kBlue);
209   h_effsteponlyHFEcutsITS1Donly->SetMarkerStyle(25);
210   h_effsteponlyHFEcutsITS1Donly->Draw("colz");
211  
212   
213   canvascomparison->cd(4);
214   TH2D* h_effsteponlyPID1Donly = (TH2D *) efficiencystepHFEcutsTRD->Project(0,3);
215   h_effsteponlyPID1Donly->SetTitle("");
216   h_effsteponlyPID1Donly->SetStats(0);
217   h_effsteponlyPID1Donly->SetLineColor(kBlue);
218   h_effsteponlyPID1Donly->SetMarkerColor(kBlue);
219   h_effsteponlyPID1Donly->SetMarkerStyle(25);
220   h_effsteponlyPID1Donly->Draw("colz");
221  
222
223   /////////////////////
224   // Different cuts
225   /////////////////////
226   TCanvas * canvascomparisonbis = new TCanvas("ITSTPCrefitStepMC","ITSTPCrefitStepMC",1000,700);
227   canvascomparisonbis->Divide(2,2);
228  
229   canvascomparisonbis->cd(1);
230   TH2D* h_effstepMCkineITSTPC1Donly =  (TH2D *)efficiencystepMCkineITSTPC->Project(0,3);
231   h_effstepMCkineITSTPC1Donly->SetTitle("");
232   h_effstepMCkineITSTPC1Donly->SetStats(0);
233   h_effstepMCkineITSTPC1Donly->SetLineColor(kBlue);
234   h_effstepMCkineITSTPC1Donly->SetMarkerColor(kBlue);
235   h_effstepMCkineITSTPC1Donly->SetMarkerStyle(25);
236   h_effstepMCkineITSTPC1Donly->Draw("colz");
237   
238
239   canvascomparisonbis->cd(2);
240   TH2D* h_effstepMCPrim1Donly = (TH2D *) efficiencystepMCPrim->Project(0,3);
241   h_effstepMCPrim1Donly->SetTitle("");
242   h_effstepMCPrim1Donly->SetStats(0);
243   h_effstepMCPrim1Donly->SetLineColor(kBlue);
244   h_effstepMCPrim1Donly->SetMarkerColor(kBlue);
245   h_effstepMCPrim1Donly->SetMarkerStyle(25);
246   h_effstepMCPrim1Donly->Draw("colz");
247   
248   canvascomparisonbis->cd(3);
249   TH2D* h_effstepMCHFEcutsITS1Donly = (TH2D *) efficiencystepMCHFEcutsITS->Project(0,3);
250   h_effstepMCHFEcutsITS1Donly->SetTitle("");
251   h_effstepMCHFEcutsITS1Donly->SetStats(0);
252   h_effstepMCHFEcutsITS1Donly->SetLineColor(kBlue);
253   h_effstepMCHFEcutsITS1Donly->SetMarkerColor(kBlue);
254   h_effstepMCHFEcutsITS1Donly->SetMarkerStyle(25);
255   h_effstepMCHFEcutsITS1Donly->Draw("colz");
256   
257   canvascomparisonbis->cd(4);
258   TH2D* h_effstepMCPID1Donly = (TH2D *) efficiencystepMCHFEcutsTRD->Project(0,3);
259   h_effstepMCPID1Donly->SetTitle("");
260   h_effstepMCPID1Donly->SetStats(0);
261   h_effstepMCPID1Donly->SetLineColor(kBlue);
262   h_effstepMCPID1Donly->SetMarkerColor(kBlue);
263   h_effstepMCPID1Donly->SetMarkerStyle(25);
264   h_effstepMCPID1Donly->Draw("colz");
265   
266
267   ///////////
268   // PID
269   ///////////
270   TCanvas * canvasc1Don = new TCanvas("PIDstepESD","PIDstepESD",1000,700);
271   canvasc1Don->Divide(3,1);
272   
273   canvasc1Don->cd(1);
274   TH2D* h_effstepPID1Donly = (TH2D *) efficiencystepPIDTOF->Project(0,3);
275   h_effstepPID1Donly->SetTitle("");
276   h_effstepPID1Donly->SetStats(0);
277   h_effstepPID1Donly->SetLineColor(kBlue);
278   h_effstepPID1Donly->SetMarkerColor(kBlue);
279   h_effstepPID1Donly->SetMarkerStyle(25);
280   gPad->SetLogz();
281   h_effstepPID1Donly->Draw("colz");
282   
283   canvasc1Don->cd(2);
284   TH2D* h_effstepPID1Donlyy = (TH2D *) efficiencystepPIDTPC->Project(0,3);
285   h_effstepPID1Donlyy->SetTitle("");
286   h_effstepPID1Donlyy->SetStats(0);
287   h_effstepPID1Donlyy->SetLineColor(kBlue);
288   h_effstepPID1Donlyy->SetMarkerColor(kBlue);
289   h_effstepPID1Donlyy->SetMarkerStyle(25);
290   gPad->SetLogz();
291   h_effstepPID1Donlyy->Draw("colz");
292   
293   canvasc1Don->cd(3);
294   TH2D* h_effstepPID1Donlyyy = (TH2D *) efficiencystepPIDall->Project(0,3);
295   h_effstepPID1Donlyyy->SetTitle("");
296   h_effstepPID1Donlyyy->SetStats(0);
297   h_effstepPID1Donlyyy->SetLineColor(kBlue);
298   h_effstepPID1Donlyyy->SetMarkerColor(kBlue);
299   h_effstepPID1Donlyyy->SetMarkerStyle(25);
300   gPad->SetLogz();
301   h_effstepPID1Donlyyy->Draw("colz");
302   
303   //
304
305   TCanvas * canvasc1DonMC = new TCanvas("PIDstepMC","PIDstepMC",1000,700);
306   canvasc1DonMC->Divide(3,1);
307   
308   canvasc1DonMC->cd(1);
309   TH2D* h_effstepPID1DMConly = (TH2D *) efficiencystepMCPIDTOF->Project(0,3);
310   h_effstepPID1DMConly->SetTitle("");
311   h_effstepPID1DMConly->SetStats(0);
312  
313   h_effstepPID1DMConly->SetLineColor(kBlue);
314  
315   h_effstepPID1DMConly->SetMarkerColor(kBlue);
316  
317   h_effstepPID1DMConly->SetMarkerStyle(25);
318   gPad->SetLogz(); 
319   h_effstepPID1DMConly->Draw("colz");
320  
321
322   canvasc1DonMC->cd(2);
323   TH2D* h_effstepPID1DMConlyy = (TH2D *) efficiencystepMCPIDTPC->Project(0,3);
324   h_effstepPID1DMConlyy->SetTitle("");
325   h_effstepPID1DMConlyy->SetStats(0);
326   h_effstepPID1DMConlyy->SetLineColor(kBlue);
327   h_effstepPID1DMConlyy->SetMarkerColor(kBlue);
328   h_effstepPID1DMConlyy->SetMarkerStyle(25);
329   gPad->SetLogz();
330   h_effstepPID1DMConlyy->Draw("colz");
331  
332
333   canvasc1DonMC->cd(3);
334   TH2D* h_effstepPID1DMConlyyy = (TH2D *) efficiencystepMCPIDall->Project(0,3);
335   h_effstepPID1DMConlyyy->SetTitle("");
336   h_effstepPID1DMConlyyy->SetStats(0);
337   h_effstepPID1DMConlyyy->SetLineColor(kBlue);
338   h_effstepPID1DMConlyyy->SetMarkerColor(kBlue);
339   h_effstepPID1DMConlyyy->SetMarkerStyle(25);
340   gPad->SetLogz();
341   h_effstepPID1DMConlyyy->Draw("colz");
342
343
344   //////////
345   // all
346   //////////
347
348   TCanvas * canvasall = new TCanvas("AllMC","AllMC",1000,700);
349   //canvasall->Divide(2,1);
350
351   canvasall->cd(1);
352   TH2D* h_effstepPID1Dall = (TH2D *) efficiencystepMCall->Project(0,3);
353   h_effstepPID1Dall->SetTitle("");
354   h_effstepPID1Dall->SetStats(0);
355   h_effstepPID1Dall->SetLineColor(kBlue);
356   h_effstepPID1Dall->SetMarkerColor(kBlue);
357   h_effstepPID1Dall->SetMarkerStyle(25);
358   gPad->SetLogz();
359   h_effstepPID1Dall->Draw("colz");
360
361   
362  
363
364
365 }
366
367 // ====================================================================
368 TList *GetResults(const Char_t *testfile,const Char_t *plus){
369   //
370   // read output
371   //
372   TFile *f = TFile::Open(testfile);
373   if(!f || f->IsZombie()){
374     printf("File not readable\n");
375     return 0x0;
376   }
377   TString name(plus);
378   name += "HFE_Results"; 
379   printf("Name of TList %s\n",(const char*)name); 
380   TList *l = dynamic_cast<TList *>(f->Get((const char*)name));
381   if(!l){
382     printf("Output container not found\n");
383     f->Close(); delete f;
384     return 0x0;
385   } 
386   TList *returnlist = dynamic_cast<TList *>(l->Clone());
387   f->Close(); delete f;
388   return returnlist;
389 }
390 //_________________________________________________________________________
391 TObject* GetSpectrum(AliCFContainer *c, Int_t step) {
392   AliCFDataGrid* data = new AliCFDataGrid("data","",*c,step);
393   //data->SetMeasured(step);
394   return data;
395 }
396 //__________________________________________________________________________
397 TObject* GetEfficiency(AliCFContainer *c, Int_t step) {
398   TString name("eff");
399   name += step;
400   AliCFEffGrid* eff = new AliCFEffGrid((const char*)name,"",*c);
401   eff->CalculateEfficiency(step,0);
402   return eff;
403 }
404 //_________________________________________________________________________
405 TObject* GetEfficiency(AliCFContainer *c, Int_t step, Int_t step0) {
406   TString name("eff");
407   name += step;
408   name+= step0;
409   AliCFEffGrid* eff = new AliCFEffGrid((const char*)name,"",*c);
410   eff->CalculateEfficiency(step,step0);
411   return eff;
412 }
413 //_______________________________________________________________________
414 TH1I *GetNEvents(const Char_t *testfile,const Char_t *plus){
415   //
416   // read output
417   //
418   TFile *f = TFile::Open(testfile);
419   if(!f || f->IsZombie()){
420     printf("File not readable\n");
421     return 0x0;
422   }
423   TString name("nEvents");
424   name += plus; 
425   printf("Name of nEvents %s\n",(const char*)name); 
426   TH1I *l = dynamic_cast<TH1I *>(f->Get((const char*)name));
427   if(!l){
428     printf("nEvents not found\n");
429     f->Close(); delete f;
430     return 0x0;
431   } 
432   TH1I *returnlist = dynamic_cast<TH1I *>(l->Clone());
433   if(!returnlist) return 0x0;
434   returnlist->SetDirectory(0);
435   f->Close(); delete f;
436   return returnlist;
437 }
438 //________________________________________________________________________
439 void CorrectFromTheWidth(TH1D *h1) {
440   //
441   // Correct from the width of the bins --> dN/dp_{T} (GeV/c)^{-1}
442   //
443
444   TAxis *axis = h1->GetXaxis();
445   Int_t nbinX = h1->GetNbinsX();
446
447   for(Int_t i = 1; i <= nbinX; i++) {
448
449     Double_t width = axis->GetBinWidth(i);
450     Double_t content = h1->GetBinContent(i);
451     Double_t error = h1->GetBinError(i); 
452     h1->SetBinContent(i,content/width); 
453     h1->SetBinError(i,error/width);
454   }
455
456 }
457 //____________________________________________________________________________
458 AliCFContainer *GetContainerCentralitySource(AliCFContainer *container,Int_t bincentrality,Int_t source) {
459
460   Int_t *vars = new Int_t[3];
461   vars[0] = 0;
462   vars[1] = 1;
463   vars[2] = 2;
464      
465   Int_t nbBinPt = container->GetNBins(0);
466   Int_t nbBinEta = container->GetNBins(1);
467   Int_t nbBinPhi = container->GetNBins(2);
468   Int_t nbBinCharge = container->GetNBins(3);
469   Int_t nbBinSource = container->GetNBins(4);
470   Int_t nbBinCentrality = container->GetNBins(5);
471   Double_t *arrayPt = new Double_t[nbBinPt + 1];
472   Double_t *arrayEta = new Double_t[nbBinEta + 1];
473   Double_t *arrayPhi = new Double_t[nbBinPhi + 1];
474   Double_t *arrayCharge = new Double_t[nbBinCharge + 1];
475   Double_t *arraySource = new Double_t[nbBinSource + 1];
476   Double_t *arrayCentrality = new Double_t[nbBinCentrality + 1];
477   container->GetBinLimits(0,arrayPt);
478   container->GetBinLimits(1,arrayEta);
479   container->GetBinLimits(2,arrayPhi);
480   container->GetBinLimits(3,arrayCharge);
481   container->GetBinLimits(4,arraySource);
482   container->GetBinLimits(4,arrayCentrality);
483
484   Int_t sourcebin = source;
485   Int_t centralitybin = bincentrality;
486  
487
488   Double_t *varMin = new Double_t[6];
489   Double_t *varMax = new Double_t[6];
490   varMin[0] = arrayPt[0];
491   varMin[1] = arrayEta[0];
492   varMin[2] = arrayPhi[0];
493   varMin[3] = arrayCharge[0];
494   if(source == -1) sourcebin = 0;
495   if(bincentrality == -1) centralitybin = 0;
496   varMin[4] = arraySource[sourcebin];
497   varMin[5] = arraySource[centralitybin];
498   varMax[0] = arrayPt[nbBinPt];
499   varMax[1] = arrayEta[nbBinEta];
500   varMax[2] = arrayPhi[nbBinPhi];
501   varMax[3] = arrayCharge[nbBinCharge];
502   if(source == -1) sourcebin = nbBinSource;
503   if(bincentrality == -1) centralitybin = nbBinCentrality;
504   varMax[4] = arraySource[sourcebin];
505   varMax[5] = arraySource[centralitybin];
506
507   //printf("Nb of bin charge %d\n",nbBinCharge);
508   //for(Int_t nb = 0; nb <= nbBinCharge; nb++) {
509   //  printf("charge %f for %d\n",arrayCharge[nb],nb);
510   //}
511
512   AliCFContainer *k = container->MakeSlice(3,vars,varMin,varMax);
513   k->SetName("Other");
514
515   return k;
516
517 }
518 //____________________________________________________________________________
519 AliCFContainer *GetContainerSourceAsFunctionOfCentrality(AliCFContainer *container,Int_t source) {
520
521   Int_t *vars = new Int_t[4];
522   vars[0] = 0;
523   vars[1] = 1;
524   vars[2] = 2;
525   vars[3] = 5;
526      
527   Int_t nbBinPt = container->GetNBins(0);
528   Int_t nbBinEta = container->GetNBins(1);
529   Int_t nbBinPhi = container->GetNBins(2);
530   Int_t nbBinCharge = container->GetNBins(3);
531   Int_t nbBinSource = container->GetNBins(4);
532   Int_t nbBinCentrality = container->GetNBins(5);
533   Double_t *arrayPt = new Double_t[nbBinPt + 1];
534   Double_t *arrayEta = new Double_t[nbBinEta + 1];
535   Double_t *arrayPhi = new Double_t[nbBinPhi + 1];
536   Double_t *arrayCharge = new Double_t[nbBinCharge + 1];
537   Double_t *arraySource = new Double_t[nbBinSource + 1];
538   Double_t *arrayCentrality = new Double_t[nbBinCentrality + 1];
539   container->GetBinLimits(0,arrayPt);
540   container->GetBinLimits(1,arrayEta);
541   container->GetBinLimits(2,arrayPhi);
542   container->GetBinLimits(3,arrayCharge);
543   container->GetBinLimits(4,arraySource);
544   container->GetBinLimits(5,arrayCentrality);
545
546   Int_t sourcebin = source;
547  
548
549   Double_t *varMin = new Double_t[6];
550   Double_t *varMax = new Double_t[6];
551   varMin[0] = arrayPt[0];
552   varMin[1] = arrayEta[0];
553   varMin[2] = arrayPhi[0];
554   varMin[3] = arrayCharge[0];
555   if(source == -1) sourcebin = 0;
556   varMin[4] = arraySource[sourcebin];
557   varMin[5] = arrayCentrality[0];
558   varMax[0] = arrayPt[nbBinPt];
559   varMax[1] = arrayEta[nbBinEta];
560   varMax[2] = arrayPhi[nbBinPhi];
561   varMax[3] = arrayCharge[nbBinCharge];
562   if(source == -1) sourcebin = nbBinSource;
563   varMax[4] = arraySource[sourcebin];
564   varMax[5] = arrayCentrality[nbBinCentrality];
565
566   //printf("Nb of bin charge %d\n",nbBinCharge);
567   //for(Int_t nb = 0; nb <= nbBinCharge; nb++) {
568   //  printf("charge %f for %d\n",arrayCharge[nb],nb);
569   //}
570
571   AliCFContainer *k = container->MakeSlice(4,vars,varMin,varMax);
572   k->SetName("OtherAsFunctionOfCentrality");
573
574   return k;
575
576 }