]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQAChecker.cxx
Implementing destructor
[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 Double_t * AliITSQAChecker::Check(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     Double_t * rv = new Double_t[AliRecoParam::kNSpecies] ; 
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     return rv ; 
251   }  // end of ESD QA
252   
253   Double_t * retval = new Double_t[AliRecoParam::kNSpecies] ; 
254   //____________________________________________________________________________
255
256   Double_t spdCheck[AliRecoParam::kNSpecies] ;
257   Double_t sddCheck[AliRecoParam::kNSpecies] ;
258   Double_t ssdCheck[AliRecoParam::kNSpecies] ;
259
260
261
262     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
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             retval[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
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               //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]) );
287               delete []stepSPD;
288               retval[specie]=spdCheck[specie];
289             }
290             //drift
291             if(fDet == 0 || fDet == 2) {
292               fSDDChecker->SetTaskOffset(fSDDOffset);
293               Double_t histoSDD=double(GetSDDHisto());
294               Double_t *stepSDD=new Double_t[AliQAv1::kNBIT];
295               CreateStepForBit(histoSDD,stepSDD);
296               fSDDChecker->SetStepBit(stepSDD);
297               sddCheck[specie] = fSDDChecker->Check(index, list[specie], recoParam);
298               if(sddCheck[specie]>fUpTestValue[AliQAv1::kFATAL]||sddCheck[specie]<0.)
299                 {
300                   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)));
301                   sddCheck[specie]=fUpTestValue[AliQAv1::kFATAL];
302                 }
303               //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]) );
304               delete []stepSDD;
305               if(sddCheck[specie]>retval[specie])retval[specie]=sddCheck[specie];  
306             }
307             //strip
308             if(fDet == 0 || fDet == 3) {
309               fSSDChecker->SetTaskOffset(fSSDOffset);
310               Double_t histoSSD=double(GetSSDHisto());
311               Double_t *stepSSD=new Double_t[AliQAv1::kNBIT];
312               CreateStepForBit(histoSSD,stepSSD);
313               fSSDChecker->SetStepBit(stepSSD);
314               ssdCheck[specie] = fSSDChecker->Check(index, list[specie], recoParam);
315               if(ssdCheck[specie]>fUpTestValue[AliQAv1::kFATAL]||ssdCheck[specie]<0.)
316                 {
317                   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)));
318                   ssdCheck[specie]=fUpTestValue[AliQAv1::kFATAL];
319                 }
320               //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]) );
321               delete [] stepSSD;
322               if(ssdCheck[specie]>retval[specie])retval[specie]=ssdCheck[specie];
323             }
324             
325             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],retval[specie]));
326             // here merging part for common ITS QA result
327             // 
328           }//end entries
329       }//end if event specie
330     }//end for
331     return retval;  
332 }
333
334
335 //____________________________________________________________________________
336 void AliITSQAChecker::SetTaskOffset(Int_t SPDOffset, Int_t SDDOffset, Int_t SSDOffset)
337 {
338   //Setting the 3 offsets for each task called
339   fSPDOffset = SPDOffset;
340   fSDDOffset = SDDOffset;
341   fSSDOffset = SSDOffset;
342 }
343
344 //____________________________________________________________________________
345 void AliITSQAChecker::SetHisto(Int_t SPDhisto, Int_t SDDhisto, Int_t SSDhisto)
346 {
347   //Setting the 3 offsets for each task called
348   fSPDHisto = SPDhisto;
349   fSDDHisto = SDDhisto;
350   fSSDHisto = SSDhisto;
351 }
352
353  //____________________________________________________________________________
354  void AliITSQAChecker::SetDetTaskOffset(Int_t subdet,Int_t offset)
355  {
356    switch(subdet){
357    case 1:
358      SetSPDTaskOffset(offset);
359      break;
360    case 2:
361      SetSDDTaskOffset(offset);
362      break;
363    case 3:
364      SetSSDTaskOffset(offset);
365      break;
366    default:
367      AliWarning("No specific (SPD,SDD or SSD) subdetector correspond to to this number!!! all offsets set to zero for all the detectors\n");
368      SetTaskOffset(0, 0, 0);
369      break;
370    }
371  }
372
373  //____________________________________________________________________________
374  void AliITSQAChecker::SetDetHisto(Int_t subdet,Int_t histo)
375  {
376    switch(subdet){
377    case 1:
378      SetSPDHisto(histo);
379      break;
380    case 2:
381      SetSDDHisto(histo);
382      break;
383    case 3:
384      SetSSDHisto(histo);
385      break;
386    default:
387      AliWarning("No specific (SPD,SDD or SSD) subdetector correspond to to this number!!! all offsets set to zero for all the detectors\n");
388      SetHisto(0, 0, 0);
389      break;
390    }
391  }
392
393 //_____________________________________________________________________________
394
395 void AliITSQAChecker::InitQACheckerLimits()
396 {
397   
398   AliInfo("Setting of tolerance values\n");
399
400   Float_t lowtolerancevalue[AliQAv1::kNBIT];
401
402   Float_t hightolerancevalue[AliQAv1::kNBIT];
403   for(Int_t bit=0;bit<AliQAv1::kNBIT;bit++)
404     {
405       lowtolerancevalue[bit]=(bit*1000.);
406       hightolerancevalue[bit]=((bit+1.)*1000.);
407     }
408   SetHiLo(hightolerancevalue,lowtolerancevalue);
409   //  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]  ));
410
411   if(fDet == 0 || fDet == 1) {
412     fSPDChecker->SetSPDLimits( lowtolerancevalue,hightolerancevalue );
413   }
414   if(fDet == 0 || fDet == 2) {
415     fSDDChecker->SetSDDLimits( lowtolerancevalue,hightolerancevalue );
416   }
417   if(fDet == 0 || fDet == 3) {
418     fSSDChecker->SetSSDLimits( lowtolerancevalue,hightolerancevalue );
419   }
420
421
422   
423 }
424
425
426 //_____________________________________________________________________________
427
428 void AliITSQAChecker::CreateStepForBit(Double_t histonumb,Double_t *steprange)
429 {
430   for(Int_t bit=0;bit < AliQAv1::kNBIT; bit++)
431     {       
432       //printf("%i\t %f \t %f \t %f \n",bit, fUpTestValue[bit],fLowTestValue[AliQAv1::kINFO],histonumb);
433       steprange[bit]=double((fUpTestValue[bit] - fLowTestValue[AliQAv1::kINFO])/histonumb);
434       //printf("%i\t %f \t %f \t %f \t %f\n",bit, fUpTestValue[bit],fLowTestValue[AliQAv1::kINFO],histonumb,steprange[bit] );
435     }
436   //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]));
437 }
438
439
440 //_____________________________________________________________________________
441 void AliITSQAChecker::SetQA(AliQAv1::ALITASK_t index, Double_t * value) const
442 {
443
444   AliQAv1 * qa = AliQAv1::Instance(index) ;
445
446
447   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
448
449     if (! qa->IsEventSpecieSet(AliRecoParam::ConvertIndex(specie)))
450       continue ;
451     if (  value == NULL ) { // No checker is implemented, set all QA to Fatal
452       qa->Set(AliQAv1::kFATAL, specie) ; 
453     } else {
454       if ( value[specie] > fLowTestValue[AliQAv1::kFATAL] && value[specie] <= fUpTestValue[AliQAv1::kFATAL] ) 
455         qa->Set(AliQAv1::kFATAL, AliRecoParam::ConvertIndex(specie)) ; 
456       else if ( value[specie] > fLowTestValue[AliQAv1::kERROR] && value[specie] <= fUpTestValue[AliQAv1::kERROR]  )
457         qa->Set(AliQAv1::kERROR, AliRecoParam::ConvertIndex(specie)) ; 
458       else if ( value[specie] > fLowTestValue[AliQAv1::kWARNING] && value[specie] <= fUpTestValue[AliQAv1::kWARNING]  )
459         qa->Set(AliQAv1::kWARNING, AliRecoParam::ConvertIndex(specie)) ;
460       else if ( value[specie] > fLowTestValue[AliQAv1::kINFO] && value[specie] <= fUpTestValue[AliQAv1::kINFO] ) 
461         qa->Set(AliQAv1::kINFO, AliRecoParam::ConvertIndex(specie)) ;   
462       //else if(value[specie]==0) qa->Set(AliQAv1::kFATAL, AliRecoParam::ConvertIndex(specie)) ; //no ckeck has been done
463     }
464     qa->ShowStatus(AliQAv1::kITS,index,AliRecoParam::ConvertIndex(specie));
465   }//end for
466
467 }
468
469