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