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