1 /**************************************************************************
2 * Copyright(c) 2007-2009, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 // *****************************************
19 // Checks the quality assurance
20 // by comparing with reference data
21 // P. Cerello Apr 2008
24 // --- ROOT system ---
30 // --- AliRoot header files ---
31 #include "AliITSQASPDChecker.h"
32 #include "AliITSQADataMakerRec.h"
35 ClassImp(AliITSQASPDChecker)
36 //__________________________________________________________________
37 AliITSQASPDChecker::AliITSQASPDChecker() :
46 //__________________________________________________________________
47 AliITSQASPDChecker& AliITSQASPDChecker::operator = (const AliITSQASPDChecker& qac )
50 this->~AliITSQASPDChecker();
51 new(this) AliITSQASPDChecker(qac);
54 //__________________________________________________________________
55 AliITSQASPDChecker::~AliITSQASPDChecker() {
56 if(fStepBitSPD) delete[] fStepBitSPD ;
57 if(fLowSPDValue)delete[]fLowSPDValue;
58 if(fHighSPDValue) delete[]fHighSPDValue;
59 if(fImage) delete[]fImage;
62 //__________________________________________________________________
63 Double_t AliITSQASPDChecker::Check(AliQAv1::ALITASK_t index, TObjArray * list, const AliDetectorRecoParam * /*recoParam*/)
66 // General methods for SPD Cheks to be used in RAWS and REC ALITASK_t
69 AliDebug(2, Form("AliITSQASPDChecker called with offset: %d\n", fSubDetOffset));
73 // Checks for ALITASK_t AliQAv1::kRAW
74 if(index == AliQAv1::kRAW) {
75 return CheckRawData(list);
77 if (list->GetEntries() == 0) {
78 test = 1.; // nothing to check
84 while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
86 TString histName = hdata->GetName();
87 if (!histName.Contains("_SPD")) continue;
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));
96 AliDebug(AliQAv1::GetQADebugLevel(), "No RecPoints in lay1");
98 else if(histName.Contains("ModPattern")) {
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++;
104 AliDebug(2, Form("%s: Entries = %d number of empty modules = %d",
105 hdata->GetName(),(Int_t)hdata->GetEntries(),ndead));
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));
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()));}
119 // else AliDebug(AliQAv1::GetQADebugLevel(), Form("%s -> %f", hdata->GetName(), rv));
124 AliError("Data type cannot be processed") ;
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
139 AliDebug(AliQAv1::GetQADebugLevel(), Form("Test Result = %f", test));
143 //__________________________________________________________________
144 Double_t AliITSQASPDChecker::CheckRawData(const TObjArray * list) {
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
154 // basic checks on input data
156 AliError("NO histogram list for RAWS");
160 if(list->GetEntries() == 0) {
161 AliWarning("No histograms in RAW list \n");
165 // loop over the raw data histograms
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;
173 while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
175 TString histName = hdata->GetName();
176 if(!histName.Contains("SPD")) continue;
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"));
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);
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;
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");
198 p->SetFillColor(kRed);
199 p->AddText("Data Format NOT OK. Please call the expert!");
201 response = fHighSPDValue[AliQAv1::kFATAL];
207 p->SetFillColor(kGreen);
212 } // data format 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"));
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);
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];
233 p->SetFillColor(kGreen);
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);
248 //__________________________________________________________________
249 void AliITSQASPDChecker::SetTaskOffset(Int_t TaskOffset)
251 // Offset for SPD within ITS QA
252 fSubDetOffset = TaskOffset;
255 //__________________________________________________________________
256 void AliITSQASPDChecker::SetStepBit(const Double_t *steprange)
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++)
262 fStepBitSPD[bit]=steprange[bit];
267 //__________________________________________________________________
268 void AliITSQASPDChecker::SetSPDLimits(const Float_t *lowvalue, const Float_t * highvalue)
270 // SPD limints for QA bit within general ITS QA
271 fLowSPDValue = new Float_t[AliQAv1::kNBIT];
272 fHighSPDValue= new Float_t[AliQAv1::kNBIT];
274 for(Int_t bit=0;bit<AliQAv1::kNBIT;bit++)
276 fLowSPDValue[bit]=lowvalue[bit];
277 fHighSPDValue[bit]= highvalue[bit];
281 //__________________________________________________________________
282 Bool_t AliITSQASPDChecker::MakeSPDImage( TObjArray ** list, AliQAv1::TASKINDEX_t task, AliQAv1::MODE_t mode)
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
290 case AliQAv1::kRAWS:{
291 val = MakeSPDRawsImage(list, task,mode);
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:;
304 //AliQAChecker::Instance()->GetDetQAChecker(0)->MakeImage(list,task,mode);
308 case AliQAv1::kNULLTASKINDEX:; case AliQAv1::kNTASKINDEX:
309 {AliWarning(Form("No histograms for these tasks ( %s ) \n", AliQAv1::GetTaskName(task).Data())); val = kFALSE;}
314 //_______________________________________________________________________
315 Bool_t AliITSQASPDChecker::MakeSPDRawsImage(TObjArray ** list, AliQAv1::TASKINDEX_t task, AliQAv1::MODE_t mode )
318 // create layout of the histograms used in the DQM
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;}
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) ;
332 fImage[esIndex]->Clear() ;
333 fImage[esIndex]->SetTitle(title) ;
334 fImage[esIndex]->cd();
336 TPaveText someText(0.015, 0.015, 0.98, 0.98);
337 someText.AddText(title);
339 fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat()), "ps") ;
340 fImage[esIndex]->Clear() ;
344 fImage[esIndex]->Divide(nx, ny) ;
348 fImage[esIndex]->cd(npad);
349 fImage[esIndex]->cd(npad)->SetBorderMode(0) ;
351 TIter next(list[esIndex]);
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);
362 hist->SetOption("colz") ;
365 if(name.Contains("MEB")) {
366 fImage[esIndex]->cd(2) ;
367 gPad->SetBorderMode(0) ;
368 gPad->SetBottomMargin(0.25);
369 hist->SetOption("colz") ;
372 if(name.Contains("SPDFastOrCorrelation")){
373 fImage[esIndex]->cd(3) ;
374 gPad->SetBorderMode(0) ;
375 hist->SetOption("") ;
379 if(name.Contains("SPDHitMapStaveChipInner")){
380 fImage[esIndex]->cd(4) ;
381 gPad->SetBorderMode(0) ;
382 gPad->SetRightMargin(0.25);
383 hist->SetOption("colz") ;
386 if(name.Contains("SPDHitMapStaveChipOuter")){
387 fImage[esIndex]->cd(5) ;
388 gPad->SetBorderMode(0) ;
389 gPad->SetRightMargin(0.25);
390 hist->SetOption("colz") ;
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") ;
405 fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat()), "ps") ;