]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/dNdEta/dNdEtaAnalysis.cxx
adding macro that applies the corrections and visualizes the result
[u/mrichter/AliRoot.git] / PWG0 / dNdEta / dNdEtaAnalysis.cxx
1 /* $Id$ */
2
3 #include "dNdEtaAnalysis.h"
4
5 #include <TFile.h>
6 #include <TH3F.h>
7 #include <TH2F.h>
8 #include <TH1F.h>
9 #include <TMath.h>
10 #include <TCanvas.h>
11 #include <TCollection.h>
12 #include <TIterator.h>
13 #include <TList.h>
14 #include <TLegend.h>
15 #include <TLine.h>
16
17 #include "AlidNdEtaCorrection.h"
18 #include <AliCorrection.h>
19 #include <AliPWG0Helper.h>
20 #include <AliCorrectionMatrix2D.h>
21 #include <AliCorrectionMatrix3D.h>
22
23 //____________________________________________________________________
24 ClassImp(dNdEtaAnalysis)
25
26 //____________________________________________________________________
27 dNdEtaAnalysis::dNdEtaAnalysis() :
28   TNamed(),
29   fData(0),
30   fMult(0),
31   fPtDist(0)
32 {
33   // default constructor
34
35   for (Int_t i=0; i<kVertexBinning; ++i)
36   {
37     fdNdEta[i] = 0;
38     fdNdEtaPtCutOffCorrected[i] = 0;
39   }
40 }
41
42 //____________________________________________________________________
43 dNdEtaAnalysis::dNdEtaAnalysis(Char_t* name, Char_t* title, AliPWG0Helper::AnalysisMode analysisMode) :
44   TNamed(name, title),
45   fData(0),
46   fMult(0),
47   fPtDist(0)
48 {
49   // constructor
50
51   // TODO this binning has to be the same than in AliCorrection, somehow passed?!
52   Float_t binLimitsPt[] = {0.0, 0.05, 0.1, 0.125, 0.15, 0.175, 0.2, 0.225, 0.25, 0.275, 0.3, 0.325, 0.35, 0.375, 0.4, 0.425, 0.45, 0.475, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0, 1.5, 2.0, 5.0, 10.0, 100.0};
53
54   fData = new AliCorrection("Analysis", Form("%s Analysis", title), analysisMode);
55
56   // do not add this hists to the directory
57   Bool_t oldStatus = TH1::AddDirectoryStatus();
58   TH1::AddDirectory(kFALSE);
59
60   fMult = new TH1F("TriggeredMultiplicity", "Triggered Events;raw multiplicity;entries", 1000, -0.5, 999.5);
61
62   fdNdEta[0] = new TH1F("dNdEta", "dN_{ch}/d#eta;#eta;dN_{ch}/d#eta", 40, -2, 2);
63
64   fdNdEtaPtCutOffCorrected[0] = dynamic_cast<TH1F*> (fdNdEta[0]->Clone(Form("%s_corrected", fdNdEta[0]->GetName())));
65
66   for (Int_t i=1; i<kVertexBinning; ++i)
67   {
68     fdNdEta[i]    = dynamic_cast<TH1F*> (fdNdEta[0]->Clone(Form("%s_%d", fdNdEta[0]->GetName(), i)));
69     fdNdEtaPtCutOffCorrected[i]    = dynamic_cast<TH1F*> (fdNdEtaPtCutOffCorrected[0]->Clone(Form("%s_%d", fdNdEtaPtCutOffCorrected[0]->GetName(), i)));
70   }
71
72   fPtDist = new TH1F("Pt", "p_{T} distribution;p_{T} [GeV/c];#frac{dN}{d#eta dp_{T}} [c/GeV]", 28, binLimitsPt);
73
74   TH1::AddDirectory(oldStatus);
75 }
76
77 //____________________________________________________________________
78 dNdEtaAnalysis::~dNdEtaAnalysis()
79 {
80   // destructor
81
82   if (fData)
83   {
84     delete fData;
85     fData = 0;
86   }
87
88   if (fMult)
89   {
90     delete fMult;
91     fMult = 0;
92   }
93
94   for (Int_t i=0; i<kVertexBinning; ++i)
95   {
96     if (fdNdEta[i])
97     {
98       delete fdNdEta[i];
99       fdNdEta[i] = 0;
100     }
101     if (fdNdEtaPtCutOffCorrected[i])
102     {
103       delete fdNdEtaPtCutOffCorrected[i];
104       fdNdEtaPtCutOffCorrected[i] = 0;
105     }
106   }
107
108   if (fPtDist)
109   {
110     delete fPtDist;
111     fPtDist = 0;
112   }
113 }
114
115 //_____________________________________________________________________________
116 dNdEtaAnalysis::dNdEtaAnalysis(const dNdEtaAnalysis &c) :
117   TNamed(c),
118   fData(0),
119   fMult(0),
120   fPtDist(0)
121 {
122   //
123   // dNdEtaAnalysis copy constructor
124   //
125
126   ((dNdEtaAnalysis &) c).Copy(*this);
127 }
128
129 //_____________________________________________________________________________
130 dNdEtaAnalysis &dNdEtaAnalysis::operator=(const dNdEtaAnalysis &c)
131 {
132   //
133   // Assignment operator
134   //
135
136   if (this != &c) ((dNdEtaAnalysis &) c).Copy(*this);
137   return *this;
138 }
139
140 //_____________________________________________________________________________
141 void dNdEtaAnalysis::Copy(TObject &c) const
142 {
143   //
144   // Copy function
145   //
146
147   dNdEtaAnalysis& target = (dNdEtaAnalysis &) c;
148
149   target.fData = dynamic_cast<AliCorrection*> (fData->Clone());
150   target.fMult = dynamic_cast<TH1F*> (fMult->Clone());
151
152   for (Int_t i=0; i<kVertexBinning; ++i)
153   {
154     target.fdNdEta[i] = dynamic_cast<TH1F*> (fdNdEta[i]->Clone());
155     target.fdNdEtaPtCutOffCorrected[i] = dynamic_cast<TH1F*> (fdNdEtaPtCutOffCorrected[i]->Clone());
156   }
157
158   target.fPtDist = dynamic_cast<TH1F*> (fPtDist->Clone());
159
160   TNamed::Copy((TNamed &) c);
161 }
162
163 //____________________________________________________________________
164 void dNdEtaAnalysis::FillTrack(Float_t vtx, Float_t eta, Float_t pt)
165 {
166   // fills a track into the histograms
167
168   fData->GetTrackCorrection()->FillMeas(vtx, eta, pt);
169 }
170
171 //____________________________________________________________________
172 void dNdEtaAnalysis::FillEvent(Float_t vtx, Float_t n)
173 {
174   // fills an event into the histograms
175
176   fData->GetEventCorrection()->FillMeas(vtx, n);
177 }
178
179 //____________________________________________________________________
180 void dNdEtaAnalysis::FillTriggeredEvent(Float_t n)
181 {
182   // fills a triggered event into the histograms
183
184   fMult->Fill(n);
185 }
186
187 //____________________________________________________________________
188 void dNdEtaAnalysis::Finish(AlidNdEtaCorrection* correction, Float_t ptCut, AlidNdEtaCorrection::CorrectionType correctionType, Int_t multCut)
189 {
190   //
191   // correct with the given correction values and calculate dNdEta and pT distribution
192   // the corrections that are applied can be steered by the flag correctionType
193   // the measured result is not used up to a multiplicity of multCut (the bin at multCut is the first that is used!)
194   //
195
196   // TODO put tag somewhere which corrections have been applied
197
198   if (multCut > 1)
199   {
200     Printf("ERROR: A bigger multiplicity cut than 1 is not possible in the current implementation");
201     return;
202   }
203
204   // set corrections to 1
205   fData->SetCorrectionToUnity();
206
207   if (correction && correctionType != AlidNdEtaCorrection::kNone)
208   {
209     TH3F* trackCorr = fData->GetTrackCorrection()->GetCorrectionHistogram();
210     TH2F* eventCorr = fData->GetEventCorrection()->GetCorrectionHistogram();
211
212     if (correctionType >= AlidNdEtaCorrection::kTrack2Particle)
213       trackCorr->Multiply(correction->GetTrack2ParticleCorrection()->GetTrackCorrection()->GetCorrectionHistogram());
214
215     if (correctionType >= AlidNdEtaCorrection::kVertexReco)
216     {
217       trackCorr->Multiply(correction->GetVertexRecoCorrection()->GetTrackCorrection()->GetCorrectionHistogram());
218       eventCorr->Multiply(correction->GetVertexRecoCorrection()->GetEventCorrection()->GetCorrectionHistogram());
219
220       // set bin with multiplicity 0 to 1 (correction has no meaning in this bin)
221       for (Int_t i=0; i<=eventCorr->GetNbinsX()+1; i++)
222         eventCorr->SetBinContent(i, 1, 1);
223     }
224
225     switch (correctionType)
226     {
227       case AlidNdEtaCorrection::kINEL :
228       {
229         trackCorr->Multiply(correction->GetTriggerBiasCorrectionINEL()->GetTrackCorrection()->GetCorrectionHistogram());
230         eventCorr->Multiply(correction->GetTriggerBiasCorrectionINEL()->GetEventCorrection()->GetCorrectionHistogram());
231         break;
232       }
233       case AlidNdEtaCorrection::kNSD :
234       {
235         trackCorr->Multiply(correction->GetTriggerBiasCorrectionNSD()->GetTrackCorrection()->GetCorrectionHistogram());
236         eventCorr->Multiply(correction->GetTriggerBiasCorrectionNSD()->GetEventCorrection()->GetCorrectionHistogram());
237         break;
238       }
239       case AlidNdEtaCorrection::kND :
240       {
241         trackCorr->Multiply(correction->GetTriggerBiasCorrectionND()->GetTrackCorrection()->GetCorrectionHistogram());
242         eventCorr->Multiply(correction->GetTriggerBiasCorrectionND()->GetEventCorrection()->GetCorrectionHistogram());
243         break;
244       }
245       default : break;
246     }
247   }
248   else
249     printf("INFO: No correction applied\n");
250
251   fData->Multiply();
252
253   if (correctionType >= AlidNdEtaCorrection::kVertexReco)
254   {
255     // There are no events with vertex that have 0 multiplicity, therefore
256     //   populate bin with 0 multiplicity with the following idea:
257     //     alpha = triggered events with vertex at a given vertex position / all triggered events with vertex
258     //     triggered events without vertex and 0 multiplicity at a given vertex position = alpha * all triggered events with 0 multiplicity
259     //   afterwards we still correct for the trigger efficiency
260
261     TH2* measuredEvents = fData->GetEventCorrection()->GetMeasuredHistogram();
262     TH2* correctedEvents = fData->GetEventCorrection()->GetGeneratedHistogram();
263
264     //new TCanvas; correctedEvents->DrawCopy("TEXT");
265
266     // start above 0 mult. bin with integration
267     TH1* vertexDist = measuredEvents->ProjectionX("vertexdist_measured", 2);
268     //new TCanvas; vertexDist->DrawCopy();
269
270     Int_t allEventsWithVertex = (Int_t) vertexDist->Integral(0, vertexDist->GetNbinsX()+1); // include under/overflow!
271     Int_t triggeredEventsWith0Mult = (Int_t) fMult->GetBinContent(1);
272
273     Printf("%d triggered events with 0 mult. -- %d triggered events with vertex", triggeredEventsWith0Mult, allEventsWithVertex);
274
275     for (Int_t i = 1; i <= measuredEvents->GetNbinsX(); i++)
276     {
277       Double_t alpha = vertexDist->GetBinContent(i) / allEventsWithVertex;
278       Double_t events = alpha * triggeredEventsWith0Mult;
279
280       // multiply with trigger correction if set above
281       events *= fData->GetEventCorrection()->GetCorrectionHistogram()->GetBinContent(i, 1);
282
283       Printf("Bin %d, alpha is %.2f, number of events with 0 mult. are %.2f", i, alpha, events);
284
285       correctedEvents->SetBinContent(i, 1, events);
286     }
287
288     //new TCanvas; correctedEvents->DrawCopy("TEXT");
289   }
290
291   fData->PrintInfo(ptCut);
292
293   TH3F* dataHist = fData->GetTrackCorrection()->GetGeneratedHistogram();
294
295   // integrate multiplicity axis out (include under/overflow bins!!!)
296   TH2F* tmp = fData->GetEventCorrection()->GetGeneratedHistogram();
297
298   TH1D* vertexHist = (TH1D*) tmp->ProjectionX("_px", 0, tmp->GetNbinsY() + 1, "e");
299
300   // create pt hist
301   {
302     // reset all ranges
303     dataHist->GetXaxis()->SetRange(0, 0);
304     dataHist->GetYaxis()->SetRange(0, 0);
305     dataHist->GetZaxis()->SetRange(0, 0);
306
307     // vtx cut
308     Int_t vertexBinBegin = dataHist->GetXaxis()->FindBin(-5);
309     Int_t vertexBinEnd = dataHist->GetXaxis()->FindBin(5);
310     dataHist->GetXaxis()->SetRange(vertexBinBegin, vertexBinEnd);
311     Float_t nEvents = vertexHist->Integral(vertexBinBegin, vertexBinEnd);
312
313     if (nEvents > 0)
314     {
315       // eta cut
316       dataHist->GetYaxis()->SetRange(dataHist->GetYaxis()->FindBin(-0.8), dataHist->GetYaxis()->FindBin(0.8));
317       Float_t etaWidth = 1.6;
318
319       TH1D* ptHist = dynamic_cast<TH1D*> (dataHist->Project3D("ze"));
320
321       for (Int_t i=1; i<=fPtDist->GetNbinsX(); ++i)
322       {
323         Float_t binSize = fPtDist->GetBinWidth(i);
324         fPtDist->SetBinContent(i, ptHist->GetBinContent(i) / binSize / nEvents / etaWidth);
325         fPtDist->SetBinError(i, ptHist->GetBinError(i) / binSize / nEvents / etaWidth);
326       }
327
328       delete ptHist;
329     }
330     else
331       printf("ERROR: nEvents is 0!\n");
332   }
333
334   // reset all ranges
335   dataHist->GetXaxis()->SetRange(0, 0);
336   dataHist->GetYaxis()->SetRange(0, 0);
337   dataHist->GetZaxis()->SetRange(0, 0);
338
339   // integrate over pt (with pt cut)
340   Int_t ptLowBin = 1;
341   if (ptCut > 0)
342     ptLowBin = dataHist->GetZaxis()->FindBin(ptCut);
343
344   dataHist->GetZaxis()->SetRange(ptLowBin, dataHist->GetZaxis()->GetNbins()+1);
345   printf("pt range %d %d\n", ptLowBin, dataHist->GetZaxis()->GetNbins()+1);
346   TH2D* vtxVsEta = dynamic_cast<TH2D*> (dataHist->Project3D("yx2e"));
347
348   dataHist->GetZaxis()->SetRange(0, 0);
349   vtxVsEta->GetXaxis()->SetTitle(dataHist->GetXaxis()->GetTitle());
350   vtxVsEta->GetYaxis()->SetTitle(dataHist->GetYaxis()->GetTitle());
351
352   if (vtxVsEta == 0)
353   {
354     printf("ERROR: pt integration failed\n");
355     return;
356   }
357
358   //new TCanvas; vtxVsEta->DrawCopy("COLZ");
359   //vtxVsEta->Rebin2D(1, 4);
360
361   const Float_t vertexRange = 9.99;
362
363   for (Int_t iEta=1; iEta<=vtxVsEta->GetNbinsY(); iEta++)
364   {
365     // do we have several histograms for different vertex positions?
366     Int_t vertexBinGlobalBegin = vertexHist->GetXaxis()->FindBin(-vertexRange);
367     Int_t vertexBinWidth = (vertexHist->GetXaxis()->FindBin(vertexRange) - vertexBinGlobalBegin + 1) / (kVertexBinning-1);
368     //printf("vertexBinGlobalBegin = %d, vertexBinWidth = %d\n", vertexBinGlobalBegin, vertexBinWidth);
369     for (Int_t vertexPos=0; vertexPos<kVertexBinning; ++vertexPos)
370     {
371       Int_t vertexBinBegin = vertexBinGlobalBegin;
372       Int_t vertexBinEnd   = vertexBinGlobalBegin + vertexBinWidth * (kVertexBinning-1);
373
374       // the first histogram is always for the whole vertex range
375       if (vertexPos > 0)
376       {
377         vertexBinBegin = vertexBinGlobalBegin + vertexBinWidth * (vertexPos-1);
378         vertexBinEnd = vertexBinBegin + vertexBinWidth;
379       }
380
381       //printf("vertexBinBegin = %d, vertexBinEnd = %d\n", vertexBinBegin, vertexBinEnd);
382
383       Float_t totalEvents = vertexHist->Integral(vertexBinBegin, vertexBinEnd - 1);
384       if (totalEvents == 0)
385       {
386         printf("WARNING: No events for hist %d %d %d\n", vertexPos, vertexBinBegin, vertexBinEnd);
387         continue;
388       }
389
390       Float_t sum = 0;
391       Float_t sumError2 = 0;
392       for (Int_t iVtx = vertexBinBegin; iVtx < vertexBinEnd; iVtx++)
393       {      
394         if (vtxVsEta->GetBinContent(iVtx, iEta) != 0)
395         {
396           sum = sum + vtxVsEta->GetBinContent(iVtx, iEta);
397
398         if (sumError2 > 10e30)
399           printf("WARNING: sum of error2 is dangerously large - be prepared for crash... ");
400
401           sumError2 = sumError2 + TMath::Power(vtxVsEta->GetBinError(iVtx, iEta),2);
402         }
403       }
404
405       Float_t ptCutOffCorrection = 1;
406       if (correction && ptCut > 0)
407         ptCutOffCorrection = correction->GetMeasuredFraction(correctionType, ptCut, vtxVsEta->GetYaxis()->GetBinCenter(iEta));
408
409       if (ptCutOffCorrection <= 0)
410       {
411         printf("UNEXPECTED: ptCutOffCorrection is %f for hist %d %d %d\n", ptCutOffCorrection, vertexPos, vertexBinBegin, vertexBinEnd);
412         continue;
413       }
414
415       //printf("Eta: %d Vertex Range: %d %d, Event Count %f, Track Sum: %f, Track Sum corrected: %f\n", iEta, vertexBinBegin, vertexBinEnd, totalEvents, sum, sum / ptCutOffCorrection);
416
417       Int_t bin = fdNdEta[vertexPos]->FindBin(vtxVsEta->GetYaxis()->GetBinCenter(iEta));
418       if (bin > 0 && bin < fdNdEta[vertexPos]->GetNbinsX())
419       {
420         Float_t dndeta = sum / totalEvents;
421         Float_t error  = TMath::Sqrt(sumError2) / totalEvents;
422
423         dndeta = dndeta / fdNdEta[vertexPos]->GetBinWidth(bin);
424         error  = error / fdNdEta[vertexPos]->GetBinWidth(bin);
425
426         fdNdEta[vertexPos]->SetBinContent(bin, dndeta);
427         fdNdEta[vertexPos]->SetBinError(bin, error);
428
429         dndeta /= ptCutOffCorrection;
430         error  /= ptCutOffCorrection;
431
432         fdNdEtaPtCutOffCorrected[vertexPos]->SetBinContent(bin, dndeta);
433         fdNdEtaPtCutOffCorrected[vertexPos]->SetBinError(bin, error);
434
435         //Printf("Bin %d has dN/deta = %f", bin, dndeta);
436       }
437     }
438   }
439
440   //new TCanvas; fdNdEta[0]->DrawCopy();
441 }
442
443 //____________________________________________________________________
444 void dNdEtaAnalysis::SaveHistograms()
445 {
446   // save the histograms to a directory with the name of this class (retrieved from TNamed)
447
448   gDirectory->mkdir(GetName());
449   gDirectory->cd(GetName());
450
451   if (fData)
452   {
453     fData->SaveHistograms();
454   }
455   else
456     printf("dNdEtaAnalysis::SaveHistograms: UNEXPECTED: fData is 0\n");
457
458   if (fMult)
459   {
460     fMult->Write();
461   }
462   else
463     printf("dNdEtaAnalysis::SaveHistograms: UNEXPECTED: fMult is 0\n");
464
465   if (fPtDist)
466     fPtDist       ->Write();
467   else
468     printf("dNdEtaAnalysis::SaveHistograms: UNEXPECTED: fPtDist is 0\n");
469
470   for (Int_t i=0; i<kVertexBinning; ++i)
471   {
472     if (fdNdEta[i])
473       fdNdEta[i]->Write();
474     else
475       printf("dNdEtaAnalysis::SaveHistograms: UNEXPECTED: fdNdEta[%d] is 0\n", i);
476
477     if (fdNdEtaPtCutOffCorrected[i])
478       fdNdEtaPtCutOffCorrected[i]->Write();
479     else
480       printf("dNdEtaAnalysis::SaveHistograms: UNEXPECTED: fdNdEtaPtCutOffCorrected[%d] is 0\n", i);
481   }
482
483   gDirectory->cd("../");
484 }
485
486 void dNdEtaAnalysis::LoadHistograms(const Char_t* dir)
487 {
488   // loads the histograms from a directory with the name of this class (retrieved from TNamed)
489
490   if (!dir)
491     dir = GetName();
492
493   gDirectory->cd(dir);
494
495   fData->LoadHistograms();
496   fMult = dynamic_cast<TH1F*> (gDirectory->Get(fMult->GetName()));
497
498   for (Int_t i=0; i<kVertexBinning; ++i)
499   {
500     fdNdEta[i] = dynamic_cast<TH1F*> (gDirectory->Get(fdNdEta[i]->GetName()));
501     fdNdEtaPtCutOffCorrected[i] = dynamic_cast<TH1F*> (gDirectory->Get(fdNdEtaPtCutOffCorrected[i]->GetName()));
502   }
503
504   fPtDist = dynamic_cast<TH1F*> (gDirectory->Get(fPtDist->GetName()));
505
506   gDirectory->cd("../");
507 }
508
509 //____________________________________________________________________
510 void dNdEtaAnalysis::DrawHistograms(Bool_t simple)
511 {
512   // draws the histograms
513
514   if (!simple)
515   {
516     if (fData)
517       fData->DrawHistograms(GetName());
518
519     TCanvas* canvas = new TCanvas(Form("%s_dNdEtaAnalysis", GetName()), Form("%s_dNdEtaAnalysis", GetName()), 800, 400);
520     canvas->Divide(2, 1);
521
522     canvas->cd(1);
523     if (fdNdEtaPtCutOffCorrected[0])
524       fdNdEtaPtCutOffCorrected[0]->DrawCopy();
525
526     if (fdNdEta[0])
527     {
528       fdNdEta[0]->SetLineColor(kRed);
529       fdNdEta[0]->DrawCopy("SAME");
530     }
531
532     canvas->cd(2);
533     if (fPtDist)
534       fPtDist->DrawCopy();
535   }
536
537     // histograms for different vertices?
538   if (kVertexBinning > 0)
539   {
540       // doesnt work, but i dont get it, giving up...
541     TCanvas* canvas2 = new TCanvas(Form("%s_dNdEtaAnalysisVtx", GetName()), Form("%s_dNdEtaAnalysisVtx", GetName()), 450, 450);
542     TCanvas* canvas3 = 0;
543     if (!simple)
544       canvas3 = new TCanvas(Form("%s_dNdEtaAnalysisVtx_noptcutoff", GetName()), Form("%s_dNdEtaAnalysisVtx_noptcutoff", GetName()), 450, 450);
545
546     //Int_t yPads = (Int_t) TMath::Ceil(((Double_t) kVertexBinning - 1) / 2);
547     //printf("%d\n", yPads);
548     //canvas2->Divide(2, yPads);
549
550     TLegend* legend = new TLegend(0.4, 0.2, 0.6, 0.4);
551
552     for (Int_t i=0; i<kVertexBinning; ++i)
553     {
554       if (fdNdEtaPtCutOffCorrected[i])
555       {
556         canvas2->cd();
557
558         fdNdEtaPtCutOffCorrected[i]->SetLineColor(i+1);
559         fdNdEtaPtCutOffCorrected[i]->DrawCopy((i == 0) ? "" : "SAME");
560         legend->AddEntry(fdNdEtaPtCutOffCorrected[i], (i == 0) ? "Vtx All" : Form("Vtx Bin %d", i-1));
561       }
562       if (canvas3 && fdNdEta[i])
563       {
564         canvas3->cd();
565
566         fdNdEta[i]->SetLineColor(i+1);
567         fdNdEta[i]->DrawCopy((i == 0) ? "" : "SAME");
568       }
569     }
570
571     canvas2->cd();
572     legend->Draw();
573     canvas2->SaveAs(Form("%s_%s.gif", canvas2->GetName(), GetName()));
574
575     if (canvas3)
576     {
577       canvas3->cd();
578       legend->Draw();
579     }
580   }
581
582   if (kVertexBinning == 3)
583   {
584      TH1* clone = dynamic_cast<TH1*> (fdNdEtaPtCutOffCorrected[1]->Clone("clone"));
585      TH1* clone2 = dynamic_cast<TH1*> (fdNdEtaPtCutOffCorrected[2]->Clone("clone2"));
586
587      if (clone && clone2)
588      {
589         TCanvas* canvas4 = new TCanvas(Form("%s_dNdEtaAnalysisVtxRatios", GetName()), Form("%s_dNdEtaAnalysisVtxRatios", GetName()), 450, 450);
590
591         clone->Divide(fdNdEtaPtCutOffCorrected[0]);
592         clone->GetYaxis()->SetRangeUser(0.95, 1.05);
593         clone->DrawCopy();
594
595         clone2->Divide(fdNdEtaPtCutOffCorrected[0]);
596         clone2->DrawCopy("SAME");
597
598         TLine* line = new TLine(-1, 1, 1, 1);
599         line->Draw();
600
601         canvas4->SaveAs(Form("%s_%s.gif", canvas4->GetName(), GetName()));
602      }
603    }
604 }
605
606 Long64_t dNdEtaAnalysis::Merge(TCollection* list)
607 {
608   // Merges a list of dNdEtaAnalysis objects with this one.
609   // This is needed for PROOF.
610   // Returns the number of merged objects (including this)
611
612   if (!list)
613     return 0;
614
615   if (list->IsEmpty())
616     return 1;
617
618   TIterator* iter = list->MakeIterator();
619   TObject* obj;
620
621   // sub collections
622   const Int_t nCollections = 2 * kVertexBinning + 3; // 3 standalone hists, 3 arrays of size kVertexBinning
623   TList* collections[nCollections];
624   for (Int_t i=0; i<nCollections; ++i)
625     collections[i] = new TList;
626
627   Int_t count = 0;
628   while ((obj = iter->Next()))
629   {
630     dNdEtaAnalysis* entry = dynamic_cast<dNdEtaAnalysis*> (obj);
631     if (entry == 0)
632       continue;
633
634     collections[0]->Add(entry->fData);
635     collections[2]->Add(entry->fMult);
636     collections[3]->Add(entry->fPtDist);
637
638     for (Int_t i=0; i<kVertexBinning; ++i)
639     {
640       collections[3+i]->Add(entry->fdNdEta[i]);
641       collections[3+kVertexBinning+i]->Add(entry->fdNdEtaPtCutOffCorrected[i]);
642     }
643
644     ++count;
645   }
646
647   fData->Merge(collections[0]);
648   fMult->Merge(collections[1]);
649   fPtDist->Merge(collections[2]);
650
651   for (Int_t i=0; i<kVertexBinning; ++i)
652   {
653     fdNdEta[i]->Merge(collections[3+i]);
654     fdNdEtaPtCutOffCorrected[i]->Merge(collections[3+kVertexBinning+i]);
655   }
656
657   for (Int_t i=0; i<nCollections; ++i)
658     delete collections[i];
659
660   return count+1;
661 }