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"AliQACheckerBase.h"
32 #include "AliITSQASPDChecker.h"
33 #include "AliITSQADataMakerRec.h"
36 ClassImp(AliITSQASPDChecker)
37 //__________________________________________________________________
38 AliITSQASPDChecker::AliITSQASPDChecker() :
47 //__________________________________________________________________
48 AliITSQASPDChecker& AliITSQASPDChecker::operator = (const AliITSQASPDChecker& qac )
51 this->~AliITSQASPDChecker();
52 new(this) AliITSQASPDChecker(qac);
55 //__________________________________________________________________
56 AliITSQASPDChecker::~AliITSQASPDChecker() {
57 if(fStepBitSPD) delete[] fStepBitSPD ;
58 if(fLowSPDValue)delete[]fLowSPDValue;
59 if(fHighSPDValue) delete[]fHighSPDValue;
60 if(fImage) delete[]fImage;
63 //__________________________________________________________________
64 Double_t AliITSQASPDChecker::Check(AliQAv1::ALITASK_t index, TObjArray * list, const AliDetectorRecoParam * /*recoParam*/)
67 // General methods for SPD Cheks to be used in RAWS and REC ALITASK_t
70 AliDebug(2, Form("AliITSQASPDChecker called with offset: %d\n", fSubDetOffset));
74 // Checks for ALITASK_t AliQAv1::kRAW
75 if(index == AliQAv1::kRAW) {
76 return CheckRawData(list);
78 if (list->GetEntries() == 0) {
79 test = 1.; // nothing to check
85 while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
87 TString histName = hdata->GetName();
88 if (!histName.Contains("_SPD")) continue;
90 if (hdata->GetEntries()>0) rv = 1;
91 if (histName.Contains("LayPattern")) {
92 if (hdata->GetBinContent(1)) {
93 Double_t ratio=hdata->GetBinContent(2)/hdata->GetBinContent(1);
94 AliDebug(2, Form("%s: ratio RecPoints lay2 / lay1 = %f", hdata->GetName(), ratio));
97 AliDebug(AliQAv1::GetQADebugLevel(), "No RecPoints in lay1");
99 else if(histName.Contains("ModPattern")) {
101 for(Int_t ibin=0;ibin<hdata->GetNbinsX();ibin++) {
102 if(histName.Contains("SPD1") && ibin<80 && hdata->GetBinContent(ibin+1)>0) ndead++;
103 if(histName.Contains("SPD2") && ibin>79 && hdata->GetBinContent(ibin+1)>0) ndead++;
105 AliDebug(2, Form("%s: Entries = %d number of empty modules = %d",
106 hdata->GetName(),(Int_t)hdata->GetEntries(),ndead));
108 else if(histName.Contains("SizeYvsZ")) {
109 Double_t meanz=hdata->GetMean(1);
110 Double_t meany=hdata->GetMean(2);
111 Double_t rmsz=hdata->GetRMS(1);
112 Double_t rmsy=hdata->GetRMS(2);
113 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s: Cluster sizeY mean = %f rms = %f", hdata->GetName(),meany,rmsy));
114 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s: Cluster sizeZ mean = %f rms = %f", hdata->GetName(),meanz,rmsz));
116 else if(histName.Contains("SPDMultiplicity")) {
117 AliDebug(2, Form("%s: Events = %d mean = %f rms = %f",
118 hdata->GetName(),(Int_t)hdata->GetEntries(),hdata->GetMean(),hdata->GetRMS()));}
120 // else AliDebug(AliQAv1::GetQADebugLevel(), Form("%s -> %f", hdata->GetName(), rv));
125 AliError("Data type cannot be processed") ;
130 if (AliITSQADataMakerRec::AreEqual(test,0)) {
131 AliWarning("Histograms are there, but they are all empty: setting flag to kWARNING");
132 test = fHighSPDValue[AliQAv1::kWARNING]; //upper limit value to set kWARNING flag for a task
140 AliDebug(AliQAv1::GetQADebugLevel(), Form("Test Result = %f", test));
144 //__________________________________________________________________
145 Double_t AliITSQASPDChecker::CheckRawData(const TObjArray * list) {
147 // Checks on the raw data histograms [ preliminary version ]
148 // The output of this method is the fraction of SPD histograms which are processed by the checker.
149 // The methods returns fHighSPDValue[AliQAv1::kFATAL] in case of data format errors or MEB errors
155 // basic checks on input data
157 AliError("NO histogram list for RAWS");
161 if(list->GetEntries() == 0) {
162 AliWarning("No histograms in RAW list \n");
166 // loop over the raw data histograms
169 Double_t totalHistos = 0;
170 Double_t goodHistos = 0; // number of histograms which passed the checks
171 Double_t response =0;
172 Bool_t fatalProblem = kFALSE;
174 while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
176 TString histName = hdata->GetName();
177 if(!histName.Contains("SPD")) continue;
180 if(histName.Contains("SPDErrorsAll")){
181 if(hdata->GetListOfFunctions()->GetEntries()<1) hdata->GetListOfFunctions()->Add(new TPaveText(0.2,0.23,0.7,0.5,"NDC"));
183 for(Int_t i=0; i<hdata->GetListOfFunctions()->GetEntries(); i++){
184 TString funcName = hdata->GetListOfFunctions()->At(i)->ClassName();
185 if(funcName.Contains("TPaveText")){
186 TPaveText *p = (TPaveText*)hdata->GetListOfFunctions()->At(i);
189 if(hdata->Integral(0,hdata->GetNbinsX())>0){
190 Bool_t isHighMult = kFALSE;
191 for(Int_t ieq=0; ieq<20; ieq++){
192 if(hdata->GetBinContent(ieq+1,17+1)>0 && hdata->GetBinContent(ieq+1,20+1)>0) isHighMult = kTRUE;
195 p->SetFillColor(kOrange);
196 p->AddText("High occupancy in a chip detected (-> errors type 17,20 and 0 are present). ");
197 p->AddText("ONLY IF OTHER error types are present CALL the expert");
199 p->SetFillColor(kRed);
200 p->AddText("Data Format NOT OK. Please call the expert!");
202 response = fHighSPDValue[AliQAv1::kFATAL];
208 p->SetFillColor(kGreen);
213 } // data format error
216 else if(histName.Contains("MEB")){
217 if(hdata->GetListOfFunctions()->GetEntries()<1) hdata->GetListOfFunctions()->Add(new TPaveText(0.2,0.23,0.7,0.5,"NDC"));
219 for(Int_t i=0; i<hdata->GetListOfFunctions()->GetEntries(); i++){
220 TString funcName = hdata->GetListOfFunctions()->At(i)->ClassName();
221 if(funcName.Contains("TPaveText")){
222 TPaveText *p = (TPaveText*)hdata->GetListOfFunctions()->At(i);
225 if(hdata->GetEntries()>0){
226 p->SetFillColor(kRed);
227 p->AddText("MEB problem could be present. Please check if SPD is in READY state.");
228 p->AddText("If SPD is in -READY- state, please notify it to the expert.");
229 response = fHighSPDValue[AliQAv1::kFATAL];
234 p->SetFillColor(kGreen);
244 if(!fatalProblem) response = goodHistos/totalHistos;
245 // printf("n histos %f - good ones %f ----> ratio %f , fatal response %i\n",totalHistos,goodHistos,goodHistos/totalHistos,(Int_t)fatalProblem);
249 //__________________________________________________________________
250 void AliITSQASPDChecker::SetTaskOffset(Int_t TaskOffset)
252 // Offset for SPD within ITS QA
253 fSubDetOffset = TaskOffset;
256 //__________________________________________________________________
257 void AliITSQASPDChecker::SetStepBit(const Double_t *steprange)
259 // Step bit for SPD within ITS QA
260 fStepBitSPD = new Double_t[AliQAv1::kNBIT];
261 for(Int_t bit=0;bit<AliQAv1::kNBIT;bit++)
263 fStepBitSPD[bit]=steprange[bit];
268 //__________________________________________________________________
269 void AliITSQASPDChecker::SetSPDLimits(const Float_t *lowvalue, const Float_t * highvalue)
271 // SPD limints for QA bit within general ITS QA
272 fLowSPDValue = new Float_t[AliQAv1::kNBIT];
273 fHighSPDValue= new Float_t[AliQAv1::kNBIT];
275 for(Int_t bit=0;bit<AliQAv1::kNBIT;bit++)
277 fLowSPDValue[bit]=lowvalue[bit];
278 fHighSPDValue[bit]= highvalue[bit];
282 //__________________________________________________________________
283 Bool_t AliITSQASPDChecker::MakeSPDImage( TObjArray ** list, AliQAv1::TASKINDEX_t task, AliQAv1::MODE_t mode)
287 fImage=(TCanvas**)AliQAChecker::Instance()->GetDetQAChecker(0)->GetImage();
288 //create the image for raws and recpoints. In the other case, the default methodof CheckerBase class will be used
291 case AliQAv1::kRAWS:{
292 val = MakeSPDRawsImage(list, task,mode);
295 case AliQAv1::kRECPOINTS:;
296 case AliQAv1::kHITS:;
297 case AliQAv1::kESDS:;
298 case AliQAv1::kDIGITS:;
299 case AliQAv1::kDIGITSR:;
300 case AliQAv1::kSDIGITS:;
301 case AliQAv1::kTRACKSEGMENTS:;
302 case AliQAv1::kRECPARTICLES:;
305 //AliQAChecker::Instance()->GetDetQAChecker(0)->MakeImage(list,task,mode);
309 case AliQAv1::kNULLTASKINDEX:; case AliQAv1::kNTASKINDEX:
310 {AliWarning(Form("No histograms for these tasks ( %s ) \n", AliQAv1::GetTaskName(task).Data())); val = kFALSE;}
315 //_______________________________________________________________________
316 Bool_t AliITSQASPDChecker::MakeSPDRawsImage(TObjArray ** list, AliQAv1::TASKINDEX_t task, AliQAv1::MODE_t mode )
319 // create layout of the histograms used in the DQM
323 for (Int_t esIndex = 0 ; esIndex < AliRecoParam::kNSpecies ; esIndex++) {
324 //printf("-------------------------> %i \n", esIndex);
325 if (! AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))->IsEventSpecieSet(AliRecoParam::ConvertIndex(esIndex)) || list[esIndex]->GetEntries() == 0)
326 {printf ("Nothing for %s \n", AliRecoParam::GetEventSpecieName(esIndex)); continue;}
328 const Char_t * title = Form("QA_%s_%s_%s", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(esIndex)) ;
329 if ( !fImage[esIndex] ) {
330 fImage[esIndex] = new TCanvas(title, title,1280,980) ;
333 fImage[esIndex]->Clear() ;
334 fImage[esIndex]->SetTitle(title) ;
335 fImage[esIndex]->cd();
337 TPaveText someText(0.015, 0.015, 0.98, 0.98);
338 someText.AddText(title);
340 fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat()), "ps") ;
341 fImage[esIndex]->Clear() ;
345 fImage[esIndex]->Divide(nx, ny) ;
349 fImage[esIndex]->cd(npad);
350 fImage[esIndex]->cd(npad)->SetBorderMode(0) ;
352 TIter next(list[esIndex]);
354 while ( (hist=static_cast<TH1*>(next())) ) {
355 //gPad=fImage[esIndex]->cd(npad)->GetPad(npad);
356 if(!hist->TestBit(AliQAv1::GetImageBit())) continue;
357 TString name(hist->GetName());
358 if(name.Contains("SPDErrorsAll")) {
359 fImage[esIndex]->cd(1) ;
360 gPad->SetBorderMode(0) ;
361 gPad->SetRightMargin(0.25);
363 hist->SetOption("text") ;
366 if(name.Contains("MEB")) {
367 fImage[esIndex]->cd(2) ;
368 gPad->SetBorderMode(0) ;
369 gPad->SetBottomMargin(0.25);
370 hist->SetOption("text") ;
373 if(name.Contains("SPDFastOrCorrelation")){
374 fImage[esIndex]->cd(3) ;
375 gPad->SetBorderMode(0) ;
376 hist->SetOption("") ;
380 if(name.Contains("SPDHitMapStaveChipInner")){
381 fImage[esIndex]->cd(4) ;
382 gPad->SetBorderMode(0) ;
383 gPad->SetRightMargin(0.25);
384 hist->SetOption("colz") ;
387 if(name.Contains("SPDHitMapStaveChipOuter")){
388 fImage[esIndex]->cd(5) ;
389 gPad->SetBorderMode(0) ;
390 gPad->SetRightMargin(0.25);
391 hist->SetOption("colz") ;
394 if(name.Contains("SPDFastOrMapStaveChip")){
395 fImage[esIndex]->cd(6) ;
396 gPad->SetBorderMode(0) ;
397 gPad->SetRightMargin(0.25);
398 gPad->SetBottomMargin(0.25);
399 hist->SetOption("colz") ;
406 fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat()), "ps") ;