]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliCentralCorrAcceptance.cxx
Fix histogram names in central corrections task.
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliCentralCorrAcceptance.cxx
1 //
2 // This class contains the acceptance correction due to dead channels 
3 // 
4 //
5 #include "AliCentralCorrAcceptance.h"
6 #include <TBrowser.h>
7 #include <TH1D.h>
8 #include <AliLog.h>
9 #include <iostream>
10
11 //____________________________________________________________________
12 AliCentralCorrAcceptance::AliCentralCorrAcceptance()
13   : fArray(), 
14     fVertexAxis(0,0,0)
15 {
16   // 
17   // Default constructor 
18   //
19   fArray.SetOwner(kTRUE);
20   fArray.SetName("acceptance");
21   fVertexAxis.SetName("vtxAxis");
22   fVertexAxis.SetTitle("v_{z} [cm]");
23   
24 }
25 //____________________________________________________________________
26 AliCentralCorrAcceptance::AliCentralCorrAcceptance(const 
27                                                AliCentralCorrAcceptance& o)
28   : TObject(o), 
29     fArray(o.fArray), 
30     fVertexAxis(o.fVertexAxis.GetNbins(), o.fVertexAxis.GetXmin(), 
31                 o.fVertexAxis.GetXmax())
32 {
33   // 
34   // Copy constructor 
35   // 
36   // Parameters:
37   //    o Object to copy from 
38   //
39   fVertexAxis.SetName("vtxAxis");
40   fVertexAxis.SetTitle("v_{z} [cm]");
41 }
42 //____________________________________________________________________
43 AliCentralCorrAcceptance::~AliCentralCorrAcceptance()
44 {
45   //
46   // Destructor 
47   // 
48   //
49   fArray.Clear();
50 }
51 //____________________________________________________________________
52 AliCentralCorrAcceptance&
53 AliCentralCorrAcceptance::operator=(const AliCentralCorrAcceptance& o)
54 {
55   // 
56   // Assignment operator 
57   // 
58   // Parameters:
59   //    o Object to assign from 
60   // 
61   // Return:
62   //    Reference to this object 
63   //
64   fArray        = o.fArray;
65   SetVertexAxis(o.fVertexAxis);
66
67   return *this;
68 }
69 //____________________________________________________________________
70 TH1D*
71 AliCentralCorrAcceptance::GetCorrection(Double_t v) const
72 {
73   // 
74   // Get the acceptance correction @f$ a_{r,v}@f$ 
75   // 
76   // Parameters:
77   //    d  Detector number (1-3)
78   //    r  Ring identifier (I or O)
79   //    v  Primary interaction point @f$z@f$ coordinate
80   // 
81   // Return:
82   //    The correction @f$ a_{r,v}@f$ 
83   //
84   Int_t b = FindVertexBin(v);
85   if (b <= 0) return 0;
86   return GetCorrection(UShort_t(b));
87 }
88 //____________________________________________________________________
89 TH1D*
90 AliCentralCorrAcceptance::GetCorrection(UShort_t b) const
91 {
92   // 
93   // Get the acceptance correction @f$ a_{r,v}@f$ 
94   // 
95   // Parameters:
96   //    d  Detector number (1-3)
97   //    r  Ring identifier (I or O)
98   //    b  Bin corresponding to the primary interaction point 
99   //           @f$z@f$ coordinate (1 based)
100   // 
101   // Return:
102   //    The correction @f$ a_{r,v}@f$ 
103   //
104   
105   // TObjArray* ringArray = GetRingArray(d, r);
106   //if (!ringArray) return 0;
107
108   if (b <= 0 || b > fArray.GetEntriesFast()) {
109     AliWarning(Form("vertex bin %d out of range [1,%d]", 
110                     b, fArray.GetEntriesFast()));
111     return 0;
112   }
113
114   TObject* o = fArray.At(b-1);
115   if (!o) { 
116     AliWarning(Form("No dead channels map found for SPD in vertex bin %d",
117                     b));
118     return 0;
119   }
120   return static_cast<TH1D*>(o);
121 }
122   
123 //____________________________________________________________________
124 Int_t
125 AliCentralCorrAcceptance::FindVertexBin(Double_t v) const
126 {
127   // 
128   // Find the vertex bin that corresponds to the passed vertex 
129   // 
130   // Parameters:
131   //    vertex The interaction points @f$z@f$-coordinate 
132   // 
133   // Return:
134   //    Vertex bin in @f$[1,N_{\mbox{vertex}}]@f$ or negative if 
135   // out of range 
136   //
137   if (fVertexAxis.GetNbins() <= 0) { 
138     AliWarning("No vertex array defined");
139     return 0;
140   }
141   Int_t bin = const_cast<TAxis&>(fVertexAxis).FindBin(v);
142   if (bin <= 0 || bin > fVertexAxis.GetNbins()) { 
143     AliWarning(Form("vertex %+8.4f out of range [%+8.4f,%+8.4f]",
144                     v, fVertexAxis.GetXmin(), fVertexAxis.GetXmax()));
145     return 0;
146   }
147   return bin;
148 }
149
150 //____________________________________________________________________
151 Bool_t
152 AliCentralCorrAcceptance::SetCorrection(UShort_t b, TH1D*  h) 
153 {
154   // 
155   // Set the acceptance correction @f$ a_{r,v}(\eta)@f$ 
156   // Note, that the object takes ownership of the passed pointer.
157   // 
158   // Parameters:
159   //    b    Bin corresponding to the primary interaction point 
160   //             @f$z@f$ coordinate  (1 based)
161   //    h    @f$ a_{r,v}(\eta)@f$ 
162   // 
163   // Return:
164   //    true if operation succeeded 
165   //
166   
167   if (b <= 0 || b > fVertexAxis.GetNbins()) { 
168     AliWarning(Form("Vertex bin %3d out of range [1,%3d]", 
169                     b, fVertexAxis.GetNbins()));
170     return false;
171   }
172   
173   h->SetName(Form("acc_vtxbin%03d",  b));
174   h->SetTitle(Form("Acceptance correction [%+5.1f<v_{z}<%+5.1f]", 
175                    fVertexAxis.GetBinLowEdge(b), 
176                    fVertexAxis.GetBinUpEdge(b)));
177   h->SetXTitle("#eta");
178   h->SetYTitle("dN_{ch}/d#eta / sum_i N_{ch,i}");
179   h->SetFillStyle(3001);
180   h->SetDirectory(0);
181   h->SetStats(0);
182   fArray.AddAtAndExpand(h, b-1);
183   return kTRUE;
184 }
185 //____________________________________________________________________
186 Bool_t
187 AliCentralCorrAcceptance::SetCorrection(Double_t v, TH1D*  h) 
188 {
189   // 
190   // Set the acceptance correction @f$ a_{r,v}(\eta)@f$.
191   // Note, that the object takes ownership of the passed pointer.
192   // 
193   // Parameters
194   //    v    Primary interaction point @f$z@f$ coordinate  
195   //    h    @f$ a_{r,v}(\eta)@f$ 
196   // 
197   // Return:
198   //    true if operation succeeded 
199   //
200   Int_t b = FindVertexBin(v);
201   if (b <= 0 || b > fVertexAxis.GetNbins()) { 
202     AliWarning(Form("Vertex %+8.4f out of range [%+8.4f,%+8.4f]", 
203                     v, fVertexAxis.GetXmin(), fVertexAxis.GetXmax()));
204     return false;
205   }
206   return SetCorrection(UShort_t(b), h);
207 }
208 //____________________________________________________________________
209 void
210 AliCentralCorrAcceptance::Browse(TBrowser* b)
211 {
212   // 
213   // Browse this object in the browser
214   // 
215   // Parameters:
216   //    b 
217   //
218   b->Add(&fArray);
219   b->Add(&fVertexAxis);
220 }
221 //____________________________________________________________________
222 void
223 AliCentralCorrAcceptance::Print(Option_t* /* option */) const
224 {
225   // 
226   // Print this object 
227   // 
228   // Parameters:
229   //    option 
230   //  
231   std::cout << "  Acceptance correction due to dead channels" 
232             << "   # of vertex bins: "  << fVertexAxis.GetNbins() << "\n"
233             << "   Vertex range:     [" << fVertexAxis.GetXmin() 
234             << "," << fVertexAxis.GetXmax() << "]\n" 
235             << "   Histograms:\n"
236             << "    ";
237   TIter next(&fArray);
238   TObject* o = 0;
239   while ((o = next())) std::cout << o->GetName() << " ";
240   std::cout << std::endl;
241 }
242     
243 //____________________________________________________________________
244 //
245 // EOF
246 //