]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/AliHMPIDQAChecker.cxx
Fix the config file
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDQAChecker.cxx
CommitLineData
e4e6ae61 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>
fdd7f404 29#include <TF1.h>
e4e6ae61 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"
4e25ac79 38#include "AliQAv1.h"
b990acb3 39#include "AliQAChecker.h"
40#include "AliHMPIDQAChecker.h"
fdd7f404 41#include "AliCDBEntry.h"
fc7e0df2 42#include "AliQAManager.h"
e4e6ae61 43
b990acb3 44ClassImp(AliHMPIDQAChecker)
e4e6ae61 45
fdd7f404 46//_________________________________________________________________
4e25ac79 47Double_t * AliHMPIDQAChecker::Check(AliQAv1::ALITASK_t index, TObjArray ** list)
fdd7f404 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
57acd2d2 53 Double_t * check = new Double_t[AliRecoParam::kNSpecies] ;
54
57acd2d2 55 AliInfo(Form("Fix needed ....."));
4e25ac79 56 char * detOCDBDir = Form("HMPID/%s/%s", AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ;
fc7e0df2 57 AliCDBEntry *QARefRec = AliQAManager::QAManager()->Get(detOCDBDir);
fdd7f404 58 if( !QARefRec){
57acd2d2 59 AliInfo("QA reference data NOT retrieved for Reconstruction check. No HMPIDChecker ...exiting");
60 return check;
fdd7f404 61 }
62
63// checking for empy histograms
57acd2d2 64 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
65 check[specie] = 1.0;
4e25ac79 66 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) )
57acd2d2 67 continue ;
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 // checking rec points
4e25ac79 74 if(index == AliQAv1::kREC) check[specie] = CheckRecPoints(list[specie],(TObjArray *)QARefRec->GetObject());
fdd7f404 75
57acd2d2 76 //default check response. It will be changed when reasonable checks will be considered
77 else check[specie] = 0.7 ; // /-> Corresponds to kINFO see AliQACheckerBase::Run
78 }
fdd7f404 79 return check;
80
81}
82//_________________________________________________________________
83Double_t AliHMPIDQAChecker::CheckEntries(TObjArray * list) const
84{
85 //
86 // check on the QA histograms on the input list:
87 //
88
89 Double_t test = 0.0 ;
90 Int_t count = 0 ;
91
92 if (list->GetEntries() == 0){
93 test = 1. ; // nothing to check
94 }
95 else {
96 TIter next(list) ;
97 TH1 * hdata ;
98 count = 0 ;
99 while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
100 if (hdata) {
101 Double_t rv = 0.;
102 //Printf("hitogram %s has entries: %f ",hdata->GetName(),hdata->GetEntries());
103 if(hdata->GetEntries()>0)rv=1;
104 count++ ;
105 test += rv ;
106 }
107 else{
108 AliError("Data type cannot be processed") ;
109 }
110
111 }
112 if (count != 0) {
113 if (test==0) {
114 AliWarning("Histograms are booked for THIS specific Task, but they are all empty: setting flag to kWARNING");
115 test = 0.; //upper limit value to set kWARNING flag for a task
116 }
117 else test = 1 ;
118 }
119 }
120
121 return test ;
122}
123//_________________________________________________________________
124
125Double_t AliHMPIDQAChecker::CheckRecPoints(TObjArray *listrec, TObjArray *listref) const
126{
0f46609e 127 //
fbe0560f 128 // check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
129 //
fdd7f404 130
131 Float_t checkresponse = 0;
132
133 Float_t counter = 0 ;
134 TIter next(listrec) ;
135 TH1* histo;
136 while ( (histo = dynamic_cast<TH1 *>(next())) ) {
3762c8fb 137 //PH The histogram should have at least 3 bins with entries
138 Int_t nbinsabove = 0;
326d0a78 139 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);
3762c8fb 140 if (histo->GetBinContent(ibin)>0) nbinsabove++;
141 }
142
143 if( nbinsabove < 3 ) counter++;
fdd7f404 144 else {
145 TString h = histo->GetTitle();
146 if(h.Contains("Zoom")){
326d0a78 147 histo->Fit("expo","LQ0","",1,50);
fdd7f404 148 if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
149 }
150 if(h.Contains("size MIP")) if(TMath::Abs(histo->GetMean()-5) > 2) counter++;
151 if(h.Contains("size Phots")) if(TMath::Abs(histo->GetMean()-2) > 2) counter++;
152 if(h.Contains("distribution")) if(histo->KolmogorovTest((TH1F *)listref->At(0))<0.8) counter++;
5379c4a3 153 AliDebug(AliQAv1::GetQADebugLevel(),Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));
fdd7f404 154 }
155 }
156 Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listrec->GetEntries())
157
158 if(response < 0.1) checkresponse = 0.7; // <10% of the check histograms show a failing check -> Corresponds to kINFO see AliQACheckerBase::Run
159 else if(response < 0.5) checkresponse = 0.4; // 50% of the check histograms show a failing check -> Corresponds to kWARNING see AliQACheckerBase::Run
160 else checkresponse = 0.001; // > 50% of the check histograms show a failing check -> Corresponds to kERROR see AliQACheckerBase::Run
161 return checkresponse;
162}
163//________________________________________________________________
e4e6ae61 164