]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSOnlineSPDscanAnalyzer.cxx
New class which is used by the SPD online detector algorithm (Henrik)
[u/mrichter/AliRoot.git] / ITS / AliITSOnlineSPDscanAnalyzer.cxx
1 ////////////////////////////////////////////////////////////
2 // Author: Henrik Tydesjo                                 //
3 // This class is used in the detector algorithm framework //
4 // to process the data stored in special container files  //
5 // (see AliITSOnlineSPDscan). For instance, minimum       //
6 // threshold values can be calculated.                    //
7 ////////////////////////////////////////////////////////////
8
9 #include "AliITSOnlineSPDscanAnalyzer.h"
10 #include "AliITSOnlineSPDscan.h"
11 #include "AliITSOnlineSPDscanSingle.h"
12 #include "AliITSOnlineSPDscanMultiple.h"
13 #include "AliITSOnlineSPDscanMeanTh.h"
14 #include "AliITSOnlineCalibrationSPDhandler.h"
15 #include "AliITSRawStreamSPD.h"
16 #include <TStyle.h>
17 #include <TMath.h>
18 #include <TF1.h>
19 #include <TGraph.h>
20 #include <TH2F.h>
21
22 Double_t itsSpdErrorf(Double_t *x, Double_t *par){
23   if (par[2]<0) par[2]=0;
24   Double_t val = par[2]+(0.12*256*32-par[2])*(0.5+0.5*TMath::Erf((x[0]-par[0])/par[1]/sqrt(2.)));
25   return val;
26 }
27 //Double_t itsSpdErrorfOrig(Double_t *x, Double_t *par){
28 //  return 0.5+0.5*TMath::Erf((x[0]-par[0])/par[1]/sqrt(2.));
29 //}
30
31
32 AliITSOnlineSPDscanAnalyzer::AliITSOnlineSPDscanAnalyzer(Char_t *fileName) :
33   fType(99),fDacId(99),fScanObj(NULL),fTriggers(NULL),
34   fOverWrite(kFALSE),fNoiseThreshold(0.01),fNoiseMinimumEvents(100),
35   fMinNrStepsBeforeIncrease(5),fMinIncreaseFromBaseLine(2),fStepDownDacSafe(2),fMaxBaseLineLevel(10)
36 {
37   // constructor
38   sprintf(fFileName,"%s",fileName);
39   for (UInt_t chipNr=0; chipNr<11; chipNr++) {
40     for (UInt_t hs=0; hs<6; hs++) {
41       fMeanMultiplicity[hs][chipNr]=NULL;
42       fHitEventEfficiency[hs][chipNr]=NULL;
43     }
44   }
45   for (Int_t module=0; module<240; module++) {
46     fHandler[module]=NULL;
47   }
48
49   Init();
50 }
51
52 AliITSOnlineSPDscanAnalyzer::AliITSOnlineSPDscanAnalyzer(const AliITSOnlineSPDscanAnalyzer& handle) :
53   fType(99),fDacId(99),fScanObj(NULL),fTriggers(NULL),
54   fOverWrite(kFALSE),fNoiseThreshold(0.01),fNoiseMinimumEvents(100),
55   fMinNrStepsBeforeIncrease(5),fMinIncreaseFromBaseLine(2),fStepDownDacSafe(2),fMaxBaseLineLevel(10)
56 {
57   // copy constructor, only copies the filename (not the processed data)
58   sprintf(fFileName,"%s",handle.fFileName);
59
60   fScanObj=NULL;
61   fType=99;
62   fDacId=99;
63   for (UInt_t chipNr=0; chipNr<11; chipNr++) {
64     for (UInt_t hs=0; hs<6; hs++) {
65       fMeanMultiplicity[hs][chipNr]=NULL;
66       fHitEventEfficiency[hs][chipNr]=NULL;
67     }
68   }
69   fTriggers=NULL;
70   for (Int_t module=0; module<240; module++) {
71     fHandler[module]=NULL;
72   }
73
74   Init();
75 }
76
77 AliITSOnlineSPDscanAnalyzer::~AliITSOnlineSPDscanAnalyzer() {
78   // destructor
79   for (UInt_t hs=0; hs<6; hs++) {
80     for (UInt_t chipNr=0; chipNr<11; chipNr++) {
81       if (fMeanMultiplicity[hs][chipNr]!=NULL) {
82         delete fMeanMultiplicity[hs][chipNr];
83       }
84       if (fHitEventEfficiency[hs][chipNr]!=NULL) {
85         delete fHitEventEfficiency[hs][chipNr];
86       }
87     }
88   }
89   if (fTriggers!=NULL) delete fTriggers;
90   if (fScanObj!=NULL) delete fScanObj;
91   for (Int_t module=0; module<240; module++) {
92     if (fHandler[module]!=NULL) {
93       delete fHandler[module];
94     }
95   }
96 }
97
98 AliITSOnlineSPDscanAnalyzer& AliITSOnlineSPDscanAnalyzer::operator=(const AliITSOnlineSPDscanAnalyzer& handle) {
99   // assignment operator, only copies the filename (not the processed data)
100   if (this!=&handle) {
101     for (UInt_t hs=0; hs<6; hs++) {
102       for (UInt_t chipNr=0; chipNr<11; chipNr++) {
103         if (fMeanMultiplicity[hs][chipNr]!=NULL) {
104           delete fMeanMultiplicity[hs][chipNr];
105         }
106         if (fHitEventEfficiency[hs][chipNr]!=NULL) {
107           delete fHitEventEfficiency[hs][chipNr];
108         }
109       }
110     }
111     if (fTriggers!=NULL) delete fTriggers;
112     if (fScanObj!=NULL) delete fScanObj;
113     for (Int_t module=0; module<240; module++) {
114       if (fHandler[module]!=NULL) {
115         delete fHandler[module];
116       }
117     }
118    
119     sprintf(fFileName,"%s",handle.fFileName);
120
121     fScanObj=NULL;
122     fType=99;
123     fDacId=99;
124     for (UInt_t chipNr=0; chipNr<11; chipNr++) {
125       for (UInt_t hs=0; hs<6; hs++) {
126         fMeanMultiplicity[hs][chipNr]=NULL;
127         fHitEventEfficiency[hs][chipNr]=NULL;
128       }
129     }
130     fTriggers=NULL;
131     for (Int_t module=0; module<240; module++) {
132       fHandler[module]=NULL;
133     }
134     
135     Init();    
136   }
137   return *this;
138 }
139
140 void AliITSOnlineSPDscanAnalyzer::Init() {
141   // first checks type of container and then initializes container obj
142   FILE* fp0 = fopen(fFileName, "r");
143   if (fp0 == NULL) {
144     return;
145   }
146   else {
147     fclose(fp0);
148   }
149   fScanObj = new AliITSOnlineSPDscan(fFileName);
150   fType = fScanObj->GetType();
151   delete fScanObj;
152
153   // init container
154   switch(fType) {
155   case kUNIMA:
156   case kNOISE:
157     fScanObj = new AliITSOnlineSPDscanSingle(fFileName);
158     break;
159   case kMINTH:
160   case kDAC:
161   case kDELAY:
162     fScanObj = new AliITSOnlineSPDscanMultiple(fFileName);
163     fDacId = ((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacId();
164     break;
165   case kMEANTH:
166     fScanObj = new AliITSOnlineSPDscanMeanTh(fFileName);
167     fDacId = ((AliITSOnlineSPDscanMeanTh*)fScanObj)->GetDacId();
168     break;
169   default:
170     printf("Type %d not defined!\n",fType);
171     fScanObj=NULL;
172     return;
173     break;
174   }
175
176   // set some default values (these should later be read from text file)
177   fOverWrite=kFALSE;
178   fNoiseThreshold=0.01;
179   fNoiseMinimumEvents=100;
180   fMinNrStepsBeforeIncrease=6;
181   fMinIncreaseFromBaseLine=2;
182   fStepDownDacSafe=2;
183   fMaxBaseLineLevel=10;
184
185 }
186
187
188 Bool_t AliITSOnlineSPDscanAnalyzer::ProcessDeadPixels(Char_t *oldcalibDir) {
189   // process dead pixel data
190   if (fScanObj==NULL) {
191     printf("No data!\n");
192     return kFALSE;
193   }
194   // should be type kUNIMA
195   if (fType!=kUNIMA) {
196     printf("Dead pixels only for scan type %d\n",kUNIMA);
197     return kFALSE;
198   }
199
200   Int_t nrDead[240];
201   for (Int_t i=0; i<240; i++) {
202     nrDead[i]=0;
203   }
204   UInt_t routerNr = fScanObj->GetRouterNr();
205   UInt_t rowStart = fScanObj->GetRowStart();
206   UInt_t rowEnd   = fScanObj->GetRowEnd();
207   for (UInt_t hs=0; hs<6; hs++) {
208     for (UInt_t chipNr=0; chipNr<10; chipNr++) {
209       if (fScanObj->GetChipPresent(hs,chipNr) && fScanObj->GetAverageMultiplicity(0,hs,chipNr)>0) { // check the status of the chippresent parameter in the mood header!!!!!!!!!!!!!!!!!!!!!!!!!!!!
210         UInt_t module = AliITSRawStreamSPD::GetModuleNumber(routerNr,hs*2+chipNr/5);
211         if (!fHandler[module]) {
212           fHandler[module] = new AliITSOnlineCalibrationSPDhandler(module);
213         }
214         fHandler[module]->SetFileLocation(oldcalibDir);
215         fHandler[module]->ReadFromFile();
216         if (fOverWrite) {fHandler[module]->ResetDead();}
217         for (UInt_t col=0; col<32; col++) {
218           for (UInt_t row=rowStart; row<=rowEnd; row++) {
219             if (col!=1 && col!=9 && col!=17 && col!=25) { //exclude test columns!!!
220               if (fScanObj->GetHits(0,hs,chipNr,col,row)==0) {
221                 Int_t newCol = 32*(chipNr%5) + col;
222                 if (fHandler[module]->SetDeadPixel(newCol,row)) {
223                   nrDead[module]++;
224                 }
225               }
226             }
227           }
228         }
229       }
230     }
231   }
232   return kTRUE;
233 }
234
235
236
237 Bool_t AliITSOnlineSPDscanAnalyzer::ProcessNoisyPixels(Char_t *oldcalibDir) {
238   // process noisy pixel data
239   if (fScanObj==NULL) {
240     printf("No data!\n");
241     return kFALSE;
242   }
243   // should be type kNOISE
244   if (fType != kNOISE) {
245     printf("Noisy pixels only for scan type %d\n",kNOISE);
246     return kFALSE;
247   }
248
249   Int_t nrNoisy[240];
250   for (Int_t i=0; i<240; i++) {
251     nrNoisy[i]=0;
252   }
253   if (fScanObj->GetTriggers(0)<fNoiseMinimumEvents) {
254     printf("Error: Process noisy: Too few events.\n"); 
255     return kFALSE;
256   }
257   UInt_t routerNr = fScanObj->GetRouterNr();
258   for (UInt_t hs=0; hs<6; hs++) {
259     for (UInt_t chipNr=0; chipNr<10; chipNr++) {
260       if (fScanObj->GetChipPresent(hs,chipNr)) { // check the status of the chippresent parameter in the mood header!!!!!!!!!!!!!!!!!!!!!!!!!!!!
261         UInt_t module = AliITSRawStreamSPD::GetModuleNumber(routerNr,hs*2+chipNr/5);
262         for (UInt_t col=0; col<32; col++) {
263           for (UInt_t row=0; row<256; row++) {
264             if (fScanObj->GetHitsEfficiency(0,hs,chipNr,col,row)>fNoiseThreshold) {
265               if (!fHandler[module]) {
266                 fHandler[module] = new AliITSOnlineCalibrationSPDhandler(module);
267               }
268               fHandler[module]->SetFileLocation(oldcalibDir);
269               fHandler[module]->ReadFromFile();
270               if (fOverWrite) {fHandler[module]->ResetNoisy();}
271               Int_t newCol = 32*(chipNr%5) + col;
272               if (fHandler[module]->SetNoisyPixel(newCol,row)) {
273                 nrNoisy[module]++;
274               }
275             }
276           }
277         }
278       }
279     }
280   }
281   return kTRUE;
282 }
283
284 Bool_t AliITSOnlineSPDscanAnalyzer::SaveDeadNoisyPixels(UInt_t module, Char_t *calibDir) {
285   // save dead and noisy pixels to file in dir calibDir
286   if (fHandler[module]!=NULL) {
287     fHandler[module]->SetFileLocation(calibDir);
288     fHandler[module]->WriteToFile();
289     return kTRUE;
290   }
291   return kFALSE;
292 }
293
294 Int_t AliITSOnlineSPDscanAnalyzer::GetDelay(UInt_t hs, UInt_t chipNr) {
295   // get delay
296   if (hs>=6 || chipNr>10) return -1;
297   if (fScanObj==NULL) {
298     printf("No data!\n");
299     return -1;
300   }
301   // should be type kDELAY or kDAC with id 42 (delay_ctrl)
302   if (fType!=kDELAY && (fType!=kDAC || fDacId!=42)) {
303     printf("Delay only for scan type %d or %d and dac_id 42\n",kDELAY,kDAC);
304     return -1;
305   }
306   if (fMeanMultiplicity[hs][chipNr]==NULL) {
307     if (!ProcessMeanMultiplicity()) {
308       return -1;
309     }
310   }
311
312   Char_t funcName[30];
313   sprintf(funcName,"Fit delay func HS%d CHIP%d",hs,chipNr);
314   Int_t minDac = ((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(0);
315   Int_t maxDac = ((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(fScanObj->GetNSteps()-1);
316   TF1* delayFunc = new TF1(funcName,"gaus",minDac,maxDac);
317   fMeanMultiplicity[hs][chipNr]->Fit(funcName,"Q0");
318   Double_t mean = delayFunc->GetParameter(1);
319   //  Double_t sigma = fDelayFunc[hs][chipNr]->GetParameter(2);
320   delete delayFunc;
321   if (mean>minDac && mean<maxDac) {
322     return (Int_t) (mean+0.5);
323   }
324   else {
325     return -1;
326   }
327 }
328
329 Int_t AliITSOnlineSPDscanAnalyzer::GetNrNoisyUnima(UInt_t hs, UInt_t chipNr) {
330   // in case of a uniformity scan, returns the nr of noisy pixels, (here > 200 hits)
331   if (hs>=6 || chipNr>10) return -1;
332   if (fScanObj==NULL) {
333     printf("No data!\n");
334     return kFALSE;
335   }
336   // should be type kUNIMA
337   if (fType != kUNIMA) {
338     printf("Noisy pixels Unima only for scan type %d\n",kUNIMA);
339     return kFALSE;
340   }
341   if (fScanObj->GetTriggers(0)!=25600) {
342     printf("Error: Process noisy unima: Incorrect number of events (!=25600.\n"); 
343     return kFALSE;
344   }
345
346   Int_t nrNoisy=0;
347   if (fScanObj->GetChipPresent(hs,chipNr)) { // check the status of the chippresent parameter in the mood header!!!!!!!!!!!!!!!!!!!!!!!!!!!!
348     for (UInt_t col=0; col<32; col++) {
349       for (UInt_t row=0; row<256; row++) {
350         if (fScanObj->GetHits(0,hs,chipNr,col,row)>200) {
351           nrNoisy++;
352         }
353       }
354     }
355   }
356   else {
357     return -1;
358   }
359   return nrNoisy;
360 }
361
362 Int_t AliITSOnlineSPDscanAnalyzer::FindLastMinThDac(UInt_t hs, UInt_t chipNr) {
363   // returns dac value where fMinIncreaseFromBaseLine reached
364   if (hs>=6 || chipNr>10) return -1;
365   if (fMeanMultiplicity[hs][chipNr]==NULL) {
366     if (!ProcessMeanMultiplicity()) {
367       return -1;
368     }
369   }
370   Double_t firstVal, dummy1;
371   fMeanMultiplicity[hs][chipNr]->GetPoint(0,dummy1,firstVal);
372   UInt_t step=0;
373   while (step<fScanObj->GetNSteps()-1) {
374     Double_t graphVal, dummy2;
375     fMeanMultiplicity[hs][chipNr]->GetPoint(step+1,dummy2,graphVal);
376     if (graphVal>firstVal+fMinIncreaseFromBaseLine) break;
377     step++;
378   }
379   if (step==fScanObj->GetNSteps()-1) return -1;
380   return ((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(step);
381 }
382
383 Int_t AliITSOnlineSPDscanAnalyzer::FindClosestLowerStep(Float_t dacValueInput) {
384   // returns step closest (lower) to a dacvalue 
385   UInt_t step=0;
386   while (step<fScanObj->GetNSteps()-1) {
387     Int_t dacVal = ((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(step+1);
388     if (dacVal>=dacValueInput) break;
389     step++;
390   }
391   return step;
392 }
393
394 Float_t AliITSOnlineSPDscanAnalyzer::GetCompareLine(UInt_t step, UInt_t hs, UInt_t chipNr, Float_t basePar2) {
395   // returns value to compare mean mult with (when finding min th)
396   if (hs>=6 || chipNr>10) return -1;
397   if (step<fMinNrStepsBeforeIncrease) return -1;
398   Float_t baseLine = basePar2;
399   if (baseLine<0) baseLine=0;
400   Float_t baseAdd;
401   Double_t baseM=0;
402   Double_t baseS=0;
403   Double_t d,m;
404   for (UInt_t st=1;st<2*step/3;st++) { // skip first point...
405     fMeanMultiplicity[hs][chipNr]->GetPoint(st,d,m);
406     baseM+=m-baseLine;
407     baseS+=(m-baseLine)*(m-baseLine);
408   }
409   baseAdd=2*sqrt( baseS/(2*step/3-1) - (baseM/(2*step/3-1))*(baseM/(2*step/3-1)) );
410   baseAdd+=0.03; // magic number
411   if (baseAdd>fMinIncreaseFromBaseLine) baseAdd=fMinIncreaseFromBaseLine;
412   return baseLine + baseAdd;
413 }
414
415 Int_t AliITSOnlineSPDscanAnalyzer::GetMinTh(UInt_t hs, UInt_t chipNr) {
416   // calculates and returns the minimum threshold
417   if (hs>=6 || chipNr>10) return -1;
418   if (fScanObj==NULL) {
419     printf("No data!\n");
420     return -1;
421   }
422   // should be type  kMINTH  or  kDAC with id 39 (pre_vth)
423   if (fType!=kMINTH && (fType!=kDAC || fDacId!=39)) {
424     printf("MinTh only for scan type %d OR %d with dac_id 39\n",kMINTH,kDAC);
425     return -1;
426   }
427   if (fMeanMultiplicity[hs][chipNr]==NULL) {
428     if (!ProcessMeanMultiplicity()) {
429       return -1;
430     }
431   }
432
433   Int_t lastDac = FindLastMinThDac(hs,chipNr);
434   if (lastDac==-1) {
435     printf("GetMinTh: HS%d, Chip%d: Increase of Mean Multiplicity by %1.2f never reached.\n",hs,chipNr,fMinIncreaseFromBaseLine);
436     return -1;
437   }
438
439   Int_t minDac = ((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(0);
440   Char_t funcName[30];
441   sprintf(funcName,"Fit minth func HS%d CHIP%d",hs,chipNr);
442   TF1 *minThFunc = new TF1(funcName,itsSpdErrorf,100,500,3);
443   minThFunc->SetParameter(0,lastDac+10);
444   minThFunc->SetParameter(1,2);
445   minThFunc->SetParameter(2,0);
446   minThFunc->SetParName(0,"Mean");
447   minThFunc->SetParName(1,"Sigma");
448   minThFunc->SetParName(2,"BaseLine");
449   minThFunc->SetLineWidth(1);
450   if (fMeanMultiplicity[hs][chipNr]==NULL) {
451     if (!ProcessMeanMultiplicity()) {
452       return -1;
453     }
454   }
455   fMeanMultiplicity[hs][chipNr]->Fit(funcName,"Q0","",minDac,lastDac);
456
457   //  Double_t mean = fMinThFunc[hs][chipNr]->GetParameter(0);
458   //  Double_t sigma = fMinThFunc[hs][chipNr]->GetParameter(1);
459   Double_t baseLine = minThFunc->GetParameter(2);
460   delete minThFunc;
461
462   if (baseLine>fMaxBaseLineLevel) {
463     printf("GetMinTh: HS%d, Chip%d: BaseLine too large (%1.2f>%1.2f)\n",hs,chipNr,baseLine,fMaxBaseLineLevel);
464     return -1;
465   }
466   UInt_t step=FindClosestLowerStep(lastDac);
467   Float_t compareLine = GetCompareLine(step,hs,chipNr,baseLine);
468   if (compareLine==-1) {
469     printf("GetMinTh: HS%d, Chip%d: Not enough steps (%d<%d) before increase to get a compare line.\n",hs,chipNr,step,fMinNrStepsBeforeIncrease);
470     return -1;
471   }
472
473   Double_t mult, dummy;
474   mult=1000;
475   while (mult > compareLine && step>0) {
476     fMeanMultiplicity[hs][chipNr]->GetPoint(step,dummy,mult);
477     step--;
478   }
479   Int_t minth = ((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(step+1)-fStepDownDacSafe;
480
481   if (step>0) {
482     return minth;
483   }
484   else {
485     printf("GetMinTh: HS%d, Chip%d: Did not find a point below the compare line (%f).\n",hs,chipNr,compareLine);
486     return -1;
487   }
488 }
489
490
491
492 Bool_t AliITSOnlineSPDscanAnalyzer::ProcessMeanMultiplicity() {
493   // process mean multiplicity data
494   if (fScanObj==NULL) {
495     printf("No data!\n");
496     return kFALSE;
497   }
498   for (UInt_t step=0; step<fScanObj->GetNSteps(); step++) {
499     for (UInt_t hs=0; hs<6; hs++) {
500       for (UInt_t chipNr=0; chipNr<11; chipNr++) {
501         //      if (fScanObj->GetChipPresent(hs,chipNr)) { // check the status of the chippresent parameter in the mood header!!!!!!!!!!!!!!!!!!!!!!!!!!!!
502         if (step==0) {
503           if (fMeanMultiplicity[hs][chipNr]!=NULL) {
504             delete fMeanMultiplicity[hs][chipNr];
505           }
506           fMeanMultiplicity[hs][chipNr] = new TGraph();
507         }
508         Float_t multiplMean=fScanObj->GetAverageMultiplicity(step,hs,chipNr);
509         if (fType==kMINTH || fType==kMEANTH || fType==kDAC || fType==kDELAY) {
510           fMeanMultiplicity[hs][chipNr]->SetPoint(step,((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(step),multiplMean);
511         }
512         else {
513           fMeanMultiplicity[hs][chipNr]->SetPoint(step,0,multiplMean);
514         }
515       }
516       //      }
517     }
518   }
519   return kTRUE;
520 }
521
522 TGraph* AliITSOnlineSPDscanAnalyzer::GetMeanMultiplicityG(UInt_t hs, UInt_t chipNr) {
523   // returns mean multiplicity graph
524   if (hs>=6 || chipNr>10) return NULL;
525   if (fMeanMultiplicity[hs][chipNr]==NULL) {
526     if (!ProcessMeanMultiplicity()) {
527       return NULL;
528     }
529   }
530   return fMeanMultiplicity[hs][chipNr];
531 }
532
533 Bool_t AliITSOnlineSPDscanAnalyzer::ProcessHitEventEfficiency() {
534   // process hit event efficiency
535   if (fScanObj==NULL) {
536     printf("No data!\n");
537     return kFALSE;
538   }
539   for (UInt_t step=0; step<fScanObj->GetNSteps(); step++) {
540     for (UInt_t hs=0; hs<6; hs++) {
541       for (UInt_t chipNr=0; chipNr<11; chipNr++) {
542         //      if (fScanObj->GetChipPresent(hs,chipNr)) { // check the status of the chippresent parameter in the mood header!!!!!!!!!!!!!!!!!!!!!!!!!!!!
543         if (step==0) {
544           if (fHitEventEfficiency[hs][chipNr]!=NULL) {
545             delete fHitEventEfficiency[hs][chipNr];
546           }
547           fHitEventEfficiency[hs][chipNr] = new TGraph();
548         }
549         Float_t efficiency=fScanObj->GetHitEventsEfficiency(step,hs,chipNr);
550         if (fType==kMINTH || fType==kMEANTH || fType==kDAC || fType==kDELAY) {
551           fHitEventEfficiency[hs][chipNr]->SetPoint(step,((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(step),efficiency);
552         }
553         else {
554           fHitEventEfficiency[hs][chipNr]->SetPoint(step,0,efficiency);
555         }
556       }
557       //      }
558     }
559   }
560   return kTRUE;
561 }
562
563 TGraph* AliITSOnlineSPDscanAnalyzer::GetHitEventEfficiencyG(UInt_t hs, UInt_t chipNr) {
564   // returns hit event efficiency graph
565   if (hs>=6 || chipNr>10) return NULL;
566   if (fHitEventEfficiency[hs][chipNr]==NULL) {
567     if (!ProcessHitEventEfficiency()) {
568       return NULL;
569     }
570   }
571   return fHitEventEfficiency[hs][chipNr];
572 }
573
574
575 Bool_t AliITSOnlineSPDscanAnalyzer::ProcessNrTriggers() {
576   // process nr of triggers data
577   if (fScanObj==NULL) {
578     printf("No data!\n");
579     return kFALSE;
580   }
581   for (UInt_t step=0; step<fScanObj->GetNSteps(); step++) {
582     if (step==0) {
583       if (fTriggers!=NULL) {
584         delete fTriggers;
585       }
586       fTriggers = new TGraph();
587     }
588     if (fType==kMINTH || fType==kMEANTH || fType==kDAC || fType==kDELAY) {
589       fTriggers->SetPoint(step,((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(step),fScanObj->GetTriggers(step));
590     }
591     else {
592       fTriggers->SetPoint(step,0,fScanObj->GetTriggers(step));
593     }
594   }
595   return kTRUE;
596 }
597
598 TGraph* AliITSOnlineSPDscanAnalyzer::GetNrTriggersG() {
599   // returns nr of triggers graph
600   if (fTriggers==NULL) {
601     if (!ProcessNrTriggers()) {
602       return NULL;
603     }
604   }
605   return fTriggers;
606 }
607
608 Bool_t AliITSOnlineSPDscanAnalyzer::GetHalfStavePresent(UInt_t hs) {
609   // returns half stave present info
610   if (hs<6 && fScanObj!=NULL) {
611     Int_t chipstatus=0;
612     for (Int_t chip=0; chip<10; chip++) {
613       chipstatus+=fScanObj->GetChipPresent(hs,chip);
614     }
615     if (chipstatus>0) return kTRUE;
616   }
617   return kFALSE;
618 }
619
620 AliITSOnlineCalibrationSPDhandler* AliITSOnlineSPDscanAnalyzer::GetOnlineCalibrationHandler(UInt_t module) {
621   // returns a pointer to the AliITSOnlineCalibrationSPDhandler
622   if (module<240) return fHandler[module]; 
623   else return NULL;
624 }
625
626 UInt_t AliITSOnlineSPDscanAnalyzer::GetRouterNr() {
627   // returns the router nr of scan obj
628   if (fScanObj!=NULL) return fScanObj->GetRouterNr(); 
629   else return 99;
630 }
631
632 TH2F* AliITSOnlineSPDscanAnalyzer::GetHitMapTot(UInt_t step) {
633   // creates and returns a pointer to a hitmap histo (half sector style a la spdmood)
634   if (fScanObj==NULL) {
635     printf("No data!\n");
636     return NULL;
637   }
638   Char_t histoname[50];
639   if (fType==kMINTH || fType==kMEANTH || fType==kDAC) {
640     sprintf(histoname,"Router %d , DAC %d",GetRouterNr(),((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(step));
641   }
642   else {
643     sprintf(histoname,"Router %d ",GetRouterNr());
644   }
645   TH2F* fHitMapTot = new TH2F(histoname,histoname,32*10,-0.5,32*10-0.5,256*6,-0.5,256*6-0.5);
646   fHitMapTot->SetNdivisions(-10,"X");
647   fHitMapTot->SetNdivisions(-006,"Y");
648   fHitMapTot->SetTickLength(0,"X");
649   fHitMapTot->SetTickLength(0,"Y");
650   fHitMapTot->GetXaxis()->SetLabelColor(gStyle->GetCanvasColor());
651   fHitMapTot->GetYaxis()->SetLabelColor(gStyle->GetCanvasColor());
652   for (UInt_t hs=0; hs<6; hs++) {
653     for (UInt_t chipNr=0; chipNr<10; chipNr++) {
654       for (UInt_t col=0; col<32; col++) {
655         for (UInt_t row=0; row<256; row++) {
656           fHitMapTot->Fill(chipNr*32+col,(5-hs)*256+row,fScanObj->GetHits(step,hs,chipNr,col,row));
657         }
658       }
659     }
660   }
661   return fHitMapTot;
662 }