]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FORWARD/analysis2/AliFMDCorrector.cxx
Various fixes to deal with centrality
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / AliFMDCorrector.cxx
1 //
2 // This class calculates the exclusive charged particle density
3 // in each for the 5 FMD rings. 
4 //
5 #include "AliFMDCorrector.h"
6 #include <AliESDFMD.h>
7 #include <TAxis.h>
8 #include <TList.h>
9 #include <TMath.h>
10 #include "AliForwardCorrectionManager.h"
11 // #include "AliFMDCorrDoubleHit.h"
12 #include "AliFMDCorrVertexBias.h"
13 #include "AliFMDCorrMergingEfficiency.h"
14 #include "AliFMDCorrAcceptance.h"
15 #include "AliLog.h"
16 #include <TH2D.h>
17 #include <TROOT.h>
18 #include <iostream>
19 #include <iomanip>
20
21 ClassImp(AliFMDCorrector)
22 #if 0
23 ; // For Emacs
24 #endif 
25
26 //____________________________________________________________________
27 AliFMDCorrector::AliFMDCorrector()
28   : TNamed(), 
29     fRingHistos(),
30     fUseSecondaryMap(true),
31     fUseVertexBias(false),
32     fUseAcceptance(true),
33     fUseMergingEfficiency(false),
34     fDebug(0)
35 {
36   // Constructor
37 }
38
39 //____________________________________________________________________
40 AliFMDCorrector::AliFMDCorrector(const char* title)
41   : TNamed("fmdCorrector", title), 
42     fRingHistos(), 
43     fUseSecondaryMap(true),
44     fUseVertexBias(false),
45     fUseAcceptance(true),
46     fUseMergingEfficiency(false),
47     fDebug(0)
48 {
49   // Constructor 
50   // 
51   // Parameters: 
52   //   title   Title
53   fRingHistos.SetName(GetName());
54   fRingHistos.Add(new RingHistos(1, 'I'));
55   fRingHistos.Add(new RingHistos(2, 'I'));
56   fRingHistos.Add(new RingHistos(2, 'O'));
57   fRingHistos.Add(new RingHistos(3, 'I'));
58   fRingHistos.Add(new RingHistos(3, 'O'));
59 }
60
61 //____________________________________________________________________
62 AliFMDCorrector::AliFMDCorrector(const AliFMDCorrector& o)
63   : TNamed(o), 
64     fRingHistos(), 
65     fUseSecondaryMap(o.fUseSecondaryMap),
66     fUseVertexBias(o.fUseVertexBias),
67     fUseAcceptance(o.fUseAcceptance),
68     fUseMergingEfficiency(o.fUseMergingEfficiency),
69     fDebug(o.fDebug)
70 {
71   // Copy constructor 
72   // 
73   // Parameters: 
74   //  o  Object to copy from 
75   TIter    next(&o.fRingHistos);
76   TObject* obj = 0;
77   while ((obj = next())) fRingHistos.Add(obj);
78 }
79
80 //____________________________________________________________________
81 AliFMDCorrector::~AliFMDCorrector()
82 {
83   // Destructor 
84   // 
85   // 
86   fRingHistos.Delete();
87 }
88
89 //____________________________________________________________________
90 AliFMDCorrector&
91 AliFMDCorrector::operator=(const AliFMDCorrector& o)
92 {
93   // Assignment operator 
94   // 
95   // Parameters:
96   //   o   Object to assign from 
97   TNamed::operator=(o);
98
99   fDebug   = o.fDebug;
100   fRingHistos.Delete();
101   fUseSecondaryMap = o.fUseSecondaryMap;
102   fUseVertexBias = o.fUseVertexBias;
103   fUseAcceptance = o.fUseAcceptance;
104   fUseMergingEfficiency = o.fUseMergingEfficiency;
105   TIter    next(&o.fRingHistos);
106   TObject* obj = 0;
107   while ((obj = next())) fRingHistos.Add(obj);
108   
109   return *this;
110 }
111
112 //____________________________________________________________________
113 void
114 AliFMDCorrector::Init(const TAxis&)
115 {
116   //
117   // Initialize this object
118   //
119   // Parameters:
120   //    etaAxis Eta axis to use
121   //
122   if (!fUseSecondaryMap)
123     AliWarning("Secondary maps not used - BE CAREFUL");
124   if (!fUseVertexBias)
125     AliWarning("Vertex bias not used");
126   if (!fUseAcceptance)
127     AliWarning("Acceptance from dead-channels not used");
128 }
129
130 //____________________________________________________________________
131 AliFMDCorrector::RingHistos*
132 AliFMDCorrector::GetRingHistos(UShort_t d, Char_t r) const
133 {
134   // 
135   // Get the ring histogram container 
136   // Parameters:
137   //    d Detector
138   //    r Ring 
139   // 
140   // Return:
141   //    Ring histogram container 
142   //
143   Int_t idx = -1;
144   switch (d) { 
145   case 1: idx = 0; break;
146   case 2: idx = 1 + (r == 'I' || r == 'i' ? 0 : 1); break;
147   case 3: idx = 3 + (r == 'I' || r == 'i' ? 0 : 1); break;
148   }
149   if (idx < 0 || idx >= fRingHistos.GetEntries()) return 0;
150   
151   return static_cast<RingHistos*>(fRingHistos.At(idx));
152 }
153     
154 //____________________________________________________________________
155 Bool_t
156 AliFMDCorrector::Correct(AliForwardUtil::Histos& hists,
157                          UShort_t                vtxbin)
158 {
159   // 
160   // Do the calculations 
161   // Parameters:
162   //    hists    Cache of histograms 
163   //    vtxBin   Vertex bin 
164   // 
165   // Return:
166   //    true on successs 
167   //
168   AliForwardCorrectionManager& fcm = AliForwardCorrectionManager::Instance();
169
170   UShort_t uvb = vtxbin;
171   for (UShort_t d=1; d<=3; d++) { 
172     UShort_t nr = (d == 1 ? 1 : 2);
173     for (UShort_t q=0; q<nr; q++) { 
174       Char_t      r  = (q == 0 ? 'I' : 'O');
175       TH2D*       h  = hists.Get(d,r);
176       RingHistos* rh = GetRingHistos(d,r);
177
178       if (fUseSecondaryMap) {
179         TH2D*  bg = fcm.GetSecondaryMap()->GetCorrection(d, r, uvb);
180         if (!bg) {
181           AliWarning(Form("No secondary correction for FMDM%d%c in vertex bin %d",
182                           d, r, uvb));
183           continue;
184         }
185         // Divide by primary/total ratio
186         h->Divide(bg);
187       }
188       if (fUseVertexBias) {
189         TH2D*  ef = fcm.GetVertexBias()->GetCorrection(r, uvb);
190         if (!ef) {
191           AliWarning(Form("No event %s vertex bias correction in vertex bin %d",
192                           (r == 'I' || r == 'i' ? "inner" : "outer"), uvb));
193           continue;
194         }
195         // Divide by the event selection efficiency
196         h->Divide(ef);
197       }
198       if (fUseAcceptance) {
199         TH2D*  ac = fcm.GetAcceptance()->GetCorrection(d, r, uvb);
200         if (!ac) {
201           AliWarning(Form("No acceptance correction for FMD%d%c in vertex bin %d",
202                         d, r, uvb));
203           continue;
204         }
205         // Divide by the acceptance correction
206         h->Divide(ac);
207       }
208
209
210
211       
212
213       if (fUseMergingEfficiency) {
214         if (!fcm.GetMergingEfficiency()) { 
215           AliWarning("No merging efficiencies");
216           continue;
217         }
218         TH1D* sf = fcm.GetMergingEfficiency()->GetCorrection(d,r,uvb);
219         if (!fcm.GetMergingEfficiency()->GetCorrection(d,r,uvb)) { 
220           AliWarning(Form("No merging efficiency for FMD%d%c at vertex bin %d",
221                           d, r, uvb));
222           continue;
223         }
224         
225       
226         for (Int_t ieta = 1; ieta <= h->GetNbinsX(); ieta++) {
227           Float_t c  = sf->GetBinContent(ieta);
228           Float_t ec = sf->GetBinError(ieta);
229           
230           if (c == 0) continue;
231           
232           for (Int_t iphi = 1; iphi <= h->GetNbinsY(); iphi++) { 
233             Double_t m  = h->GetBinContent(ieta, iphi) / c;
234             Double_t em = h->GetBinError(ieta, iphi);
235           
236             Double_t e  = TMath::Sqrt(em * em + (m * ec) * (m * ec));
237             
238             h->SetBinContent(ieta,iphi,m);
239             h->SetBinError(ieta,iphi,e);
240           }
241         }
242       }
243       //HHD
244       /*
245       TH2D*  bg = fcm.GetSecondaryMap()->GetCorrection(d, r, uvb);
246       TH2D   hRing("hring","hring",bg->GetNbinsX(),
247                    bg->GetXaxis()->GetXmin(),
248                    bg->GetXaxis()->GetXmax(),
249                    bg->GetNbinsY(),
250                    bg->GetYaxis()->GetXmin(),
251                    bg->GetYaxis()->GetXmax());
252       
253       Int_t edgebin[4] = {0,0,0,0};
254       for(Int_t ii = 1; ii <=bg->GetNbinsX(); ii++) {
255         for(Int_t jj = 1; jj <=bg->GetNbinsY(); jj++) {
256           Float_t bgcor = bg->GetBinContent(ii,jj);
257           if(bgcor<0.1) continue;
258           if(edgebin[0] == 0) edgebin[0] = ii;
259           if(edgebin[0] == ii) continue;
260           if(edgebin[0] > 0 && edgebin[1] == 0) edgebin[1] = ii;
261           if(edgebin[0]>0 && edgebin[1]>0) break; 
262         }
263       }
264       for(Int_t ii = bg->GetNbinsX(); ii >= 1;  ii--) {
265         for(Int_t jj = 1; jj <=bg->GetNbinsY(); jj++) {
266           Float_t bgcor = bg->GetBinContent(ii,jj);
267           if(bgcor<0.1) continue;
268           if(edgebin[2] == 0) edgebin[2] = ii;
269           if(edgebin[2] == ii) continue;
270           if(edgebin[2] > 0 && edgebin[3] == 0) edgebin[3] = ii;
271           if(edgebin[2]>0 && edgebin[3]>0) break; 
272         }
273       }
274       for(Int_t ii = 1; ii <=bg->GetNbinsX(); ii++) {
275         for(Int_t jj = 1; jj <=bg->GetNbinsY(); jj++) {
276           Float_t data = h->GetBinContent(ii,jj);
277           if(data <0.000001) continue;
278           if(edgebin[0] == ii || edgebin[1] == ii || edgebin[2] == ii || edgebin[3] == ii) continue;
279           hRing.SetBinContent(ii,jj,data);
280           hRing.SetBinError(ii,jj,h->GetBinError(ii,jj));
281         }
282       }
283       
284       //std::cout<<edgebin[0]<<"  "<<edgebin[1]<<"  "<<edgebin[2]<<"   "<<edgebin[3]<<std::endl;
285       */
286       
287       rh->fDensity->Add(h);
288     }
289   }
290   
291   return kTRUE;
292 }
293
294 //____________________________________________________________________
295 void
296 AliFMDCorrector::ScaleHistograms(const TList* dir, Int_t nEvents)
297 {
298   // 
299   // Scale the histograms to the total number of events 
300   // Parameters:
301   //    dir     Where the output is stored
302   //    nEvents Number of events 
303   //
304   if (nEvents <= 0) return;
305   TList* d = static_cast<TList*>(dir->FindObject(GetName()));
306   if (!d) return;
307
308   TIter    next(&fRingHistos);
309   RingHistos* o = 0;
310   while ((o = static_cast<RingHistos*>(next())))
311     o->ScaleHistograms(d, nEvents);
312 }
313 //____________________________________________________________________
314 void
315 AliFMDCorrector::DefineOutput(TList* dir)
316 {
317   // 
318   // Output diagnostic histograms to directory 
319   // 
320   // Parameters:
321   //    dir List to write in
322   //  
323   TList* d = new TList;
324   d->SetName(GetName());
325   dir->Add(d);
326   TIter    next(&fRingHistos);
327   RingHistos* o = 0;
328   while ((o = static_cast<RingHistos*>(next()))) {
329     o->Output(d);
330   }
331 }
332
333 //____________________________________________________________________
334 void
335 AliFMDCorrector::Print(Option_t* /* option */) const
336 {
337   // 
338   // Print information
339   // Parameters:
340   //    option Not used 
341   //
342   char ind[gROOT->GetDirLevel()+1];
343   for (Int_t i = 0; i < gROOT->GetDirLevel(); i++) ind[i] = ' ';
344   ind[gROOT->GetDirLevel()] = '\0';
345   std::cout << ind << "AliFMDCorrector: " << GetName() <<  "\n"
346             << std::boolalpha
347             << ind << " Use secondary maps:     " << fUseSecondaryMap << "\n"
348             << ind << " Use vertex bias:        " << fUseVertexBias << "\n"
349             << ind << " Use acceptance:         " << fUseAcceptance << "\n"
350             << ind << " Use merging efficiency: " << fUseMergingEfficiency
351             << std::noboolalpha << std::endl;
352 }
353
354 //====================================================================
355 AliFMDCorrector::RingHistos::RingHistos()
356   : AliForwardUtil::RingHistos(), 
357     fDensity(0)
358 {
359   // Constructor 
360   //
361   // 
362 }
363
364 //____________________________________________________________________
365 AliFMDCorrector::RingHistos::RingHistos(UShort_t d, Char_t r)
366   : AliForwardUtil::RingHistos(d,r), 
367     fDensity(0)
368 {
369   // 
370   // Constructor
371   // Parameters:
372   //    d detector
373   //    r ring 
374   //
375   fDensity = new TH2D(Form("FMD%d%c_Primary_Density", d, r), 
376                       Form("in FMD%d%c", d, r), 
377                       200, -4, 6, (r == 'I' || r == 'i' ? 20 : 40), 
378                       0, 2*TMath::Pi());
379   fDensity->SetDirectory(0);
380   fDensity->SetXTitle("#eta");
381   fDensity->SetYTitle("#phi [radians]");
382   fDensity->SetZTitle("Primary N_{ch} density");
383 }
384 //____________________________________________________________________
385 AliFMDCorrector::RingHistos::RingHistos(const RingHistos& o)
386   : AliForwardUtil::RingHistos(o), 
387     fDensity(o.fDensity)
388 {
389   // 
390   // Copy constructor 
391   // Parameters:
392   //    o Object to copy from 
393   //
394 }
395
396 //____________________________________________________________________
397 AliFMDCorrector::RingHistos&
398 AliFMDCorrector::RingHistos::operator=(const RingHistos& o)
399 {
400   // 
401   // Assignment operator 
402   // Parameters:
403   //    o Object to assign from 
404   // 
405   // Return:
406   //    Reference to this 
407   //
408   AliForwardUtil::RingHistos::operator=(o);
409   
410   if (fDensity) delete fDensity;
411   
412   fDensity = static_cast<TH2D*>(o.fDensity->Clone());
413   
414   return *this;
415 }
416 //____________________________________________________________________
417 AliFMDCorrector::RingHistos::~RingHistos()
418 {
419   // 
420   // Destructor 
421   //
422   if (fDensity) delete fDensity;
423 }
424
425 //____________________________________________________________________
426 void
427 AliFMDCorrector::RingHistos::Output(TList* dir)
428 {
429   // 
430   // Make output 
431   // Parameters:
432   //    dir Where to put it 
433   //
434   TList* d = DefineOutputList(dir);
435   d->Add(fDensity);
436 }
437
438 //____________________________________________________________________
439 void
440 AliFMDCorrector::RingHistos::ScaleHistograms(TList* dir, Int_t /*nEvents*/)
441
442   // 
443   // Scale the histograms to the total number of events 
444   // Parameters:
445   //    dir     where the output is stored
446   //    nEvents Number of events 
447   //
448   TList* l = GetOutputList(dir);
449   if (!l) return; 
450
451   //TH1* density = GetOutputHist(l,Form("%s_Primary_Density", fName.Data()));
452   //if (density) density->Scale(1./nEvents);
453 }
454
455 //____________________________________________________________________
456 //
457 // EOF
458 //
459           
460
461