]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliFMDCorrAcceptance.cxx
Added possibility to get phi acceptance as a 1D histogram
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliFMDCorrAcceptance.cxx
1 //
2 // This class contains the acceptance correction due to dead channels 
3 // 
4 //
5 #include "AliFMDCorrAcceptance.h"
6 #include <TBrowser.h>
7 #include <TH2D.h>
8 #include <AliLog.h>
9 #include <iostream>
10
11 //____________________________________________________________________
12 AliFMDCorrAcceptance::AliFMDCorrAcceptance()
13   : fRingArray(), 
14     fCache(0),
15     fVertexAxis(0,0,0),
16     fHasOverflow(false)
17 {
18   // 
19   // Default constructor 
20   //
21   fRingArray.SetOwner(kTRUE);
22   fRingArray.SetName("rings");
23   fVertexAxis.SetName("vtxAxis");
24   fVertexAxis.SetTitle("v_{z} [cm]");
25   
26 }
27 //____________________________________________________________________
28 AliFMDCorrAcceptance::AliFMDCorrAcceptance(const 
29                                                AliFMDCorrAcceptance& o)
30   : TObject(o), 
31     fRingArray(o.fRingArray), 
32     fCache(o.fCache),
33     fVertexAxis(o.fVertexAxis.GetNbins(), o.fVertexAxis.GetXmin(), 
34                 o.fVertexAxis.GetXmax()),
35     fHasOverflow(o.fHasOverflow)
36 {
37   // 
38   // Copy constructor 
39   // 
40   // Parameters:
41   //    o Object to copy from 
42   //
43   fVertexAxis.SetName("vtxAxis");
44   fVertexAxis.SetTitle("v_{z} [cm]");
45 }
46 //____________________________________________________________________
47 AliFMDCorrAcceptance::~AliFMDCorrAcceptance()
48 {
49   //
50   // Destructor 
51   // 
52   //
53   fRingArray.Clear();
54   if (fCache) fCache->Clear();
55 }
56 //____________________________________________________________________
57 AliFMDCorrAcceptance&
58 AliFMDCorrAcceptance::operator=(const AliFMDCorrAcceptance& o)
59 {
60   // 
61   // Assignment operator 
62   // 
63   // Parameters:
64   //    o Object to assign from 
65   // 
66   // Return:
67   //    Reference to this object 
68   //
69   fRingArray        = o.fRingArray;
70   fCache            = o.fCache;
71   fHasOverflow      = o.fHasOverflow;
72   SetVertexAxis(o.fVertexAxis);
73
74   return *this;
75 }
76 //____________________________________________________________________
77 TH2D*
78 AliFMDCorrAcceptance::GetCorrection(UShort_t d, Char_t r, Double_t v) const
79 {
80   // 
81   // Get the acceptance correction @f$ a_{r,v}@f$ 
82   // 
83   // Parameters:
84   //    d  Detector number (1-3)
85   //    r  Ring identifier (I or O)
86   //    v  Primary interaction point @f$z@f$ coordinate
87   // 
88   // Return:
89   //    The correction @f$ a_{r,v}@f$ 
90   //
91   Int_t b = FindVertexBin(v);
92   if (b <= 0) return 0;
93   return GetCorrection(d, r, UShort_t(b));
94 }
95 //____________________________________________________________________
96 TH2D*
97 AliFMDCorrAcceptance::GetCorrection(UShort_t d, Char_t r, UShort_t b) const
98 {
99   // 
100   // Get the acceptance correction @f$ a_{r,v}@f$ 
101   // 
102   // Parameters:
103   //    d  Detector number (1-3)
104   //    r  Ring identifier (I or O)
105   //    b  Bin corresponding to the primary interaction point 
106   //           @f$z@f$ coordinate (1 based)
107   // 
108   // Return:
109   //    The correction @f$ a_{r,v}@f$ 
110   //
111   return static_cast<TH2D*>(GetObject(fRingArray, d, r, b));
112 }
113 //____________________________________________________________________
114 TH1D*
115 AliFMDCorrAcceptance::GetPhiAcceptance(UShort_t d, Char_t r, Double_t v) const
116 {
117   // 
118   // Get the acceptance correction @f$ a_{r,v}@f$ 
119   // 
120   // Parameters:
121   //    d  Detector number (1-3)
122   //    r  Ring identifier (I or O)
123   //    v  Primary interaction point @f$z@f$ coordinate
124   // 
125   // Return:
126   //    The correction @f$ a_{r,v}@f$ 
127   //
128   Int_t b = FindVertexBin(v);
129   if (b <= 0) return 0;
130   return GetPhiAcceptance(d, r, UShort_t(b));
131 }
132 //____________________________________________________________________
133 TH1D*
134 AliFMDCorrAcceptance::GetPhiAcceptance(UShort_t d, Char_t r, UShort_t b) const
135 {
136   // 
137   // Get the acceptance correction @f$ a_{r,v}@f$ 
138   // 
139   // Parameters:
140   //    d  Detector number (1-3)
141   //    r  Ring identifier (I or O)
142   //    b  Bin corresponding to the primary interaction point 
143   //           @f$z@f$ coordinate (1 based)
144   // 
145   // Return:
146   //    The correction @f$ a_{r,v}@f$ 
147   //
148   if (!fHasOverflow) return 0;
149   if (!fCache) FillCache();
150   return static_cast<TH1D*>(GetObject(*fCache, d, r, b));
151 }
152   
153 //____________________________________________________________________
154 Int_t
155 AliFMDCorrAcceptance::FindVertexBin(Double_t v) const
156 {
157   // 
158   // Find the vertex bin that corresponds to the passed vertex 
159   // 
160   // Parameters:
161   //    vertex The interaction points @f$z@f$-coordinate 
162   // 
163   // Return:
164   //    Vertex bin in @f$[1,N_{\mbox{vertex}}]@f$ or negative if 
165   // out of range 
166   //
167   if (fVertexAxis.GetNbins() <= 0) { 
168     AliWarning("No vertex array defined");
169     return 0;
170   }
171   Int_t bin = const_cast<TAxis&>(fVertexAxis).FindBin(v);
172   if (bin <= 0 || bin > fVertexAxis.GetNbins()) { 
173     AliWarning(Form("vertex %+8.4f out of range [%+8.4f,%+8.4f]",
174                     v, fVertexAxis.GetXmin(), fVertexAxis.GetXmax()));
175     return 0;
176   }
177   return bin;
178 }
179 //____________________________________________________________________
180 Int_t
181 AliFMDCorrAcceptance::GetRingIndex(UShort_t d, Char_t r) const
182 {
183   // Get the index corresponding to the given ring 
184   // 
185   // Parameters:
186   //    d Detector
187   //    r Ring 
188   // 
189   // Return:
190   //    Index (0 based) or negative in case of errors
191   //
192   switch (d) {
193   case 1:  return 0;
194   case 2:  return (r == 'I' || r == 'i' ? 1 : 2); break;  
195   case 3:  return (r == 'I' || r == 'i' ? 3 : 4); break;  
196
197   }
198   AliWarning(Form("Index for FMD%d%c not found", d, r));
199   return -1;
200 }
201 //____________________________________________________________________
202 TObject*
203 AliFMDCorrAcceptance::GetObject(const TObjArray& m, UShort_t d, 
204                                 Char_t r, UShort_t b) const
205 {
206   // 
207   // Get the object @f$ a_{r,v}@f$ 
208   // 
209   // Parameters:
210   //    m  Mother list
211   //    d  Detector number (1-3)
212   //    r  Ring identifier (I or O)
213   //    b  Bin corresponding to the primary interaction point 
214   //           @f$z@f$ coordinate (1 based)
215   // 
216   // Return:
217   //    The correction @f$ a_{r,v}@f$ 
218   //
219   TObjArray* ringArray = GetRingArray(m, d, r);
220   if (!ringArray) return 0;
221
222   if (b <= 0 || b > ringArray->GetEntriesFast()) {
223     AliWarning(Form("vertex bin %d out of range [1,%d]", 
224                     b, ringArray->GetEntriesFast()));
225     return 0;
226   }
227
228   TObject* o = ringArray->At(b-1);
229   if (o) return o;
230
231   AliWarning(Form("No dead channels map found for FMD%d%c in vertex bin %d",
232                     d,r,b));
233   return 0;
234 }
235 //____________________________________________________________________
236 TObjArray*
237 AliFMDCorrAcceptance::GetRingArray(const TObjArray& m, 
238                                    UShort_t d, Char_t r) const
239 {
240   // 
241   // Get the ring array corresponding to the specified ring
242   // 
243   // Parameters:
244   //    d Detector 
245   //    r Ring 
246   // 
247   // Return:
248   //    Pointer to ring array, or null in case of problems
249   //
250   Int_t idx = GetRingIndex(d,r);
251   if (idx < 0) return 0;
252   
253   TObject* o = m.At(idx);
254   if (!o) { 
255     AliWarning(Form("No array found for FMD%d%c", d, r));
256     return 0;
257   }
258
259   return static_cast<TObjArray*>(o);
260 }
261 //____________________________________________________________________
262 TObjArray*
263 AliFMDCorrAcceptance::GetOrMakeRingArray(TObjArray& m, 
264                                          UShort_t d, Char_t r) const
265 {
266   // 
267   // Get the ring array corresponding to the specified ring
268   // 
269   // Parameters:
270   //    d Detector 
271   //    r Ring 
272   // 
273   // Return:
274   //    Pointer to ring array, or newly created container 
275   //
276   Int_t idx = GetRingIndex(d,r);
277   if (idx < 0) return 0;
278   
279   TObject* o = m.At(idx);
280   if (!o) { 
281     TObjArray* a = new TObjArray(fVertexAxis.GetNbins());
282     a->SetName(Form("FMD%d%c", d, r));
283     a->SetOwner(kTRUE);
284     m.AddAtAndExpand(a, idx);
285     return a;
286   }
287
288   return static_cast<TObjArray*>(m.At(idx));
289 }
290
291 //____________________________________________________________________
292 Bool_t
293 AliFMDCorrAcceptance::SetCorrection(UShort_t d, Char_t r, 
294                                     UShort_t b, TH2D*  h) 
295 {
296   // 
297   // Set the acceptance correction @f$ a_{r,v}(\eta)@f$ 
298   // Note, that the object takes ownership of the passed pointer.
299   // 
300   // Parameters:
301   //    d    Detector number (1-3)
302   //    r    Ring identifier (I or O)
303   //    b    Bin corresponding to the primary interaction point 
304   //             @f$z@f$ coordinate  (1 based)
305   //    h    @f$ a_{r,v}(\eta)@f$ 
306   // 
307   // Return:
308   //    true if operation succeeded 
309   //
310   TObjArray* ringArray = GetOrMakeRingArray(fRingArray, d, r);
311   if (!ringArray) return false;
312   
313   if (b <= 0 || b > fVertexAxis.GetNbins()) { 
314     AliWarning(Form("Vertex bin %3d out of range [1,%3d]", 
315                     b, fVertexAxis.GetNbins()));
316     return false;
317   }
318   h->SetName(Form("FMD%d%c_vtxbin%03d", d, r, b));
319   h->SetTitle(Form("Acceptance correction for FMD%d%c "
320                    "in vertex bin %d [%+8.4f,%+8.4f]", 
321                    d, r, b, fVertexAxis.GetBinLowEdge(b), 
322                    fVertexAxis.GetBinUpEdge(b)));
323   h->SetXTitle("#eta");
324   h->SetYTitle("N_{strips,OK}/N_{strips}");
325   h->SetFillStyle(3001);
326   h->SetDirectory(0);
327   h->SetStats(0);
328   ringArray->AddAtAndExpand(h, b-1);
329   return kTRUE;
330 }
331 //____________________________________________________________________
332 Bool_t
333 AliFMDCorrAcceptance::SetCorrection(UShort_t d, Char_t r, 
334                                     Double_t v, TH2D*  h) 
335 {
336   // 
337   // Set the acceptance correction @f$ a_{r,v}(\eta)@f$.
338   // Note, that the object takes ownership of the passed pointer.
339   // 
340   // Parameters:
341   //    d    Detector number (1-3)
342   //    r    Ring identifier (I or O)
343   //    v    Primary interaction point @f$z@f$ coordinate  
344   //    h    @f$ a_{r,v}(\eta)@f$ 
345   // 
346   // Return:
347   //    true if operation succeeded 
348   //
349   Int_t b = FindVertexBin(v);
350   if (b <= 0 || b > fVertexAxis.GetNbins()) { 
351     AliWarning(Form("Vertex %+8.4f out of range [%+8.4f,%+8.4f]", 
352                     v, fVertexAxis.GetXmin(), fVertexAxis.GetXmax()));
353     return false;
354   }
355   return SetCorrection(d, r, UShort_t(b), h);
356 }
357
358 //____________________________________________________________________
359 void
360 AliFMDCorrAcceptance::FillCache() const
361 {
362   if (fCache) return;
363
364   fCache = new TObjArray;
365   fCache->SetOwner(kTRUE);
366   fCache->SetName("cache");
367
368   Int_t nV = fVertexAxis.GetNbins();
369   for (UShort_t v = 1; v <= nV; v++) {
370     for(UShort_t d = 1; d <= 3;d++) { 
371       UShort_t nR = (d == 1 ? 1 : 2);
372       for (UShort_t q = 0; q < nR; q++) { 
373         Char_t   r  = (q == 0 ? 'I' : 'O');
374
375         TObjArray* a = GetOrMakeRingArray(*fCache, d, r);
376
377         TH2D* corr = GetCorrection(d, r, v);
378         if (!corr) continue;
379
380         Int_t nY = corr->GetNbinsY();
381         TH1D* h  = corr->ProjectionX("tmp", nY+1, nY+1, "");
382         h->SetName(Form("FMD%d%c_vtxbin%03d", d, r, v));
383         h->SetTitle(Form("#phi acceptance correction for FMD%d%c "
384                            "in vertex bin %d [%+8.4f,%+8.4f]", 
385                            d, r, v, fVertexAxis.GetBinLowEdge(v), 
386                            fVertexAxis.GetBinUpEdge(v)));
387         h->SetXTitle("#eta");
388         h->SetYTitle("N_{strips}/N_{strips,OK}");
389         h->SetFillStyle(3001);
390         h->SetDirectory(0);
391         h->SetStats(0);
392         a->AddAtAndExpand(h,v-1);
393
394         if (fHasOverflow) continue;
395
396         // Construct the overflow bin from 
397         Int_t nX = corr->GetNbinsX();
398         for (Int_t eta = 1; eta <= nX; eta++) { 
399           Double_t sum = 0;
400           for (Int_t phi = 1; phi <= nY; phi++) 
401             sum += corr->GetBinContent(eta, phi);
402           if (nY <= 0) continue;
403           h->SetBinContent(eta, nY/sum);
404         } // for eta 
405       } // for q 
406     } // for d 
407   } // for v 
408 }
409 //____________________________________________________________________
410 void
411 AliFMDCorrAcceptance::Browse(TBrowser* b)
412 {
413   // 
414   // Browse this object in the browser
415   // 
416   // Parameters:
417   //    b 
418   //
419   b->Add(&fRingArray);
420   if (fCache) b->Add(fCache);
421   b->Add(&fVertexAxis);
422 }
423 //____________________________________________________________________
424 void
425 AliFMDCorrAcceptance::Print(Option_t* option) const
426 {
427   // 
428   // Print this object 
429   // 
430   // Parameters:
431   //    option 
432   //  
433   std::cout << "Acceptance correction due to dead channels" << std::endl;
434   fRingArray.Print(option);
435   fVertexAxis.Print(option);
436 }
437     
438 //____________________________________________________________________
439 //
440 // EOF
441 //