]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/AliCorrection.cxx
fixing wrong EOL
[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) : TNamed(name, title),
34   fEventCorr(0),
35   fTrackCorr(0)
36 {
37   // constructor initializing tnamed
38
39   Float_t binLimitsPt[] = {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};
40   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};
41   //Float_t binLimitsVtx[] = {-20,-15,-10,-6,-3,0,3,6,10,15,20};
42   //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};
43   Float_t binLimitsVtx[] = {-20,-15,-10,-8,-6,-4,-2,0,2,4,6,8,10,15,20};
44   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,
45                             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};
46
47   TH3F* dummyBinning = new TH3F("","",14, binLimitsVtx, 40, binLimitsEta , 28, binLimitsPt);
48
49   fEventCorr = new AliCorrectionMatrix2D("EventCorrection", Form("%s EventCorrection", title), 14, binLimitsVtx, 22, binLimitsN);
50   fTrackCorr = new AliCorrectionMatrix3D("TrackCorrection", Form("%s TrackCorrection", title), dummyBinning);
51
52   delete dummyBinning;
53
54   fEventCorr->SetAxisTitles("vtx z [cm]", "Ntracks");
55   fTrackCorr->SetAxisTitles("vtx z [cm]", "#eta", "p_{T} [GeV/c]");
56 }
57
58 //____________________________________________________________________
59 AliCorrection::AliCorrection(const AliCorrection& c) : TNamed(c),
60   fEventCorr(0),
61   fTrackCorr(0)
62 {
63   // copy constructor
64   ((AliCorrection &)c).Copy(*this);
65 }
66
67 //____________________________________________________________________
68 AliCorrection::~AliCorrection()
69 {
70   //
71   // destructor
72   //
73
74   if (fEventCorr)
75   {
76     delete fEventCorr;
77     fEventCorr = 0;
78   }
79
80   if (fTrackCorr)
81   {
82     delete fTrackCorr;
83     fTrackCorr = 0;
84   }
85 }
86
87 //____________________________________________________________________
88 AliCorrection &AliCorrection::operator=(const AliCorrection &c)
89 {
90   // assigment operator
91
92   if (this != &c)
93     ((AliCorrection &) c).Copy(*this);
94
95   return *this;
96 }
97
98 //____________________________________________________________________
99 void AliCorrection::Copy(TObject& c) const
100 {
101   // copy function
102
103   AliCorrection& target = (AliCorrection &) c;
104
105   if (fEventCorr)
106     target.fEventCorr = dynamic_cast<AliCorrectionMatrix2D*> (fEventCorr->Clone());
107
108   if (fTrackCorr)
109     target.fTrackCorr = dynamic_cast<AliCorrectionMatrix3D*> (fTrackCorr->Clone());
110 }
111
112 //____________________________________________________________________
113 Long64_t AliCorrection::Merge(TCollection* list)
114 {
115   // Merge a list of AliCorrection objects with this (needed for
116   // PROOF). 
117   // Returns the number of merged objects (including this).
118
119   if (!list)
120     return 0;
121   
122   if (list->IsEmpty())
123     return 1;
124
125   TIterator* iter = list->MakeIterator();
126   TObject* obj;
127
128   // collections of measured and generated histograms
129   TList* collectionEvent = new TList;
130   TList* collectionTrack = new TList;
131   
132   Int_t count = 0;
133   while ((obj = iter->Next())) {
134     
135     AliCorrection* entry = dynamic_cast<AliCorrection*> (obj);
136     if (entry == 0) 
137       continue;
138
139     collectionEvent->Add(entry->fEventCorr);
140     collectionTrack->Add(entry->fTrackCorr);
141
142     count++;
143   }
144   fEventCorr->Merge(collectionEvent);
145   fTrackCorr->Merge(collectionTrack);
146
147   delete collectionEvent;
148   delete collectionTrack;
149
150   return count+1;
151 }
152
153 //____________________________________________________________________
154 void AliCorrection::Divide()
155 {
156   //
157   // divide the histograms to get the correction
158   //
159   
160   if (!fEventCorr || !fTrackCorr)
161     return;
162     
163   fEventCorr->Divide();
164   fTrackCorr->Divide();
165
166   Int_t emptyBins = fTrackCorr->CheckEmptyBins(-9.99, 9.99, -0.79, 0.79, 0.3, 9.9);
167   printf("INFO: In the central region the track correction of %s has %d empty bins\n", GetName(), emptyBins);
168 }
169
170 //____________________________________________________________________
171 void AliCorrection::Add(AliCorrection* aCorrectionToAdd, Float_t c)
172 {
173   //
174   // add to measured and generated the measured and generated of aCorrectionToAdd
175   // with the weight c
176
177   fEventCorr->Add(aCorrectionToAdd->GetEventCorrection(),c);
178   fTrackCorr->Add(aCorrectionToAdd->GetTrackCorrection(),c);
179 }
180
181
182 //____________________________________________________________________
183 Bool_t AliCorrection::LoadHistograms(const Char_t* dir)
184 {
185   //
186   // loads the histograms from a file
187   // if dir is empty a directory with the name of this object is taken (like in SaveHistogram)
188   //
189
190   if (!fEventCorr || !fTrackCorr)
191     return kFALSE;
192
193   if (!dir)
194     dir = GetName();
195
196   if (!gDirectory->cd(dir))
197     return kFALSE;
198
199   Bool_t success = fEventCorr->LoadHistograms();
200   success &= fTrackCorr->LoadHistograms();
201
202   gDirectory->cd("..");
203
204   return success;
205 }
206
207 //____________________________________________________________________
208 void AliCorrection::SaveHistograms()
209 {
210   //
211   // saves the histograms in a directory with the name of this object (GetName)
212   //
213   
214   gDirectory->mkdir(GetName());
215   gDirectory->cd(GetName());
216
217   if (fEventCorr)
218     fEventCorr->SaveHistograms();
219
220   if (fTrackCorr)
221     fTrackCorr->SaveHistograms();
222     
223   gDirectory->cd("..");
224 }
225
226 //____________________________________________________________________
227 void AliCorrection::ReduceInformation()
228 {
229   // this function deletes the measured and generated histograms to reduce the amount of data
230   // in memory
231
232   if (!fEventCorr || !fTrackCorr)
233     return;
234
235   fEventCorr->ReduceInformation();
236   fTrackCorr->ReduceInformation();
237 }
238
239 //____________________________________________________________________
240 void AliCorrection::Reset(Option_t* option)
241 {
242   // resets the histograms
243
244   if (fEventCorr)
245     fEventCorr->Reset(option);
246
247   if (fTrackCorr)
248     fTrackCorr->Reset(option);
249 }
250
251 //____________________________________________________________________
252 void AliCorrection::DrawHistograms(const Char_t* name)
253 {
254   // draws the corrections
255
256   if (!name)
257     name = GetName();
258
259   if (fEventCorr)
260     fEventCorr->DrawHistograms(Form("%s event", name));
261
262   if (fTrackCorr)
263     fTrackCorr->DrawHistograms(Form("%s track", name));
264 }
265
266 //____________________________________________________________________
267 void AliCorrection::SetCorrectionToUnity()
268 {
269   // set the corrections to unity
270
271   if (fEventCorr)
272     fEventCorr->SetCorrectionToUnity();
273
274   if (fTrackCorr)
275     fTrackCorr->SetCorrectionToUnity();
276 }
277
278 //____________________________________________________________________
279 void AliCorrection::Multiply()
280 {
281   // call Multiply
282
283   if (fEventCorr)
284   {
285     fEventCorr->Multiply();
286     // now we manually copy the overflow bin of the y axis (multiplicity) over. This is important to get the event count correct
287     TH2F* hist = fEventCorr->GetMeasuredHistogram();
288     for (Int_t x = 1; x <= hist->GetNbinsX(); ++x)
289       fEventCorr->GetGeneratedHistogram()->SetBinContent(x, hist->GetNbinsY() + 1, hist->GetBinContent(x, hist->GetNbinsY() + 1));
290   }
291
292   if (fTrackCorr)
293     fTrackCorr->Multiply();
294 }