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