]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDCalibDrawer.cxx
monitor dQdl (PID) at TRD entrance
[u/mrichter/AliRoot.git] / FMD / AliFMDCalibDrawer.cxx
1 #include <AliFMDCalibDrawer.h>
2 #include <TH1.h>
3 #include <TH2.h>
4 #include <TPad.h>
5 #include <TCanvas.h>
6 #include <AliFMDParameters.h>
7 #include <AliCDBManager.h>
8 #include <AliLog.h>
9 // #include <AliFMDCalibPedestal.h>
10 // #include <AliFMDCalibSampleRate.h>
11 // #include <AliFMDCalibStripRange.h>
12
13 //____________________________________________________________________
14 void
15 AliFMDCalibDrawer::Init(Int_t runNo, const char* ocdb)
16 {
17   AliCDBManager* cdb = AliCDBManager::Instance();
18   cdb->SetRun(runNo);
19   if (ocdb && ocdb[0] != '\0') cdb->SetDefaultStorage(ocdb);
20
21   AliFMDParameters* pars = AliFMDParameters::Instance();
22   pars->Init(kTRUE);
23 }
24
25 //____________________________________________________________________
26 Double_t
27 AliFMDCalibDrawer::GetHistMax(EWhat what) const
28 {
29   switch (what) {
30   case kPedestal:  return 150;
31   case kNoise:     return 5;
32   case kGain:      return 4;
33   case kDead:      return 1.5;
34   case kRate:      return 5;
35   case kRange:     return 128;
36   case kZeroSuppression:  return 10;
37   }
38   return 1024;
39 }
40 //____________________________________________________________________
41 Double_t
42 AliFMDCalibDrawer::GetHistMin(EWhat what) const
43 {
44   switch (what) {
45   case kPedestal:  return 0;
46   case kNoise:     return 0;
47   case kGain:      return 0;
48   case kDead:      return -.5;
49   case kRate:      return 0;
50   case kRange:     return -1;
51   case kZeroSuppression:  return -1;
52   }
53   return -1;
54 }
55   
56 //____________________________________________________________________
57 const char*
58 AliFMDCalibDrawer::GetHistName(EWhat what) const
59 {
60   switch (what) {
61   case kPedestal:  return "peds"; 
62   case kNoise:     return "noise";
63   case kGain:      return "gain";
64   case kDead:      return "dead";
65   case kRate:      return "rate";
66   case kRange:     return "range";
67   case kZeroSuppression:  return "thrs";
68   }
69   return "unknown";
70 }
71 //____________________________________________________________________
72 const char*
73 AliFMDCalibDrawer::GetHistTitle(EWhat what) const
74 {
75   switch (what) {
76   case kPedestal:  return "Pedestal & noise"; 
77   case kNoise:     return "Noise";
78   case kGain:      return "Gain";
79   case kDead:      return "Dead map";
80   case kRate:      return "Sample rate";
81   case kRange:     return "Strip range";
82   case kZeroSuppression:  return "ZS threshold";
83   }
84   return "Unknown";
85 }
86   
87 //____________________________________________________________________
88 Int_t
89 AliFMDCalibDrawer::GetRingColor(UShort_t d, Char_t r) const
90 {
91   return ((d == 1 ? kRed : (d == 2 ? kGreen : kBlue))
92           + ((r == 'I' || r == 'i') ? 2 : -3));
93 }
94 //____________________________________________________________________
95 void
96 AliFMDCalibDrawer::SetAttributes(TH1* ret, EWhat what, 
97                                  UShort_t d, Char_t r) const
98 {
99   ret->SetFillColor(GetRingColor(d,r));
100   ret->SetMarkerColor(GetRingColor(d,r));
101   ret->SetLineColor(GetRingColor(d,r));
102   ret->SetFillStyle(3001);
103   ret->SetMarkerStyle(20);
104   ret->SetLineStyle(1);
105   ret->SetStats(0);
106   ret->SetDirectory(0);
107   ret->SetMinimum(GetHistMin(what));
108   ret->SetMaximum(GetHistMax(what));
109 }
110 //____________________________________________________________________
111 TH2D*
112 AliFMDCalibDrawer::Make2D(EWhat what, UShort_t d, Char_t r) const
113 {
114   UShort_t q  = (r == 'I' || r == 'i') ? 0 : 1;
115   UShort_t nY = q == 0 ?  20 :  40;
116   UShort_t nX = q == 0 ? 512 : 256;
117   TString  n(Form("%s_FMD%d%c", GetHistName(what), d, r));
118   TString  t(Form("%s for FMD%d%c", GetHistTitle(what),  d, r));
119   TH2D*    ret = new TH2D(n, t, nX, 0, nX, nY, 0, nY);
120   ret->SetXTitle("Strip #");
121   ret->SetYTitle("Sector #");
122   ret->SetZTitle(GetHistTitle(what));
123   ret->GetXaxis()->SetNdivisions(1600 + (q == 0 ? 4 : 2),false);
124   ret->GetYaxis()->SetNdivisions(nY, false);
125   SetAttributes(ret, what, d, r);
126   return ret;
127
128 }
129
130 //____________________________________________________________________
131 TH1D*
132 AliFMDCalibDrawer::Make1D(EWhat what, UShort_t d, Char_t r, UShort_t s) const  
133 {
134   UShort_t q    = (r == 'I' || r == 'i') ? 0 : 1;
135   UShort_t nStr = (r == 'I' || r == 'i') ? 512 : 256;
136   TString  n(Form("%s_FMD%d%c_%02d", GetHistName(what), d, r, s));
137   TString  t(Form("%s for FMD%d%c[%2d]", GetHistTitle(what), d, r, s));
138   TH1D* ret = new TH1D(n, t, nStr, 0, nStr);
139   ret->SetXTitle("Strip #");
140   ret->SetYTitle(GetHistTitle(what));
141   ret->GetXaxis()->SetNdivisions(1600 + (q == 0 ? 4 : 2),false);
142   SetAttributes(ret, what, d, r);
143   return ret;
144 }
145
146 //____________________________________________________________________
147 void
148 AliFMDCalibDrawer::GetNumber(EWhat what, UShort_t d, Char_t r, UShort_t s, 
149                              UShort_t t, Double_t& val, Double_t& err) const
150 {
151   AliFMDParameters* pars = AliFMDParameters::Instance();
152   err = 0;
153   switch (what) {
154   case kPedestal:
155     val = pars->GetPedestal(d, r, s, t);
156     err = pars->GetPedestalWidth(d, r, s, t);
157     break;
158   case kNoise: val = pars->GetPedestalWidth(d, r, s, t); break;
159   case kGain:  val = pars->GetPulseGain(d, r, s, t); break;
160   case kDead:  val = pars->IsDead(d, r, s, t) ? 0 : 1; break;
161   case kRate:  val = pars->GetSampleRate(d, r, s, t); break;
162   case kRange:
163     val = pars->GetMaxStrip(d,r,s,t);
164     err = pars->GetMinStrip(d,r,s,t);
165     break;
166   case kZeroSuppression:  val = pars->GetZeroSuppression(d,r,s,t); break;
167   default: 
168     AliError(Form("Unknown quantity: %d - nothing returned", what));
169     break;
170   }
171 }
172
173 //____________________________________________________________________
174 TH1*
175 AliFMDCalibDrawer::FillRing(EWhat what, UShort_t d, Char_t r) const
176 {
177   if (what == kRate || what == kRange) { 
178     TH1D* ret = new TH1D(Form("%s_FMD%d%c", GetHistName(what), d, r), 
179                          Form("%s for FMD%d%c", GetHistTitle(what), d, r), 
180                          2, -.5, 1.5);
181     UShort_t nSec = (r == 'I' || r == 'i') ? 20 : 40;
182     ret->GetXaxis()->SetBinLabel(1, Form("Top [%02d-%02d]", 0,      nSec/2-1));
183     ret->GetXaxis()->SetBinLabel(2, Form("Bottom [%02d-%02d]",nSec/2-1,nSec-1));
184     ret->SetYTitle(GetHistTitle(what));
185     SetAttributes(ret, what, d, r);
186     if (what == kRate) { 
187       ret->SetBarWidth(0.8);
188       ret->SetBarOffset(0.1);
189     }
190     else 
191       ret->SetMarkerSize(0);
192     for (UShort_t s = 0; s < ret->GetNbinsX(); s++) { 
193       Double_t val; 
194       Double_t err; 
195       GetNumber(what, d, r, s*(nSec/2), 0, val, err);
196       if (what == kRange) { 
197         Double_t min = err;
198         err          = (val - err) / 2;
199         val          = min + err;
200       }
201       ret->SetBinContent(s+1, val);
202       ret->SetBinError(s+1, err);
203     }
204     return ret;
205   }
206   if (what == kZeroSuppression) { 
207     UShort_t nY  = (r == 'I' || r == 'i') ?  20 :  40;
208     UShort_t nX  = ((r == 'I' || r == 'i') ? 512 : 256) / 128;
209     TH2D*    ret = new TH2D(Form("%s_FMD%d%c", GetHistName(what), d, r), 
210                             Form("%s for FMD%d%c", GetHistTitle(what), d, r),
211                             nX, 0, nX, nY, 0, nY);
212     ret->SetXTitle("Channel #");
213     ret->SetYTitle("Sector #");
214     ret->SetZTitle(GetHistTitle(what));
215     ret->GetXaxis()->SetNdivisions(nX, false);
216     ret->GetYaxis()->SetNdivisions(nY, false);
217     SetAttributes(ret, what, d, r);
218     for (UShort_t s = 0; s < ret->GetNbinsY(); s++) { 
219       for (UShort_t t = 0; t < ret->GetNbinsX(); t++) { 
220         Double_t val; 
221         Double_t err; 
222         GetNumber(what, d, r, s, t*128, val, err);
223         ret->SetBinContent(t+1, s+1, val);
224         ret->SetBinError(t+1, s+1, err);
225       }
226     }
227     return ret;
228   }
229     
230   TH2D* ret = Make2D(what, d, r);
231   for (UShort_t s = 0; s < ret->GetNbinsY(); s++) { 
232     for (UShort_t t = 0; t < ret->GetNbinsX(); t++) { 
233       Double_t val; 
234       Double_t err; 
235       GetNumber(what, d, r, s, t, val, err);
236       ret->SetBinContent(t+1, s+1, val);
237       ret->SetBinError(t+1, s+1, err);
238     }
239   }
240   return ret;
241 }
242 //____________________________________________________________________
243 TH1*
244 AliFMDCalibDrawer::FillSector(EWhat what,UShort_t d, Char_t r, UShort_t s) const
245 {
246   if (what == kZeroSuppression) { 
247     UShort_t nX  = ((r == 'I' || r == 'i') ? 512 : 256) / 128;
248     TH1D*    ret = new TH1D(Form("%s_FMD%d%c[%02d]", 
249                                  GetHistName(what), d, r, s), 
250                             Form("%s for FMD%d%c[%02d]", 
251                                  GetHistTitle(what), d, r, s),
252                             nX, 0, nX);
253     ret->SetXTitle("Channel #");
254     ret->SetYTitle(GetHistTitle(what));
255     SetAttributes(ret, what, d, r);
256     ret->GetXaxis()->SetNdivisions(nX, false);
257     for (UShort_t t = 0; t < ret->GetNbinsX(); t++) { 
258       Double_t val; 
259       Double_t err; 
260       GetNumber(what, d, r, s, t*128, val, err);
261       ret->SetBinContent(t+1, s+1, val);
262       if (err >= 0) ret->SetBinError(t+1, s+1, err);
263     }
264     return ret;
265   }
266
267   TH1D* ret = Make1D(what, d, r, s);
268   for (UShort_t t = 0; t < ret->GetNbinsX(); t++) { 
269     Double_t val; 
270     Double_t err; 
271     GetNumber(what, d, r, s, t, val, err);
272     ret->SetBinContent(t+1, s+1, val);
273     if (err >= 0) ret->SetBinError(t+1, s+1, err);
274   }
275   return ret;
276 }
277
278 //____________________________________________________________________
279 void
280 AliFMDCalibDrawer::DrawOne(EWhat what, Short_t d, Char_t r, 
281                            Short_t s, Short_t /*t*/) const
282 {
283   Info("DrawOne", "Will draw %s for D=%d, R=%c, S=%d",
284        GetHistTitle(what), d, (r == '\0' ? '-' : r), s);
285   TVirtualPad* tp = gPad;
286   if (!tp) tp = TCanvas::MakeDefCanvas();
287   if (d <= 0 || d > 3) { 
288     // All detectors - need to split current pad
289     tp->Divide(1,3,0,0);
290     for (UShort_t dd = 1; dd <= 3; dd++) { 
291       tp->cd(dd);
292       // Call our selves 
293       DrawOne(what, dd, '\0', -1, -1);
294     }
295   }
296   else if ((r != 'I' && r != 'i') && (r != 'O' && r != 'o')) { 
297     // All rings in a detector.  Split current pad in two 
298     tp->SetName(Form("FMD%d", d));
299     tp->Divide(2,1,0,0);
300     for (UShort_t q = 0; q < 2; q++) { 
301       tp->cd(q+1);
302       Char_t       rr = q == 0 ? 'I' : 'O';
303       // Call ourselves 
304       DrawOne(what, d, rr, -1, -1);
305     }
306   }
307   else if (what == kRate || what == kRange || s < 0 || s > 39) {
308     // A single ring - and special case for sample rate and strip range
309     tp->SetName(Form("FMD%d%c", d, r));
310     tp->SetFillColor(0);
311     tp->SetFillStyle(0);
312     tp->SetRightMargin(0.10);
313     if (d == 1 && (r == 'O' || r == 'o')) return;
314     TH1* h = FillRing(what, d, r);
315     switch (what) { 
316     case kRate:  h->Draw("bar"); break;
317     case kRange: h->Draw("e3"); break;
318     default: h->Draw("colz"); break;
319     }
320   }
321   else {
322     // Finally, we're down to a single sector.  In this case, we just do
323     // a single 1D histogram 
324     TH1* h = FillSector(what, d, r, s);
325     h->Draw("hist e");
326   }
327   tp->Update();
328   tp->Modified();
329   tp->cd();
330 }
331
332
333
334 //____________________________________________________________________
335 //
336 // EOF
337 //