]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQASPDChecker.cxx
Coverity corrections.
[u/mrichter/AliRoot.git] / ITS / AliITSQASPDChecker.cxx
1 /**************************************************************************
2  * Copyright(c) 2007-2009, 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 /* $Id $ */
17
18 // *****************************************
19 //  Checks the quality assurance 
20 //  by comparing with reference data
21 //  P. Cerello Apr 2008
22 //  INFN Torino
23
24 // --- ROOT system ---
25 #include "TH1.h"
26 #include "TString.h"
27 #include "TList.h"
28 #include "TCanvas.h"
29
30 // --- AliRoot header files ---
31 #include "AliITSQASPDChecker.h"
32 #include "AliITSQADataMakerRec.h"
33 #include "AliLog.h"
34
35 ClassImp(AliITSQASPDChecker)
36 //__________________________________________________________________
37 AliITSQASPDChecker::AliITSQASPDChecker() : 
38  TObject(),
39  fSubDetOffset(0), 
40  fStepBitSPD(NULL),
41  fLowSPDValue(NULL),
42  fHighSPDValue(NULL),
43  fImage(NULL) 
44  {
45  }
46 //__________________________________________________________________
47 AliITSQASPDChecker& AliITSQASPDChecker::operator = (const AliITSQASPDChecker& qac ) 
48 {
49   // Equal operator.
50   this->~AliITSQASPDChecker();
51   new(this) AliITSQASPDChecker(qac);
52   return *this;
53 }
54 //__________________________________________________________________
55 AliITSQASPDChecker::~AliITSQASPDChecker() {
56 if(fStepBitSPD) delete[] fStepBitSPD ;
57 if(fLowSPDValue)delete[]fLowSPDValue;
58 if(fHighSPDValue) delete[]fHighSPDValue;
59 if(fImage) delete[]fImage;
60
61
62 //__________________________________________________________________
63 Double_t AliITSQASPDChecker::Check(AliQAv1::ALITASK_t index, TObjArray * list, const AliDetectorRecoParam * /*recoParam*/)
64 {
65 //
66 // General methods for SPD Cheks to be used in RAWS and REC ALITASK_t
67 //
68
69   AliDebug(2, Form("AliITSQASPDChecker called with offset: %d\n", fSubDetOffset));
70
71   Double_t test = 0.0;
72   Int_t count = 0;
73   // Checks for ALITASK_t AliQAv1::kRAW
74   if(index == AliQAv1::kRAW) {
75   return CheckRawData(list);
76   } else {
77   if (list->GetEntries() == 0) {
78     test = 1.; // nothing to check
79   }
80   else {
81     TIter next(list);
82     TH1 * hdata;
83     count = 0;
84     while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
85       if (hdata) {
86         TString histName = hdata->GetName();
87         if (!histName.Contains("_SPD")) continue;
88         Double_t rv = 0.;
89         if (hdata->GetEntries()>0) rv = 1;
90         if (histName.Contains("LayPattern")) {
91          if (hdata->GetBinContent(1)) {
92            Double_t ratio=hdata->GetBinContent(2)/hdata->GetBinContent(1);
93            AliDebug(2, Form("%s: ratio RecPoints lay2 / lay1 = %f", hdata->GetName(), ratio));
94          }
95          else
96            AliDebug(AliQAv1::GetQADebugLevel(), "No RecPoints in lay1");
97        }
98         else if(histName.Contains("ModPattern")) {
99            Int_t ndead=0;
100            for(Int_t ibin=0;ibin<hdata->GetNbinsX();ibin++) {
101              if(histName.Contains("SPD1") && ibin<80 && hdata->GetBinContent(ibin+1)>0) ndead++;
102              if(histName.Contains("SPD2") && ibin>79 && hdata->GetBinContent(ibin+1)>0) ndead++;
103            }
104            AliDebug(2, Form("%s: Entries = %d  number of empty modules = %d", 
105                         hdata->GetName(),(Int_t)hdata->GetEntries(),ndead));
106         }
107         else if(histName.Contains("SizeYvsZ")) {
108            Double_t meanz=hdata->GetMean(1);
109            Double_t meany=hdata->GetMean(2);
110            Double_t rmsz=hdata->GetRMS(1);
111            Double_t rmsy=hdata->GetRMS(2);
112            AliDebug(AliQAv1::GetQADebugLevel(), Form("%s: Cluster sizeY mean = %f  rms = %f", hdata->GetName(),meany,rmsy));
113            AliDebug(AliQAv1::GetQADebugLevel(), Form("%s: Cluster sizeZ mean = %f  rms = %f", hdata->GetName(),meanz,rmsz));
114         }
115         else if(histName.Contains("SPDMultiplicity")) {
116            AliDebug(2, Form("%s: Events = %d  mean = %f  rms = %f",
117                         hdata->GetName(),(Int_t)hdata->GetEntries(),hdata->GetMean(),hdata->GetRMS()));}
118
119        // else AliDebug(AliQAv1::GetQADebugLevel(), Form("%s -> %f", hdata->GetName(), rv));
120         count++;
121         test += rv;
122       }
123       else {
124         AliError("Data type cannot be processed") ;
125       }
126     }
127
128     if (count != 0) {
129       if (AliITSQADataMakerRec::AreEqual(test,0)) {
130         AliWarning("Histograms are there, but they are all empty: setting flag to kWARNING");
131         test = fHighSPDValue[AliQAv1::kWARNING];  //upper limit value to set kWARNING flag for a task
132       }
133       else {
134         test /= count;
135       }
136     }
137   }
138 }
139   AliDebug(AliQAv1::GetQADebugLevel(), Form("Test Result = %f", test));
140   return test ;
141
142 }
143 //__________________________________________________________________
144 Double_t AliITSQASPDChecker::CheckRawData(const TObjArray * list) {
145 //
146 // Checks on the raw data histograms [ preliminary version ]
147 // The output of this method is the fraction of SPD histograms which are processed by the checker. 
148 // The methods returns fHighSPDValue[AliQAv1::kFATAL] in case of data format errors or MEB errors
149 // 
150 // A. Mastroserio
151
152 Double_t test =0;
153
154 // basic checks on input data
155 if(!list) {
156  AliError("NO histogram list for RAWS");
157  return test;
158  }
159
160  if(list->GetEntries() == 0) {
161  AliWarning("No histograms in RAW list \n");
162  return test;
163  }
164
165 // loop over the raw data histograms
166 TIter next(list);
167 TH1 * hdata;
168 Double_t totalHistos = 0;
169 Double_t goodHistos = 0; // number of histograms which passed the checks
170 Double_t response =0;
171 Bool_t fatalProblem = kFALSE;
172
173 while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
174  if (hdata) {
175         TString histName = hdata->GetName();
176         if(!histName.Contains("SPD")) continue;
177         totalHistos++;
178         // data format error
179        if(histName.Contains("SPDErrorsAll")){
180         if(hdata->GetListOfFunctions()->GetEntries()<1) hdata->GetListOfFunctions()->Add(new TPaveText(0.2,0.23,0.7,0.5,"NDC"));
181
182            for(Int_t i=0; i<hdata->GetListOfFunctions()->GetEntries(); i++){
183            TString funcName = hdata->GetListOfFunctions()->At(i)->ClassName();
184            if(funcName.Contains("TPaveText")){
185             TPaveText *p = (TPaveText*)hdata->GetListOfFunctions()->At(i);
186              p->Clear();
187
188             if(hdata->Integral(0,hdata->GetNbinsX())>0){
189              Bool_t isHighMult = kFALSE;
190              for(Int_t ieq=0; ieq<20; ieq++){
191               if(hdata->GetBinContent(ieq+1,17+1)>0 && hdata->GetBinContent(ieq+1,20+1)>0) isHighMult = kTRUE;
192              }
193              if(isHighMult) {
194               p->SetFillColor(kOrange);
195               p->AddText("High occupancy in a chip detected (-> errors type 17,20 and 0 are present). ");
196               p->AddText("ONLY IF OTHER error types are present CALL the expert");
197               } else {
198              p->SetFillColor(kRed);
199              p->AddText("Data Format NOT OK. Please call the expert!");
200              }
201             response = fHighSPDValue[AliQAv1::kFATAL];
202             fatalProblem=kTRUE;
203             continue;
204             } // if errors 
205            else {
206              p->Clear();
207              p->SetFillColor(kGreen);
208              p->AddText("OK");
209             }
210            } // TPaveText
211          } // list entries   
212         } // data format error
213
214         // MEB error
215         else if(histName.Contains("MEB")){
216           if(hdata->GetListOfFunctions()->GetEntries()<1) hdata->GetListOfFunctions()->Add(new TPaveText(0.2,0.23,0.7,0.5,"NDC"));
217
218            for(Int_t i=0; i<hdata->GetListOfFunctions()->GetEntries(); i++){
219            TString funcName = hdata->GetListOfFunctions()->At(i)->ClassName();
220            if(funcName.Contains("TPaveText")){
221              TPaveText *p = (TPaveText*)hdata->GetListOfFunctions()->At(i);
222              p->Clear();
223
224           if(hdata->GetEntries()>0){
225              p->SetFillColor(kRed);
226              p->AddText("MEB problem could be present. Please check if SPD is in READY state.");
227              p->AddText("If SPD is in -READY- state, please notify it to the expert."); 
228             response = fHighSPDValue[AliQAv1::kFATAL];
229             fatalProblem=kTRUE;
230             continue;
231
232            } else {
233              p->SetFillColor(kGreen);
234              p->AddText("OK");
235               }   
236
237             } // pave text
238            } // list 
239          }
240         goodHistos++;
241         }
242      }
243     if(!fatalProblem) response = goodHistos/totalHistos;
244    // printf("n histos %f - good ones %f ----> ratio %f , fatal response %i\n",totalHistos,goodHistos,goodHistos/totalHistos,(Int_t)fatalProblem);
245     return response;
246 }
247
248 //__________________________________________________________________
249 void AliITSQASPDChecker::SetTaskOffset(Int_t TaskOffset)
250 {
251 // Offset for SPD within ITS QA
252   fSubDetOffset = TaskOffset;
253 }
254
255 //__________________________________________________________________
256 void AliITSQASPDChecker::SetStepBit(const Double_t *steprange) 
257 {
258 // Step bit for SPD within ITS QA
259   fStepBitSPD = new Double_t[AliQAv1::kNBIT];
260   for(Int_t bit=0;bit<AliQAv1::kNBIT;bit++)
261     {
262       fStepBitSPD[bit]=steprange[bit];
263     }
264 }
265
266
267 //__________________________________________________________________
268 void  AliITSQASPDChecker::SetSPDLimits(const Float_t *lowvalue, const Float_t * highvalue)
269 {
270 // SPD limints for QA bit within general ITS QA
271   fLowSPDValue = new Float_t[AliQAv1::kNBIT];
272   fHighSPDValue= new Float_t[AliQAv1::kNBIT];
273
274   for(Int_t bit=0;bit<AliQAv1::kNBIT;bit++)
275     {
276       fLowSPDValue[bit]=lowvalue[bit];
277       fHighSPDValue[bit]= highvalue[bit];
278     }
279
280 }
281 //__________________________________________________________________
282 Bool_t  AliITSQASPDChecker::MakeSPDImage( TObjArray ** list, AliQAv1::TASKINDEX_t task, AliQAv1::MODE_t mode)
283 {
284   Bool_t val=kFALSE;
285
286   fImage=(TCanvas**)AliQAChecker::Instance()->GetDetQAChecker(0)->GetImage();
287   //create the image for raws and recpoints. In the other case, the default methodof CheckerBase class will be used
288   switch(task)
289     {
290     case AliQAv1::kRAWS:{
291       val = MakeSPDRawsImage(list, task,mode);
292     }
293       break;
294     case AliQAv1::kRECPOINTS:;
295     case AliQAv1::kHITS:; 
296     case AliQAv1::kESDS:; 
297     case AliQAv1::kDIGITS:;
298     case AliQAv1::kDIGITSR:;
299     case AliQAv1::kSDIGITS:;
300     case AliQAv1::kTRACKSEGMENTS:;
301     case AliQAv1::kRECPARTICLES:; 
302     default:
303     {
304        //AliQAChecker::Instance()->GetDetQAChecker(0)->MakeImage(list,task,mode);
305       val = kFALSE;
306     }
307     break;
308     case AliQAv1::kNULLTASKINDEX:; case  AliQAv1::kNTASKINDEX: 
309       {AliWarning(Form("No histograms for these tasks ( %s ) \n", AliQAv1::GetTaskName(task).Data())); val = kFALSE;}
310       break;
311     }
312  return val; 
313 }
314 //_______________________________________________________________________
315 Bool_t AliITSQASPDChecker::MakeSPDRawsImage(TObjArray ** list, AliQAv1::TASKINDEX_t task, AliQAv1::MODE_t mode )
316 {
317     //
318     // create layout of the histograms used in the DQM
319     //
320   
321   
322     for (Int_t esIndex = 0 ; esIndex < AliRecoParam::kNSpecies ; esIndex++) {
323       //printf("-------------------------> %i \n", esIndex);
324       if (! AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))->IsEventSpecieSet(AliRecoParam::ConvertIndex(esIndex)) || list[esIndex]->GetEntries() == 0) 
325           {printf ("Nothing for %s \n", AliRecoParam::GetEventSpecieName(esIndex)); continue;}
326       else{
327         const Char_t * title = Form("QA_%s_%s_%s", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(esIndex)) ; 
328         if ( !fImage[esIndex] ) {
329           fImage[esIndex] = new TCanvas(title, title,1280,980) ;
330         }
331         
332         fImage[esIndex]->Clear() ; 
333         fImage[esIndex]->SetTitle(title) ; 
334         fImage[esIndex]->cd();
335  
336         TPaveText someText(0.015, 0.015, 0.98, 0.98);
337         someText.AddText(title);
338         someText.Draw(); 
339         fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat()), "ps") ; 
340         fImage[esIndex]->Clear() ; 
341         Int_t nx =3; 
342         Int_t ny =2; 
343         
344         fImage[esIndex]->Divide(nx, ny) ; 
345          
346         TH1* hist = NULL ;
347         Int_t npad = 1 ; 
348         fImage[esIndex]->cd(npad); 
349         fImage[esIndex]->cd(npad)->SetBorderMode(0) ;
350         
351         TIter next(list[esIndex]);
352
353         while ( (hist=static_cast<TH1*>(next())) ) {
354           //gPad=fImage[esIndex]->cd(npad)->GetPad(npad);
355           if(!hist->TestBit(AliQAv1::GetImageBit())) continue;
356           TString name(hist->GetName());
357            if(name.Contains("SPDErrorsAll")) {
358             fImage[esIndex]->cd(1) ; 
359             gPad->SetBorderMode(0) ;  
360             gPad->SetRightMargin(0.25);
361             hist->SetStats(0);
362             hist->SetOption("colz") ;
363             hist->DrawCopy();  
364            }     
365            if(name.Contains("MEB")) {
366             fImage[esIndex]->cd(2) ; 
367             gPad->SetBorderMode(0) ;  
368             gPad->SetBottomMargin(0.25);
369             hist->SetOption("colz") ;
370             hist->DrawCopy();  
371            }     
372            if(name.Contains("SPDFastOrCorrelation")){
373             fImage[esIndex]->cd(3) ; 
374             gPad->SetBorderMode(0) ;  
375             hist->SetOption("") ;
376             hist->DrawCopy();               
377            }
378            
379             if(name.Contains("SPDHitMapStaveChipInner")){
380             fImage[esIndex]->cd(4) ; 
381             gPad->SetBorderMode(0) ;  
382             gPad->SetRightMargin(0.25);
383             hist->SetOption("colz") ;
384             hist->DrawCopy();               
385            }
386             if(name.Contains("SPDHitMapStaveChipOuter")){
387             fImage[esIndex]->cd(5) ; 
388             gPad->SetBorderMode(0) ;  
389             gPad->SetRightMargin(0.25);
390             hist->SetOption("colz") ;
391             hist->DrawCopy();               
392            }
393             if(name.Contains("SPDFastOrMapStaveChip")){
394             fImage[esIndex]->cd(6) ; 
395             gPad->SetBorderMode(0) ;  
396             gPad->SetRightMargin(0.25);
397             gPad->SetBottomMargin(0.25);
398             hist->SetOption("colz") ;
399             hist->DrawCopy();               
400            }
401           
402           
403         }
404         
405         fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat()), "ps") ; 
406       }
407     }
408   return kTRUE;
409 }