]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQAChecker.cxx
Two bug fixes: mismatch between subdetector names and fix for a return statement...
[u/mrichter/AliRoot.git] / ITS / AliITSQAChecker.cxx
1 /**************************************************************************
2  * Copyright(c) 2007-2009, 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 // *****************************************
19 //  Checks the quality assurance 
20 //  by comparing with reference data
21 //  W.Ferrarese  P.Cerello  Mag 2008
22 //  INFN Torino
23
24 // --- ROOT system ---
25 #include "TH1.h"
26 #include <Riostream.h>
27
28 // --- AliRoot header files ---
29 #include "AliITSQAChecker.h"
30 #include "AliITSQASPDChecker.h"
31 #include "AliITSQASDDChecker.h"
32 #include "AliITSQASSDChecker.h"
33 #include "AliITSQADataMakerRec.h"
34
35 ClassImp(AliITSQAChecker)
36
37 //____________________________________________________________________________
38 AliITSQAChecker::AliITSQAChecker(Bool_t kMode, Short_t subDet, Short_t ldc) :
39 AliQACheckerBase("ITS","SDD Quality Assurance Checker"),
40 fkOnline(0),
41 fDet(0),  
42 fLDC(0),
43 fSPDOffset(0), 
44 fSDDOffset(0), 
45 fSSDOffset(0),
46 fSPDHisto(0),
47 fSDDHisto(0),
48 fSSDHisto(0),
49 fSPDChecker(0),  // SPD Checker
50 fSDDChecker(0),  // SDD Checker
51 fSSDChecker(0)  // SSD Checker
52
53 {
54   // Standard constructor
55   fkOnline = kMode; fDet = subDet; fLDC = ldc;
56   if(fDet == 0 || fDet == 1) {
57     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQAChecker::Create SPD Checker\n");
58     fSPDChecker = new AliITSQASPDChecker();
59   }
60   if(fDet == 0 || fDet == 2) {
61     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQAChecker::Create SDD Checker\n");
62     fSDDChecker = new AliITSQASDDChecker();
63   }
64   if(fDet == 0 || fDet == 3) {
65     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQAChecker::Create SSD Checker\n");
66     fSSDChecker = new AliITSQASSDChecker();
67   }
68   InitQACheckerLimits();
69 }
70
71 //____________________________________________________________________________
72 AliITSQAChecker::AliITSQAChecker(const AliITSQAChecker& qac):
73 AliQACheckerBase(qac.GetName(), qac.GetTitle()), 
74 fkOnline(qac.fkOnline), 
75 fDet(qac.fDet), 
76 fLDC(qac.fLDC), 
77 fSPDOffset(qac.fSPDOffset), 
78 fSDDOffset(qac.fSDDOffset), 
79 fSSDOffset(qac.fSSDOffset), 
80 fSPDHisto(qac.fSPDHisto),
81 fSDDHisto(qac.fSDDHisto),
82 fSSDHisto(qac.fSSDHisto),
83 fSPDChecker(qac.fSPDChecker), 
84 fSDDChecker(qac.fSDDChecker), 
85 fSSDChecker(qac.fSSDChecker)
86 {
87   // copy constructor
88   AliError("Copy should not be used with this class\n");
89 }
90 //____________________________________________________________________________
91 AliITSQAChecker& AliITSQAChecker::operator=(const AliITSQAChecker& qac){
92   // assignment operator
93   this->~AliITSQAChecker();
94   new(this)AliITSQAChecker(qac);
95   return *this;
96 }
97
98
99 //____________________________________________________________________________
100 AliITSQAChecker::~AliITSQAChecker(){
101   // destructor
102   if(fSPDChecker)delete fSPDChecker;
103   if(fSDDChecker)delete fSDDChecker;
104   if(fSSDChecker)delete fSSDChecker;
105
106 }
107 //____________________________________________________________________________
108 void AliITSQAChecker::Check(Double_t * rv, AliQAv1::ALITASK_t index, TObjArray ** list, const AliDetectorRecoParam * recoParam)
109 {
110
111
112   // basic checks on the QA histograms on the input list
113   //for the ITS subdetectorQA (Raws Digits Hits RecPoints SDigits) return the worst value of the three result
114   if(index == AliQAv1::kESD){
115     
116     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
117       rv[specie] = 0.0 ; 
118       if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
119         continue ; 
120       AliDebug(AliQAv1::GetQADebugLevel(),"Checker for ESD");
121       Int_t tested = 0;
122       Int_t empty = 0;
123       // The following flags are set to kTRUE if the corresponding
124       // QA histograms exceed a given quality threshold
125       Bool_t cluMapSA = kFALSE;
126       Bool_t cluMapMI = kFALSE;
127       Bool_t cluMI = kFALSE;
128       Bool_t cluSA = kFALSE;
129       Bool_t verSPDZ = kFALSE;
130       if (list[specie]->GetEntries() == 0) {
131         rv[specie] = 0.; // nothing to check
132       }
133       else {
134         Double_t *stepbit=new Double_t[AliQAv1::kNBIT];
135         Double_t histonumb= list[specie]->GetEntries();
136         CreateStepForBit(histonumb,stepbit); 
137         TIter next1(list[specie]);
138         TH1 * hdata;
139         Int_t nskipped=0;
140         Bool_t skipped[6]={kFALSE,kFALSE,kFALSE,kFALSE,kFALSE,kFALSE};
141         // look for layers that we wanted to skip
142         while ( (hdata = dynamic_cast<TH1 *>(next1())) ) {
143           if(!hdata) continue;
144           TString hname = hdata->GetName();
145           if(!hname.Contains("hESDSkippedLayers")) continue;
146           for(Int_t k=1; k<7; k++) {
147             if(hdata->GetBinContent(k)>0) { 
148               nskipped++; 
149               skipped[k-1]=kTRUE; 
150             } 
151           } 
152         }
153         TIter next(list[specie]);
154         while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
155           if(hdata){
156             TString hname = hdata->GetName();
157             Double_t entries = hdata->GetEntries();
158             ++tested;
159             if(!(entries>0.))++empty;
160             AliDebug(AliQAv1::GetQADebugLevel(),Form("ESD hist name %s - entries %12.1g",hname.Data(),entries));
161             if(hname.Contains("hESDClusterMapSA") && entries>0.){
162               cluMapSA = kTRUE;
163               AliDebug(AliQAv1::GetQADebugLevel(),Form("Processing histogram %s",hname.Data()));
164               // Check if there are layers with anomalously low 
165               // contributing points to SA reconstructed tracks
166               for(Int_t k=1;k<7;k++){
167                 // check if the layer was skipped
168                 if(skipped[k-1]) continue;
169                 if(hdata->GetBinContent(k)<0.5*(entries/6.)){
170                   cluMapSA = kFALSE;
171                   AliDebug(AliQAv1::GetQADebugLevel(),Form("SA tracks have few points on layer %d - look at histogram hESDClustersSA",k));
172                 }
173               }  
174             }
175
176             else if(hname.Contains("hESDClusterMapMI") && entries>0.){
177               // Check if there are layers with anomalously low 
178               // contributing points to MI reconstructed tracks
179               AliDebug(AliQAv1::GetQADebugLevel(),Form("Processing histogram %s",hname.Data()));
180               cluMapMI = kTRUE;
181               for(Int_t k=1;k<7;k++){
182                 // check if the layer was skipped
183                 if(skipped[k-1]) continue;
184                 if(hdata->GetBinContent(k)<0.5*(entries/6.)){
185                   cluMapMI = kFALSE;
186                   AliDebug(AliQAv1::GetQADebugLevel(),Form("MI tracks have few points on layer %d - look at histogram hESDClustersMI",k));
187                 }
188               }  
189             }
190
191             else if(hname.Contains("hESDClustersMI") && entries>0.){
192               // Check if 6 clusters MI tracks are the majority
193               AliDebug(AliQAv1::GetQADebugLevel(),Form("Processing histogram %s",hname.Data()));
194               cluMI = kTRUE;
195               Double_t maxlaytracks = hdata->GetBinContent(7-nskipped);
196               for(Int_t k=2; k<7-nskipped; k++){
197                 if(hdata->GetBinContent(k)>maxlaytracks){
198                   cluMI = kFALSE;
199                   AliDebug(AliQAv1::GetQADebugLevel(),Form("MI Tracks with %d clusters are more than tracks with %d clusters. Look at histogram hESDClustersMI",k-1,6-nskipped));
200                 }
201               }
202             }
203
204             else if(hname.Contains("hESDClustersSA") && entries>0.){
205               // Check if 6 clusters SA tracks are the majority
206               AliDebug(AliQAv1::GetQADebugLevel(),Form("Processing histogram %s",hname.Data()));
207               cluSA = kTRUE;
208               Double_t maxlaytracks = hdata->GetBinContent(7-nskipped);
209               for(Int_t k=2; k<7-nskipped; k++){
210                 if(hdata->GetBinContent(k)>maxlaytracks){
211                   cluSA = kFALSE;
212                   AliDebug(AliQAv1::GetQADebugLevel(), Form("SA Tracks with %d clusters are more than tracks with %d clusters. Look at histogram hESDClustersSA",k-1,6-nskipped));
213                 }
214               }
215             }
216
217             else if(hname.Contains("hSPDVertexZ") && entries>0.){
218               // Check if average Z vertex coordinate is -5 < z < 5 cm
219               AliDebug(AliQAv1::GetQADebugLevel(),Form("Processing histogram %s",hname.Data()));
220               verSPDZ = kTRUE;
221               if(hdata->GetMean()<-5. && hdata->GetMean()>5.){
222                 verSPDZ = kFALSE;
223                 AliDebug(AliQAv1::GetQADebugLevel(),Form("Average z vertex coordinate is at z= %10.4g cm",hdata->GetMean()));
224               }
225             }
226           }
227           else{
228             AliError("ESD Checker - invalid data type");
229           }
230         }
231         rv[specie] = 0.;
232         if(tested>0){
233           if(tested == empty){
234             rv[specie] = 2500.; // set to error
235             AliWarning(Form("All ESD histograms are empty - specie=%d",specie));
236           }
237           else {
238             rv[specie] = 2500.-1500.*(static_cast<Double_t>(tested-empty)/static_cast<Double_t>(tested)); // INFO if all histos are filled
239             if(cluMapSA)rv[specie]-=200.;
240             if(cluMapMI)rv[specie]-=200.;
241             if(cluMI)rv[specie]-=200.;
242             if(cluSA)rv[specie]-=200.;
243             if(verSPDZ)rv[specie]-=199.;  // down to 1 if everything is OK
244           }
245         }
246       }
247       //     AliDebug(AliQAv1::GetQADebugLevel(), Form("ESD - Tested %d histograms, Return value %f \n",tested,rv[specie]));
248       AliInfo(Form("ESD - Tested %d histograms, Return value %f \n",tested,rv[specie]));
249     }
250   }  // end of ESD QA
251   else{
252     
253     //____________________________________________________________________________
254
255     Double_t spdCheck[AliRecoParam::kNSpecies] ;
256     Double_t sddCheck[AliRecoParam::kNSpecies] ;
257     Double_t ssdCheck[AliRecoParam::kNSpecies] ;
258
259
260
261     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
262       if ( !AliQAv1::Instance()->IsEventSpecieSet(specie)) continue; 
263       if ( AliQAv1::Instance()->IsEventSpecieSet(specie) ) {
264         Double_t histotot=list[specie]->GetEntries();
265         if(histotot!=0)
266           {
267             spdCheck[specie]=0.;
268             sddCheck[specie]=0.;
269             ssdCheck[specie]=0.;
270             rv[specie] = 0.0 ;// 
271             //pixel
272             if(fDet == 0 || fDet == 1) {
273               fSPDChecker->SetTaskOffset(fSPDOffset);
274               //printf("spdoffset = %i \n",fSPDOffset );
275               Double_t histoSPD=double(GetSPDHisto());
276               if(AliITSQADataMakerRec::AreEqual(histoSPD,0)==kFALSE){
277                 Double_t *stepSPD=new Double_t[AliQAv1::kNBIT];
278                 CreateStepForBit(histoSPD,stepSPD);
279                 fSPDChecker->SetStepBit(stepSPD);
280                 spdCheck[specie] = fSPDChecker->Check(index, list[specie], recoParam);
281                 if(spdCheck[specie]>fUpTestValue[AliQAv1::kFATAL]||spdCheck[specie]<0.)
282                   {
283                     AliInfo(Form("SPD check result for %s  is out of range (%f)!!! Retval of specie %s is sit to -1\n ",AliQAv1::GetAliTaskName(index),spdCheck[specie],AliRecoParam::GetEventSpecieName(specie)));
284                     spdCheck[specie]=fUpTestValue[AliQAv1::kFATAL];
285                   }
286                 delete []stepSPD;
287               }//end check SPD entries
288               else{spdCheck[specie]=fUpTestValue[AliQAv1::kFATAL];}
289               rv[specie]=spdCheck[specie];
290             }//end SPD check
291             //drift
292             if(fDet == 0 || fDet == 2) {
293               fSDDChecker->SetTaskOffset(fSDDOffset);
294               Double_t histoSDD=double(GetSDDHisto());
295               if(AliITSQADataMakerRec::AreEqual(histoSDD,0)==kFALSE){
296                 Double_t *stepSDD=new Double_t[AliQAv1::kNBIT];
297                 CreateStepForBit(histoSDD,stepSDD);
298                 fSDDChecker->SetStepBit(stepSDD);
299                 sddCheck[specie] = fSDDChecker->Check(index, list[specie], recoParam);  
300                 if(sddCheck[specie]>fUpTestValue[AliQAv1::kFATAL]||sddCheck[specie]<0.)
301                   {
302                     AliInfo(Form("SDD check result for %s  is out of range (%f)!!! Retval of specie %s is sit to -1\n ",AliQAv1::GetAliTaskName(index),sddCheck[specie],AliRecoParam::GetEventSpecieName(specie)));
303                     sddCheck[specie]=fUpTestValue[AliQAv1::kFATAL];
304                   }
305                 delete []stepSDD;
306               }//end check SDD entries
307               else{ssdCheck[specie]=fUpTestValue[AliQAv1::kFATAL];}
308               if(sddCheck[specie]>rv[specie])rv[specie]=sddCheck[specie];  
309             }//end SDD
310             //strip
311             if(fDet == 0 || fDet == 3) {
312               fSSDChecker->SetTaskOffset(fSSDOffset);
313               Double_t histoSSD=double(GetSSDHisto());
314               if(AliITSQADataMakerRec::AreEqual(histoSSD,0)==kFALSE){
315               Double_t *stepSSD=new Double_t[AliQAv1::kNBIT];
316               CreateStepForBit(histoSSD,stepSSD);
317               fSSDChecker->SetStepBit(stepSSD);
318               ssdCheck[specie] = fSSDChecker->Check(index, list[specie], recoParam);
319               if(ssdCheck[specie]>fUpTestValue[AliQAv1::kFATAL]||ssdCheck[specie]<0.)
320                 {
321                   AliInfo(Form("SSD check result for %s is out of range (%f)!!! Retval of specie %s is sit to -1\n ",AliQAv1::GetAliTaskName(index),ssdCheck[specie],AliRecoParam::GetEventSpecieName(specie)));
322                   ssdCheck[specie]=fUpTestValue[AliQAv1::kFATAL];
323                 }
324               delete [] stepSSD;
325               }//end check SSD entries
326               else{ssdCheck[specie]=fUpTestValue[AliQAv1::kFATAL];}
327               if(ssdCheck[specie]>rv[specie])rv[specie]=ssdCheck[specie];
328             }//end SSD
329             
330             AliInfo(Form("Check result for %s: \n\t  SPD %f \n\t  SDD %f \n\t  SSD %f \n Check result %f \n ",AliQAv1::GetAliTaskName(index),spdCheck[specie],sddCheck[specie],ssdCheck[specie],rv[specie]));
331             // here merging part for common ITS QA result
332             // 
333           }//end entries
334       }//end if event specie
335     }//end for
336   }
337 }
338
339 //____________________________________________________________________________
340 void AliITSQAChecker::SetTaskOffset(Int_t SPDOffset, Int_t SDDOffset, Int_t SSDOffset)
341 {
342   //Setting the 3 offsets for each task called
343   fSPDOffset = SPDOffset;
344   fSDDOffset = SDDOffset;
345   fSSDOffset = SSDOffset;
346 }
347
348 //____________________________________________________________________________
349 void AliITSQAChecker::SetHisto(Int_t SPDhisto, Int_t SDDhisto, Int_t SSDhisto)
350 {
351   //Setting the 3 offsets for each task called
352   fSPDHisto = SPDhisto;
353   fSDDHisto = SDDhisto;
354   fSSDHisto = SSDhisto;
355 }
356
357  //____________________________________________________________________________
358  void AliITSQAChecker::SetDetTaskOffset(Int_t subdet,Int_t offset)
359  {
360    switch(subdet){
361    case 1:
362      SetSPDTaskOffset(offset);
363      break;
364    case 2:
365      SetSDDTaskOffset(offset);
366      break;
367    case 3:
368      SetSSDTaskOffset(offset);
369      break;
370    default:
371      AliWarning("No specific (SPD,SDD or SSD) subdetector correspond to to this number!!! all offsets set to zero for all the detectors\n");
372      SetTaskOffset(0, 0, 0);
373      break;
374    }
375  }
376
377  //____________________________________________________________________________
378  void AliITSQAChecker::SetDetHisto(Int_t subdet,Int_t histo)
379  {
380    switch(subdet){
381    case 1:
382      SetSPDHisto(histo);
383      break;
384    case 2:
385      SetSDDHisto(histo);
386      break;
387    case 3:
388      SetSSDHisto(histo);
389      break;
390    default:
391      AliWarning("No specific (SPD,SDD or SSD) subdetector correspond to to this number!!! all offsets set to zero for all the detectors\n");
392      SetHisto(0, 0, 0);
393      break;
394    }
395  }
396
397 //_____________________________________________________________________________
398
399 void AliITSQAChecker::InitQACheckerLimits()
400 {
401   
402   AliInfo("Setting of tolerance values\n");
403
404   Float_t lowtolerancevalue[AliQAv1::kNBIT];
405
406   Float_t hightolerancevalue[AliQAv1::kNBIT];
407   for(Int_t bit=0;bit<AliQAv1::kNBIT;bit++)
408     {
409       lowtolerancevalue[bit]=(bit*1000.);
410       hightolerancevalue[bit]=((bit+1.)*1000.);
411     }
412   SetHiLo(hightolerancevalue,lowtolerancevalue);
413   //  AliInfo(Form("Range Value  \n INFO    -> %f <  value <  %f \n WARNING -> %f <  value <= %f \n ERROR   -> %f <  value <= %f \n FATAL   -> %f <= value <  %f \n", fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO], fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING], fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR], fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]  ));
414
415   if(fDet == 0 || fDet == 1) {
416     fSPDChecker->SetSPDLimits( lowtolerancevalue,hightolerancevalue );
417   }
418   if(fDet == 0 || fDet == 2) {
419     fSDDChecker->SetSDDLimits( lowtolerancevalue,hightolerancevalue );
420   }
421   if(fDet == 0 || fDet == 3) {
422     fSSDChecker->SetSSDLimits( lowtolerancevalue,hightolerancevalue );
423   }
424
425
426   
427 }
428
429
430 //_____________________________________________________________________________
431
432 void AliITSQAChecker::CreateStepForBit(Double_t histonumb,Double_t *steprange)
433 {
434   for(Int_t bit=0;bit < AliQAv1::kNBIT; bit++)
435     {       
436       //printf("%i\t %f \t %f \t %f \n",bit, fUpTestValue[bit],fLowTestValue[AliQAv1::kINFO],histonumb);
437       steprange[bit]=double((fUpTestValue[bit] - fLowTestValue[AliQAv1::kINFO])/histonumb);
438       //printf("%i\t %f \t %f \t %f \t %f\n",bit, fUpTestValue[bit],fLowTestValue[AliQAv1::kINFO],histonumb,steprange[bit] );
439     }
440   //AliInfo(Form("StepBitValue:numner of histo %f\n\t INFO %f \t WARNING %f \t ERROR %f \t FATAL %f \n",histonumb, steprange[AliQAv1::kINFO],steprange[AliQAv1::kWARNING],steprange[AliQAv1::kERROR],steprange[AliQAv1::kFATAL]));
441 }
442
443
444 //_____________________________________________________________________________
445 void AliITSQAChecker::SetQA(AliQAv1::ALITASK_t index, Double_t * value) const
446 {
447
448   AliQAv1 * qa = AliQAv1::Instance(index) ;
449
450
451   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
452
453     if (! qa->IsEventSpecieSet(AliRecoParam::ConvertIndex(specie)))
454       continue ;
455     if (  value == NULL ) { // No checker is implemented, set all QA to Fatal
456       qa->Set(AliQAv1::kFATAL, specie) ; 
457     } else {
458       if ( value[specie] > fLowTestValue[AliQAv1::kFATAL] && value[specie] <= fUpTestValue[AliQAv1::kFATAL] ) 
459         qa->Set(AliQAv1::kFATAL, AliRecoParam::ConvertIndex(specie)) ; 
460       else if ( value[specie] > fLowTestValue[AliQAv1::kERROR] && value[specie] <= fUpTestValue[AliQAv1::kERROR]  )
461         qa->Set(AliQAv1::kERROR, AliRecoParam::ConvertIndex(specie)) ; 
462       else if ( value[specie] > fLowTestValue[AliQAv1::kWARNING] && value[specie] <= fUpTestValue[AliQAv1::kWARNING]  )
463         qa->Set(AliQAv1::kWARNING, AliRecoParam::ConvertIndex(specie)) ;
464       else if ( value[specie] > fLowTestValue[AliQAv1::kINFO] && value[specie] <= fUpTestValue[AliQAv1::kINFO] ) 
465         qa->Set(AliQAv1::kINFO, AliRecoParam::ConvertIndex(specie)) ;   
466       //else if(value[specie]==0) qa->Set(AliQAv1::kFATAL, AliRecoParam::ConvertIndex(specie)) ; //no ckeck has been done
467     }
468     qa->ShowStatus(AliQAv1::kITS,index,AliRecoParam::ConvertIndex(specie));
469   }//end for
470
471 }
472
473