changed binning: accepts events in z +- 30cm and eta +- 3
[u/mrichter/AliRoot.git] / PWG0 / AliCorrection.cxx
1 /* $Id$ */
2
3 // ------------------------------------------------------
4 //
5 // Class to handle corrections.
6 //
7 // ------------------------------------------------------
8 //
9
10 #include <TFile.h>
11 #include <TCanvas.h>
12 #include <TH3F.h>
13 #include <TH2F.h>
14
15 #include <AliLog.h>
16 #include "AliCorrectionMatrix2D.h"
17 #include "AliCorrectionMatrix3D.h"
18
19 #include "AliCorrection.h"
20
21 //____________________________________________________________________
22 ClassImp(AliCorrection)
23
24 //____________________________________________________________________
25 AliCorrection::AliCorrection() : TNamed(),
26   fEventCorr(0),
27   fTrackCorr(0)
28 {
29   // default constructor
30 }
31
32 //____________________________________________________________________
33 AliCorrection::AliCorrection(const Char_t* name, const Char_t* title, AliPWG0Helper::AnalysisMode analysisMode) : TNamed(name, title),
34   fEventCorr(0),
35   fTrackCorr(0)
36 {
37   // constructor initializing tnamed
38
39   Float_t* binLimitsPt = 0;
40   Int_t nBinsPt = 0;
41
42   // different binnings, better solution could be anticipated...
43   if (analysisMode == AliPWG0Helper::kTPC || analysisMode == AliPWG0Helper::kTPCITS)
44   {
45     static Float_t binLimitsPtTmp[] = {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};
46     binLimitsPt = (Float_t*) binLimitsPtTmp;
47     nBinsPt = 28;
48   }
49   else if (analysisMode == AliPWG0Helper::kSPD)
50   {
51     static Float_t binLimitsPtTmp[] = {-0.5, 100.0};
52     binLimitsPt = (Float_t*) binLimitsPtTmp;
53     nBinsPt = 1;
54   }
55
56   if (!binLimitsPt)
57   {
58     Printf("FATAL: Invalid binning");
59     return;
60   }
61
62   Float_t binLimitsN[]   = {-0.5, 0.5, 1.5, 2.5, 3.5, 4.5, 5.5, 6.5, 7.5, 8.5, 9.5, 10.5, 12.5, 14.5, 16.5, 18.5, 20.5, 25.5, 30.5, 40.5, 50.5, 100.5, 300.5};
63   //Float_t binLimitsVtx[] = {-20,-15,-10,-6,-3,0,3,6,10,15,20};
64   //Float_t binLimitsVtx[] = {-20,-19,-18,-17,-16,-15,-14,-13,-12,-11,-10,-9,-8,-7,-6,-5,-4,-3,-2,-1,0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16,17,18,19,20};
65   Float_t binLimitsVtx[] = {-30,-25,-20,-15,-10,-8,-6,-4,-2,0,2,4,6,8,10,15,20,25,30};
66   Float_t binLimitsEta[] = {-3.0,-2.8,-2.6,-2.4,-2.2,
67                             -2.0,-1.8,-1.6,-1.4,-1.2,
68                             -1.0,-0.8,-0.6,-0.4,-0.2, 0.0,
69                              0.2, 0.4, 0.6, 0.8, 1.0,
70                              1.2, 1.4, 1.6, 1.8, 2.0,
71                              2.2, 2.4, 2.6, 2.8, 3.0};
72 //  Float_t binLimitsEta[] = {-2,-1.8,-1.6,-1.4,-1.2,-1.0,-0.8,-0.6,-0.4,-0.2,0.0,0.2,0.4,0.6,0.8,1.0,1.2,1.4,1.6,1.8,2.0};
73 //  Float_t binLimitsEta[] = {-2,-1.9,-1.8,-1.7,-1.6,-1.5,-1.4,-1.3,-1.2,-1.1,-1.0,-0.9,-0.8,-0.7,-0.6,-0.5,-0.4,-0.3,-0.2,-0.1,0.0,0.1,0.2,0.3,0.4,0.5,0.6,0.7,0.8,0.9,1.0,1.1,1.2,1.3,1.4,1.5,1.6,1.7,1.8,1.9,2.0};
74 //  Float_t binLimitsEta[] = { -7.0, -6.9, -6.8, -6.7, -6.6, -6.5, -6.4, -6.3, -6.2, -6.1, -6.0, -5.9, -5.8, -5.7, -5.6, -5.5, -5.4, -5.3, -5.2, -5.1, -5.0, -4.9, -4.8, -4.7, -4.6, -4.5, -4.4, -4.3, -4.2, -4.1, -4.0, -3.9, -3.8, -3.7, -3.6, -3.5, -3.4, -3.3, -3.2, -3.1, -3.0, -2.9, -2.8, -2.7, -2.6, -2.5, -2.4, -2.3, -2.2, -2.1, -2.0, -1.9, -1.8, -1.7, -1.6, -1.5, -1.4, -1.3, -1.2, -1.1, -1.0, -0.9, -0.8, -0.7, -0.6, -0.5, -0.4, -0.3, -0.2, -0.1, -0.0, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0, 1.1, 1.2, 1.3, 1.4, 1.5, 1.6, 1.7, 1.8, 1.9, 2.0, 2.1, 2.2, 2.3, 2.4, 2.5, 2.6, 2.7, 2.8, 2.9, 3.0, 3.1, 3.2, 3.3, 3.4, 3.5, 3.6, 3.7, 3.8, 3.9, 4.0, 4.1, 4.2, 4.3, 4.4, 4.5, 4.6, 4.7, 4.8, 4.9, 5.0, 5.1, 5.2, 5.3, 5.4, 5.5, 5.6, 5.7, 5.8, 5.9, 6.0, 6.1, 6.2, 6.3, 6.4, 6.5, 6.6, 6.7, 6.8, 6.9, 7.0 };
75
76   TH3F* dummyBinning = new TH3F("dummyBinning","dummyBinning",18, binLimitsVtx, 30, binLimitsEta , nBinsPt, binLimitsPt);
77
78   fEventCorr = new AliCorrectionMatrix2D("EventCorrection", Form("%s EventCorrection", fTitle.Data()), 18, binLimitsVtx, 22, binLimitsN);
79   fTrackCorr = new AliCorrectionMatrix3D("TrackCorrection", Form("%s TrackCorrection", fTitle.Data()), dummyBinning);
80
81   delete dummyBinning;
82
83   fEventCorr->SetAxisTitles("vtx z [cm]", "Ntracks");
84   fTrackCorr->SetAxisTitles("vtx z [cm]", "#eta", "p_{T} [GeV/c]");
85 }
86
87 //____________________________________________________________________
88 AliCorrection::AliCorrection(const AliCorrection& c) : TNamed(c),
89   fEventCorr(0),
90   fTrackCorr(0)
91 {
92   // copy constructor
93   ((AliCorrection &)c).Copy(*this);
94 }
95
96 //____________________________________________________________________
97 AliCorrection::~AliCorrection()
98 {
99   //
100   // destructor
101   //
102
103   if (fEventCorr)
104   {
105     delete fEventCorr;
106     fEventCorr = 0;
107   }
108
109   if (fTrackCorr)
110   {
111     delete fTrackCorr;
112     fTrackCorr = 0;
113   }
114 }
115
116 //____________________________________________________________________
117 AliCorrection &AliCorrection::operator=(const AliCorrection &c)
118 {
119   // assigment operator
120
121   if (this != &c)
122     ((AliCorrection &) c).Copy(*this);
123
124   return *this;
125 }
126
127 //____________________________________________________________________
128 void AliCorrection::Copy(TObject& c) const
129 {
130   // copy function
131
132   AliCorrection& target = (AliCorrection &) c;
133
134   if (fEventCorr)
135     target.fEventCorr = dynamic_cast<AliCorrectionMatrix2D*> (fEventCorr->Clone());
136
137   if (fTrackCorr)
138     target.fTrackCorr = dynamic_cast<AliCorrectionMatrix3D*> (fTrackCorr->Clone());
139 }
140
141 //____________________________________________________________________
142 Long64_t AliCorrection::Merge(TCollection* list)
143 {
144   // Merge a list of AliCorrection objects with this (needed for
145   // PROOF). 
146   // Returns the number of merged objects (including this).
147
148   if (!list)
149     return 0;
150   
151   if (list->IsEmpty())
152     return 1;
153
154   TIterator* iter = list->MakeIterator();
155   TObject* obj;
156
157   // collections of measured and generated histograms
158   TList* collectionEvent = new TList;
159   TList* collectionTrack = new TList;
160   
161   Int_t count = 0;
162   while ((obj = iter->Next())) {
163     
164     AliCorrection* entry = dynamic_cast<AliCorrection*> (obj);
165     if (entry == 0) 
166       continue;
167
168     collectionEvent->Add(entry->fEventCorr);
169     collectionTrack->Add(entry->fTrackCorr);
170
171     count++;
172   }
173   fEventCorr->Merge(collectionEvent);
174   fTrackCorr->Merge(collectionTrack);
175
176   delete collectionEvent;
177   delete collectionTrack;
178
179   return count+1;
180 }
181
182 //____________________________________________________________________
183 void AliCorrection::Divide()
184 {
185   //
186   // divide the histograms to get the correction
187   //
188   
189   if (!fEventCorr || !fTrackCorr)
190     return;
191     
192   fEventCorr->Divide();
193   fTrackCorr->Divide();
194
195   Int_t emptyBins = fTrackCorr->CheckEmptyBins(-9.99, 9.99, -0.79, 0.79, 0.3, 4.9);
196   printf("INFO: In the central region the track correction of %s has %d empty bins\n", GetName(), emptyBins);
197 }
198
199 //____________________________________________________________________
200 void AliCorrection::Add(AliCorrection* aCorrectionToAdd, Float_t c)
201 {
202   //
203   // add to measured and generated the measured and generated of aCorrectionToAdd
204   // with the weight c
205
206   fEventCorr->Add(aCorrectionToAdd->GetEventCorrection(),c);
207   fTrackCorr->Add(aCorrectionToAdd->GetTrackCorrection(),c);
208 }
209
210
211 //____________________________________________________________________
212 Bool_t AliCorrection::LoadHistograms(const Char_t* dir)
213 {
214   //
215   // loads the histograms from a file
216   // if dir is empty a directory with the name of this object is taken (like in SaveHistogram)
217   //
218
219   if (!fEventCorr || !fTrackCorr)
220     return kFALSE;
221
222   if (!dir)
223     dir = GetName();
224
225   if (!gDirectory->cd(dir))
226     return kFALSE;
227
228   Bool_t success = fEventCorr->LoadHistograms();
229   success &= fTrackCorr->LoadHistograms();
230
231   gDirectory->cd("..");
232
233   return success;
234 }
235
236 //____________________________________________________________________
237 void AliCorrection::SaveHistograms()
238 {
239   //
240   // saves the histograms in a directory with the name of this object (GetName)
241   //
242   
243   gDirectory->mkdir(GetName());
244   gDirectory->cd(GetName());
245
246   if (fEventCorr)
247     fEventCorr->SaveHistograms();
248
249   if (fTrackCorr)
250     fTrackCorr->SaveHistograms();
251     
252   gDirectory->cd("..");
253 }
254
255 //____________________________________________________________________
256 void AliCorrection::ReduceInformation()
257 {
258   // this function deletes the measured and generated histograms to reduce the amount of data
259   // in memory
260
261   if (!fEventCorr || !fTrackCorr)
262     return;
263
264   fEventCorr->ReduceInformation();
265   fTrackCorr->ReduceInformation();
266 }
267
268 //____________________________________________________________________
269 void AliCorrection::Reset(Option_t* option)
270 {
271   // resets the histograms
272
273   if (fEventCorr)
274     fEventCorr->Reset(option);
275
276   if (fTrackCorr)
277     fTrackCorr->Reset(option);
278 }
279
280 //____________________________________________________________________
281 void AliCorrection::DrawHistograms(const Char_t* name)
282 {
283   // draws the corrections
284
285   if (!name)
286     name = GetName();
287
288   if (fEventCorr)
289     fEventCorr->DrawHistograms(Form("%s event", name));
290
291   if (fTrackCorr)
292     fTrackCorr->DrawHistograms(Form("%s track", name));
293 }
294
295 void AliCorrection::DrawOverview(const char* canvasName)
296 {
297   // draw projection of the corrections
298   //   to the 3 axis of fTrackCorr and 2 axis of fEventCorr
299
300   TString canvasNameTmp(GetName());
301   if (canvasName)
302     canvasNameTmp = canvasName;
303
304   TCanvas* canvas = new TCanvas(canvasNameTmp, canvasNameTmp, 1200, 800);
305   canvas->Divide(3, 2);
306
307   if (fTrackCorr) {
308     canvas->cd(1);
309     fTrackCorr->Get1DCorrectionHistogram("x", 0.3, 5, -1, 1)->DrawCopy()->GetYaxis()->SetRangeUser(0, 10);
310
311     canvas->cd(2);
312     fTrackCorr->Get1DCorrectionHistogram("y", 0.3, 5, 0, 0)->DrawCopy()->GetYaxis()->SetRangeUser(0, 10);
313
314     canvas->cd(3);
315     fTrackCorr->Get1DCorrectionHistogram("z", 0, -1, -1, 1)->DrawCopy()->GetYaxis()->SetRangeUser(0, 10);
316   }
317
318   if (fEventCorr)
319   {
320     canvas->cd(4);
321     fEventCorr->Get1DCorrectionHistogram("x")->DrawCopy();
322
323     canvas->cd(5);
324     fEventCorr->Get1DCorrectionHistogram("y")->DrawCopy()->GetXaxis()->SetRangeUser(0, 30);
325   }
326 }
327
328 //____________________________________________________________________
329 void AliCorrection::SetCorrectionToUnity()
330 {
331   // set the corrections to unity
332
333   if (fEventCorr)
334     fEventCorr->SetCorrectionToUnity();
335
336   if (fTrackCorr)
337     fTrackCorr->SetCorrectionToUnity();
338 }
339
340 //____________________________________________________________________
341 void AliCorrection::Multiply()
342 {
343   // call Multiply
344
345   if (fEventCorr)
346   {
347     fEventCorr->Multiply();
348     // now we manually copy the overflow bin of the y axis (multiplicity) over. This is important to get the event count correct
349     TH2F* hist = fEventCorr->GetMeasuredHistogram();
350     for (Int_t x = 1; x <= hist->GetNbinsX(); ++x)
351       fEventCorr->GetGeneratedHistogram()->SetBinContent(x, hist->GetNbinsY() + 1, hist->GetBinContent(x, hist->GetNbinsY() + 1));
352   }
353
354   if (fTrackCorr)
355     fTrackCorr->Multiply();
356 }
357
358 //____________________________________________________________________
359 void AliCorrection::Scale(Double_t factor)
360 {
361   // scales the two contained corrections
362
363   fEventCorr->Scale(factor);
364   fTrackCorr->Scale(factor);
365 }
366
367 //____________________________________________________________________
368 void AliCorrection::PrintInfo(Float_t ptCut)
369 {
370   // prints some stats
371
372   TH3* measured = GetTrackCorrection()->GetMeasuredHistogram();
373   TH3* generated = GetTrackCorrection()->GetGeneratedHistogram();
374
375   TH2* measuredEvents = GetEventCorrection()->GetMeasuredHistogram();
376   TH2* generatedEvents = GetEventCorrection()->GetGeneratedHistogram();
377
378   Printf("AliCorrection::PrintInfo: Whole phasespace:");
379
380   Printf("tracks measured: %.1f tracks generated: %.1f, events measured: %.1f, events generated %.1f", measured->Integral(), generated->Integral(), measuredEvents->Integral(), generatedEvents->Integral());
381
382   Printf("AliCorrection::PrintInfo: Values in |eta| < 0.8, |vtx-z| < 10 and pt > %.2f:", ptCut);
383
384   Float_t tracksM = measured->Integral(measured->GetXaxis()->FindBin(-9.9), measured->GetXaxis()->FindBin(9.9), measured->GetYaxis()->FindBin(-0.79), measured->GetYaxis()->FindBin(0.79), measured->GetZaxis()->FindBin(ptCut), measured->GetZaxis()->GetNbins());
385   Float_t tracksG = generated->Integral(generated->GetXaxis()->FindBin(-9.9), generated->GetXaxis()->FindBin(9.9), generated->GetYaxis()->FindBin(-0.79), generated->GetYaxis()->FindBin(0.79), generated->GetZaxis()->FindBin(ptCut), generated->GetZaxis()->GetNbins());
386
387   Float_t eventsM = measuredEvents->Integral(measuredEvents->GetXaxis()->FindBin(-9.9), measuredEvents->GetXaxis()->FindBin(9.9), 1, measuredEvents->GetNbinsY());
388   Float_t eventsG = generatedEvents->Integral(generatedEvents->GetXaxis()->FindBin(-9.9), generatedEvents->GetXaxis()->FindBin(9.9), 1, generatedEvents->GetNbinsY());
389
390   Printf("tracks measured: %.1f tracks generated: %.1f, events measured: %.1f, events generated %.1f", tracksM, tracksG, eventsM, eventsG);
391
392   if (tracksM > 0)
393     Printf("Effective average correction factor for TRACKS: %.3f", tracksG / tracksM);
394   if (eventsM > 0)
395   Printf("Effective average correction factor for EVENTS: %.3f", eventsG / eventsM);
396
397   if (eventsM > 0 && eventsG > 0)
398   {
399     // normalize to number of events;
400     tracksM /= eventsM;
401     tracksG /= eventsG;
402
403     Printf("%.2f tracks/event measured, %.2f tracks/event after correction --> effective average correction factor is %.3f (pt cut %.2f GeV/c)", tracksM, tracksG, tracksG / tracksM, ptCut);
404   }
405 }