]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/AliCorrection.cxx
3b8d455b9e359702a2e0210c319bfff936d9005d
[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   Int_t nBinsN = 22;
63   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};
64   //Int_t nBinsN = 52;
65   //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, 11.5, 12.5, 13.5, 14.5, 15.5, 16.5, 17.5, 18.5, 19.5, 20.5, 21.5, 22.5, 23.5, 24.5, 25.5, 26.5, 27.5, 28.5, 29.5, 30.5, 31.5, 32.5, 33.5, 34.5, 35.5, 36.5, 37.5, 38.5, 39.5, 40.5, 41.5, 42.5, 43.5, 44.5, 45.5, 46.5, 47.5, 48.5, 49.5, 50.5, 100.5, 300.5};
66
67   //Float_t binLimitsVtx[] = {-20,-15,-10,-6,-3,0,3,6,10,15,20};
68   //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};
69   Float_t binLimitsVtx[] = {-30,-25,-20,-15,-10,-8,-6,-4,-2,0,2,4,6,8,10,15,20,25,30};
70   Float_t binLimitsEta[] = {-3.0,-2.8,-2.6,-2.4,-2.2,
71                             -2.0,-1.8,-1.6,-1.4,-1.2,
72                             -1.0,-0.8,-0.6,-0.4,-0.2, 0.0,
73                              0.2, 0.4, 0.6, 0.8, 1.0,
74                              1.2, 1.4, 1.6, 1.8, 2.0,
75                              2.2, 2.4, 2.6, 2.8, 3.0};
76 //  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};
77 //  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};
78 //  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 };
79
80   TH3F* dummyBinning = new TH3F("dummyBinning","dummyBinning",18, binLimitsVtx, 30, binLimitsEta , nBinsPt, binLimitsPt);
81
82   fEventCorr = new AliCorrectionMatrix2D("EventCorrection", Form("%s EventCorrection", fTitle.Data()), 18, binLimitsVtx, nBinsN, binLimitsN);
83   fTrackCorr = new AliCorrectionMatrix3D("TrackCorrection", Form("%s TrackCorrection", fTitle.Data()), dummyBinning);
84
85   delete dummyBinning;
86
87   fEventCorr->SetAxisTitles("vtx z [cm]", "Ntracks");
88   fTrackCorr->SetAxisTitles("vtx z [cm]", "#eta", "p_{T} [GeV/c]");
89 }
90
91 //____________________________________________________________________
92 AliCorrection::AliCorrection(const AliCorrection& c) : TNamed(c),
93   fEventCorr(0),
94   fTrackCorr(0)
95 {
96   // copy constructor
97   ((AliCorrection &)c).Copy(*this);
98 }
99
100 //____________________________________________________________________
101 AliCorrection::~AliCorrection()
102 {
103   //
104   // destructor
105   //
106
107   if (fEventCorr)
108   {
109     delete fEventCorr;
110     fEventCorr = 0;
111   }
112
113   if (fTrackCorr)
114   {
115     delete fTrackCorr;
116     fTrackCorr = 0;
117   }
118 }
119
120 //____________________________________________________________________
121 AliCorrection &AliCorrection::operator=(const AliCorrection &c)
122 {
123   // assigment operator
124
125   if (this != &c)
126     ((AliCorrection &) c).Copy(*this);
127
128   return *this;
129 }
130
131 //____________________________________________________________________
132 void AliCorrection::Copy(TObject& c) const
133 {
134   // copy function
135
136   AliCorrection& target = (AliCorrection &) c;
137
138   if (fEventCorr)
139     target.fEventCorr = dynamic_cast<AliCorrectionMatrix2D*> (fEventCorr->Clone());
140
141   if (fTrackCorr)
142     target.fTrackCorr = dynamic_cast<AliCorrectionMatrix3D*> (fTrackCorr->Clone());
143 }
144
145 //____________________________________________________________________
146 Long64_t AliCorrection::Merge(TCollection* list)
147 {
148   // Merge a list of AliCorrection objects with this (needed for
149   // PROOF). 
150   // Returns the number of merged objects (including this).
151
152   if (!list)
153     return 0;
154   
155   if (list->IsEmpty())
156     return 1;
157
158   TIterator* iter = list->MakeIterator();
159   TObject* obj;
160
161   // collections of measured and generated histograms
162   TList* collectionEvent = new TList;
163   TList* collectionTrack = new TList;
164   
165   Int_t count = 0;
166   while ((obj = iter->Next())) {
167     
168     AliCorrection* entry = dynamic_cast<AliCorrection*> (obj);
169     if (entry == 0) 
170       continue;
171
172     collectionEvent->Add(entry->fEventCorr);
173     collectionTrack->Add(entry->fTrackCorr);
174
175     count++;
176   }
177   fEventCorr->Merge(collectionEvent);
178   fTrackCorr->Merge(collectionTrack);
179
180   delete collectionEvent;
181   delete collectionTrack;
182
183   return count+1;
184 }
185
186 //____________________________________________________________________
187 void AliCorrection::Divide()
188 {
189   //
190   // divide the histograms to get the correction
191   //
192   
193   if (!fEventCorr || !fTrackCorr)
194     return;
195     
196   fEventCorr->Divide();
197   fTrackCorr->Divide();
198
199   Int_t emptyBins = fTrackCorr->CheckEmptyBins(-9.99, 9.99, -0.79, 0.79, 0.3, 4.9);
200   printf("INFO: In the central region the track correction of %s has %d empty bins\n", GetName(), emptyBins);
201 }
202
203 //____________________________________________________________________
204 void AliCorrection::Add(AliCorrection* aCorrectionToAdd, Float_t c)
205 {
206   //
207   // add to measured and generated the measured and generated of aCorrectionToAdd
208   // with the weight c
209
210   fEventCorr->Add(aCorrectionToAdd->GetEventCorrection(),c);
211   fTrackCorr->Add(aCorrectionToAdd->GetTrackCorrection(),c);
212 }
213
214
215 //____________________________________________________________________
216 Bool_t AliCorrection::LoadHistograms(const Char_t* dir)
217 {
218   //
219   // loads the histograms from a file
220   // if dir is empty a directory with the name of this object is taken (like in SaveHistogram)
221   //
222
223   if (!fEventCorr || !fTrackCorr)
224     return kFALSE;
225
226   if (!dir)
227     dir = GetName();
228
229   if (!gDirectory->cd(dir))
230     return kFALSE;
231
232   Bool_t success = fEventCorr->LoadHistograms();
233   success &= fTrackCorr->LoadHistograms();
234
235   gDirectory->cd("..");
236
237   return success;
238 }
239
240 //____________________________________________________________________
241 void AliCorrection::SaveHistograms()
242 {
243   //
244   // saves the histograms in a directory with the name of this object (GetName)
245   //
246   
247   gDirectory->mkdir(GetName());
248   gDirectory->cd(GetName());
249
250   if (fEventCorr)
251     fEventCorr->SaveHistograms();
252
253   if (fTrackCorr)
254     fTrackCorr->SaveHistograms();
255     
256   gDirectory->cd("..");
257 }
258
259 //____________________________________________________________________
260 void AliCorrection::ReduceInformation()
261 {
262   // this function deletes the measured and generated histograms to reduce the amount of data
263   // in memory
264
265   if (!fEventCorr || !fTrackCorr)
266     return;
267
268   fEventCorr->ReduceInformation();
269   fTrackCorr->ReduceInformation();
270 }
271
272 //____________________________________________________________________
273 void AliCorrection::Reset(Option_t* option)
274 {
275   // resets the histograms
276
277   if (fEventCorr)
278     fEventCorr->Reset(option);
279
280   if (fTrackCorr)
281     fTrackCorr->Reset(option);
282 }
283
284 //____________________________________________________________________
285 void AliCorrection::DrawHistograms(const Char_t* name)
286 {
287   // draws the corrections
288
289   if (!name)
290     name = GetName();
291
292   if (fEventCorr)
293     fEventCorr->DrawHistograms(Form("%s event", name));
294
295   if (fTrackCorr)
296     fTrackCorr->DrawHistograms(Form("%s track", name));
297 }
298
299 void AliCorrection::DrawOverview(const char* canvasName)
300 {
301   // draw projection of the corrections
302   //   to the 3 axis of fTrackCorr and 2 axis of fEventCorr
303
304   TString canvasNameTmp(GetName());
305   if (canvasName)
306     canvasNameTmp = canvasName;
307
308   TCanvas* canvas = new TCanvas(canvasNameTmp, canvasNameTmp, 1200, 800);
309   canvas->Divide(3, 2);
310
311   if (fTrackCorr) {
312     canvas->cd(1);
313     fTrackCorr->Get1DCorrectionHistogram("x", 0.3, 5, -1, 1)->DrawCopy()->GetYaxis()->SetRangeUser(0, 10);
314
315     canvas->cd(2);
316     fTrackCorr->Get1DCorrectionHistogram("y", 0.3, 5, 0, 0)->DrawCopy()->GetYaxis()->SetRangeUser(0, 10);
317
318     canvas->cd(3);
319     fTrackCorr->Get1DCorrectionHistogram("z", 0, -1, -1, 1)->DrawCopy()->GetYaxis()->SetRangeUser(0, 10);
320   }
321
322   if (fEventCorr)
323   {
324     canvas->cd(4);
325     fEventCorr->Get1DCorrectionHistogram("x")->DrawCopy();
326
327     canvas->cd(5);
328     fEventCorr->Get1DCorrectionHistogram("y")->DrawCopy()->GetXaxis()->SetRangeUser(0, 30);
329   }
330 }
331
332 //____________________________________________________________________
333 void AliCorrection::SetCorrectionToUnity()
334 {
335   // set the corrections to unity
336
337   if (fEventCorr)
338     fEventCorr->SetCorrectionToUnity();
339
340   if (fTrackCorr)
341     fTrackCorr->SetCorrectionToUnity();
342 }
343
344 //____________________________________________________________________
345 void AliCorrection::Multiply()
346 {
347   // call Multiply
348
349   if (fEventCorr)
350   {
351     fEventCorr->Multiply();
352     // now we manually copy the overflow bin of the y axis (multiplicity) over. This is important to get the event count correct
353     TH2F* hist = fEventCorr->GetMeasuredHistogram();
354     for (Int_t x = 1; x <= hist->GetNbinsX(); ++x)
355       fEventCorr->GetGeneratedHistogram()->SetBinContent(x, hist->GetNbinsY() + 1, hist->GetBinContent(x, hist->GetNbinsY() + 1));
356   }
357
358   if (fTrackCorr)
359     fTrackCorr->Multiply();
360 }
361
362 //____________________________________________________________________
363 void AliCorrection::Scale(Double_t factor)
364 {
365   // scales the two contained corrections
366
367   fEventCorr->Scale(factor);
368   fTrackCorr->Scale(factor);
369 }
370
371 //____________________________________________________________________
372 void AliCorrection::PrintInfo(Float_t ptCut)
373 {
374   // prints some stats
375
376   TH3* measured = GetTrackCorrection()->GetMeasuredHistogram();
377   TH3* generated = GetTrackCorrection()->GetGeneratedHistogram();
378
379   TH2* measuredEvents = GetEventCorrection()->GetMeasuredHistogram();
380   TH2* generatedEvents = GetEventCorrection()->GetGeneratedHistogram();
381
382   Printf("AliCorrection::PrintInfo: Whole phasespace:");
383
384   Printf("tracks measured: %.1f tracks generated: %.1f, events measured: %.1f, events generated %.1f", measured->Integral(), generated->Integral(), measuredEvents->Integral(), generatedEvents->Integral());
385
386   Printf("AliCorrection::PrintInfo: Values in |eta| < 0.8, |vtx-z| < 10 and pt > %.2f:", ptCut);
387
388   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());
389   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());
390
391   Float_t eventsM = measuredEvents->Integral(measuredEvents->GetXaxis()->FindBin(-9.9), measuredEvents->GetXaxis()->FindBin(9.9), 1, measuredEvents->GetNbinsY());
392   Float_t eventsG = generatedEvents->Integral(generatedEvents->GetXaxis()->FindBin(-9.9), generatedEvents->GetXaxis()->FindBin(9.9), 1, generatedEvents->GetNbinsY());
393
394   Printf("tracks measured: %.1f tracks generated: %.1f, events measured: %.1f, events generated %.1f", tracksM, tracksG, eventsM, eventsG);
395
396   if (tracksM > 0)
397     Printf("Effective average correction factor for TRACKS: %.3f", tracksG / tracksM);
398   if (eventsM > 0)
399   Printf("Effective average correction factor for EVENTS: %.3f", eventsG / eventsM);
400
401   if (eventsM > 0 && eventsG > 0)
402   {
403     // normalize to number of events;
404     tracksM /= eventsM;
405     tracksG /= eventsG;
406
407     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);
408   }
409 }