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