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