]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQASDDChecker.cxx
Updated version of ITS QA Checker and related modifications (Melinda)
[u/mrichter/AliRoot.git] / ITS / AliITSQASDDChecker.cxx
1
2 /**************************************************************************
3  * Copyright(c) 2007-2009, ALICE Experiment at CERN, All rights reserved. *
4  *                                                                        *
5  * Author: The ALICE Off-line Project.                                    *
6  * Contributors are mentioned in the code where appropriate.              *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 /* $Id$ */
18
19 // *****************************************
20 //  Checks the quality assurance 
21 //  by comparing with reference data
22 //  P. Cerello Apr 2008
23 //  INFN Torino
24
25 // --- ROOT system ---
26 #include "TH1.h"
27 #include <TCanvas.h>
28 #include <TH1D.h>
29 #include <TH2.h>
30 #include <TF1.h>
31
32
33
34 // --- AliRoot header files ---
35 #include "AliITSQASDDChecker.h"
36 #include "AliLog.h"
37 #include "AliCDBEntry.h"
38 #include "AliQAManager.h"
39 #include "AliQACheckerBase.h"
40 #include "TSystem.h"
41 #include "AliITSCalibrationSDD.h"
42 #include "AliITSgeomTGeo.h"
43
44 ClassImp(AliITSQASDDChecker)
45 //__________________________________________________________________
46 AliITSQASDDChecker& AliITSQASDDChecker::operator = (const AliITSQASDDChecker& qac ) 
47 {
48   // Equal operator.
49   this->~AliITSQASDDChecker();
50   new(this) AliITSQASDDChecker(qac);
51   return *this;
52 }
53
54 AliITSQASDDChecker::~AliITSQASDDChecker() 
55 {
56   if(fStepBitSDD) 
57     {
58       delete[] fStepBitSDD ;
59       fStepBitSDD = NULL;
60     }
61   if(fLowSDDValue)
62     {
63       delete[]fLowSDDValue;
64       fLowSDDValue=NULL;
65     }
66   if(fHighSDDValue)
67     { 
68       delete[]fHighSDDValue;
69       fHighSDDValue=NULL;
70     }
71   if(fCalibration)
72     {
73       delete fCalibration;
74       fCalibration=NULL;
75     }
76 } // dtor
77
78 //__________________________________________________________________
79 Double_t AliITSQASDDChecker::Check(AliQAv1::ALITASK_t index, TObjArray * list/*, AliDetectorRecoParam *recoparam*/) 
80 {  
81   AliInfo(Form("AliITSQASDDChecker called with offset: %d\n", fSubDetOffset) );
82
83   AliDebug(1,Form("AliITSQASDDChecker called with offset: %d\n", fSubDetOffset));
84
85   Double_t test = 0.;
86   TH1 *hdata=NULL;
87   Double_t entries=0.;
88   Double_t entries2[2];
89   for(Int_t i=0;i<2;i++)entries2[i]=0.;
90   switch(index)
91     {
92
93     case AliQAv1::kRAW:
94       AliInfo(Form("Check on %s\n",AliQAv1::GetAliTaskName(index)));
95       if(!fCalibration){
96         AliCDBEntry *calibSDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD");
97         Bool_t cacheStatus = AliCDBManager::Instance()->GetCacheFlag();
98         if(!calibSDD)
99           {
100             AliError("Calibration object retrieval failed! SDD will not be processed");
101             fCalibration = NULL;
102             test= fHighSDDValue[AliQAv1::kWARNING];
103           }
104         fCalibration = (TObjArray *)calibSDD->GetObject();
105         
106         if(!cacheStatus)calibSDD->SetObject(NULL);
107         calibSDD->SetOwner(kTRUE);
108         if(!cacheStatus)
109           {
110             delete calibSDD;
111           }
112       }
113       AliInfo("Calib SDD Created\n ");
114       
115       if (list->GetEntries() == 0.){ //check if the list is empty
116         //printf("test = %f \t value %f\n",test,fHighSDDValue[AliQAv1::kFATAL]);
117         test=test+fHighSDDValue[AliQAv1::kFATAL];
118         break;
119       }//end if getentries
120       else{
121         TIter next(list);
122         Int_t offset = 0;
123         for(offset =0;offset < fSubDetOffset; offset++){hdata = dynamic_cast<TH1*>(next());}//end for
124         Int_t emptymodules[2];
125         Int_t filledmodules[2];
126         Int_t emptyladders[2];
127         Int_t filledladders[2];
128         for(Int_t i=0;i<2;i++){
129           emptymodules[i]=0;
130           filledmodules[i]=0;
131           emptyladders[i]=0;
132           filledladders[i]=0;
133         }
134         TH1 *hmodule=NULL;
135         TH2 *hlayer[2];
136         for(Int_t i=0;i<2;i++)hlayer[i]=NULL;    
137         while( hdata = dynamic_cast<TH1* >(next()) ){
138           if (hdata){
139             TString hname=hdata->GetName();
140             if(hname.Contains("SDDchargeMap"))continue;
141             if(hname.Contains("SDDModPattern")){
142               hmodule=hdata;
143               entries= hdata->GetEntries();
144               if(entries==0){
145                 AliWarning(Form("===================>>>>>> No entries in  %s \n",hname.Data()));
146                 //printf("test = %f \t value %f\n",test,fHighSDDValue[AliQAv1::kFATAL]);
147                 test=test+fStepBitSDD[AliQAv1::kFATAL];
148               }//endif entries
149               else{
150                 int modmax=hdata->GetNbinsX();
151                 Int_t empty=0;
152                 Int_t filled=0;
153                 Double_t content=0;
154                 for(Int_t i=1;i<=modmax;i++){
155                   content=hdata->GetBinContent(i);
156                   if(content==0.){empty++;}
157                   else if(content!=0.){filled++;}
158                 }//end for
159                 AliInfo(Form(" %s : empty modules %i \t filled modules %i",hname.Data(), empty, filled));
160               }//end else pattern entries !=0         
161             }//end if modpattern
162             else if(hname.Contains("SDDphizL3")||hname.Contains("SDDphizL4")){
163               Int_t layer=0;
164               if(hname.Contains("3"))layer=0;
165               else  if(hname.Contains("4"))layer=1;
166               entries2[layer]=hdata->GetEntries();
167               if(entries2[layer]==0){
168                 AliWarning(Form("===================>>>>>> No entries in  %s \n",hname.Data()));
169                 //printf("test = %f \t value %f\n",test,fStepBitSDD[AliQAv1::kFATAL]);
170                 test=test+fStepBitSDD[AliQAv1::kFATAL];
171                 if(entries==0){ 
172                   //return test; 
173                   //break;
174                 }
175               }//end if getentries
176               else{
177                 Int_t layer=0;
178                 if(hname.Contains("3"))layer=0;
179                 else  if(hname.Contains("4"))layer=1;
180                 hlayer[layer]=dynamic_cast<TH2*>(hdata); 
181                 hlayer[layer]->RebinX(2);
182                 int modmay=hlayer[layer]->GetNbinsY();
183                 TH1D* hproj= hlayer[layer]->ProjectionY();
184                 Double_t ladcontent=0;
185                 for(Int_t i=1;i<=modmay;i++) {//loop on the ladders
186                   ladcontent=hproj->GetBinContent(i);
187                   if(ladcontent==0){emptyladders[layer]++;}
188                   else if(ladcontent!=0){filledladders[layer]++;} 
189                 }//end for
190                 AliInfo(Form(" %s : empty ladders %i \t filled ladders %i\n",hname.Data(), emptyladders[layer], filledladders[layer]));//end else layer 3
191                 delete hproj;
192                 hproj=NULL;     
193               }//end else entries !=0         
194             }//end if layer 3
195           }//end if hdata       
196         }//end while
197         if(entries==0.&&entries2[0]==0.&&entries2[1]==0.) break;
198         else{
199           if(hmodule||(hlayer[0]&&hlayer[1])){
200             Int_t excluded=0;
201             Int_t active=0;
202             Int_t exactive=0;//excluded but taking data
203             //AliITSCalibrationSDD *cal;
204             for(Int_t imod=0;imod<fgknSDDmodules;imod++){
205               Int_t lay=0;
206               Int_t lad=0;
207               Int_t det=0;
208               Int_t module=0;
209               module=imod+fgkmodoffset;
210               AliITSCalibrationSDD * cal=(AliITSCalibrationSDD*)fCalibration->At(imod);
211               if(cal==0) { delete cal; continue;}
212               AliITSgeomTGeo::GetModuleId(module,lay,lad,det);
213               if (cal->IsBad()){
214                 excluded++;
215                 Double_t content=0.;
216                 Double_t contentlayer[2];
217                 for(Int_t i=0;i<2;i++)contentlayer[i]=0.;
218                 if(hmodule)content=hmodule->GetBinContent(imod+1);//if expert bit is active the histogram has been created 
219                 contentlayer[lay-3]=hlayer[lay-3]->GetBinContent(det,lad);
220                 if(content!=0.||contentlayer[lay-3]!=0.)
221                   {
222                     filledmodules[lay-3]++;
223                     AliWarning(Form("The module %d (layer %i, ladder %i det %i ) excluded from the acquisition, took data \n ",module,lay,lad,det));
224                     exactive++;
225                   }
226                 else if(content==0.&&contentlayer[lay-3]==0.)emptymodules[lay-3]++;
227                 //AliInfo(Form("The module %d (layer %i, ladder %i det %i ) is bad, content %f content layer %f  filled modules position %d ",module,lay,lad,det,contentlayer[lay-3],content,lay-3) );
228               }//end if bad
229               else
230                 {
231                   Double_t contentgood=0.;
232                   active++;
233                   //printf("lay: %i\t det %i \t lad %i \n",lay,det,lad );
234                   contentgood=hlayer[lay-3]->GetBinContent(det,lad);
235                   if(contentgood==0.){emptymodules[lay-3]++;}
236                   else if(contentgood!=0.){filledmodules[lay-3]++;}
237                 }
238               
239               //delete cal;
240               //cal=NULL;
241             }//end for
242             for(Int_t i=0;i<2;i++){AliInfo(Form(" %s :Layer %i \tempty modules %i \t filled modules %i\n",hlayer[i]->GetName(), i+3,emptymodules[i], filledmodules[i]));}//end else layers
243             if(exactive==0){
244               AliInfo(Form("All the active modules (%i) are in acquisition. The number of excluded modules are %i \n",active,excluded));
245               test=fHighSDDValue[AliQAv1::kINFO];}
246             if(exactive!=0){
247               AliWarning(Form("%i modules excluded from the acquisition took data. Active modules%i \n ",exactive,active));
248               test=fHighSDDValue[AliQAv1::kWARNING];
249             }
250             if(excluded==exactive){
251               AliWarning(Form("All the modules exluded from the acquisition (%d) took data!  Active modules %i\n",excluded,active));
252               test=fHighSDDValue[AliQAv1::kWARNING];
253             }
254             if(active==0){
255               AliWarning(Form("No modules took data: excluded %i \t exactive %i \n", excluded, exactive)); 
256               test=fHighSDDValue[AliQAv1::kFATAL];
257             }
258           } 
259         }
260       }//end getentries !=0
261       //delete calSDD;
262       
263       
264       
265       //delete calibSDD;
266       //delete calSDD;
267    
268       break;
269     case AliQAv1::kNULLTASK:
270       AliInfo(Form("No Check on %s\n",AliQAv1::GetAliTaskName(index))); 
271       test=1.;
272       break;
273     case AliQAv1::kSIM:
274       AliInfo(Form("Check on %s\n",AliQAv1::GetAliTaskName(index))); 
275       test=1.;
276       break;
277     case AliQAv1::kREC:
278       AliInfo(Form("Check on %s\n",AliQAv1::GetAliTaskName(index))); 
279       test=1.;
280       break;
281     case AliQAv1::kANA:
282       AliInfo(Form("Check on %s\n",AliQAv1::GetAliTaskName(index)));
283       test=1.; 
284       break;
285     case AliQAv1::kESD:
286       AliInfo(Form("Check on %s\n",AliQAv1::GetAliTaskName(index)));
287       test=1.; 
288       break;
289     case AliQAv1::kNTASK:
290       AliInfo(Form("No Check on %s\n",AliQAv1::GetAliTaskName(index))); 
291       test=1.;
292       break;
293       
294     }//end switch
295   delete hdata;
296   return test;  
297 }
298  
299 //__________________________________________________________________
300 void AliITSQASDDChecker::SetTaskOffset(Int_t taskoffset)
301 {
302   fSubDetOffset = taskoffset;
303 }
304
305
306 //__________________________________________________________________
307 void AliITSQASDDChecker::SetStepBit(Double_t *steprange)
308 {
309
310   fStepBitSDD = new Double_t[AliQAv1::kNBIT];
311   for(Int_t bit=0;bit<AliQAv1::kNBIT;bit++)
312     {
313       fStepBitSDD[bit]=steprange[bit];
314     }
315 }
316
317 //__________________________________________________________________
318 void  AliITSQASDDChecker::SetSDDLimits(Float_t *lowvalue, Float_t * highvalue)
319 {
320
321   fLowSDDValue = new Float_t[AliQAv1::kNBIT];
322   fHighSDDValue= new Float_t[AliQAv1::kNBIT];
323
324   for(Int_t bit=0;bit<AliQAv1::kNBIT;bit++)
325     {
326       fLowSDDValue[bit]=lowvalue[bit];
327       fHighSDDValue[bit]= highvalue[bit];
328     }
329
330 }
331
332