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