]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliFMDCorrector.cxx
Merge branch 'master_patch'
[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 "AliFMDCorrSecondaryMap.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(false),
34     fUseMergingEfficiency(false),
35     fDebug(0)
36 {
37   // Constructor
38   DGUARD(fDebug, 3, "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(false),
48     fUseMergingEfficiency(false),
49     fDebug(0)
50 {
51   // Constructor 
52   // 
53   // Parameters: 
54   //   title   Title
55   DGUARD(fDebug, 3, "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, 3, "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::SetupForData(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 void
164 AliFMDCorrector::DivideMap(TH2* num, const TH2* denom,
165                            Bool_t alsoUnderOver) const
166 {
167   // 
168   // Implement TH1::Divide but 
169   // - Assume compatible histograms 
170   // - Unless third argument is true, do not divide over/under flow bins
171   // 
172   if (!num || !denom) return;
173
174   Int_t first = (alsoUnderOver ? 0 : 1);
175   Int_t lastX = num->GetNbinsX() + (alsoUnderOver ? 1 : 0);
176   Int_t lastY = num->GetNbinsY() + (alsoUnderOver ? 1 : 0);
177   
178   for (Int_t ix = first; ix <= lastX; ix++) {
179     for (Int_t iy = first; iy <= lastY; iy++) { 
180       Int_t    bin = num->GetBin(ix,iy);
181       Double_t c0  = num->GetBinContent(bin);
182       Double_t c1  = denom->GetBinContent(bin);
183       if (!c1) { 
184         num->SetBinContent(bin,0);
185         num->SetBinError(bin, 0);
186         continue;
187       }
188       Double_t w   = c0 / c1;
189       Double_t e0  = num->GetBinError(bin);
190       Double_t e1  = denom->GetBinError(bin);
191       Double_t c12 = c1*c1;
192       Double_t e2  = (e0*e0*c1*c1 + e1*e1*c0*c0)/(c12*c12);
193       
194       num->SetBinContent(bin, w);
195       num->SetBinError(bin, TMath::Sqrt(e2));
196     }
197   }
198 }
199 //____________________________________________________________________
200 Bool_t
201 AliFMDCorrector::Correct(AliForwardUtil::Histos& hists,
202                          UShort_t                vtxbin)
203 {
204   // 
205   // Do the calculations 
206   // Parameters:
207   //    hists    Cache of histograms 
208   //    vtxBin   Vertex bin 
209   // 
210   // Return:
211   //    true on successs 
212   //
213   DGUARD(fDebug, 1, "Correct histograms of AliFMDCorrector");
214   AliForwardCorrectionManager& fcm = AliForwardCorrectionManager::Instance();
215
216   UShort_t uvb = vtxbin;
217   for (UShort_t d=1; d<=3; d++) { 
218     UShort_t nr = (d == 1 ? 1 : 2);
219     for (UShort_t q=0; q<nr; q++) { 
220       Char_t      r  = (q == 0 ? 'I' : 'O');
221       TH2D*       h  = hists.Get(d,r);
222       RingHistos* rh = GetRingHistos(d,r);
223
224       if (fUseSecondaryMap) {
225         TH2D*  bg = fcm.GetSecondaryMap()->GetCorrection(d, r, uvb);
226         if (!bg) {
227           AliWarning(Form("No secondary correction for FMDM%d%c "
228                           "in vertex bin %d", d, r, uvb));
229           continue;
230         }
231         // Divide by primary/total ratio
232         DivideMap(h, bg, false);
233       }
234       if (fUseVertexBias) {
235         TH2D*  ef = fcm.GetVertexBias()->GetCorrection(r, uvb);
236         if (!ef) {
237           AliWarning(Form("No event %s vertex bias correction in vertex bin %d",
238                           (r == 'I' || r == 'i' ? "inner" : "outer"), uvb));
239           continue;
240         }
241         // Divide by the event selection efficiency
242         DivideMap(h, ef, false);
243       }
244       if (fUseAcceptance) {
245         TH2D*  ac = fcm.GetAcceptance()->GetCorrection(d, r, uvb);
246         if (!ac) {
247           AliWarning(Form("No acceptance correction for FMD%d%c in "
248                           "vertex bin %d", d, r, uvb));
249           continue;
250         }
251         // Fill overflow bin with ones 
252         for (Int_t i = 1; i <= h->GetNbinsX(); i++) 
253           h->SetBinContent(i, h->GetNbinsY()+1, 1);
254
255         // Divide by the acceptance correction - 
256         DivideMap(h, ac, fcm.GetAcceptance()->HasOverflow());
257       }
258
259       if (fUseMergingEfficiency) {
260         if (!fcm.GetMergingEfficiency()) { 
261           AliWarning("No merging efficiencies");
262           continue;
263         }
264         TH1D* sf = fcm.GetMergingEfficiency()->GetCorrection(d,r,uvb);
265         if (!fcm.GetMergingEfficiency()->GetCorrection(d,r,uvb)) { 
266           AliWarning(Form("No merging efficiency for FMD%d%c at vertex bin %d",
267                           d, r, uvb));
268           continue;
269         }
270         
271       
272         for (Int_t ieta = 1; ieta <= h->GetNbinsX(); ieta++) {
273           Float_t c  = sf->GetBinContent(ieta);
274           Float_t ec = sf->GetBinError(ieta);
275                   
276           for (Int_t iphi = 1; iphi <= h->GetNbinsY(); iphi++) { 
277             if (c == 0) {
278               h->SetBinContent(ieta,iphi,0);
279               h->SetBinError(ieta,iphi,0);
280               continue;
281             }
282
283             Double_t m  = h->GetBinContent(ieta, iphi) / c;
284             Double_t em = h->GetBinError(ieta, iphi);
285           
286             Double_t e  = TMath::Sqrt(em * em + (m * ec) * (m * ec));
287             
288             h->SetBinContent(ieta,iphi,m);
289             h->SetBinError(ieta,iphi,e);
290           }
291         }
292       }
293       
294       rh->fDensity->Add(h);
295     }
296   }
297   
298   return kTRUE;
299 }
300
301 //____________________________________________________________________
302 void
303 AliFMDCorrector::Terminate(const TList* dir, TList* output, Int_t nEvents)
304 {
305   // 
306   // Scale the histograms to the total number of events 
307   // Parameters:
308   //    dir     Where the output is stored
309   //    nEvents Number of events 
310   //
311   DGUARD(fDebug, 1, "Scale histograms of AliFMDCorrector");
312   if (nEvents <= 0) return;
313   TList* d = static_cast<TList*>(dir->FindObject(GetName()));
314   if (!d) return;
315
316   TList* out = new TList;
317   out->SetName(d->GetName());
318   out->SetOwner();
319
320   TIter    next(&fRingHistos);
321   RingHistos* o = 0;
322   THStack* sums = new THStack("sums", "Sums of ring results");
323   while ((o = static_cast<RingHistos*>(next()))) {
324     o->Terminate(d, nEvents);
325     if (!o->fDensity) { 
326       Warning("Terminate", "No density from %s", o->GetName());
327       continue;
328     }
329     TH1D* sum = o->fDensity->ProjectionX(o->GetName(), 1, 
330                                          o->fDensity->GetNbinsY(),"e");
331     sum->Scale(1., "width");
332     sum->SetTitle(o->GetName());
333     sum->SetDirectory(0);
334     sum->SetYTitle("#sum N_{ch,primary}");
335     sums->Add(sum);
336   }
337   out->Add(sums);
338   output->Add(out);
339 }
340 //____________________________________________________________________
341 void
342 AliFMDCorrector::CreateOutputObjects(TList* dir)
343 {
344   // 
345   // Output diagnostic histograms to directory 
346   // 
347   // Parameters:
348   //    dir List to write in
349   //  
350   DGUARD(fDebug, 1, "Define output of AliFMDCorrector");
351   TList* d = new TList;
352   d->SetName(GetName());
353   dir->Add(d);
354
355   d->Add(AliForwardUtil::MakeParameter("secondary", fUseSecondaryMap));
356   d->Add(AliForwardUtil::MakeParameter("vertexBias", fUseVertexBias));
357   d->Add(AliForwardUtil::MakeParameter("acceptance", fUseAcceptance));
358   d->Add(AliForwardUtil::MakeParameter("merging", fUseMergingEfficiency));
359   
360
361   TIter    next(&fRingHistos);
362   RingHistos* o = 0;
363   while ((o = static_cast<RingHistos*>(next()))) {
364     o->CreateOutputObjects(d);
365   }
366 }
367
368 #define PF(N,V,...)                                     \
369   AliForwardUtil::PrintField(N,V, ## __VA_ARGS__)
370 #define PFB(N,FLAG)                             \
371   do {                                                                  \
372     AliForwardUtil::PrintName(N);                                       \
373     std::cout << std::boolalpha << (FLAG) << std::noboolalpha << std::endl; \
374   } while(false)
375 #define PFV(N,VALUE)                                    \
376   do {                                                  \
377     AliForwardUtil::PrintName(N);                       \
378     std::cout << (VALUE) << std::endl; } while(false)
379
380 //____________________________________________________________________
381 void
382 AliFMDCorrector::Print(Option_t* /* option */) const
383 {
384   // 
385   // Print information
386   // Parameters:
387   //    option Not used 
388   //
389   AliForwardUtil::PrintTask(*this);
390   gROOT->IncreaseDirLevel();
391   PFB("Use secondary maps",             fUseSecondaryMap );
392   PFB("Use vertex bias",                fUseVertexBias );
393   PFB("Use acceptance",                 fUseAcceptance );
394   PFB("Use merging efficiency",         fUseMergingEfficiency);
395   gROOT->DecreaseDirLevel();  
396 }
397
398 //====================================================================
399 AliFMDCorrector::RingHistos::RingHistos()
400   : AliForwardUtil::RingHistos(), 
401     fDensity(0)
402 {
403   // Constructor 
404   //
405   // 
406 }
407
408 //____________________________________________________________________
409 AliFMDCorrector::RingHistos::RingHistos(UShort_t d, Char_t r)
410   : AliForwardUtil::RingHistos(d,r), 
411     fDensity(0)
412 {
413   // 
414   // Constructor
415   // Parameters:
416   //    d detector
417   //    r ring 
418   //
419   fDensity = new TH2D("primaryDensity", 
420                       "#sum N_{ch,primary}/(#Delta#eta#Delta#phi)", 
421                       200, -4, 6, (r == 'I' || r == 'i' ? 20 : 40), 
422                       0, 2*TMath::Pi());
423   fDensity->SetDirectory(0);
424   fDensity->SetMarkerColor(Color());
425   fDensity->Sumw2();
426   fDensity->SetXTitle("#eta");
427   fDensity->SetYTitle("#phi [radians]");
428   fDensity->SetZTitle("Primary N_{ch} density");
429 }
430 //____________________________________________________________________
431 AliFMDCorrector::RingHistos::RingHistos(const RingHistos& o)
432   : AliForwardUtil::RingHistos(o), 
433     fDensity(o.fDensity)
434 {
435   // 
436   // Copy constructor 
437   // Parameters:
438   //    o Object to copy from 
439   //
440 }
441
442 //____________________________________________________________________
443 AliFMDCorrector::RingHistos&
444 AliFMDCorrector::RingHistos::operator=(const RingHistos& o)
445 {
446   // 
447   // Assignment operator 
448   // Parameters:
449   //    o Object to assign from 
450   // 
451   // Return:
452   //    Reference to this 
453   //
454   if (&o == this) return *this; 
455   AliForwardUtil::RingHistos::operator=(o);
456   
457   if (fDensity) delete fDensity;
458   
459   fDensity = static_cast<TH2D*>(o.fDensity->Clone());
460   
461   return *this;
462 }
463 //____________________________________________________________________
464 AliFMDCorrector::RingHistos::~RingHistos()
465 {
466   // 
467   // Destructor 
468   //
469   // if (fDensity) delete fDensity;
470 }
471
472 //____________________________________________________________________
473 void
474 AliFMDCorrector::RingHistos::CreateOutputObjects(TList* dir)
475 {
476   // 
477   // Make output 
478   // Parameters:
479   //    dir Where to put it 
480   //
481   TList* d = DefineOutputList(dir);
482   d->Add(fDensity);
483 }
484
485 //____________________________________________________________________
486 void
487 AliFMDCorrector::RingHistos::Terminate(TList* dir, Int_t nEvents)
488
489   // 
490   // Scale the histograms to the total number of events 
491   // Parameters:
492   //    dir     where the output is stored
493   //    nEvents Number of events 
494   //
495   TList* l = GetOutputList(dir);
496   if (!l) return; 
497
498   TH2D* density = static_cast<TH2D*>(GetOutputHist(l,"primaryDensity"));
499   if (density) density->Scale(1./nEvents);
500   fDensity = density;
501 }
502
503 //____________________________________________________________________
504 //
505 // EOF
506 //
507           
508
509