]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/AliCorrectionMatrix.cxx
AliHMPIDDigitN no longer needed
[u/mrichter/AliRoot.git] / PWG0 / AliCorrectionMatrix.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 <TH2F.h>
13
14 #include <AliLog.h>
15
16 #include "AliCorrectionMatrix.h"
17
18 //____________________________________________________________________
19 ClassImp(AliCorrectionMatrix)
20
21 //____________________________________________________________________
22 AliCorrectionMatrix::AliCorrectionMatrix() : TNamed(),
23   fhMeas(0),
24   fhGene(0),
25   fhCorr(0)
26 {
27   // default constructor
28 }
29
30 //____________________________________________________________________
31 AliCorrectionMatrix::AliCorrectionMatrix(const Char_t* name, const Char_t* title) : TNamed(name, title),
32   fhMeas(0),
33   fhGene(0),
34   fhCorr(0)
35 {
36   // constructor initializing tnamed
37 }
38
39 //____________________________________________________________________
40 AliCorrectionMatrix::AliCorrectionMatrix(const AliCorrectionMatrix& c) : TNamed(c),
41   fhMeas(0),
42   fhGene(0),
43   fhCorr(0)
44 {
45   // copy constructor
46   ((AliCorrectionMatrix &)c).Copy(*this);
47 }
48
49 //____________________________________________________________________
50 AliCorrectionMatrix::~AliCorrectionMatrix()
51 {
52   //
53   // destructor
54   //
55
56   if (fhMeas)
57   {
58     delete fhMeas;
59     fhMeas = 0;
60   }
61
62   if (fhGene)
63   {
64     delete fhGene;
65     fhGene = 0;
66   }
67
68   if (fhCorr)
69   {
70     delete fhCorr;
71     fhCorr = 0;
72   }
73 }
74
75 //____________________________________________________________________
76 AliCorrectionMatrix &AliCorrectionMatrix::operator=(const AliCorrectionMatrix &c)
77 {
78   // assigment operator
79
80   if (this != &c)
81     ((AliCorrectionMatrix &) c).Copy(*this);
82
83   return *this;
84 }
85
86 //____________________________________________________________________
87 void AliCorrectionMatrix::Copy(TObject& c) const
88 {
89   // copy function
90
91   AliCorrectionMatrix& target = (AliCorrectionMatrix &) c;
92
93   if (fhMeas)
94     target.fhMeas = dynamic_cast<TH1*> (fhMeas->Clone());
95
96   if (fhGene)
97     target.fhGene = dynamic_cast<TH1*> (fhGene->Clone());
98
99   if (fhCorr)
100     target.fhCorr = dynamic_cast<TH1*> (fhCorr->Clone());
101 }
102
103 //________________________________________________________________________
104 void AliCorrectionMatrix::SetAxisTitles(const Char_t* titleX, const Char_t* titleY, const Char_t* titleZ)
105 {
106   //
107   // method for setting the axis titles of the histograms
108   //
109
110   fhMeas ->SetXTitle(titleX);  fhMeas ->SetYTitle(titleY);  fhMeas ->SetZTitle(titleZ);
111   fhGene ->SetXTitle(titleX);  fhGene ->SetYTitle(titleY);  fhGene ->SetZTitle(titleZ);
112   fhCorr ->SetXTitle(titleX);  fhCorr ->SetYTitle(titleY);  fhCorr ->SetZTitle(titleZ);
113 }
114
115 //____________________________________________________________________
116 Long64_t AliCorrectionMatrix::Merge(TCollection* list)
117 {
118   // Merge a list of AliCorrectionMatrix objects with this (needed for
119   // PROOF). 
120   // Returns the number of merged objects (including this).
121
122   if (!list)
123     return 0;
124   
125   if (list->IsEmpty())
126     return 1;
127
128   TIterator* iter = list->MakeIterator();
129   TObject* obj;
130
131   // collections of measured and generated histograms
132   TList* collectionMeas = new TList;
133   TList* collectionGene = new TList;
134   
135   Int_t count = 0;
136   while ((obj = iter->Next())) {
137     
138     AliCorrectionMatrix* entry = dynamic_cast<AliCorrectionMatrix*> (obj);
139     if (entry == 0) 
140       continue;
141
142     collectionMeas->Add(entry->GetMeasuredHistogram());
143     collectionGene->Add(entry->GetGeneratedHistogram());
144
145     count++;
146   }
147   fhMeas->Merge(collectionMeas);
148   fhGene->Merge(collectionGene);
149
150   delete collectionMeas;
151   delete collectionGene;
152
153   return count+1;
154 }
155
156 //____________________________________________________________________
157 void AliCorrectionMatrix::Divide()
158 {
159   //
160   // divides generated by measured to get the correction
161   //
162
163   if (!fhMeas || !fhGene || !fhCorr)
164     return;
165
166   fhCorr->Divide(fhGene, fhMeas, 1, 1, "B");
167
168   Int_t emptyBins = 0;
169   for (Int_t x=1; x<=fhCorr->GetNbinsX(); ++x)
170     for (Int_t y=1; y<=fhCorr->GetNbinsY(); ++y)
171       for (Int_t z=1; z<=fhCorr->GetNbinsZ(); ++z)
172         if (fhCorr->GetBinContent(x, y, z) == 0)
173           ++emptyBins;
174
175   if (emptyBins > 0)
176     printf("INFO: In %s we have %d empty bins (of %d) in the correction map\n", GetTitle(), emptyBins, fhCorr->GetNbinsX() * fhCorr->GetNbinsY() * fhCorr->GetNbinsZ());
177 }
178
179 //____________________________________________________________________
180 void AliCorrectionMatrix::Multiply()
181 {
182   //
183   // multiplies measured with correction to get the generated
184   //
185
186   if (!fhMeas || !fhGene || !fhCorr)
187     return;
188
189   fhGene->Multiply(fhMeas, fhCorr, 1, 1, "B");
190 }
191
192 //____________________________________________________________________
193 Bool_t AliCorrectionMatrix::LoadHistograms(const Char_t* dir)
194 {
195   //
196   // loads the histograms from a file
197   // if dir is empty a directory with the name of this object is taken (like in SaveHistogram)
198   //
199
200   if (!dir)
201     dir = GetName();
202
203   if (!gDirectory->cd(dir))
204     return kFALSE;
205
206   if (fhGene)
207   {
208     delete fhGene;
209     fhGene=0;
210   }
211
212   if (fhCorr)
213   {
214     delete fhCorr;
215     fhCorr=0;
216   }
217
218   if (fhMeas)
219   {
220     delete fhMeas;
221     fhMeas=0;
222   }
223
224   fhMeas  = dynamic_cast<TH1*> (gDirectory->Get("measured"));
225   if (!fhMeas)
226     Info("LoadHistograms", "No measured hist available");
227
228   fhGene  = dynamic_cast<TH1*> (gDirectory->Get("generated"));
229   if (!fhGene)
230     Info("LoadHistograms", "No generated hist available");
231
232   fhCorr  = dynamic_cast<TH1*> (gDirectory->Get("correction"));
233
234   Bool_t success = kTRUE;
235   if (!fhCorr)
236   {
237     Info("LoadHistograms", "No correction hist available");
238     success = kFALSE;
239   }
240
241   gDirectory->cd("..");
242
243   return success;
244 }
245
246 //____________________________________________________________________
247 void AliCorrectionMatrix::SaveHistograms()
248 {
249   //
250   // saves the histograms
251   //
252
253   gDirectory->mkdir(GetName());
254   gDirectory->cd(GetName());
255
256   if (fhMeas)
257     fhMeas ->Write();
258
259   if (fhGene)
260     fhGene ->Write();
261
262   if (fhCorr)
263     fhCorr->Write();
264
265   gDirectory->cd("..");
266 }
267
268 //____________________________________________________________________
269 void AliCorrectionMatrix::DrawHistograms(const Char_t* canvasName)
270 {
271   //
272   // draws all histograms on one TCanvas
273   // if canvasName is 0 the name of this object is taken
274   //
275
276   if (!canvasName)
277     canvasName = Form("%s_canvas", GetName());
278
279   TCanvas* canvas = new TCanvas(canvasName, GetTitle(), 1200, 400);
280   canvas->Divide(3, 1);
281
282   canvas->cd(1);
283   if (fhMeas)
284     fhMeas->Draw("COLZ");
285
286   canvas->cd(2);
287   if (fhGene)
288   {
289     // work around ROOT bug #22011
290     if (fhGene->GetEntries() == 0)
291       fhGene->SetEntries(1);
292     fhGene->Draw("COLZ");
293   }
294
295   canvas->cd(3);
296   if (fhCorr)
297     fhCorr->Draw("COLZ");
298 }
299
300 //____________________________________________________________________
301 void AliCorrectionMatrix::ReduceInformation()
302 {
303   // this function deletes the measured and generated histograms to reduce the amount of data
304   // in memory
305
306   if (fhMeas)
307   {
308     delete fhMeas;
309     fhMeas = 0;
310   }
311
312   if (fhGene)
313   {
314     delete fhGene;
315     fhGene = 0;
316   }
317 }
318
319 //____________________________________________________________________
320 void AliCorrectionMatrix::Reset(Option_t* option)
321 {
322   // resets the histograms
323
324   if (fhGene)
325     fhGene->Reset(option);
326
327   if (fhMeas)
328     fhMeas->Reset(option);
329
330   if (fhCorr)
331     fhCorr->Reset(option);
332 }
333
334 //____________________________________________________________________
335 void AliCorrectionMatrix::SetCorrectionToUnity()
336 {
337   // sets the correction matrix to unity
338
339   if (!fhCorr)
340     return;
341
342   for (Int_t x=1; x<=fhCorr->GetNbinsX(); ++x)
343     for (Int_t y=1; y<=fhCorr->GetNbinsY(); ++y)
344       for (Int_t z=1; z<=fhCorr->GetNbinsZ(); ++z)
345         fhCorr->SetBinContent(x, y, z, 1);
346 }