]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/AliCorrection.cxx
-- new class to store track tree
[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 #include <TMath.h>
15
16 #include <AliLog.h>
17 #include "AliCorrectionMatrix2D.h"
18 #include "AliCorrectionMatrix3D.h"
19
20 #include "AliCorrection.h"
21
22 //____________________________________________________________________
23 ClassImp(AliCorrection)
24
25 //____________________________________________________________________
26 AliCorrection::AliCorrection() : TNamed(),
27   fEventCorr(0),
28   fTrackCorr(0)
29 {
30   // default constructor
31 }
32
33 //____________________________________________________________________
34 AliCorrection::AliCorrection(const Char_t* name, const Char_t* title, AliPWG0Helper::AnalysisMode analysisMode) : TNamed(name, title),
35   fEventCorr(0),
36   fTrackCorr(0)
37 {
38   // constructor initializing tnamed
39
40   Float_t* binLimitsPt = 0;
41   Int_t nBinsPt = 0;
42
43   // different binnings, better solution could be anticipated...
44   if ((analysisMode & AliPWG0Helper::kTPC || analysisMode & AliPWG0Helper::kTPCITS) && (analysisMode & AliPWG0Helper::kFieldOn))
45   {
46     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};
47     //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.1, 1.2, 1.3, 1.4, 1.5, 2.0, 5.0, 10.0, 100.0};
48     binLimitsPt = (Float_t*) binLimitsPtTmp;
49     nBinsPt = 28;
50   }
51   else if (analysisMode & AliPWG0Helper::kSPD || !(analysisMode & AliPWG0Helper::kFieldOn))
52   {
53     static Float_t binLimitsPtTmp[] = {-0.5, 100.0};
54     binLimitsPt = (Float_t*) binLimitsPtTmp;
55     nBinsPt = 1;
56   }
57
58   if (!binLimitsPt)
59   {
60     Printf("FATAL: Invalid binning");
61     return;
62   }
63
64   // third axis for track histogram
65   Int_t nBinsN2 = 1;
66   Float_t binLimitsN2[]   = {-0.5, 1000};
67   const char* title3 = "Ntracks";
68   
69   //Int_t nBinsN2 = 21;
70   //Float_t binLimitsN2[]   = {-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, 1000.5};
71   
72   // phi
73   //Int_t nBinsN2 = 36;
74   //Float_t binLimitsN2[]   = {0.000, 0.175, 0.349, 0.524, 0.698, 0.873, 1.047, 1.222, 1.396, 1.571, 1.745, 1.920, 2.094, 2.269, 2.443, 2.618, 2.793, 2.967, 3.142, 3.316, 3.491, 3.665, 3.840, 4.014, 4.189, 4.363, 4.538, 4.712, 4.887, 5.061, 5.236, 5.411, 5.585, 5.760, 5.934, 6.109, TMath::Pi() * 2};
75   //const char* title3 = "#phi";
76
77   // mult axis for event histogram
78   Int_t nBinsN = 22;
79   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};
80   //Int_t nBinsN = 8;
81   //Float_t binLimitsN[]   = {-0.5, 0.5, 1.5, 2.5, 3.5, 4.5, 5.5, 6.5, 300.5};
82   //Int_t nBinsN = 52;
83   //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};
84
85   //Float_t binLimitsVtx[] = {-20,-15,-10,-6,-3,0,3,6,10,15,20};
86   //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};
87   //Float_t binLimitsVtx[] = {-30,-25,-20,-15,-10,-8,-6,-4,-2,0,2,4,6,8,10,15,20,25,30};
88   Float_t binLimitsVtx[] = {-30,-25,-20,-15,-10,-7,-5.5,-4,-3,-2,-1,0,1,2,3,4,5.5,7,10,15,20,25,30};
89   /*Float_t binLimitsEta[] = {-3.0,-2.8,-2.6,-2.4,-2.2,
90                             -2.0,-1.8,-1.6,-1.4,-1.2,
91                             -1.0,-0.8,-0.6,-0.4,-0.2, 0.0,
92                              0.2, 0.4, 0.6, 0.8, 1.0,
93                              1.2, 1.4, 1.6, 1.8, 2.0,
94                              2.2, 2.4, 2.6, 2.8, 3.0};*/
95   //Float_t binLimitsEta[] = { -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 };
96   /*Float_t binLimitsEta[] = {-3,-2.75,-2.5,-2.25,-2.0,
97                             -1.75,-1.5,-1.25,
98                             -1,-0.75,-0.5,-0.25,
99                             0,0.25,0.5,0.75,
100                             1.0,1.25,1.5,1.75,
101                             2.0,2.25,2.5,2.75,3.0}; */
102   //Float_t binLimitsEta[] = {-2.8,-2.4,-2,-1.6,-1.2,-0.8,-0.4,0,0.4,0.8,1.2,1.6,2.0,2.4,2.8}; 
103   //Float_t binLimitsEta[] = {-2.8,-2.4,-2,-1.6,-1.2,-0.8,-0.5,0,0.5,0.8,1.2,1.6,2.0,2.4,2.8}; 
104
105   //Float_t binLimitsEta[] = {-3,-2.6,-2.2,-1.8,-1.4,-1,-0.6,-0.2,0.2,0.6,1,1.4,1.8,2.2,2.6,3.0};
106   //1. standard
107   //  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};
108   // 2. MB Working Group definition
109   Float_t binLimitsEta[] = {-2, -1.8, -1.6, -1.4, -1.2, -1., -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};
110 //  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 };
111
112   fEventCorr = new AliCorrectionMatrix2D("EventCorrection", Form("%s EventCorrection", fTitle.Data()), 22, binLimitsVtx, nBinsN, binLimitsN);
113
114   if (analysisMode & AliPWG0Helper::kSPD)
115   {
116     // 1. standard
117     //TH3F* dummyBinning = new TH3F("dummyBinning","dummyBinning",22, binLimitsVtx, 40, binLimitsEta, nBinsN2, binLimitsN2);
118     // 2. MB Working Group definition
119     TH3F* dummyBinning = new TH3F("dummyBinning","dummyBinning",22, binLimitsVtx, 20, binLimitsEta, nBinsN2, binLimitsN2);
120     fTrackCorr = new AliCorrectionMatrix3D("TrackCorrection", Form("%s TrackCorrection", fTitle.Data()), dummyBinning);
121     fTrackCorr->SetAxisTitles("vtx-z (cm)", "#eta", title3);
122     delete dummyBinning;
123   }
124   else
125   {
126     // 1. standard
127     //TH3F* dummyBinning = new TH3F("dummyBinning","dummyBinning",22, binLimitsVtx, 40, binLimitsEta , nBinsPt, binLimitsPt);
128     // 2. MB Working Group definition
129     TH3F* dummyBinning = new TH3F("dummyBinning","dummyBinning",22, binLimitsVtx, 20, binLimitsEta , nBinsPt, binLimitsPt);
130     fTrackCorr = new AliCorrectionMatrix3D("TrackCorrection", Form("%s TrackCorrection", fTitle.Data()), dummyBinning);
131     fTrackCorr->SetAxisTitles("vtx-z (cm)", "#eta", "p_{T} (GeV/c)");
132     delete dummyBinning;
133   }
134
135
136   fEventCorr->SetAxisTitles("vtx-z (cm)", "Ntracks");
137 }
138
139 //____________________________________________________________________
140 AliCorrection::AliCorrection(const AliCorrection& c) : TNamed(c),
141   fEventCorr(0),
142   fTrackCorr(0)
143 {
144   // copy constructor
145   ((AliCorrection &)c).Copy(*this);
146 }
147
148 //____________________________________________________________________
149 AliCorrection::~AliCorrection()
150 {
151   //
152   // destructor
153   //
154
155   if (fEventCorr)
156   {
157     delete fEventCorr;
158     fEventCorr = 0;
159   }
160
161   if (fTrackCorr)
162   {
163     delete fTrackCorr;
164     fTrackCorr = 0;
165   }
166 }
167
168 //____________________________________________________________________
169 AliCorrection &AliCorrection::operator=(const AliCorrection &c)
170 {
171   // assigment operator
172
173   if (this != &c)
174     ((AliCorrection &) c).Copy(*this);
175
176   return *this;
177 }
178
179 //____________________________________________________________________
180 void AliCorrection::Copy(TObject& c) const
181 {
182   // copy function
183
184   AliCorrection& target = (AliCorrection &) c;
185
186   if (fEventCorr)
187     target.fEventCorr = dynamic_cast<AliCorrectionMatrix2D*> (fEventCorr->Clone());
188
189   if (fTrackCorr)
190     target.fTrackCorr = dynamic_cast<AliCorrectionMatrix3D*> (fTrackCorr->Clone());
191 }
192
193 //____________________________________________________________________
194 Long64_t AliCorrection::Merge(TCollection* list)
195 {
196   // Merge a list of AliCorrection objects with this (needed for
197   // PROOF). 
198   // Returns the number of merged objects (including this).
199
200   if (!list)
201     return 0;
202   
203   if (list->IsEmpty())
204     return 1;
205
206   TIterator* iter = list->MakeIterator();
207   TObject* obj;
208
209   // collections of measured and generated histograms
210   TList* collectionEvent = new TList;
211   TList* collectionTrack = new TList;
212   
213   Int_t count = 0;
214   while ((obj = iter->Next())) {
215     
216     AliCorrection* entry = dynamic_cast<AliCorrection*> (obj);
217     if (entry == 0) 
218       continue;
219
220     collectionEvent->Add(entry->fEventCorr);
221     collectionTrack->Add(entry->fTrackCorr);
222
223     count++;
224   }
225   fEventCorr->Merge(collectionEvent);
226   fTrackCorr->Merge(collectionTrack);
227
228   delete collectionEvent;
229   delete collectionTrack;
230
231   return count+1;
232 }
233
234 //____________________________________________________________________
235 void AliCorrection::Divide()
236 {
237   //
238   // divide the histograms to get the correction
239   //
240   
241   if (!fEventCorr || !fTrackCorr)
242     return;
243     
244   fEventCorr->Divide();
245   fTrackCorr->Divide();
246
247   Int_t emptyBins = fTrackCorr->CheckEmptyBins(-9.99, 9.99, -0.79, 0.79, 0.2, 4.9);
248   printf("INFO: In the central region the track correction of %s has %d empty bins\n", GetName(), emptyBins);
249 }
250
251 //____________________________________________________________________
252 void AliCorrection::Add(AliCorrection* aCorrectionToAdd, Float_t c)
253 {
254   //
255   // add to measured and generated the measured and generated of aCorrectionToAdd
256   // with the weight c
257
258   fEventCorr->Add(aCorrectionToAdd->GetEventCorrection(),c);
259   fTrackCorr->Add(aCorrectionToAdd->GetTrackCorrection(),c);
260 }
261
262
263 //____________________________________________________________________
264 Bool_t AliCorrection::LoadHistograms(const Char_t* dir)
265 {
266   //
267   // loads the histograms from a file
268   // if dir is empty a directory with the name of this object is taken (like in SaveHistogram)
269   //
270
271   if (!fEventCorr || !fTrackCorr)
272     return kFALSE;
273
274   if (!dir)
275     dir = GetName();
276
277   if (!gDirectory->cd(dir))
278     return kFALSE;
279
280   Bool_t success = fEventCorr->LoadHistograms();
281   success &= fTrackCorr->LoadHistograms();
282
283   gDirectory->cd("..");
284
285   return success;
286 }
287
288 //____________________________________________________________________
289 void AliCorrection::SaveHistograms()
290 {
291   //
292   // saves the histograms in a directory with the name of this object (GetName)
293   //
294   
295   gDirectory->mkdir(GetName());
296   gDirectory->cd(GetName());
297
298   if (fEventCorr)
299     fEventCorr->SaveHistograms();
300
301   if (fTrackCorr)
302     fTrackCorr->SaveHistograms();
303     
304   gDirectory->cd("..");
305 }
306
307 //____________________________________________________________________
308 void AliCorrection::ReduceInformation()
309 {
310   // this function deletes the measured and generated histograms to reduce the amount of data
311   // in memory
312
313   if (!fEventCorr || !fTrackCorr)
314     return;
315
316   fEventCorr->ReduceInformation();
317   fTrackCorr->ReduceInformation();
318 }
319
320 //____________________________________________________________________
321 void AliCorrection::Reset(Option_t* option)
322 {
323   // resets the histograms
324
325   if (fEventCorr)
326     fEventCorr->Reset(option);
327
328   if (fTrackCorr)
329     fTrackCorr->Reset(option);
330 }
331
332 //____________________________________________________________________
333 void AliCorrection::DrawHistograms(const Char_t* name)
334 {
335   // draws the corrections
336
337   if (!name)
338     name = GetName();
339
340   if (fEventCorr)
341     fEventCorr->DrawHistograms(Form("%s event", name));
342
343   if (fTrackCorr)
344     fTrackCorr->DrawHistograms(Form("%s track", name));
345 }
346
347 void AliCorrection::DrawOverview(const char* canvasName)
348 {
349   // draw projection of the corrections
350   //   to the 3 axis of fTrackCorr and 2 axis of fEventCorr
351
352   TString canvasNameTmp(GetName());
353   if (canvasName)
354     canvasNameTmp = canvasName;
355
356   TCanvas* canvas = new TCanvas(canvasNameTmp, canvasNameTmp, 1200, 1200);
357   canvas->Divide(3, 3);
358
359   if (fTrackCorr) {
360     canvas->cd(1);
361     fTrackCorr->Get2DCorrectionHistogram("xy", 0, 0)->DrawCopy("COLZ")->GetZaxis()->SetRangeUser(0, 10);
362     
363     canvas->cd(2);
364     fTrackCorr->Get2DCorrectionHistogram("yz", 0, 0)->DrawCopy("COLZ")->GetZaxis()->SetRangeUser(0, 10);
365     
366     canvas->cd(3);
367     fTrackCorr->Get2DCorrectionHistogram("zx", 0, 0)->DrawCopy("COLZ")->GetZaxis()->SetRangeUser(0, 10);
368     
369     canvas->cd(4);
370     fTrackCorr->Get1DCorrectionHistogram("x", 0.3, 5, -1, 1)->DrawCopy()->GetYaxis()->SetRangeUser(0, 10);
371
372     canvas->cd(5);
373     fTrackCorr->Get1DCorrectionHistogram("y", 0.3, 5, 0, 0)->DrawCopy()->GetYaxis()->SetRangeUser(0, 10);
374
375     canvas->cd(6);
376     fTrackCorr->Get1DCorrectionHistogram("z", 0, -1, -1, 1)->DrawCopy()->GetYaxis()->SetRangeUser(0, 10);
377   }
378
379   if (fEventCorr)
380   {
381     canvas->cd(7);
382     fEventCorr->GetCorrectionHistogram()->DrawCopy("COLZ")->GetYaxis()->SetRangeUser(0, 30);
383     
384     canvas->cd(8);
385     fEventCorr->Get1DCorrectionHistogram("x")->DrawCopy();
386
387     canvas->cd(9);
388     fEventCorr->Get1DCorrectionHistogram("y")->DrawCopy()->GetXaxis()->SetRangeUser(0, 30);
389   }
390 }
391
392 //____________________________________________________________________
393 void AliCorrection::SetCorrectionToUnity()
394 {
395   // set the corrections to unity
396
397   if (fEventCorr)
398     fEventCorr->SetCorrectionToUnity();
399
400   if (fTrackCorr)
401     fTrackCorr->SetCorrectionToUnity();
402 }
403
404 //____________________________________________________________________
405 void AliCorrection::Multiply()
406 {
407   // call Multiply
408
409   if (fEventCorr)
410   {
411     fEventCorr->Multiply();
412     // now we manually copy the overflow bin of the y axis (multiplicity) over. This is important to get the event count correct
413     TH2* hist = fEventCorr->GetMeasuredHistogram();
414     for (Int_t x = 1; x <= hist->GetNbinsX(); ++x)
415       fEventCorr->GetGeneratedHistogram()->SetBinContent(x, hist->GetNbinsY() + 1, hist->GetBinContent(x, hist->GetNbinsY() + 1));
416   }
417
418   if (fTrackCorr)
419     fTrackCorr->Multiply();
420 }
421
422 //____________________________________________________________________
423 void AliCorrection::Scale(Double_t factor)
424 {
425   // scales the two contained corrections
426
427   fEventCorr->Scale(factor);
428   fTrackCorr->Scale(factor);
429 }
430
431 //____________________________________________________________________
432 void AliCorrection::PrintStats(Float_t zRange, Float_t etaRange, Float_t ptCut)
433 {
434   // prints statistics and effective correction factors
435
436   Printf("AliCorrection::PrintInfo: Values in |eta| < %.2f, |vtx-z| < %.2f and pt > %.2f:", etaRange, zRange, ptCut);
437
438   // prevent to be on bin border
439   zRange -= 0.1;
440   etaRange -= 0.1;
441
442   TH3* measured = GetTrackCorrection()->GetMeasuredHistogram();
443   TH3* generated = GetTrackCorrection()->GetGeneratedHistogram();
444
445   TH2* measuredEvents = GetEventCorrection()->GetMeasuredHistogram();
446   TH2* generatedEvents = GetEventCorrection()->GetGeneratedHistogram();
447
448   Float_t tracksM = measured->Integral(measured->GetXaxis()->FindBin(-zRange), measured->GetXaxis()->FindBin(zRange), measured->GetYaxis()->FindBin(-etaRange), measured->GetYaxis()->FindBin(etaRange), measured->GetZaxis()->FindBin(ptCut), measured->GetZaxis()->GetNbins());
449   Float_t tracksG = generated->Integral(generated->GetXaxis()->FindBin(-zRange), generated->GetXaxis()->FindBin(zRange), generated->GetYaxis()->FindBin(-etaRange), generated->GetYaxis()->FindBin(etaRange), generated->GetZaxis()->FindBin(ptCut), generated->GetZaxis()->GetNbins());
450
451   Float_t eventsM = measuredEvents->Integral(measuredEvents->GetXaxis()->FindBin(-zRange), measuredEvents->GetXaxis()->FindBin(zRange), 1, measuredEvents->GetNbinsY());
452   Float_t eventsG = generatedEvents->Integral(generatedEvents->GetXaxis()->FindBin(-zRange), generatedEvents->GetXaxis()->FindBin(zRange), 1, generatedEvents->GetNbinsY());
453
454   Float_t eventsM1 = measuredEvents->Integral(measuredEvents->GetXaxis()->FindBin(-zRange), measuredEvents->GetXaxis()->FindBin(zRange), 2, measuredEvents->GetNbinsY());
455   Float_t eventsG1 = generatedEvents->Integral(generatedEvents->GetXaxis()->FindBin(-zRange), generatedEvents->GetXaxis()->FindBin(zRange), 2, generatedEvents->GetNbinsY());
456   
457   Printf("tracks measured: %.1f tracks generated: %.1f, events measured: %.1f, events generated %.1f, events (N>0) measured: %.1f, events (N>0) generated %.1f", tracksM, tracksG, eventsM, eventsG, eventsM1, eventsG1);
458
459   if (tracksM > 0)
460     Printf("Effective average correction factor for TRACKS: %.3f", tracksG / tracksM);
461   if (eventsM > 0)
462     Printf("Effective average correction factor for EVENTS: %.3f", eventsG / eventsM);
463
464   if (eventsM > 0 && eventsG > 0)
465   {
466     // normalize to number of events;
467     tracksM /= eventsM;
468     tracksG /= eventsG;
469
470     Printf("%.2f tracks/event measured, %.2f tracks/event after correction --> effective average correction factor is %.3f (pt cut %.2f GeV/c)", tracksM, tracksG, (tracksM > 0) ? (tracksG / tracksM) : -1, ptCut);
471   }
472 }
473
474 //____________________________________________________________________
475 void AliCorrection::PrintInfo(Float_t ptCut)
476 {
477   // prints some stats
478
479   TH3* measured = GetTrackCorrection()->GetMeasuredHistogram();
480   TH3* generated = GetTrackCorrection()->GetGeneratedHistogram();
481
482   TH2* measuredEvents = GetEventCorrection()->GetMeasuredHistogram();
483   TH2* generatedEvents = GetEventCorrection()->GetGeneratedHistogram();
484
485   Printf("AliCorrection::PrintInfo: Whole phasespace:");
486
487   Printf("tracks measured: %.1f tracks generated: %.1f, events measured: %.1f, events generated %.1f", measured->Integral(), generated->Integral(), measuredEvents->Integral(), generatedEvents->Integral());
488
489   Printf("Example centered bin: tracks measured: %.1f tracks generated: %.1f, events measured: %.1f, events generated %.1f", measured->GetBinContent(measured->GetNbinsX() / 2, measured->GetNbinsY() / 2, measured->GetNbinsZ() / 2), generated->GetBinContent(measured->GetNbinsX() / 2, measured->GetNbinsY() / 2, measured->GetNbinsZ() / 2), measuredEvents->GetBinContent(measuredEvents->GetNbinsX() / 2, measuredEvents->GetNbinsY() / 2), generatedEvents->GetBinContent(measuredEvents->GetNbinsX() / 2, measuredEvents->GetNbinsY() / 2));
490
491   PrintStats(10, 0.6, ptCut);
492   PrintStats(10, 0.8, ptCut);
493   PrintStats(10, 1.5, ptCut);
494 }
495
496 //____________________________________________________________________
497 void AliCorrection::ResetErrorsOnCorrections()
498 {
499   // resets the errors in the correction matrix
500   
501   fEventCorr->ResetErrorsOnCorrections();
502   fTrackCorr->ResetErrorsOnCorrections();  
503 }