]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0QAChecker.cxx
bug fix in reconstruction of simulated data
[u/mrichter/AliRoot.git] / T0 / AliT0QAChecker.cxx
1 /**************************************************************************
2  * Coyright(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 //  Checks the quality assurance. 
17 //  By comparing with reference data
18 //  Skeleton for T0
19 //---------------------------------------------
20 //checkig without reference data:
21 //for RAW QA all histograms should have approximatly the same 
22 //number of entries as RefPoint
23 //for Rec Points checks 
24 //  - amplitude measured by 2 methos
25 //  - online and offline T0 measurements
26 // for ESD quality of reconstruction ( and measurements):
27 // RMS of vertex and T0 less than 75ps
28 //
29 // Alla.Maevskaya@cern.ch   
30 //...
31
32 // --- ROOT system ---
33 #include <Riostream.h>
34 #include <TClass.h>
35 #include <TH1F.h> 
36 #include <TF1.h> 
37 #include <TFitResultPtr.h>
38 #include <TH2.h> 
39 #include <TIterator.h> 
40 #include <TKey.h> 
41 #include <TFile.h> 
42 #include <TMath.h>
43 #include <TString.h>
44 #include <TPaveText.h>
45 #include <TLegend.h>
46
47 // --- Standard library ---
48
49 // --- AliRoot header files ---
50 #include "AliLog.h"
51 #include "AliQAv1.h"
52 #include "AliQAChecker.h"
53 #include "AliCDBEntry.h"
54 #include "AliQAManager.h"
55 #include "AliT0QAChecker.h"
56 #include "AliQAThresholds.h"
57 #include "AliDAQ.h"
58
59 ClassImp(AliT0QAChecker)
60 //____________________________________________________________________________
61 AliT0QAChecker::AliT0QAChecker() :
62   AliQACheckerBase("T0","T0 Quality Assurance Checker"),
63   fCFDErrorThreshold(0),
64   fLEDErrorThreshold(0),
65   fQTCErrorThreshold(0),
66   fRatioCFDEffLEDEffErrorThreshold(1),
67   fQTCEfficiencyErrorThreshold(0),
68   fBCIDPeriodParam(3564),
69   fBCIDOffsetParam(37),
70   fBCIDBandWidthParam(10),
71   fTZeroAPlusCErrorThreshold(2000.0),
72   fTZeroAMinusCErrorThreshold(2000.0)
73 {
74   // Standard constructor
75   for(Int_t i=0; i<24; i++){ 
76     fMeanCFDFromGoodRunParam[i]=0; 
77     fMeanLEDFromGoodRunParam[i]=0; 
78     fMeanQTCFromGoodRunParam[i]=0; 
79   }
80
81 }
82
83 //____________________________________________________________________________
84 AliT0QAChecker::AliT0QAChecker(const AliT0QAChecker& qac):
85   AliQACheckerBase(qac.GetName(), qac.GetTitle()),
86   fCFDErrorThreshold(qac.fCFDErrorThreshold),
87   fLEDErrorThreshold(qac.fLEDErrorThreshold),
88   fQTCErrorThreshold(qac.fQTCErrorThreshold),
89   fRatioCFDEffLEDEffErrorThreshold(qac.fRatioCFDEffLEDEffErrorThreshold),
90   fQTCEfficiencyErrorThreshold(qac.fQTCEfficiencyErrorThreshold),
91   fBCIDPeriodParam(qac.fBCIDPeriodParam),
92   fBCIDOffsetParam(qac.fBCIDOffsetParam),
93   fBCIDBandWidthParam(qac.fBCIDBandWidthParam),
94   fTZeroAPlusCErrorThreshold(qac.fTZeroAPlusCErrorThreshold),
95   fTZeroAMinusCErrorThreshold(qac.fTZeroAMinusCErrorThreshold)
96 {
97   // copy constructor
98   AliError("Copy should not be used with this class\n");
99   for(Int_t i=0; i<24; i++){ 
100     fMeanCFDFromGoodRunParam[i]=qac.fMeanCFDFromGoodRunParam[i]; 
101     fMeanLEDFromGoodRunParam[i]=qac.fMeanLEDFromGoodRunParam[i]; 
102     fMeanQTCFromGoodRunParam[i]=qac.fMeanQTCFromGoodRunParam[i]; 
103   }
104
105 }
106 //____________________________________________________________________________
107 AliT0QAChecker& AliT0QAChecker::operator=(const AliT0QAChecker& qac){
108   // assignment operator
109   this->~AliT0QAChecker();
110   new(this)AliT0QAChecker(qac);
111   return *this;
112 }
113
114
115 //____________________________________________________________________________
116 AliT0QAChecker::~AliT0QAChecker(){
117   // destructor
118
119 }
120
121 //__________________________________________________________________
122 void AliT0QAChecker::Check(Double_t *  test, AliQAv1::ALITASK_t index, TObjArray ** list, const AliDetectorRecoParam * /*recoParam*/)
123 {
124
125   AliCDBManager* man = AliCDBManager::Instance();
126   //man->SetDefaultStorage(gSystem->Getenv("AMORE_CDB_URI"));
127   if(!man) return; 
128   AliCDBEntry* entry = man->Get("GRP/Calib/QAThresholds");
129   if(!entry) return; 
130   TObjArray* t0branch = (TObjArray*) entry->GetObject();
131   if(!list) return;
132   AliQAThresholds*  thresholds = (AliQAThresholds*) t0branch->FindObject("T00");
133   // here you should test that you got a non-null pointer
134
135   if(!thresholds) return;
136   if(AliDAQ::DetectorID("T0")!= thresholds->GetDetectorId()){
137     AliInfo(Form("DETECTOR ID %d DOES NOT MATCH TO TZERO",thresholds->GetDetectorId()));
138     return;
139   }
140
141   int iparam; 
142   for(int ipmt=0; ipmt<24;ipmt++){ 
143     iparam = ipmt + 1; //current consecutive number of parameter
144     if((TParameter<float>*) thresholds->GetThreshold(iparam)){ // mean CFD from a good run 
145       fMeanCFDFromGoodRunParam[ipmt] = ((TParameter<float>*) thresholds->GetThreshold(iparam))->GetVal();
146     }
147
148     iparam = ipmt + 25;
149     if((TParameter<float>*) thresholds->GetThreshold(iparam)){ // mean LED from a good run 
150       fMeanLEDFromGoodRunParam[ipmt] = ((TParameter<float>*) thresholds->GetThreshold(iparam))->GetVal();
151     } 
152     iparam = ipmt + 49;
153     if((TParameter<float>*) thresholds->GetThreshold(iparam)){ // mean QTC from a good run 
154       fMeanQTCFromGoodRunParam[ipmt] = ((TParameter<float>*) thresholds->GetThreshold(iparam))->GetVal();
155     } 
156   }
157   iparam = 73; //CFD threshold
158   if((TParameter<float>*) thresholds->GetThreshold(iparam)){ 
159     fCFDErrorThreshold = ((TParameter<float>*) thresholds->GetThreshold(iparam))->GetVal();
160   }
161   iparam = 74; //LED threshold
162   if((TParameter<float>*) thresholds->GetThreshold(iparam)){ 
163     fLEDErrorThreshold = ((TParameter<float>*) thresholds->GetThreshold(iparam))->GetVal();
164   }
165   iparam = 75; //QTC threshold
166   if((TParameter<float>*) thresholds->GetThreshold(iparam)){ 
167     fQTCErrorThreshold = ((TParameter<float>*) thresholds->GetThreshold(iparam))->GetVal();
168   }
169  
170   iparam = 82; //Error level threshold on CFD efficiency/LED efficiency ratio
171   if((TParameter<float>*) thresholds->GetThreshold(iparam)){ 
172     fRatioCFDEffLEDEffErrorThreshold = ((TParameter<float>*) thresholds->GetThreshold(iparam))->GetVal();
173   }
174   // Super-basic check on the QA histograms on the input list:
175   // look whether they are empty!
176
177   iparam = 83; //Error level threshold on QTC efficiency 
178   if((TParameter<float>*) thresholds->GetThreshold(iparam)){ 
179     fQTCEfficiencyErrorThreshold = ((TParameter<float>*) thresholds->GetThreshold(iparam))->GetVal();
180   }
181
182   iparam = 84; 
183   if((TParameter<int>*) thresholds->GetThreshold(iparam)){ 
184     fBCIDPeriodParam = ((TParameter<int>*) thresholds->GetThreshold(iparam))->GetVal();
185   }
186
187   iparam = 85; 
188   if((TParameter<int>*) thresholds->GetThreshold(iparam)){ 
189     fBCIDOffsetParam = ((TParameter<int>*) thresholds->GetThreshold(iparam))->GetVal();
190   }
191
192   iparam = 86; 
193   if((TParameter<int>*) thresholds->GetThreshold(iparam)){ 
194     fBCIDBandWidthParam = ((TParameter<int>*) thresholds->GetThreshold(iparam))->GetVal();
195   }
196
197   iparam = 87; 
198   if((TParameter<float>*) thresholds->GetThreshold(iparam)){ 
199     fTZeroAPlusCErrorThreshold = ((TParameter<float>*) thresholds->GetThreshold(iparam))->GetVal();
200   }
201
202   iparam = 88; 
203   if((TParameter<float>*) thresholds->GetThreshold(iparam)){ 
204     fTZeroAMinusCErrorThreshold = ((TParameter<float>*) thresholds->GetThreshold(iparam))->GetVal();
205   }
206
207
208   char * detOCDBDir = Form("T0/%s/%s", AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ; 
209
210   AliCDBEntry *QARefRec = AliQAManager::QAManager()->Get(detOCDBDir);
211   //  QARefRec->Dump();
212   if( !QARefRec){
213     AliInfo("QA reference data NOT retrieved for Reconstruction check. No T0 reference distribution");
214   }
215
216     
217   for(Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++){ 
218     test[specie]    = 1.0; //FK//  initiate qa flag for the whole set of histograms as good 
219   }
220
221
222   for(Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
223     if(!(AliQAv1::Instance()->IsEventSpecieSet(specie) && list[specie]) || list[specie]->GetEntries() == 0) {
224       continue;
225     }
226
227     if(index == AliQAv1::kRAW){
228
229       //if(AliRecoParam::ConvertIndex(specie) == AliRecoParam::kCalib){//      if (index == AliQAv1::kRAW )
230         //check laser data efficiencies   
231       //  Double_t qaFlag = CheckLaser(list[specie]);
232       //  if(qaFlag < test[specie]) test[specie] = qaFlag;
233       //}
234
235       //if(//AliRecoParam::ConvertIndex(specie) == AliRecoParam::kCalib   ||
236        //  AliRecoParam::ConvertIndex(specie) == AliRecoParam::kHighMult ||
237       //   AliRecoParam::ConvertIndex(specie) == AliRecoParam::kLowMult){ 
238          //AliRecoParam::ConvertIndex(specie) == AliRecoParam::kDefault ||
239          
240         //check BCID   
241       //  Double_t qaFlag = CheckBCID(list[specie]);
242       //  if(qaFlag < test[specie]) test[specie] = qaFlag;
243       //}
244
245       if(AliRecoParam::ConvertIndex(specie) == AliRecoParam::kHighMult ||
246         AliRecoParam::ConvertIndex(specie) == AliRecoParam::kLowMult){ 
247         //AliRecoParam::ConvertIndex(specie) == AliRecoParam::kDefault ||
248         //check physics 
249         Double_t qaFlag = CheckRaw(list[specie]);
250         if(qaFlag < test[specie]) test[specie] = qaFlag;
251       }
252     }
253
254     if(index == AliQAv1::kESD && AliRecoParam::Convert(specie) != AliRecoParam::kCalib){
255       test[specie] = CheckESD(list[specie]);
256     } 
257   }
258 }
259
260 //--------------------------------------------------------------------------
261 //Double_t AliT0QAChecker::CheckLaser(TObjArray *listrec) const {
262 //   
263 //  return 1.0; 
264 //}
265
266 //--------------------------------------------------------------------------
267 Double_t AliT0QAChecker::CheckRaw(TObjArray *listrec) const {
268    
269   //Fk Set drawing options for LED and CFD efficiencies from the raw data
270   TH1F *hCFDEffData = (TH1F*) listrec->UncheckedAt(207);//hRawTrigger 
271   TH1F *hLEDEffData = (TH1F*) listrec->UncheckedAt(208);//hRawTrigger 
272
273   //clean objects added at previous checks
274   EraseOldMessages((TH1*) hCFDEffData); 
275   hCFDEffData->GetListOfFunctions()->Add((TH1D*) hLEDEffData->Clone());       
276
277   TLegend leg(0.12,0.76,0.9,0.92," ","brNDC");
278   leg.SetFillStyle(0); leg.SetBorderSize(0); leg.SetTextSize(0.04); leg.SetNColumns(2);
279   leg.AddEntry((TH1D*) hCFDEffData,"CFD","p");
280   leg.AddEntry((TH1D*) hLEDEffData,"LED","p");
281   hCFDEffData->GetListOfFunctions()->Add((TLegend*) leg.Clone());             
282
283
284   //Fk Draw CFD-mean for each PMT
285   TH2F* fhCFD = (TH2F*) listrec->UncheckedAt(210);
286   TH1F* fhCFDSubtrMean = (TH1F*) listrec->UncheckedAt(231);
287
288   EraseOldMessages((TH1*) fhCFDSubtrMean); 
289   for(int ipmt=0; ipmt<24;ipmt++){ 
290     TH1F*  hProjDummy = (TH1F*) fhCFD->ProjectionY("dummy",ipmt+1,ipmt+1);
291     Float_t mean=0.0, rms=0.0;
292     GetMeanAndRmsAroundMainMaximum(mean, rms,  hProjDummy,0);
293
294     Float_t deviation = mean - fMeanCFDFromGoodRunParam[ipmt]; 
295
296     fhCFDSubtrMean->SetBinContent(ipmt+1,deviation);
297     fhCFDSubtrMean->SetBinError(ipmt+1,rms);
298       
299     delete hProjDummy;
300   }
301   TLine linelowredCFD(0, fCFDErrorThreshold, 24, fCFDErrorThreshold);    
302   linelowredCFD.SetLineColor(2);
303   linelowredCFD.SetLineStyle(3);
304   linelowredCFD.SetLineWidth(4);
305   TLine linehighredCFD(0, -fCFDErrorThreshold, 24, -fCFDErrorThreshold);    
306   linehighredCFD.SetLineColor(2);
307   linehighredCFD.SetLineStyle(3);
308   linehighredCFD.SetLineWidth(4);
309   fhCFDSubtrMean->GetListOfFunctions()->Add((TLine*) linelowredCFD.Clone());          
310   fhCFDSubtrMean->GetListOfFunctions()->Add((TLine*) linehighredCFD.Clone());         
311
312  
313  
314   //Fk Draw LED-mean for each PMT
315   TH2F* fhLED = (TH2F*) listrec->UncheckedAt(211);
316   TH1F* fhLEDSubtrMean = (TH1F*) listrec->UncheckedAt(232);
317  
318   EraseOldMessages((TH1*) fhLEDSubtrMean); 
319   for(int ipmt=0; ipmt<24;ipmt++){ 
320     TH1F*  hProjDummy = (TH1F*) fhLED->ProjectionY("dummy",ipmt+1,ipmt+1);
321     Float_t mean=0.0, rms=0.0;
322     GetMeanAndRmsAroundMainMaximum(mean, rms,  hProjDummy,1);
323     Float_t deviation = mean - fMeanLEDFromGoodRunParam[ipmt]; 
324
325     fhLEDSubtrMean->SetBinContent(ipmt+1,deviation);
326     fhLEDSubtrMean->SetBinError(ipmt+1,rms);
327       
328     delete hProjDummy;
329   }
330   TLine linelowredLED(0, fLEDErrorThreshold, 24, fLEDErrorThreshold);    
331   linelowredLED.SetLineColor(2);
332   linelowredLED.SetLineStyle(3);
333   linelowredLED.SetLineWidth(4);
334   TLine linehighredLED(0, -fLEDErrorThreshold, 24, -fLEDErrorThreshold);    
335   linehighredLED.SetLineColor(2);
336   linehighredLED.SetLineStyle(3);
337   linehighredLED.SetLineWidth(4);
338   fhLEDSubtrMean->GetListOfFunctions()->Add((TLine*) linelowredLED.Clone());          
339   fhLEDSubtrMean->GetListOfFunctions()->Add((TLine*) linehighredLED.Clone());         
340
341      
342   //Fk Draw QTC-mean for each PMT
343   TH2F* fhQTC = (TH2F*) listrec->UncheckedAt(212);
344   TH1F* fhQTCSubtrMean = (TH1F*) listrec->UncheckedAt(233);
345    
346   EraseOldMessages((TH1*) fhQTCSubtrMean); 
347   for(int ipmt=0; ipmt<24;ipmt++){ 
348     TH1F*  hProjDummy = (TH1F*) fhQTC->ProjectionY("dummy",ipmt+1,ipmt+1);
349     Float_t mean=0.0, rms=0.0;
350     GetMeanAndRmsAroundMainMaximum(mean, rms,  hProjDummy,2);
351     Float_t deviation = mean - fMeanQTCFromGoodRunParam[ipmt]; 
352
353     fhQTCSubtrMean->SetBinContent(ipmt+1,deviation);
354     fhQTCSubtrMean->SetBinError(ipmt+1,rms);
355       
356     delete hProjDummy;
357   }
358   TLine linelowredQTC(0, fQTCErrorThreshold, 24, fQTCErrorThreshold);    
359   linelowredQTC.SetLineColor(2);
360   linelowredQTC.SetLineStyle(3);
361   linelowredQTC.SetLineWidth(4);
362   TLine linehighredQTC(0, -fQTCErrorThreshold, 24, -fQTCErrorThreshold);    
363   linehighredQTC.SetLineColor(2);
364   linehighredQTC.SetLineStyle(3);
365   linehighredQTC.SetLineWidth(4);
366   fhQTCSubtrMean->GetListOfFunctions()->Add((TLine*) linelowredQTC.Clone());          
367   fhQTCSubtrMean->GetListOfFunctions()->Add((TLine*) linehighredQTC.Clone());         
368
369   //CFD and LED efficiency in range ~2000- ~3000 
370   TH1F* hCFDeffSubRange = (TH1F*) listrec->UncheckedAt(237);
371   TH1F* hEffLEDSubRange = (TH1F*) listrec->UncheckedAt(238);
372   // ratio CDF eff /LEF eff in subragne 
373   TH1F* hRatioCFDLEDeff = (TH1F*) listrec->UncheckedAt(239);//FK   
374   EraseOldMessages((TH1*) hRatioCFDLEDeff); 
375   int npmt = hRatioCFDLEDeff->GetNbinsX();
376   for(int ipmt=1;ipmt<=npmt;ipmt++){
377     Float_t c0 = hCFDeffSubRange->GetBinContent(ipmt); 
378     Float_t c1 = hEffLEDSubRange->GetBinContent(ipmt);
379     if(c1){
380       hRatioCFDLEDeff->SetBinContent(ipmt,c0/c1);  
381     }else{
382       hRatioCFDLEDeff->SetBinContent(ipmt,0);  
383     }  
384   }
385
386   TLine linelowredRatioCFDLEDeff(0, 1+fRatioCFDEffLEDEffErrorThreshold, 24, 1+fRatioCFDEffLEDEffErrorThreshold);    
387   linelowredRatioCFDLEDeff.SetLineColor(2);
388   linelowredRatioCFDLEDeff.SetLineStyle(3);
389   linelowredRatioCFDLEDeff.SetLineWidth(4);
390   TLine linehighredRatioCFDLEDeff(0, 1-fRatioCFDEffLEDEffErrorThreshold, 24, 1-fRatioCFDEffLEDEffErrorThreshold);    
391   linehighredRatioCFDLEDeff.SetLineColor(2);
392   linehighredRatioCFDLEDeff.SetLineStyle(3);
393   linehighredRatioCFDLEDeff.SetLineWidth(4);
394   hRatioCFDLEDeff->GetListOfFunctions()->Add((TLine*) linelowredRatioCFDLEDeff.Clone());              
395   hRatioCFDLEDeff->GetListOfFunctions()->Add((TLine*) linehighredRatioCFDLEDeff.Clone());             
396
397   //        PERFROM CHECKS on HISTOGRAMS
398
399   //-------- triggers -----------
400   Int_t qualityFlagTrigger = kT0Info; //init quality flag for a given histogram; 
401
402   TH1F *hTrigger = (TH1F*) listrec->UncheckedAt(169);//hRawTrigger 
403
404   // clean objects added at previous checks
405   EraseOldMessages((TH1*) hTrigger); 
406
407   if(hTrigger->Integral()>0){
408     //trigger plot does have some counts in it
409     //are Mean, ORA and ORC not empty?  
410     if( hTrigger->GetBinContent(1)<0.001 || hTrigger->GetBinContent(3)<0.001 || hTrigger->GetBinContent(4)<0.001){
411       qualityFlagTrigger = kT0Error; //no entries on diagonal
412       AliDebug(AliQAv1::GetQADebugLevel(), Form("T0: too little ORA and ORC in  %s", hTrigger->GetName() ));
413
414       TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
415       text.AddText(Form("Check ORA and ORC")); 
416       text.AddText(Form("Report problem to the T0 on-call expert")); 
417       text.SetBorderSize(0);
418       text.SetFillStyle(0);
419       hTrigger->GetListOfFunctions()->Add((TPaveText*)text.Clone());          
420     }
421
422   }else{ //Trigger histo empty
423
424     qualityFlagTrigger = kT0Error;
425     AliDebug(AliQAv1::GetQADebugLevel(), Form("T0 histogram  %s has NO entries", hTrigger->GetName() ));
426
427     TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
428     text.AddText(Form("NO ENTRIES!!!")); 
429     text.AddText(Form("If T0 is READY report")); 
430     text.AddText(Form("readout problem to the T0 on-call expert")); 
431     text.SetBorderSize(0);
432     text.SetFillStyle(0);
433     hTrigger->GetListOfFunctions()->Add((TPaveText*)text.Clone());            
434
435   }
436   //---------- CFD eff/LED eff within subrange  -----------
437   Int_t qualityFlagRatioCFDeffLEDeff = kT0Info; //init quality flag for a given histogram;
438   int nPMTs = hRatioCFDLEDeff->GetNbinsX(); 
439
440   for(int ipmt=1; ipmt<=nPMTs; ipmt++){
441     if(TMath::Abs( hRatioCFDLEDeff->GetBinContent(ipmt) -1) > fRatioCFDEffLEDEffErrorThreshold){ //mean is expected to be around 1
442       qualityFlagRatioCFDeffLEDeff = kT0Error;
443     }
444   }
445   if(qualityFlagRatioCFDeffLEDeff == kT0Error){
446     TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
447     text.AddText(Form("Problem with efficiency ratio CFD/LED !!!")); 
448     text.AddText(Form("If T0 is READY and beam is on report")); 
449     text.AddText(Form("the problem to the T0 on-call expert"));
450     text.SetBorderSize(0);
451     text.SetFillStyle(0);
452     hRatioCFDLEDeff->GetListOfFunctions()->Add((TPaveText*)text.Clone());             
453   }
454   
455
456   //---------- CFD -  mean CFD  -----------
457   Int_t qualityFlagCFDSubtr = kT0Info; //init quality flag for a given histogram;
458   nPMTs = fhCFDSubtrMean->GetNbinsX(); 
459
460   for(int ipmt=1; ipmt<=nPMTs; ipmt++){
461     if(TMath::Abs( fhCFDSubtrMean->GetBinContent(ipmt)) > fCFDErrorThreshold){
462       qualityFlagCFDSubtr = kT0Error;
463
464    }
465   }
466   if(qualityFlagCFDSubtr == kT0Error){
467     TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
468     text.AddText(Form("Problem with CFD timing  !!!")); 
469     text.AddText(Form("If T0 is READY and beam is on report")); 
470     text.AddText(Form("the problem to the T0 on-call expert"));
471     text.SetBorderSize(0);
472     text.SetFillStyle(0);
473     fhCFDSubtrMean->GetListOfFunctions()->Add((TPaveText*)text.Clone());              
474   }
475   //--------- QTC efficiency ---------------
476
477   Int_t qualityFlagEffQTC = kT0Info; //init quality flag for a given histogram;
478   TH1F *hEffQTC = (TH1F*) listrec->UncheckedAt(209);//hRawTrigger 
479   
480   EraseOldMessages((TH1*) hEffQTC); // clean objects added at previous checks
481
482   nPMTs = hEffQTC->GetNbinsX(); 
483   for(int ipmt=1; ipmt<=nPMTs; ipmt++){
484     if(TMath::Abs( hEffQTC->GetBinContent(ipmt)) < fQTCEfficiencyErrorThreshold){
485       qualityFlagEffQTC = kT0Error;
486     }
487   }
488   if( qualityFlagEffQTC == kT0Error){
489     TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
490     text.AddText(Form("Problem with QTC efficiency !!!")); 
491     text.AddText(Form("If T0 is READY and beam is on report")); 
492     text.AddText(Form("the problem to the T0 on-call expert"));
493     text.SetBorderSize(0);
494     text.SetFillStyle(0);
495     hEffQTC->GetListOfFunctions()->Add((TPaveText*)text.Clone());             
496   }
497   TLine linehighredQTCeff(0, fQTCEfficiencyErrorThreshold, 24, fQTCEfficiencyErrorThreshold);    
498   linehighredQTCeff.SetLineColor(2);
499   linehighredQTCeff.SetLineStyle(3);
500   linehighredQTCeff.SetLineWidth(4);
501   hEffQTC->GetListOfFunctions()->Add((TLine*) linehighredQTCeff.Clone());             
502  
503   //---------- BCID --------------------
504   Int_t qualityFlagBCID = kT0Info; //init quality flag for a given histogram; 
505   TH2F *hBCID = (TH2F*) listrec->UncheckedAt(236); //BCID versus TRM  BCID
506      
507   // clean objects added at previous checks
508   EraseOldMessages((TH1*)hBCID);
509
510   if(hBCID->Integral()>0){
511     //BCID does have some counts in it
512
513     Int_t nbinsX = hBCID->GetNbinsX();
514     Int_t startX = hBCID->GetXaxis()->FindBin(100); //skip region close to orbit period
515     Int_t nbinsY = hBCID->GetNbinsY();
516     Int_t binWidthX = (Int_t) hBCID->GetXaxis()->GetBinWidth(1);
517     Int_t binWidthY = (Int_t) hBCID->GetYaxis()->GetBinWidth(1);
518     double entriesOnDiagonal  = 0; //count diagonal and off diagonal entries
519     double entriesOffDiagonal = 0;
520
521     for(Int_t itrm=startX; itrm<=nbinsX; itrm++){ //BCID TRM
522       for(Int_t ibcid=1; ibcid<=nbinsY; ibcid++){ //BCID
523         if(TMath::Abs( (itrm*binWidthX - fBCIDOffsetParam) % fBCIDPeriodParam - binWidthY*ibcid) < fBCIDBandWidthParam ){ 
524           entriesOnDiagonal  += hBCID->GetBinContent(itrm,ibcid); //On  Diagonal
525           //hBCID->Fill(itrm*binWidthX,ibcid*binWidthY,0.001); // visualize the diagonal belt
526         }else{
527           entriesOffDiagonal += hBCID->GetBinContent(itrm,ibcid); //Off Diagonal
528         }
529       }
530     }
531     if(entriesOnDiagonal<1 || entriesOffDiagonal>20){
532       qualityFlagBCID = kT0Error; //no entries on diagonal
533       AliDebug(AliQAv1::GetQADebugLevel(), Form("T0   %s is not diagonal", hBCID->GetName() ));
534
535       TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
536       text.AddText(Form("Check if entries are on diagonal")); 
537       text.AddText(Form("Report readout problem to the T0 on-call expert")); 
538       text.SetBorderSize(0);
539       text.SetFillStyle(0);
540       hBCID->GetListOfFunctions()->Add((TPaveText*)text.Clone());             
541     }
542   }else{ //BCID empty
543
544     qualityFlagBCID = kT0Error;
545     AliDebug(AliQAv1::GetQADebugLevel(), Form("T0 :  %s has NO entries", hBCID->GetName() ));
546
547     TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
548     text.AddText(Form("NO ENTRIES!!!")); 
549     text.AddText(Form("If T0 is READY report")); 
550     text.AddText(Form("readout problem to the T0 on-call expert"));
551     text.SetBorderSize(0);
552     text.SetFillStyle(0);
553     hBCID->GetListOfFunctions()->Add((TPaveText*)text.Clone());       
554   }
555
556   //--------- mean versus vertex 1st ---------
557   Int_t qualityFlagMeanVersusVertex1st = kT0Info; //init quality flag for a given histogram;
558   TH2F *hMeanVersusVertex1st  = (TH2F*) listrec->UncheckedAt(220);  //fhMeanBest  time
559   EraseOldMessages((TH1*) hMeanVersusVertex1st); 
560  
561   if(hMeanVersusVertex1st->Integral()<1){
562     qualityFlagMeanVersusVertex1st = kT0Error;
563     AliDebug(AliQAv1::GetQADebugLevel(), Form("T0 histogram  %s has NO entries", hMeanVersusVertex1st->GetName() ));
564
565     TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
566     text.AddText(Form("NO ENTRIES!!!")); 
567     text.AddText(Form("If T0 is READY and beam is on report")); 
568     text.AddText(Form("the problem to the T0 on-call expert"));
569     text.SetBorderSize(0);
570     text.SetFillStyle(0);
571     hMeanVersusVertex1st->GetListOfFunctions()->Add((TPaveText*)text.Clone());        
572   }
573
574
575   //--------- vertex TVDC on ---------
576   Int_t qualityFlagVertex1stTVDCon = kT0Info; //init quality flag for a given histogram;
577   TH1F *hVertex1stTVDCon  = (TH1F*) listrec->UncheckedAt(223);  //fhMeanBest  time
578   EraseOldMessages((TH1*) hVertex1stTVDCon); 
579  
580   if(hVertex1stTVDCon->Integral()<1){
581     qualityFlagVertex1stTVDCon = kT0Error;
582     AliDebug(AliQAv1::GetQADebugLevel(), Form("T0 histogram  %s has NO entries", hVertex1stTVDCon->GetName() ));
583
584     TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
585     text.AddText(Form("NO ENTRIES!!!")); 
586     text.AddText(Form("If T0 is READY and beam is on report")); 
587     text.AddText(Form("the problem to the T0 on-call expert"));
588     text.SetBorderSize(0);
589     text.SetFillStyle(0);
590     hVertex1stTVDCon->GetListOfFunctions()->Add((TPaveText*)text.Clone());            
591   }else{
592     if(TMath::Abs(hVertex1stTVDCon->GetMean()) > fTZeroAMinusCErrorThreshold){ 
593       qualityFlagVertex1stTVDCon = kT0Error;
594
595       TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
596       text.AddText(Form("Displaced vertex")); 
597       text.AddText(Form("If T0 is READY and beam is on report")); 
598       text.AddText(Form("the problem to the T0 on-call expert"));
599       text.SetBorderSize(0);
600       text.SetFillStyle(0);
601       hVertex1stTVDCon->GetListOfFunctions()->Add((TPaveText*)text.Clone());          
602     }
603   }
604
605   //--------- vertex TVDC off ---------
606   Int_t qualityFlagVertex1stTVDCoff = kT0Info; //init quality flag for a given histogram;
607   TH1F *hVertex1stTVDCoff  = (TH1F*) listrec->UncheckedAt(225);  //fhMeanBest  time
608   EraseOldMessages((TH1*) hVertex1stTVDCoff); 
609  
610   if(hVertex1stTVDCoff->Integral()<1){
611     qualityFlagVertex1stTVDCoff = kT0Warning;
612     AliDebug(AliQAv1::GetQADebugLevel(), Form("T0 histogram  %s has NO entries", hVertex1stTVDCoff->GetName() ));
613
614     TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
615     text.AddText(Form("Warning: NO ENTRIES")); 
616     text.SetBorderSize(0);
617     text.SetFillStyle(0);
618     hVertex1stTVDCoff->GetListOfFunctions()->Add((TPaveText*)text.Clone());           
619   }
620
621
622
623   //----------- time TVDC on ---------
624     
625   Int_t qualityFlagMean1stTVDCon = kT0Info; //init quality flag for a given histogram;
626   TH1F *hMean1stTVDCon  = (TH1F*) listrec->UncheckedAt(226);  //fhMeanBest  time
627   EraseOldMessages((TH1*) hMean1stTVDCon); 
628  
629   if(hMean1stTVDCon->Integral()<1){
630     qualityFlagMean1stTVDCon = kT0Error;
631     AliDebug(AliQAv1::GetQADebugLevel(), Form("T0 histogram  %s has NO entries", hMean1stTVDCon->GetName() ));
632
633     TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
634     text.AddText(Form("NO ENTRIES!!!")); 
635     text.AddText(Form("If T0 is READY and beam is on report")); 
636     text.AddText(Form("the problem to the T0 on-call expert"));
637     text.SetBorderSize(0);
638     text.SetFillStyle(0);
639     hMean1stTVDCon->GetListOfFunctions()->Add((TPaveText*)text.Clone());              
640   }else{
641     //cout<<"Mean: "<<TMath::Abs(hMean1stTVDCon->GetMean())<<" threshold "<<fTZeroAPlusCErrorThreshold<<endl;
642     if(TMath::Abs(hMean1stTVDCon->GetMean()) > fTZeroAPlusCErrorThreshold){ 
643       qualityFlagMean1stTVDCon = kT0Error;
644
645       TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
646       text.AddText(Form("Shift of mean time")); 
647       text.AddText(Form("If T0 is READY and beam is on report")); 
648       text.AddText(Form("the problem to the T0 on-call expert"));
649       text.SetBorderSize(0);
650       text.SetFillStyle(0);
651       hMean1stTVDCon->GetListOfFunctions()->Add((TPaveText*)text.Clone());            
652     }
653   }
654
655   //----------- time TVDC off ---------
656     
657   Int_t qualityFlagMean1stTVDCoff = kT0Info; //init quality flag for a given histogram;
658   TH1F *hMean1stTVDCoff  = (TH1F*) listrec->UncheckedAt(227);  //fhMeanBest  time
659   EraseOldMessages((TH1*) hMean1stTVDCoff); 
660  
661   if(hMean1stTVDCoff->Integral()<1){
662     qualityFlagMean1stTVDCoff = kT0Warning;
663     AliDebug(AliQAv1::GetQADebugLevel(), Form("T0 histogram  %s has NO entries", hMean1stTVDCoff->GetName() ));
664
665     TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
666     text.AddText(Form("Warning: NO ENTRIES")); 
667     text.SetBorderSize(0);
668     text.SetFillStyle(0);
669     hMean1stTVDCoff->GetListOfFunctions()->Add((TPaveText*)text.Clone());             
670   }
671  
672
673
674   //----------------------- executive summary ---------------------
675   int lowestQualityFlag = (int) qualityFlagTrigger;
676   if(qualityFlagRatioCFDeffLEDeff   < lowestQualityFlag)  lowestQualityFlag = qualityFlagRatioCFDeffLEDeff;
677   if(qualityFlagCFDSubtr            < lowestQualityFlag)  lowestQualityFlag = qualityFlagCFDSubtr;
678   if(qualityFlagEffQTC              < lowestQualityFlag)  lowestQualityFlag = qualityFlagEffQTC;
679   if(qualityFlagMeanVersusVertex1st < lowestQualityFlag)  lowestQualityFlag = qualityFlagMeanVersusVertex1st;
680   if(qualityFlagBCID                < lowestQualityFlag)  lowestQualityFlag = qualityFlagBCID;
681   if(qualityFlagVertex1stTVDCon     < lowestQualityFlag)  lowestQualityFlag = qualityFlagVertex1stTVDCon;
682   if(qualityFlagVertex1stTVDCoff    < lowestQualityFlag)  lowestQualityFlag = qualityFlagVertex1stTVDCoff;
683   if(qualityFlagMean1stTVDCon       < lowestQualityFlag)  lowestQualityFlag = qualityFlagMean1stTVDCon;
684   if(qualityFlagMean1stTVDCoff      < lowestQualityFlag)  lowestQualityFlag = qualityFlagMean1stTVDCoff;
685   
686
687   return ConvertQualityFlagToDouble(lowestQualityFlag); 
688   
689 }
690
691 //--------------------------------------------------------------------------
692 Double_t AliT0QAChecker::CheckESD(TObjArray *listrec ) const
693 {
694   Float_t checkr = 0;
695   TH1 *fhESD;
696  
697   fhESD  = (TH1*) listrec->UncheckedAt(2);
698   if(fhESD){
699     AliDebug(AliQAv1::GetQADebugLevel(), Form("count %s ", fhESD->GetName()) );
700      TF1 *f1 = new TF1("f1","gaus",-1,1);
701     fhESD->Fit("f1","R","Q", -1,1);
702     Double_t par[3];
703     f1->GetParameters(&par[0]);
704     
705     TPaveText text(0.30,0.50,0.99,0.99,"NDC");
706     
707     text.AddText(Form("T0 RUN %d ",AliCDBManager::Instance()->GetRun()));
708     
709     AliDebug(AliQAv1::GetQADebugLevel(), Form("numentries %d mean %f  #sigma %f", (int)fhESD->GetEntries(),par[1], par[2]));
710     
711     
712     if (par[2] > 0.07 && par[2] < 1.) {
713       checkr=0.5;
714        text.AddText(Form("not good resolution :\n %f ns\n", par[2] ));
715        text.SetFillColor(5);
716        printf("T0 detector resolution is not good enouph: %f ns\n",par[2] );
717     }
718     if(TMath::Abs(par[1])>0.05) {
719       checkr = 0.5;
720       text.AddText(Form(" Check clock shift on %f ns", par[1]));
721       text.SetFillColor(5);
722     }
723     if (par[2] >  1. || TMath::Abs(par[1])>0.1) {
724       checkr = 0.25;
725       text.AddText(Form(" Bad resolution:\n mean %f ns sigma %f ns", par[1], par[2]));
726       text.SetFillColor(2);
727       { // RS Clean previous additions
728         TList* lstF = fhESD->GetListOfFunctions();
729         if (lstF) {
730           TObject *stats = lstF->FindObject("stats");
731           lstF->Remove(stats);
732           TObject *obj;
733           while ((obj = lstF->First())) {
734             while(lstF->Remove(obj)) { }
735             delete obj;
736           }
737           if (stats) lstF->Add(stats);
738         } 
739       }
740       fhESD->GetListOfFunctions()->Add((TPaveText*)text.Clone());       
741       AliDebug(AliQAv1::GetQADebugLevel(),
742                Form("Please, check calibration: shift= %f resolution %f test=%f\n",
743                     par[1], par[2], checkr) ) ; 
744     }
745   }
746   else
747     {
748       AliDebug(AliQAv1::GetQADebugLevel(),
749                Form("No ESD QA histogram found, nothing to check"));
750         checkr=0;
751     }
752   
753   
754   return checkr;
755 }
756
757
758 //--------------------------------------------------------------------------
759 void AliT0QAChecker::EraseOldMessages(TH1* h) const 
760 {
761   //erase the old captions 
762   TList* lstF = h->GetListOfFunctions();
763   if(lstF){  
764      TObject *stats = lstF->FindObject("stats");
765      lstF->Remove(stats);
766      TObject *obj;
767      while ((obj = lstF->First())) {
768        while(lstF->Remove(obj)) { }
769          delete obj;
770     }
771     if (stats) lstF->Add(stats);
772   }
773 }
774 //--------------------------------------------------------------------------
775 Double_t AliT0QAChecker::ConvertQualityFlagToDouble(int qualityFlag) const 
776 {
777   //covert quality flag to double
778   Double_t checkr=1.0;
779
780   switch ( qualityFlag ){
781     case kT0Info:
782         checkr = 1.0; break;
783     case kT0Warning:
784         checkr = 0.75; break;
785     case kT0Error:
786           checkr = 0.25; break;
787     case kT0Fatal:
788         checkr = -1.0; break;
789     default:
790          AliError("Invalid ecc value. FIXME !");
791          checkr = 0.25; break;
792   };
793
794   return checkr; 
795 }
796  
797 //--------------------------------------------------------------------------
798 Float_t AliT0QAChecker::GetMeanAboveThreshold(TH1F* hV, Float_t thr) const{
799   //caculate mean value of histo bins above threshold
800   Int_t nBins = hV->GetNbinsX();
801   Int_t nBinsAboveThr = 0; 
802   Float_t sum = 0;
803
804   for(Int_t ib=1;ib<=nBins;ib++){
805     Float_t val = hV->GetBinContent(ib);
806     if(val<thr){
807       sum+=val;
808       nBinsAboveThr++; 
809     }
810   }
811
812   if(nBinsAboveThr>0) return sum/nBinsAboveThr;
813   else return hV->GetMean();
814 }
815     
816 //--------------------------------------------------------------------------
817 void AliT0QAChecker::GetMeanAndRmsAroundMainMaximum(Float_t &meanHisto,Float_t &rmsHisto, TH1F *histo, int type) const{
818
819   if(!histo){
820     meanHisto=0.0;
821     rmsHisto=0.0;
822     return;
823   }
824   if(histo->Integral()<0.00001){
825     meanHisto=0.0;
826     rmsHisto=0.0;
827     return;
828   }
829
830   double nSigma      = 3.0; //n sigma window around mean and main maximum 
831   double expectedRMS = 13.0; // expected rms of the main peak in case of CFD
832   if(type == 1) expectedRMS = 34.0; //LED
833   if(type == 2) expectedRMS = 34.0; //QTC
834
835   //0) approx of mean is global maximum
836   Int_t   nb     =  histo->GetNbinsX();
837   Int_t   bmax   =  histo->GetMaximumBin();
838   Float_t xmax   =  histo->GetBinCenter(bmax);
839
840   double window = expectedRMS * nSigma;
841   //1) estimate a mean in the nSigma window around main max
842   int nlow = histo->FindBin( xmax - window);
843   int nhigh = histo->FindBin( xmax + window);
844   if(nlow<1)   nlow = 1; 
845   if(nhigh>nb) nhigh = nb;
846
847   Float_t sum=0.0, sumWeight=0.0, sumWeightSqr=0.0;
848   for(int ii=nlow;ii<=nhigh;ii++){
849     sum          += histo->GetBinContent(ii);  
850     sumWeight    += histo->GetBinContent(ii)*histo->GetBinCenter(ii);
851   }  
852   if(sum>0.0){
853     meanHisto = sumWeight/sum; //mean in 1st itteration
854   }else{
855     meanHisto = 0.0;
856     rmsHisto=0.0;
857     return;
858   }
859   
860   //2) recalculte mean and rms in the nSigma window around mean 1)
861   nlow  = histo->FindBin( meanHisto - window);
862   nhigh = histo->FindBin( meanHisto + window);
863   if(nlow<1)   nlow = 1; 
864   if(nhigh>nb) nhigh = nb;
865
866   sum=0.0; sumWeight=0.0; sumWeightSqr=0.0;
867   for(int ii=nlow;ii<=nhigh;ii++){
868     sum          += histo->GetBinContent(ii);  
869     sumWeight    += histo->GetBinContent(ii)*histo->GetBinCenter(ii);
870     sumWeightSqr += histo->GetBinContent(ii)*histo->GetBinCenter(ii)*histo->GetBinCenter(ii);
871   }
872   if(sum>0.0){
873     meanHisto = sumWeight/sum; //mean in 2nd itteration
874     rmsHisto  = sumWeightSqr/sum - meanHisto*meanHisto;//rms square
875     if(rmsHisto > 0.0) rmsHisto = sqrt(rmsHisto); //rms
876     else  rmsHisto = 0.0;
877     return;
878   }else{
879     meanHisto = 0.0;
880     rmsHisto  = 0.0;
881     return;
882   }
883 }