]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDQAChecker.cxx
Introducing event specie in QA (Yves)
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDQAChecker.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16
17 /* $Id$ */
18
19 //...
20 //  Checks the quality assurance. 
21 //  By comparing with reference data
22 //  Skeleton for HMPID
23 //...
24
25 // --- ROOT system ---
26 #include <TClass.h>
27 #include <TH1F.h> 
28 #include <TH1I.h> 
29 #include <TF1.h> 
30 #include <TIterator.h> 
31 #include <TKey.h> 
32 #include <TFile.h> 
33
34 // --- Standard library ---
35
36 // --- AliRoot header files ---
37 #include "AliLog.h"
38 #include "AliQA.h"
39 #include "AliQAChecker.h"
40 #include "AliHMPIDQAChecker.h"
41 #include "AliCDBEntry.h"
42 #include "AliCDBManager.h"
43
44 ClassImp(AliHMPIDQAChecker)
45
46 //____________________________________________________________________________
47 Double_t * AliHMPIDQAChecker::Check(AliQA::ALITASK_t /*index*/)
48 {
49   Double_t * rv = new Double_t[AliRecoParam::kNSpecies] ; 
50   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) 
51     rv[specie] = 0.0 ; 
52   return rv ;  
53 }
54
55 //_________________________________________________________________
56 Double_t * AliHMPIDQAChecker::Check(AliQA::ALITASK_t index, TObjArray ** list) 
57 {
58 //
59 // Main check function: Depending on the TASK, different checks are applied
60 // At the moment:       check for empty histograms and checks for RecPoints
61
62   Double_t * check = new Double_t[AliRecoParam::kNSpecies] ; 
63   
64 //YS THIS IS NOT CORRECT
65   AliInfo(Form("Fix needed ....."));
66   AliCDBEntry *QARefRec = AliCDBManager::Instance()->Get("HMPID/QARef/Rec");
67   if( !QARefRec){
68     AliInfo("QA reference data NOT retrieved for Reconstruction check. No HMPIDChecker  ...exiting");
69     return check;
70   }
71 //YS THIS IS NOT CORRECT
72
73 // checking for empy histograms
74   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
75     check[specie] = 1.0;
76     if ( !AliQA::Instance()->IsEventSpecieSet(specie) ) 
77       continue ; 
78     if(CheckEntries(list[specie]) == 0)  {
79       AliWarning("histograms are empty");
80       check[specie] = 0.4;//-> Corresponds to kWARNING see AliQACheckerBase::Run
81     }
82   
83     // checking rec points
84     if(index == AliQA::kREC) check[specie] = CheckRecPoints(list[specie],(TObjArray *)QARefRec->GetObject());
85
86     //default check response. It will be changed when reasonable checks will be considered
87     else check[specie] = 0.7 ; // /-> Corresponds to kINFO see AliQACheckerBase::Run 
88   }
89   return check;
90
91 }
92 //_________________________________________________________________
93 Double_t AliHMPIDQAChecker::CheckEntries(TObjArray * list) const
94 {
95   //
96   //  check on the QA histograms on the input list: 
97   // 
98
99   Double_t test = 0.0  ;
100   Int_t count = 0 ; 
101   
102   if (list->GetEntries() == 0){  
103     test = 1. ; // nothing to check
104   }
105   else {
106     TIter next(list) ; 
107     TH1 * hdata ;
108     count = 0 ; 
109     while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
110       if (hdata) { 
111         Double_t rv = 0.;
112         //Printf("hitogram %s     has entries: %f ",hdata->GetName(),hdata->GetEntries());
113         if(hdata->GetEntries()>0)rv=1; 
114         count++ ; 
115         test += rv ; 
116       }
117       else{
118         AliError("Data type cannot be processed") ;
119       }
120       
121     }
122     if (count != 0) { 
123       if (test==0) {
124         AliWarning("Histograms are booked for THIS specific Task, but they are all empty: setting flag to kWARNING");
125         test = 0.;  //upper limit value to set kWARNING flag for a task
126       }
127       else test = 1 ;
128     }
129   }
130
131   return test ; 
132 }  
133 //_________________________________________________________________
134
135 Double_t AliHMPIDQAChecker::CheckRecPoints(TObjArray *listrec, TObjArray *listref) const
136 {
137   //
138   //  check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
139   //
140
141    Float_t checkresponse = 0;
142
143    Float_t counter = 0 ;
144    TIter next(listrec) ;
145    TH1* histo;
146    while ( (histo = dynamic_cast<TH1 *>(next())) ) {
147    if( histo->GetEntries() < 3 ) counter++;
148    else {
149     TString h = histo->GetTitle();
150     if(h.Contains("Zoom")){
151     histo->Fit("expo","Q0","",1,50);
152     if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
153     }
154     if(h.Contains("size  MIP"))   if(TMath::Abs(histo->GetMean()-5) > 2) counter++;
155     if(h.Contains("size  Phots")) if(TMath::Abs(histo->GetMean()-2) > 2) counter++;
156     if(h.Contains("distribution")) if(histo->KolmogorovTest((TH1F *)listref->At(0))<0.8) counter++;
157     AliDebug(1,Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));  
158    }
159   }
160  Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listrec->GetEntries())
161  
162  if(response < 0.1) checkresponse = 0.7;      // <10% of the check histograms show a failing check -> Corresponds to kINFO see AliQACheckerBase::Run
163  else if(response < 0.5) checkresponse = 0.4; //  50%  of the check histograms show a failing check -> Corresponds to kWARNING see AliQACheckerBase::Run
164  else checkresponse = 0.001;                  // > 50% of the check histograms show a failing check -> Corresponds to kERROR see AliQACheckerBase::Run
165  return checkresponse;
166 }
167 //________________________________________________________________
168