]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDQAChecker.cxx
4ad7d09e257a90b8845b991d04ae13c3966f9698
[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 "AliQAv1.h"
39 #include "AliQAChecker.h"
40 #include "AliHMPIDQAChecker.h"
41 #include "AliCDBEntry.h"
42 #include "AliQAManager.h"
43
44 ClassImp(AliHMPIDQAChecker)
45
46 //____________________________________________________________________________
47 Double_t * AliHMPIDQAChecker::Check(AliQAv1::ALITASK_t index)
48 {
49   
50  TObjArray **list;
51
52   Double_t * rv = new Double_t[AliRecoParam::kNSpecies] ;
53  Int_t count[AliRecoParam::kNSpecies]   = { 0 }; 
54  for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++){ // species loop 
55     rv[specie] = 1.0 ; 
56  
57    TObjArray *dataList;  Int_t iList=0;  
58   
59    if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
60    continue ; 
61    if (!fDataSubDir) {
62      rv[specie] = 0. ; // nothing to check
63    } 
64     else if (!fRefSubDir && !fRefOCDBSubDir) {
65         rv[specie] = -1 ; // no reference data
66     } 
67     else {
68       TList * keyList = fDataSubDir->GetListOfKeys() ; 
69       TIter next(keyList) ; 
70       TKey * key ;
71       count[specie] = 0 ; 
72       while ( (key = static_cast<TKey *>(next())) ) {
73         iList++;
74         TObject * odata = fDataSubDir->Get(key->GetName()) ;
75         dataList->AddAt(odata,iList);
76     }// while
77   }// else 
78   list[specie] = dataList;
79  } //species loop
80  Check(index,list);
81  return rv;            
82 }
83 //_________________________________________________________________
84 Double_t * AliHMPIDQAChecker::Check(AliQAv1::ALITASK_t index, TObjArray ** list) 
85 {
86 //
87 // Main check function: Depending on the TASK, different checks are applied
88 // At the moment:       check for empty histograms and checks for RecPoints
89
90   Double_t * check = new Double_t[AliRecoParam::kNSpecies] ; 
91   
92   AliInfo(Form("Fix needed ....."));
93   char * detOCDBDir = Form("HMPID/%s/%s", AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ; 
94   AliCDBEntry *QARefRec = AliQAManager::QAManager()->Get(detOCDBDir);
95   if(!QARefRec){
96     AliInfo("QA reference data NOT retrieved for Reconstruction check. No HMPIDChecker  ...exiting");
97     return check;
98   }
99
100   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
101     check[specie] = 1.0;
102     if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
103       continue ; 
104     // checking for empy histograms
105     if(CheckEntries(list[specie]) == 0)  {
106       AliWarning("histograms are empty");
107       check[specie] = 0.4;//-> Corresponds to kWARNING see AliQACheckerBase::Run
108     }
109   
110     //check sim
111     if(index == AliQAv1::kSIM) check[specie] = CheckSim(list[specie],(TObjArray *)QARefRec->GetObject());
112
113     // checking rec points
114     if(index == AliQAv1::kREC) check[specie] = CheckRec(list[specie],(TObjArray *)QARefRec->GetObject());
115
116     //default check response. It will be changed when reasonable checks will be considered
117     else check[specie] = 0.7 ; // /-> Corresponds to kINFO see AliQACheckerBase::Run 
118   } // species loop
119
120   return check;
121
122 }
123 //_________________________________________________________________
124 Double_t AliHMPIDQAChecker::CheckEntries(TObjArray * list) const
125 {
126   //
127   //  check on the QA histograms on the input list: 
128   // 
129
130   Double_t test = 0.0  ;
131   Int_t count = 0 ; 
132   
133   if (list->GetEntries() == 0){  
134     test = 1. ; // nothing to check
135   }
136   else {
137     TIter next(list) ; 
138     TH1 * hdata ;
139     count = 0 ; 
140     while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
141       if (hdata) { 
142         Double_t rv = 0.;
143         //Printf("hitogram %s     has entries: %f ",hdata->GetName(),hdata->GetEntries());
144         if(hdata->GetEntries()>0)rv=1; 
145         count++ ; 
146         test += rv ; 
147       }
148       else{
149         AliError("Data type cannot be processed") ;
150       }
151       
152     }
153     if (count != 0) { 
154       if (test==0) {
155         AliWarning("Histograms are booked for THIS specific Task, but they are all empty: setting flag to kWARNING");
156         test = 0.;  //upper limit value to set kWARNING flag for a task
157       }
158       else test = 1 ;
159     }
160   }
161
162   return test ; 
163 }  
164 //_________________________________________________________________
165 Double_t AliHMPIDQAChecker::CheckSim(TObjArray *listsim, TObjArray *listref) const
166 {
167   //
168   //  check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
169   //
170
171    Float_t checkresponse = 0;
172
173    Float_t counter = 0 ;
174    TIter next(listsim) ;
175    TH1* histo;
176    while ( (histo = dynamic_cast<TH1 *>(next())) ) {
177      //PH The histogram should have at least 3 bins with entries
178      Int_t nbinsabove = 0;
179      for (Int_t ibin=histo->FindBin(1); ibin<=histo->FindBin(50); ibin++) { //1,50 is the fit region, see histo->Fit("expo","Q0","",1,50);
180        if (histo->GetBinContent(ibin)>0) nbinsabove++;
181      }
182
183    if( nbinsabove < 3 ) counter++;
184    else {
185     TString h = histo->GetTitle();
186     if(h.Contains("Zoom")){
187     histo->Fit("expo","LQ0","",1,50);
188     if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
189     }
190     if(h.Contains("size  MIP"))   if(TMath::Abs(histo->GetMean()-5) > 2) counter++;
191     if(h.Contains("size  Phots")) if(TMath::Abs(histo->GetMean()-2) > 2) counter++;
192     if(h.Contains("distribution")) if(histo->KolmogorovTest((TH1F *)listref->At(0))<0.8) counter++;
193     AliDebug(AliQAv1::GetQADebugLevel(),Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));  
194    }
195   }
196  Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listsim->GetEntries())
197  
198  if(response < 0.1) checkresponse = 0.7;      // <10% of the check histograms show a failing check -> Corresponds to kINFO see AliQACheckerBase::Run
199  else if(response < 0.5) checkresponse = 0.4; //  50%  of the check histograms show a failing check -> Corresponds to kWARNING see AliQACheckerBase::Run
200  else checkresponse = 0.001;                  // > 50% of the check histograms show a failing check -> Corresponds to kERROR see AliQACheckerBase::Run
201  return checkresponse;
202 }
203
204 //___________________________________________________________________________________________________
205 Double_t AliHMPIDQAChecker::CheckRec(TObjArray *listrec, TObjArray *listref) const
206 {
207   //
208   //  check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
209   //
210
211    Float_t checkresponse = 0;
212
213    Float_t counter = 0 ;
214    TIter next(listrec) ;
215    TH1* histo;
216    while ( (histo = dynamic_cast<TH1 *>(next())) ) {
217      //PH The histogram should have at least 3 bins with entries
218      Int_t nbinsabove = 0;
219      for (Int_t ibin=histo->FindBin(1); ibin<=histo->FindBin(50); ibin++) { //1,50 is the fit region, see histo->Fit("expo","Q0","",1,50);
220        if (histo->GetBinContent(ibin)>0) nbinsabove++;
221      }
222
223    if( nbinsabove < 3 ) counter++;
224    else {
225     TString h = histo->GetTitle();
226     if(h.Contains("Zoom")){
227     histo->Fit("expo","LQ0","",1,50);
228     if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
229     }
230     if(h.Contains("size  MIP"))   if(TMath::Abs(histo->GetMean()-5) > 2) counter++;
231     if(h.Contains("size  Phots")) if(TMath::Abs(histo->GetMean()-2) > 2) counter++;
232     if(h.Contains("distribution")) if(histo->KolmogorovTest((TH1F *)listref->At(0))<0.8) counter++;
233     AliDebug(AliQAv1::GetQADebugLevel(),Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));  
234    }
235   }
236  Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listrec->GetEntries())
237  
238  if(response < 0.1) checkresponse = 0.7;      // <10% of the check histograms show a failing check -> Corresponds to kINFO see AliQACheckerBase::Run
239  else if(response < 0.5) checkresponse = 0.4; //  50%  of the check histograms show a failing check -> Corresponds to kWARNING see AliQACheckerBase::Run
240  else checkresponse = 0.001;                  // > 50% of the check histograms show a failing check -> Corresponds to kERROR see AliQACheckerBase::Run
241  return checkresponse;
242 }
243