]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQAChecker.cxx
Added classes for geometry - reconstruction interface, minor fixes in other classes
[u/mrichter/AliRoot.git] / ITS / AliITSQAChecker.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 //  W.Ferrarese  P.Cerello  Mag 2008
22 //  INFN Torino
23
24 // --- ROOT system ---
25 #include "TH1.h"
26 #include <Riostream.h>
27 #include "TStyle.h"
28
29 // --- AliRoot header files ---
30 #include "AliITSQAChecker.h"
31 #include "AliITSQASPDChecker.h"
32 #include "AliITSQASDDChecker.h"
33 #include "AliITSQASSDChecker.h"
34 #include "AliITSQADataMakerRec.h"
35
36 ClassImp(AliITSQAChecker)
37
38 //____________________________________________________________________________
39 AliITSQAChecker::AliITSQAChecker(Bool_t kMode, Short_t subDet, Short_t ldc) :
40 AliQACheckerBase("ITS","SDD Quality Assurance Checker"),
41 fkOnline(0),
42 fDet(0),  
43 fLDC(0),
44 fSPDOffset(0), 
45 fSDDOffset(0), 
46 fSSDOffset(0),
47 fSPDHisto(0),
48 fSDDHisto(0),
49 fSSDHisto(0),
50 fSPDChecker(0),  // SPD Checker
51 fSDDChecker(0),  // SDD Checker
52 fSSDChecker(0)  // SSD Checker
53
54 {
55   // Standard constructor
56   fkOnline = kMode; fDet = subDet; fLDC = ldc;
57   if(fDet == 0 || fDet == 1) {
58     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQAChecker::Create SPD Checker\n");
59     fSPDChecker = new AliITSQASPDChecker();
60   }
61   if(fDet == 0 || fDet == 2) {
62     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQAChecker::Create SDD Checker\n");
63     fSDDChecker = new AliITSQASDDChecker();
64   }
65   if(fDet == 0 || fDet == 3) {
66     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQAChecker::Create SSD Checker\n");
67     fSSDChecker = new AliITSQASSDChecker();
68   }
69   InitQACheckerLimits();
70 }
71
72 //____________________________________________________________________________
73 AliITSQAChecker::~AliITSQAChecker(){
74   // destructor
75   if(fSPDChecker)delete fSPDChecker;
76   if(fSDDChecker)delete fSDDChecker;
77   if(fSSDChecker)delete fSSDChecker;
78
79 }
80 //____________________________________________________________________________
81 void AliITSQAChecker::Check(Double_t * rv, AliQAv1::ALITASK_t index, TObjArray ** list, const AliDetectorRecoParam * recoParam)
82 {
83
84
85   // basic checks on the QA histograms on the input list
86   //for the ITS subdetectorQA (Raws Digits Hits RecPoints SDigits) return the worst value of the three result
87   if(index == AliQAv1::kESD){
88     
89     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
90       rv[specie] = 0.0 ; 
91       if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
92         continue ; 
93       AliDebug(AliQAv1::GetQADebugLevel(),"Checker for ESD");
94       Int_t tested = 0;
95       Int_t empty = 0;
96       // The following flags are set to kTRUE if the corresponding
97       // QA histograms exceed a given quality threshold
98       Bool_t cluMapSA = kFALSE;
99       Bool_t cluMapMI = kFALSE;
100       Bool_t cluMI = kFALSE;
101       Bool_t cluSA = kFALSE;
102       Bool_t verSPDZ = kFALSE;
103       if (list[specie]->GetEntries() == 0) {
104         rv[specie] = 0.; // nothing to check
105       }
106       else {
107         Double_t stepbit[AliQAv1::kNBIT];
108         Double_t histonumb= list[specie]->GetEntries();
109         CreateStepForBit(histonumb,stepbit); 
110         TIter next1(list[specie]);
111         TH1 * hdata;
112         Int_t nskipped=0;
113         Bool_t skipped[6]={kFALSE,kFALSE,kFALSE,kFALSE,kFALSE,kFALSE};
114         // look for layers that we wanted to skip
115         while ( (hdata = dynamic_cast<TH1 *>(next1())) ) {
116           if(hdata){
117             TString hname = hdata->GetName();
118             if(!hname.Contains("hESDSkippedLayers")) continue;
119             for(Int_t k=1; k<7; k++) {
120               if(hdata->GetBinContent(k)>0) { 
121                 nskipped++; 
122                 skipped[k-1]=kTRUE; 
123               } 
124             } 
125           }
126         }
127         TIter next(list[specie]);
128         while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
129           if(hdata){
130             TString hname = hdata->GetName();
131             Double_t entries = hdata->GetEntries();
132             ++tested;
133             if(!(entries>0.))++empty;
134             AliDebug(AliQAv1::GetQADebugLevel(),Form("ESD hist name %s - entries %12.1g",hname.Data(),entries));
135             if(hname.Contains("hESDClusterMapSA") && entries>0.){
136               cluMapSA = kTRUE;
137               AliDebug(AliQAv1::GetQADebugLevel(),Form("Processing histogram %s",hname.Data()));
138               // Check if there are layers with anomalously low 
139               // contributing points to SA reconstructed tracks
140               for(Int_t k=1;k<7;k++){
141                 // check if the layer was skipped
142                 if(skipped[k-1]) continue;
143                 if(hdata->GetBinContent(k)<0.5*(entries/6.)){
144                   cluMapSA = kFALSE;
145                   AliDebug(AliQAv1::GetQADebugLevel(),Form("SA tracks have few points on layer %d - look at histogram hESDClustersSA",k));
146                 }
147               }  
148             }//end clustermapsa 
149             
150             else if(hname.Contains("hESDClusterMapMI") && entries>0.){
151               // Check if there are layers with anomalously low 
152               // contributing points to MI reconstructed tracks
153               AliDebug(AliQAv1::GetQADebugLevel(),Form("Processing histogram %s",hname.Data()));
154               cluMapMI = kTRUE;
155               for(Int_t k=1;k<7;k++){
156                 // check if the layer was skipped
157                 if(skipped[k-1]) continue;
158                 if(hdata->GetBinContent(k)<0.5*(entries/6.)){
159                   cluMapMI = kFALSE;
160                   AliDebug(AliQAv1::GetQADebugLevel(),Form("MI tracks have few points on layer %d - look at histogram hESDClustersMI",k));
161                 }
162               }  
163             }//end clustermapmi
164             
165             else if(hname.Contains("hESDClustersMI") && entries>0.){
166               // Check if 6 clusters MI tracks are the majority
167               AliDebug(AliQAv1::GetQADebugLevel(),Form("Processing histogram %s",hname.Data()));
168               cluMI = kTRUE;
169               Double_t maxlaytracks = hdata->GetBinContent(7-nskipped);
170               for(Int_t k=2; k<7-nskipped; k++){
171                 if(hdata->GetBinContent(k)>maxlaytracks){
172                   cluMI = kFALSE;
173                   AliDebug(AliQAv1::GetQADebugLevel(),Form("MI Tracks with %d clusters are more than tracks with %d clusters. Look at histogram hESDClustersMI",k-1,6-nskipped));
174                 }
175               }
176             }//end clustersmi
177             
178             else if(hname.Contains("hESDClustersSA") && entries>0.){
179               // Check if 6 clusters SA tracks are the majority
180               AliDebug(AliQAv1::GetQADebugLevel(),Form("Processing histogram %s",hname.Data()));
181               cluSA = kTRUE;
182               Double_t maxlaytracks = hdata->GetBinContent(7-nskipped);
183               for(Int_t k=2; k<7-nskipped; k++){
184                 if(hdata->GetBinContent(k)>maxlaytracks){
185                   cluSA = kFALSE;
186                   AliDebug(AliQAv1::GetQADebugLevel(), Form("SA Tracks with %d clusters are more than tracks with %d clusters. Look at histogram hESDClustersSA",k-1,6-nskipped));
187                 }
188               }
189             }//end clusterssa
190             
191             else if(hname.Contains("hSPDVertexZ") && entries>0.){
192               // Check if average Z vertex coordinate is -5 < z < 5 cm
193               AliDebug(AliQAv1::GetQADebugLevel(),Form("Processing histogram %s",hname.Data()));
194               verSPDZ = kTRUE;
195               if(hdata->GetMean()<-5. && hdata->GetMean()>5.){
196                 verSPDZ = kFALSE;
197                 AliDebug(AliQAv1::GetQADebugLevel(),Form("Average z vertex coordinate is at z= %10.4g cm",hdata->GetMean()));
198               }
199             }//end spdvertexz
200             
201             else{ AliError("ESD Checker - invalid data type");}//end else
202             
203             rv[specie] = 0.;
204             if(tested>0){
205               if(tested == empty){
206                 rv[specie] = 2500.; // set to error
207                 AliWarning(Form("All ESD histograms are empty - specie=%d",specie));
208               }
209               else {
210                 rv[specie] = 2500.-1500.*(static_cast<Double_t>(tested-empty)/static_cast<Double_t>(tested)); // INFO if all histos are filled
211                 if(cluMapSA)rv[specie]-=200.;
212                 if(cluMapMI)rv[specie]-=200.;
213                 if(cluMI)rv[specie]-=200.;
214                 if(cluSA)rv[specie]-=200.;
215                 if(verSPDZ)rv[specie]-=199.;  // down to 1 if everything is OK
216               }
217             }//end tested
218           }//end hdata
219         }//end while
220         //     AliDebug(AliQAv1::GetQADebugLevel(), Form("ESD - Tested %d histograms, Return value %f \n",tested,rv[specie]));
221         AliInfo(Form("ESD - Tested %d histograms, Return value %f \n",tested,rv[specie]));
222       }
223     }
224   } // end of ESD QA
225   else{
226     
227     //____________________________________________________________________________
228
229     Double_t spdCheck[AliRecoParam::kNSpecies] ;
230     Double_t sddCheck[AliRecoParam::kNSpecies] ;
231     Double_t ssdCheck[AliRecoParam::kNSpecies] ;
232
233
234
235     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
236       if ( !AliQAv1::Instance()->IsEventSpecieSet(specie)) continue; 
237       if ( AliQAv1::Instance()->IsEventSpecieSet(specie) ) {
238         Double_t histotot=list[specie]->GetEntries();
239         if(histotot!=0)
240           {
241             spdCheck[specie]=0.;
242             sddCheck[specie]=0.;
243             ssdCheck[specie]=0.;
244             rv[specie] = 0.0 ;// 
245             //pixel
246             if(fDet == 0 || fDet == 1) {
247               fSPDChecker->SetTaskOffset(fSPDOffset);
248               //printf("spdoffset = %i \n",fSPDOffset );
249               Double_t histoSPD=double(GetSPDHisto());
250               if(AliITSQADataMakerRec::AreEqual(histoSPD,0)==kFALSE){
251                 Double_t *stepSPD=new Double_t[AliQAv1::kNBIT];
252                 CreateStepForBit(histoSPD,stepSPD);
253                 fSPDChecker->SetStepBit(stepSPD);
254                 spdCheck[specie] = fSPDChecker->Check(index, list[specie], recoParam);
255                 if(spdCheck[specie]>fUpTestValue[AliQAv1::kFATAL]||spdCheck[specie]<0.)
256                   {
257                     AliInfo(Form("SPD check result for %s  is out of range (%f)!!! Retval of specie %s is sit to -1\n ",AliQAv1::GetAliTaskName(index),spdCheck[specie],AliRecoParam::GetEventSpecieName(specie)));
258                     spdCheck[specie]=fUpTestValue[AliQAv1::kFATAL];
259                   }
260                 delete []stepSPD;
261               }//end check SPD entries
262               else{spdCheck[specie]=fUpTestValue[AliQAv1::kFATAL];}
263               rv[specie]=spdCheck[specie];
264             }//end SPD check
265             //drift
266             if(fDet == 0 || fDet == 2) {
267               fSDDChecker->SetTaskOffset(fSDDOffset);
268               fSDDChecker->SetEventSpecieForCheck(specie);
269               Double_t histoSDD=double(GetSDDHisto());
270               if(AliITSQADataMakerRec::AreEqual(histoSDD,0)==kFALSE){
271                 Double_t *stepSDD=new Double_t[AliQAv1::kNBIT];
272                 CreateStepForBit(histoSDD,stepSDD);
273                 fSDDChecker->SetStepBit(stepSDD);
274                 sddCheck[specie] = fSDDChecker->Check(index, list[specie], recoParam);  
275                 if(sddCheck[specie]>fUpTestValue[AliQAv1::kFATAL]||sddCheck[specie]<0.)
276                   {
277                     AliInfo(Form("SDD check result for %s  is out of range (%f)!!! Retval of specie %s is sit to -1\n ",AliQAv1::GetAliTaskName(index),sddCheck[specie],AliRecoParam::GetEventSpecieName(specie)));
278                     sddCheck[specie]=fUpTestValue[AliQAv1::kFATAL];
279                   }
280                 delete []stepSDD;
281               }//end check SDD entries
282               else{ssdCheck[specie]=fUpTestValue[AliQAv1::kFATAL];}
283               if(sddCheck[specie]>rv[specie])rv[specie]=sddCheck[specie];  
284             }//end SDD
285             //strip
286             if(fDet == 0 || fDet == 3) {
287               fSSDChecker->SetTaskOffset(fSSDOffset);
288               Double_t histoSSD=double(GetSSDHisto());
289               if(AliITSQADataMakerRec::AreEqual(histoSSD,0)==kFALSE){
290               Double_t *stepSSD=new Double_t[AliQAv1::kNBIT];
291               CreateStepForBit(histoSSD,stepSSD);
292               fSSDChecker->SetStepBit(stepSSD);
293               ssdCheck[specie] = fSSDChecker->Check(index, list[specie], recoParam);
294               if(ssdCheck[specie]>fUpTestValue[AliQAv1::kFATAL]||ssdCheck[specie]<0.)
295                 {
296                   AliInfo(Form("SSD check result for %s is out of range (%f)!!! Retval of specie %s is sit to -1\n ",AliQAv1::GetAliTaskName(index),ssdCheck[specie],AliRecoParam::GetEventSpecieName(specie)));
297                   ssdCheck[specie]=fUpTestValue[AliQAv1::kFATAL];
298                 }
299               delete [] stepSSD;
300               }//end check SSD entries
301               else{ssdCheck[specie]=fUpTestValue[AliQAv1::kFATAL];}
302               if(ssdCheck[specie]>rv[specie])rv[specie]=ssdCheck[specie];
303             }//end SSD
304             
305             AliInfo(Form("Check result for %s: \n\t  SPD %f \n\t  SDD %f \n\t  SSD %f \n Check result %f \n ",AliQAv1::GetAliTaskName(index),spdCheck[specie],sddCheck[specie],ssdCheck[specie],rv[specie]));
306             // here merging part for common ITS QA result
307             // 
308           }//end entries
309       }//end if event specie
310     }//end for
311   }
312 }
313
314 //____________________________________________________________________________
315 void AliITSQAChecker::SetTaskOffset(Int_t SPDOffset, Int_t SDDOffset, Int_t SSDOffset)
316 {
317   //Setting the 3 offsets for each task called
318   fSPDOffset = SPDOffset;
319   fSDDOffset = SDDOffset;
320   fSSDOffset = SSDOffset;
321 }
322
323 //____________________________________________________________________________
324 void AliITSQAChecker::SetHisto(Int_t SPDhisto, Int_t SDDhisto, Int_t SSDhisto)
325 {
326   //Setting the 3 offsets for each task called
327   fSPDHisto = SPDhisto;
328   fSDDHisto = SDDhisto;
329   fSSDHisto = SSDhisto;
330 }
331
332  //____________________________________________________________________________
333  void AliITSQAChecker::SetDetTaskOffset(Int_t subdet,Int_t offset)
334  {
335    //returns the offset for each task and for each subdetector
336    switch(subdet){
337    case 1:
338      SetSPDTaskOffset(offset);
339      break;
340    case 2:
341      SetSDDTaskOffset(offset);
342      break;
343    case 3:
344      SetSSDTaskOffset(offset);
345      break;
346    default:
347      AliWarning("No specific (SPD,SDD or SSD) subdetector correspond to to this number!!! all offsets set to zero for all the detectors\n");
348      SetTaskOffset(0, 0, 0);
349      break;
350    }
351  }
352
353  //____________________________________________________________________________
354  void AliITSQAChecker::SetDetHisto(Int_t subdet,Int_t histo)
355  {
356    //set the number od histograms for the subdetector
357    switch(subdet){
358    case 1:
359      SetSPDHisto(histo);
360      break;
361    case 2:
362      SetSDDHisto(histo);
363      break;
364    case 3:
365      SetSSDHisto(histo);
366      break;
367    default:
368      AliWarning("No specific (SPD,SDD or SSD) subdetector correspond to to this number!!! all offsets set to zero for all the detectors\n");
369      SetHisto(0, 0, 0);
370      break;
371    }
372  }
373
374 //_____________________________________________________________________________
375
376 void AliITSQAChecker::InitQACheckerLimits()
377 {
378   //init the tolerance range for each QA bit 
379   AliInfo("Setting of tolerance values\n");
380
381   Float_t lowtolerancevalue[AliQAv1::kNBIT];
382
383   Float_t hightolerancevalue[AliQAv1::kNBIT];
384   for(Int_t bit=0;bit<AliQAv1::kNBIT;bit++)
385     {
386       lowtolerancevalue[bit]=(bit*1000.);
387       hightolerancevalue[bit]=((bit+1.)*1000.);
388     }
389   SetHiLo(hightolerancevalue,lowtolerancevalue);
390   //  AliInfo(Form("Range Value  \n INFO    -> %f <  value <  %f \n WARNING -> %f <  value <= %f \n ERROR   -> %f <  value <= %f \n FATAL   -> %f <= value <  %f \n", fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO], fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING], fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR], fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]  ));
391
392   if(fDet == 0 || fDet == 1) {
393     fSPDChecker->SetSPDLimits( lowtolerancevalue,hightolerancevalue );
394   }
395   if(fDet == 0 || fDet == 2) {
396     fSDDChecker->SetSDDLimits( lowtolerancevalue,hightolerancevalue );
397   }
398   if(fDet == 0 || fDet == 3) {
399     fSSDChecker->SetSSDLimits( lowtolerancevalue,hightolerancevalue );
400   }
401
402
403   
404 }
405
406
407 //_____________________________________________________________________________
408
409 void AliITSQAChecker::CreateStepForBit(Double_t histonumb,Double_t *steprange)
410 {
411   //creation of the step bit for each QA bit 
412   for(Int_t bit=0;bit < AliQAv1::kNBIT; bit++)
413     {       
414       //printf("%i\t %f \t %f \t %f \n",bit, fUpTestValue[bit],fLowTestValue[AliQAv1::kINFO],histonumb);
415       steprange[bit]=double((fUpTestValue[bit] - fLowTestValue[AliQAv1::kINFO])/histonumb);
416       //printf("%i\t %f \t %f \t %f \t %f\n",bit, fUpTestValue[bit],fLowTestValue[AliQAv1::kINFO],histonumb,steprange[bit] );
417     }
418   //AliInfo(Form("StepBitValue:numner of histo %f\n\t INFO %f \t WARNING %f \t ERROR %f \t FATAL %f \n",histonumb, steprange[AliQAv1::kINFO],steprange[AliQAv1::kWARNING],steprange[AliQAv1::kERROR],steprange[AliQAv1::kFATAL]));
419 }
420
421
422 //_____________________________________________________________________________
423 void AliITSQAChecker::SetQA(AliQAv1::ALITASK_t index, Double_t * value) const
424 {
425   //Setting of the QA tolerance values
426   AliQAv1 * qa = AliQAv1::Instance(index) ;
427
428
429   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
430
431     if (! qa->IsEventSpecieSet(AliRecoParam::ConvertIndex(specie)))
432       continue ;
433     if (  value == NULL ) { // No checker is implemented, set all QA to Fatal
434       qa->Set(AliQAv1::kFATAL, specie) ; 
435     } else {
436       if ( value[specie] > fLowTestValue[AliQAv1::kFATAL] && value[specie] <= fUpTestValue[AliQAv1::kFATAL] ) 
437         qa->Set(AliQAv1::kFATAL, AliRecoParam::ConvertIndex(specie)) ; 
438       else if ( value[specie] > fLowTestValue[AliQAv1::kERROR] && value[specie] <= fUpTestValue[AliQAv1::kERROR]  )
439         qa->Set(AliQAv1::kERROR, AliRecoParam::ConvertIndex(specie)) ; 
440       else if ( value[specie] > fLowTestValue[AliQAv1::kWARNING] && value[specie] <= fUpTestValue[AliQAv1::kWARNING]  )
441         qa->Set(AliQAv1::kWARNING, AliRecoParam::ConvertIndex(specie)) ;
442       else if ( value[specie] > fLowTestValue[AliQAv1::kINFO] && value[specie] <= fUpTestValue[AliQAv1::kINFO] ) 
443         qa->Set(AliQAv1::kINFO, AliRecoParam::ConvertIndex(specie)) ;   
444       //else if(value[specie]==0) qa->Set(AliQAv1::kFATAL, AliRecoParam::ConvertIndex(specie)) ; //no ckeck has been done
445     }
446     qa->ShowStatus(AliQAv1::kITS,index,AliRecoParam::ConvertIndex(specie));
447   }//end for
448
449 }
450
451
452 //__________________________________________________________________
453 void  AliITSQAChecker::MakeImage( TObjArray ** list, AliQAv1::TASKINDEX_t task, AliQAv1::MODE_t mode)
454 {
455   //make a summary image
456   //gStyle->SetPalette(1);
457
458   //Int_t nImages = 0 ;
459   //Int_t imageindex=0;
460   for (Int_t esIndex = 0 ; esIndex < AliRecoParam::kNSpecies ; esIndex++) {
461     if (! AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))->IsEventSpecieSet(AliRecoParam::ConvertIndex(esIndex)) ) 
462       continue ;
463     //else imageindex=esIndex;
464
465     TIter next(list[esIndex]) ;  
466     TH1 * hdata = NULL ; 
467     while ( (hdata=static_cast<TH1 *>(next())) ) {
468       TString cln(hdata->ClassName()) ; 
469       if ( ! cln.Contains("TH") )
470         continue ; 
471       if(cln.Contains("TH2")) 
472         {
473           Float_t min=hdata->GetMinimum();
474           Float_t max=hdata->GetMaximum();
475           if(max>min) hdata->SetOption("colz");
476         }
477     }
478     break ; 
479   }
480
481   Bool_t retvalue=kFALSE;
482
483   if(GetSubDet()==0) MakeITSImage(list,task, mode);
484   else if(GetSubDet()==1) 
485     {
486       retvalue=fSPDChecker->MakeSPDImage(list,task, mode) ; 
487       if(retvalue==kFALSE)AliQACheckerBase::MakeImage(list,task, mode);
488     }
489   else if(GetSubDet()==2){ retvalue=fSDDChecker->MakeSDDImage(list,task, mode) ;if(retvalue==kFALSE)AliQACheckerBase::MakeImage(list,task,mode); }
490   else if(GetSubDet()==3) 
491     {
492       retvalue=fSSDChecker->MakeSSDImage(list,task, mode) ;
493       if(retvalue==kFALSE)AliQACheckerBase::MakeImage(list,task, mode); 
494     }
495
496 }
497