]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONTrackerQAChecker.cxx
Added new static method AliAnalysisManager::GetRunFromAlienPath() that extracts the...
[u/mrichter/AliRoot.git] / MUON / AliMUONTrackerQAChecker.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 // $Id$
17
18 #include "AliMUONTrackerQAChecker.h"
19
20 /// \class AliMUONTrackerQAChecker
21 ///
22 /// Implementation of AliQACheckerBase for MCH and MTR
23 ///
24 /// For the moment we only implement the checking of raw data QA for the tracker
25 /// by looking at the occupancy at the bus patch level.
26 ///
27 /// \author Laurent Aphecetche, Subatech
28
29 #include "AliCDBManager.h"
30 #include "AliCodeTimer.h"
31 #include "AliLog.h"
32 #include "AliMUONQAIndices.h"
33 #include "AliMUONRecoParam.h"
34 #include "AliMpBusPatch.h"
35 #include "AliMpDDLStore.h"
36 #include "AliQAv1.h"
37 #include "Riostream.h"
38 #include "TAxis.h"
39 #include "TDirectory.h"
40 #include "TGaxis.h"
41 #include "TH1.h"
42 #include "TLine.h"
43 #include "TMath.h"
44 #include "TPaveText.h"
45 #include "TVirtualPad.h"
46
47 /// \cond CLASSIMP
48 ClassImp(AliMUONTrackerQAChecker)
49 /// \endcond
50
51 namespace {
52   
53   //___________________________________________________________________________  
54   int trim(Int_t n, 
55            Double_t* x,
56            Double_t alpha,
57            Double_t& tmean,
58            Double_t& tvar,
59            Double_t& min,
60            Double_t& max)
61   {
62     //
63     // Calculates the trimmed (tmean) mean
64     // of a sample (x) and estimates the variance (tvar)
65     // of that mean.
66     //
67     
68     // First check input parameters
69     
70     // number of observations
71     if ( n < 2 )
72     {
73       return -1;
74     }
75     
76     if ( alpha < 0 || alpha >= 0.5 )
77       // proportion of observations
78       // to be trimmed at each end of the sorted sample
79     {
80       return -2;
81     }
82     
83     // Input parameters are good. Let's move on.
84     
85     // Insure we use a sample sorted into ascending order.
86     
87     Int_t* indices = new Int_t[n];
88     
89     TMath::Sort(n,x,indices,kFALSE);
90     
91     Double_t* sx = new Double_t[n];
92     
93     for ( Int_t i = 0; i < n; ++i )
94     {
95       sx[i] = x[indices[i]];
96     }
97     delete[] indices;
98     
99     
100     // Number of observations trimmed at each end.
101     
102     Int_t k = TMath::FloorNint(alpha * n);
103     
104     double sum = 0.0;
105     
106     for ( Int_t i = k; i < n - k ; ++i )
107     {
108       sum += sx[i];
109     }
110     
111     tmean = sum / ( n - 2 * k );
112     
113     double t2 = 0.0;
114     
115     for ( Int_t i = k; i < n - k; ++i )
116     {
117       t2 += (sx[i] - tmean) * (sx[i] - tmean);
118     }
119     
120     tvar = (
121             t2 +
122             k * (sx[k] - tmean) * (sx[k] - tmean) +
123             k * (sx[n - k - 1] - tmean) * (sx[n - k - 1] - tmean)
124             ) / (n * n);
125     
126     // get the min and max for the non-rejected values
127     min = DBL_MAX;
128     max = 0.0;
129     
130     for ( Int_t i = k; i < n-k; ++i ) 
131     {
132       min = TMath::Min(min,sx[i]);
133       max = TMath::Max(max,sx[i]);
134     }
135     
136     delete[] sx;
137     
138     return 0;
139   }
140
141   //___________________________________________________________________________  
142   Int_t GetColorFromCheckCode(AliMUONVQAChecker::ECheckCode code)
143   {
144     const Int_t INFOCOLOR(kGreen);  // green = INFO
145     const Int_t WARNINGCOLOR(kYellow); // yellow = WARNING
146     const Int_t ERRORCOLOR(kOrange); // orange = ERROR
147     const Int_t FATALCOLOR(kRed); // red = FATAL
148     
149     if ( code == AliMUONVQAChecker::kInfo ) return INFOCOLOR;
150     else if ( code == AliMUONVQAChecker::kWarning ) return WARNINGCOLOR;
151     else if ( code ==  AliMUONVQAChecker::kFatal) return FATALCOLOR;
152     else return ERRORCOLOR;
153   }
154   
155   const char* NOTENOUGHEVENTMESSAGE = "Not enough event to judge. Please wait a bit";
156   const char* NOTIFYEXPERTMESSAGE = "PLEASE NOTIFY EXPERT !";
157   const char* ALLISFINEMESSAGE = "All is fine. Just enjoy.";
158   
159   const int NOTENOUGHEVENTLIMIT = 50;
160   
161   //___________________________________________________________________________  
162   void SetupHisto(Int_t neventsseen, Int_t neventsused, const TObjArray& messages, TH1& histo, AliMUONVQAChecker::ECheckCode code, const char* extraopt="")
163   {
164     Bool_t allIsFine(kFALSE);
165     
166     if ( code == AliMUONVQAChecker::kInfo ) 
167     {
168       allIsFine = kTRUE;
169     }    
170     
171     Double_t y1 = 0.99 - (messages.GetLast()+(allIsFine?1:0)+2)*0.075;
172     
173     TPaveText* text = new TPaveText(0.5,y1,0.99,0.99,"NDC");
174     
175     text->AddText(Form("MCH RUN %d - %d events seen - %d events used",AliCDBManager::Instance()->GetRun(),neventsseen,neventsused));
176     
177     TIter next(&messages);
178     TObjString* str;
179     
180     while ( ( str = static_cast<TObjString*>(next()) ) )
181     {
182       text->AddText(str->String());
183     }
184     
185     if ( allIsFine )
186     {
187       text->AddText(ALLISFINEMESSAGE);
188     }
189                     
190     text->SetFillColor(GetColorFromCheckCode(code));
191                       
192     Int_t color = GetColorFromCheckCode(code);
193     
194     histo.SetFillStyle(1001);
195     histo.SetFillColor(color);    
196     TString sopt("BAR");
197     sopt += extraopt;
198     histo.SetOption(sopt.Data());
199     
200     histo.SetTitle(kFALSE);
201     
202     TObject* title = histo.GetListOfFunctions()->FindObject("title");
203     if (title) title->Delete();
204
205     histo.GetListOfFunctions()->Add(text);
206   }
207   
208 }
209
210 //__________________________________________________________________
211 AliMUONTrackerQAChecker::AliMUONTrackerQAChecker() : AliMUONVQAChecker()
212 {
213         /// ctor
214 }          
215
216 //__________________________________________________________________
217 AliMUONTrackerQAChecker::~AliMUONTrackerQAChecker() 
218 {
219         /// dtor
220 }
221
222 //__________________________________________________________________
223 AliMUONTrackerQAChecker::AliMUONTrackerQAChecker(const AliMUONTrackerQAChecker& qac) : 
224     AliMUONVQAChecker(qac) 
225 {
226         /// copy ctor 
227 }   
228
229 //______________________________________________________________________________
230 AliMUONVQAChecker::ECheckCode*
231 AliMUONTrackerQAChecker::CheckRecPoints(TObjArray ** list, const AliMUONRecoParam* /*recoParam*/)
232 {
233   /// Check rec points
234   /// Very binary check for the moment. 
235   
236   AliCodeTimerAuto("",0);
237   
238   AliMUONVQAChecker::ECheckCode * rv = new AliMUONVQAChecker::ECheckCode[AliRecoParam::kNSpecies] ; 
239   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) 
240     rv[specie] =  AliMUONVQAChecker::kInfo; 
241   
242   
243   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) 
244   {
245     TH1* h = AliQAv1::GetData(list,AliMUONQAIndices::kTrackerNumberOfClustersPerDE,AliRecoParam::ConvertIndex(specie));
246
247     if ( !h ) rv[specie] = AliMUONVQAChecker::kWarning; // only a warning if histo not found, in order not to kill anything because QA failed...
248   
249     else if ( h->GetMean() == 0.0 ) rv[specie] =  MarkHisto(*h,AliMUONVQAChecker::kFatal);
250   }
251   return rv;
252 }
253
254 //______________________________________________________________________________
255 AliMUONVQAChecker::ECheckCode 
256 AliMUONTrackerQAChecker::MarkHisto(TH1& histo, AliMUONVQAChecker::ECheckCode value) const
257 {
258   /// Mark histo as originator of some QA error/warning
259   
260   if ( value != AliMUONVQAChecker::kInfo )
261   {
262     histo.SetBit(AliQAv1::GetQABit());
263   }
264   
265   return value;
266 }
267
268 //______________________________________________________________________________
269 AliMUONVQAChecker::ECheckCode*
270 AliMUONTrackerQAChecker::CheckESD(TObjArray ** list, const AliMUONRecoParam* /*recoParam*/)
271 {
272   /// Check ESD
273   
274   AliCodeTimerAuto("",0);
275   
276   AliMUONVQAChecker::ECheckCode * rv = new AliMUONVQAChecker::ECheckCode[AliRecoParam::kNSpecies] ; 
277   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) 
278     rv[specie] = AliMUONVQAChecker::kInfo; 
279   
280   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
281     
282     TH1* h = AliQAv1::GetData(list,AliMUONQAIndices::kESDnTracks,AliRecoParam::ConvertIndex(specie));
283   
284     if (!h) rv[specie] = AliMUONVQAChecker::kWarning;
285   
286     else if ( h->GetMean() == 0.0 ) rv[specie] =  MarkHisto(*h,AliMUONVQAChecker::kFatal); // no track -> fatal
287   
288     h = AliQAv1::GetData(list,AliMUONQAIndices::kESDMatchTrig,AliRecoParam::ConvertIndex(specie));
289   
290     if (!h) rv[specie] = AliMUONVQAChecker::kWarning;
291   
292     else if (h->GetMean() == 0.0 ) rv[specie] = MarkHisto(*h,AliMUONVQAChecker::kError); // no trigger matching -> error
293   }
294   return rv;
295 }
296
297 //______________________________________________________________________________
298 AliMUONVQAChecker::ECheckCode*
299 AliMUONTrackerQAChecker::CheckRaws(TObjArray ** list, const AliMUONRecoParam* recoParam)
300 {
301   /// Check raws
302
303   AliCodeTimerAuto("",0);
304   
305   if (!recoParam) return 0x0;
306   
307   AliMUONVQAChecker::ECheckCode * rv = new AliMUONVQAChecker::ECheckCode[AliRecoParam::kNSpecies] ; 
308
309   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) 
310   {
311     rv[specie] = AliMUONVQAChecker::kInfo; 
312   }
313   
314   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) 
315   {
316     TH1* hneventsseen = AliQAv1::GetData(list,AliMUONQAIndices::kTrackerNofPhysicsEventsSeen,AliRecoParam::ConvertIndex(specie));
317     TH1* hneventsused = AliQAv1::GetData(list,AliMUONQAIndices::kTrackerNofGoodPhysicsEventsUsed,AliRecoParam::ConvertIndex(specie));
318
319     if (!hneventsseen || !hneventsused ) continue;
320     
321     Int_t neventsseen = TMath::Nint(hneventsseen->GetBinContent(1));
322     Int_t neventsused = TMath::Nint(hneventsused->GetBinContent(1));
323         
324     if ( !hneventsseen || !hneventsused )
325     {
326       continue;
327     }
328     
329     AliMUONVQAChecker::ECheckCode c1 = AliMUONVQAChecker::kInfo;
330     AliMUONVQAChecker::ECheckCode c2 = AliMUONVQAChecker::kInfo;
331     AliMUONVQAChecker::ECheckCode c3 = AliMUONVQAChecker::kInfo;
332     
333     TH1* hbp = AliQAv1::GetData(list,AliMUONQAIndices::kTrackerBusPatchOccupancy,AliRecoParam::ConvertIndex(specie));    
334     TH1* hbpconfig = AliQAv1::GetData(list,AliMUONQAIndices::kTrackerBusPatchConfig,AliRecoParam::ConvertIndex(specie));
335     TH1* hddl = AliQAv1::GetData(list,AliMUONQAIndices::kTrackerDDLOccupancy,AliRecoParam::ConvertIndex(specie));    
336
337     if ( hbp && hbpconfig && hddl ) 
338     {
339       c1 = BeautifyOccupancyHistograms(*hddl,*hbp,hbpconfig,neventsseen,neventsused,*recoParam);  
340     }
341     else
342     {
343       AliError(Form("Could not BeautifyOccupancyHistograms : hddl=%p hbpconfig=%p hbp=%p",hddl,hbpconfig,hbp));
344     }
345     
346     TH1* hbuspatchtokenerrors = AliQAv1::GetData(list,AliMUONQAIndices::kTrackerBusPatchTokenLostErrors,AliRecoParam::ConvertIndex(specie));
347     TH1* hrostatus = AliQAv1::GetData(list,AliMUONQAIndices::kTrackerReadoutStatus,AliRecoParam::ConvertIndex(specie));    
348     TH1* hrostatusnorm = AliQAv1::GetData(list,AliMUONQAIndices::kTrackerReadoutStatusPerEvent,AliRecoParam::ConvertIndex(specie));
349     
350     if ( hrostatus && hrostatusnorm && hbuspatchtokenerrors )
351     {
352       c2 = BeautifyReadoutHistograms(*hrostatus,*hrostatusnorm,*hbuspatchtokenerrors,
353                                      neventsseen,neventsused,*recoParam);
354     }
355     else
356     {
357       AliError(Form("Could not BeautifyReadoutHistograms : hrostatus=%p hrostatusnorm=%p hbuspatchtokenerrors=%p",
358                     hrostatus,hrostatusnorm,hbuspatchtokenerrors));
359     }
360     
361     
362     TH1* heventsize = AliQAv1::GetData(list,AliMUONQAIndices::kTrackerDDLEventSize,AliRecoParam::ConvertIndex(specie));    
363     TH1* heventsizeperevent = AliQAv1::GetData(list,AliMUONQAIndices::kTrackerDDLEventSizePerEvent,AliRecoParam::ConvertIndex(specie));
364     
365     if ( heventsize && heventsizeperevent ) 
366     {
367       c3 = BeautifyEventsizeHistograms(*heventsize,*heventsizeperevent,
368                                        neventsseen,neventsused,*recoParam);
369     }
370     else
371     {
372       AliError(Form("Could not BeautifyEventsizeHistograms heventsize=%p heventsizeperevent=%p",heventsize,heventsizeperevent));
373     }
374     
375     if ( c1 == AliMUONVQAChecker::kFatal || c2 == AliMUONVQAChecker::kFatal || c3 == AliMUONVQAChecker::kFatal )
376     {
377       rv[specie] = AliMUONVQAChecker::kFatal;
378     }
379     else if ( c1 == AliMUONVQAChecker::kError || c2 == AliMUONVQAChecker::kError || c3 == AliMUONVQAChecker::kError )
380     {
381       rv[specie] = AliMUONVQAChecker::kError;
382     }
383     else if ( c1 == AliMUONVQAChecker::kWarning || c2 == AliMUONVQAChecker::kWarning || c3 == AliMUONVQAChecker::kWarning )
384     {
385       rv[specie] = AliMUONVQAChecker::kWarning;
386     }
387     else
388     {
389       rv[specie] = AliMUONVQAChecker::kInfo;
390     }
391   }
392   
393   return rv;
394 }
395
396 //____________________________________________________________________________ 
397 AliMUONVQAChecker::ECheckCode
398 AliMUONTrackerQAChecker::BeautifyOccupancyHistograms(TH1& hddl,
399                                                      TH1& hbp, 
400                                                      const TH1* hbuspatchconfig, 
401                                                      Int_t neventsseen,
402                                                      Int_t neventsused,
403                                                      const AliMUONRecoParam& recoParam)
404 {
405   /// Put labels, limits and so on on the TrackerBusPatchOccupancy histograms
406   /// hbuspatchconfig and hbp must have the same bin definitions
407   
408   if ( hbuspatchconfig )
409   {
410     if ( hbp.GetNbinsX() != hbuspatchconfig->GetNbinsX() ||
411         hbp.GetXaxis()->GetXmin() != hbuspatchconfig->GetXaxis()->GetXmin() ||
412         hbp.GetXaxis()->GetXmax() != hbuspatchconfig->GetXaxis()->GetXmax() )
413     {
414       AliError("hbp and hbuspatchconfig histograms are not compatible !");
415       return AliMUONVQAChecker::kFatal;
416     }
417   }
418   
419   hddl.SetStats(kFALSE);
420   hbp.SetXTitle("Absolute Bus Patch Id");
421   hbp.SetYTitle("Occupancy (percent)");
422   hbp.SetStats(kFALSE);
423   
424   Double_t xmin = hbp.GetXaxis()->GetXmin();
425   Double_t xmax = hbp.GetXaxis()->GetXmax();
426   
427   Double_t occMax(0.1); // 0.1% y-limit for the plot
428   Double_t maxToleratedOccupancy(recoParam.BuspatchOccupancyHighLimit()*100.0); 
429   Double_t minToleratedOccupancy(recoParam.BuspatchOccupancyLowLimit()*100.0);   
430   TLine* line1 = new TLine(xmin,maxToleratedOccupancy,xmax,maxToleratedOccupancy);
431   line1->SetLineColor(1);
432   line1->SetLineWidth(1);
433
434   TLine* line2 = new TLine(xmin,minToleratedOccupancy,xmax,minToleratedOccupancy);
435   line2->SetLineColor(1);
436   line2->SetLineWidth(1);
437   
438   hbp.GetListOfFunctions()->Add(line1);
439   hbp.GetListOfFunctions()->Add(line2);
440     
441   TIter next(AliMpDDLStore::Instance()->CreateBusPatchIterator());
442   AliMpBusPatch* bp(0x0);
443
444   Int_t nBusPatches(0);
445   Int_t nMissingBusPatches(0);
446   
447   while ( ( bp = static_cast<AliMpBusPatch*>(next())) )
448   {
449     ++nBusPatches;
450     
451     Int_t bin = hbp.FindBin(bp->GetId());
452     
453     if ( hbp.GetBinContent(bin) <= 0.0 ) 
454     {
455       ++nMissingBusPatches;
456     }
457   }
458   
459   next.Reset();
460   
461   Int_t ok(-1);
462   Int_t n(0);
463   Int_t nBusPatchesAboveLimit(0);
464   Int_t nBusPatchesBelowLimit(0);
465   Double_t alpha(0.1); // trim 10% of data
466   Double_t tmean(0.0),tvar(0.0);
467   Double_t ymin(0.0),ymax(0.0);
468   
469   if ( nBusPatches ) 
470   {
471     Double_t* x = new Double_t[nBusPatches];
472     
473     while ( ( bp = static_cast<AliMpBusPatch*>(next())) )
474     {
475       Int_t bin = hbp.FindBin(bp->GetId());
476       if ( hbp.GetBinContent(bin) > 0 )
477       {
478         x[n] = hbp.GetBinContent(bin);
479         ++n;
480       }
481       if ( hbp.GetBinContent(bin) > maxToleratedOccupancy )
482       {
483         ++nBusPatchesAboveLimit;
484       }
485       if ( hbp.GetBinContent(bin) < minToleratedOccupancy )
486       {
487         // check whether this buspatch has a reason to be absent (only valid
488         // if we got the config, otherwise we cannot do the test)
489         if ( hbuspatchconfig && hbuspatchconfig->GetBinContent(bin) > 0 )
490         {
491           // should be there, so it's an error
492           ++nBusPatchesBelowLimit;
493         }
494       }
495     }
496     
497     // computed the truncated mean of the occupancy values, in order to get a 
498     // reasonable y-range for the histogram (without giant peaks to the roof 
499     // for misbehaving buspatches).
500     ok = trim(n,x,alpha,tmean,tvar,ymin,ymax);
501     
502     delete[] x;
503   }
504   
505   if ( ok < 0 ) 
506   {
507     ymax = occMax;
508   }
509   else
510   {
511     ymax = TMath::Max(ymax,occMax);
512   }
513   
514   hbp.SetMaximum(ymax*1.4);
515   
516   AliMUONVQAChecker::ECheckCode rv(AliMUONVQAChecker::kInfo);
517
518   TObjArray messages;
519   messages.SetOwner(kTRUE);
520   
521   if ( neventsseen < NOTENOUGHEVENTLIMIT ) 
522   {
523     messages.Add(new TObjString(NOTENOUGHEVENTMESSAGE));
524   }
525   else
526   {
527     if ( ok < 0 ) 
528     {
529       messages.Add(new TObjString("Could not compute truncated mean. Not enough events ?"));
530       
531       if ( neventsused < TMath::Nint(0.1*neventsseen) )
532       {
533         messages.Add(new TObjString(Form("We've actually seen %d events, but used only %d !",neventsseen,neventsused)));
534         messages.Add(new TObjString("For a normal physics run, this is highly suspect !"));
535         messages.Add(new TObjString(NOTIFYEXPERTMESSAGE));
536         rv = AliMUONVQAChecker::kFatal;
537       }
538     }
539     else if (!nBusPatches)
540     {
541       messages.Add(new TObjString(Form("Could not get the total number of buspatches (%d). ERROR !!!",
542                          nBusPatches)));
543       messages.Add(new TObjString("Please check with expert !"));
544       rv = AliMUONVQAChecker::kFatal;
545     }
546     else
547     {
548       Float_t missingBusPatchFraction = nMissingBusPatches*100.0/nBusPatches;
549       Float_t aboveLimitFraction = nBusPatchesAboveLimit*100.0/nBusPatches;
550       Float_t belowLimitFraction = nBusPatchesBelowLimit*100.0/nBusPatches;
551       
552       messages.Add(new TObjString(Form("%5.2f %% of missing buspatches (%d out of %d)",missingBusPatchFraction,nMissingBusPatches,nBusPatches)));
553       messages.Add(new TObjString(Form("%5.2f %% bus patches above the %5.2f %% limit",aboveLimitFraction,maxToleratedOccupancy)));
554       messages.Add(new TObjString(Form("%5.2f %% bus patches below the %e %% limit",belowLimitFraction,minToleratedOccupancy)));
555       messages.Add(new TObjString(Form("Bus patch mean occupancy (truncated at %2d %%) is %7.2f %%",(Int_t)(alpha*100),tmean)));
556       
557       if ( aboveLimitFraction > recoParam.FractionOfBuspatchOutsideOccupancyLimit()*100.0 || 
558           belowLimitFraction > recoParam.FractionOfBuspatchOutsideOccupancyLimit()*100.0 ) 
559       {
560         rv = AliMUONVQAChecker::kError;
561       }
562       else
563       {
564         rv = AliMUONVQAChecker::kInfo;
565       }
566     }
567   }
568   
569   SetupHisto(neventsseen,neventsused,messages,hbp,rv);
570   
571   /// Make as well a version for DDL occupancy, that'll be used by the shifter
572   SetupHisto(neventsseen,neventsused,messages,hddl,rv);
573   
574   
575   hddl.SetStats(kFALSE);
576   
577   return rv;
578 }
579
580 //______________________________________________________________________________
581 AliMUONVQAChecker::ECheckCode 
582 AliMUONTrackerQAChecker::BeautifyReadoutHistograms(TH1& hrostatus,
583                                                    TH1& hrostatusnorm,
584                                                    const TH1& hbuspatchtokenerrors,
585                                                    Int_t neventsseen,
586                                                    Int_t neventsused,
587                                                    const AliMUONRecoParam& recoParam)
588 {
589   /// Normalize and put some text on the readout error histogram
590   /// Note in particular the treatment of tokenlost errors !
591   
592   hrostatusnorm.Reset();
593   
594   AliMUONVQAChecker::ECheckCode rv(AliMUONVQAChecker::kInfo);
595   
596   TObjArray messages;
597   messages.SetOwner(kTRUE);
598   
599   if ( neventsseen < NOTENOUGHEVENTLIMIT )
600   {
601     messages.Add(new TObjString(NOTENOUGHEVENTMESSAGE));
602   }
603   else
604   {
605     hrostatusnorm.Add(&hrostatus,100.0/neventsseen);
606     hrostatusnorm.SetOption("TEXT45"); // note : cannot use HIST option, otherwise the associated function (i.e. the tpavetext) is not drawn...
607     hrostatusnorm.SetMinimum(0.0);
608     
609     Double_t tokenLost = hrostatusnorm.GetBinContent(hrostatusnorm.FindBin(1.0*AliMUONQAIndices::kTrackerRawNofTokenLostErrors));
610     
611     if ( tokenLost > recoParam.TokenLostLimit() )
612     {
613       rv = AliMUONVQAChecker::kError;
614       
615       messages.Add(new TObjString("There are some token lost errors !"));
616       messages.Add(new TObjString("PLEASE CHECK THE BUSY TIME FOR MUON !"));
617       messages.Add(new TObjString("If above 5 ms please have the MUON expert"));
618       messages.Add(new TObjString("check the following bus patches :"));
619       
620       for ( Int_t i = 1; i <= hbuspatchtokenerrors.GetNbinsX(); ++i ) 
621       {
622         if ( hbuspatchtokenerrors.GetBinContent(i) > 0 ) 
623         {
624           messages.Add(new TObjString(Form("BP %4d",i)));
625         }
626       }
627     }
628     
629     if ( hrostatusnorm.GetBinContent(hrostatusnorm.FindBin(AliMUONQAIndices::kTrackerRawNofEmptyEvents)) > 25.0 )
630     {
631       messages.Add(new TObjString("Too many empty events !"));
632       messages.Add(new TObjString(NOTIFYEXPERTMESSAGE));
633       rv = AliMUONVQAChecker::kFatal;
634       
635     }
636   }
637    
638   SetupHisto(neventsseen,neventsused,messages,hrostatusnorm,rv,"text45");
639   
640   return rv;
641 }
642
643 //______________________________________________________________________________
644 AliMUONVQAChecker::ECheckCode 
645 AliMUONTrackerQAChecker::BeautifyEventsizeHistograms(TH1& heventsize,
646                                                      TH1& heventsizeperevent,
647                                                      Int_t neventsseen,
648                                                      Int_t neventsused,
649                                                      const AliMUONRecoParam& recoParam)
650 {
651   /// Normalize and put some text on the event size histogram
652
653   AliMUONVQAChecker::ECheckCode rv(AliMUONVQAChecker::kInfo);
654
655   heventsizeperevent.Reset();
656
657   TObjArray messages;
658   messages.SetOwner(kTRUE);
659   
660   if ( neventsseen < NOTENOUGHEVENTLIMIT )
661   {
662     messages.Add(new TObjString(NOTENOUGHEVENTMESSAGE));
663   }
664   else
665   {
666     heventsizeperevent.Add(&heventsize,1.0/neventsseen/1024.0);
667     heventsizeperevent.SetMinimum(0);
668     
669     Double_t totalEventSizePerEvent = heventsizeperevent.Integral();
670     
671     TString msg;
672     TString action;
673     
674     if ( totalEventSizePerEvent > recoParam.EventSizeHardLimit() ) 
675     {
676       rv = AliMUONVQAChecker::kFatal;
677       msg = "That is really too high.";
678       action = NOTIFYEXPERTMESSAGE;
679     }
680     else if ( totalEventSizePerEvent > recoParam.EventSizeSoftLimit() ) 
681     {
682       msg = "That is a bit high.";
683       action = "Please keep an eye on it.";
684       rv = AliMUONVQAChecker::kError;
685     }
686     else 
687     {
688       rv = AliMUONVQAChecker::kInfo;
689     }
690     
691     messages.Add(new TObjString(Form("<MCH event size> %5.1f KB/event\n",totalEventSizePerEvent)));
692     if (msg.Length()>0)
693     {
694       messages.Add(new TObjString(msg));
695       messages.Add(new TObjString(action));
696     }
697   }
698   
699   SetupHisto(neventsseen,neventsused,messages,heventsizeperevent,rv);
700   
701   return rv;
702 }
703
704
705