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