]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSOnlineSPDscanAnalyzer.cxx
Fix for Savannah bug 67210 (Peter Hristov)
[u/mrichter/AliRoot.git] / ITS / AliITSOnlineSPDscanAnalyzer.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 // Author: Henrik Tydesjo                                 //
20 // This class is used in the detector algorithm framework //
21 // to process the data stored in special container files  //
22 // (see AliITSOnlineSPDscan). For instance, minimum       //
23 // threshold values can be extracted.                     //
24 ////////////////////////////////////////////////////////////
25
26 #include "AliITSOnlineSPDscanAnalyzer.h"
27 #include "AliITSOnlineSPDscan.h"
28 #include "AliITSOnlineSPDscanSingle.h"
29 #include "AliITSOnlineSPDscanMultiple.h"
30 #include "AliITSOnlineSPDscanMeanTh.h"
31 #include "AliITSOnlineCalibrationSPDhandler.h"
32 #include "AliITSRawStreamSPD.h"
33 #include <TStyle.h>
34 #include <TMath.h>
35 #include <TLine.h>
36 #include <TF1.h>
37 #include <TGraph.h>
38 #include <TH2F.h>
39 #include <TError.h>
40 #include <iostream>
41 #include <fstream>
42
43 Double_t itsSpdErrorf(Double_t *x, Double_t *par){
44   if (par[2]<0) par[2]=0;
45   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.)));
46   return val;
47 }
48 //Double_t itsSpdErrorfOrig(Double_t *x, Double_t *par){
49 //  return 0.5+0.5*TMath::Erf((x[0]-par[0])/par[1]/sqrt(2.));
50 //}
51 //_________________________________________________________________________
52 Double_t itsSpdScurveForMeanTh(Double_t *x, Double_t *par){
53   if (par[2]<0) par[2]=0;
54   Double_t val = 1.- par[2]*(1.-TMath::Erf((x[0]-par[0])/par[1]/sqrt(2.)));
55 //  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.)));
56   return val;
57 }
58
59 //_________________________________________________________________________
60 AliITSOnlineSPDscanAnalyzer::AliITSOnlineSPDscanAnalyzer(const Char_t *fileName, AliITSOnlineCalibrationSPDhandler *handler, Bool_t readFromGridFile) :
61   fType(99),fDacId(99),fFileName(fileName),fScanObj(NULL),fHandler(handler),fTriggers(NULL),fTPeff(0),fTPeffHS(NULL),fDeadPixel(0),fDeadPixelHS(NULL),fNoisyPixel(0),fNoisyPixelHS(NULL),
62   fOverWrite(kFALSE),fNoiseThreshold(0.01),fNoiseMinimumEvents(100),
63   fMinNrStepsBeforeIncrease(5),fMinIncreaseFromBaseLine(2),fStepDownDacSafe(5),fMaxBaseLineLevel(10)
64 {
65   // constructor
66   for (UInt_t chipNr=0; chipNr<11; chipNr++) {
67     for (UInt_t hs=0; hs<6; hs++) {
68       fMeanMultiplicity[hs][chipNr]=NULL;
69       fHitEventEfficiency[hs][chipNr]=NULL;
70     }
71   }
72   for (UInt_t hs=0; hs<6; hs++) {
73     fTPeffChip[hs]=NULL;
74     fDeadPixelChip[hs]=NULL;
75     fNoisyPixelChip[hs]=NULL;
76   }
77
78   for (UInt_t mod=0; mod<240; mod++) {
79     fbModuleScanned[mod]=kFALSE;
80   }
81
82   Init(readFromGridFile);
83 }
84 //_________________________________________________________________________
85 AliITSOnlineSPDscanAnalyzer::AliITSOnlineSPDscanAnalyzer(const AliITSOnlineSPDscanAnalyzer& handle) :
86   fType(99),fDacId(99),fFileName("."),fScanObj(NULL),fHandler(NULL),fTriggers(NULL),fTPeff(0),fTPeffHS(NULL),fDeadPixel(0),fDeadPixelHS(NULL),fNoisyPixel(0),fNoisyPixelHS(NULL),
87   fOverWrite(kFALSE),fNoiseThreshold(0.01),fNoiseMinimumEvents(100),
88   fMinNrStepsBeforeIncrease(5),fMinIncreaseFromBaseLine(2),fStepDownDacSafe(5),fMaxBaseLineLevel(10)
89 {
90   // copy constructor, only copies the filename and params (not the processed data)
91   fFileName=handle.fFileName;
92   fOverWrite=handle.fOverWrite;
93   fNoiseThreshold=handle.fNoiseThreshold;
94   fNoiseMinimumEvents=handle.fNoiseMinimumEvents;
95   fMinNrStepsBeforeIncrease=handle.fMinNrStepsBeforeIncrease;
96   fMinIncreaseFromBaseLine=handle.fMinIncreaseFromBaseLine;
97   fStepDownDacSafe=handle.fStepDownDacSafe;
98   fMaxBaseLineLevel=handle.fMaxBaseLineLevel;
99
100   for (UInt_t chipNr=0; chipNr<11; chipNr++) {
101     for (UInt_t hs=0; hs<6; hs++) {
102       fMeanMultiplicity[hs][chipNr]=NULL;
103       fHitEventEfficiency[hs][chipNr]=NULL;
104     }
105   }
106   for (UInt_t hs=0; hs<6; hs++) {
107     fTPeffChip[hs]=NULL;
108     fDeadPixelChip[hs]=NULL;
109     fNoisyPixelChip[hs]=NULL;
110   }
111
112   for (UInt_t mod=0; mod<240; mod++) {
113     fbModuleScanned[mod]=kFALSE;
114   }
115
116   Init();
117 }
118 //_________________________________________________________________________
119 AliITSOnlineSPDscanAnalyzer::~AliITSOnlineSPDscanAnalyzer() {
120   // destructor
121   for (UInt_t hs=0; hs<6; hs++) {
122     for (UInt_t chipNr=0; chipNr<11; chipNr++) {
123       if (fMeanMultiplicity[hs][chipNr]!=NULL) {
124         delete fMeanMultiplicity[hs][chipNr];
125         fMeanMultiplicity[hs][chipNr]=NULL;
126       }
127       if (fHitEventEfficiency[hs][chipNr]!=NULL) {
128         delete fHitEventEfficiency[hs][chipNr];
129         fHitEventEfficiency[hs][chipNr]=NULL;
130       }
131     }
132   }
133
134   if (fTriggers!=NULL) {
135     delete fTriggers;
136     fTriggers=NULL;
137   }
138
139   DeleteUniformityHistograms();
140
141   if (fScanObj!=NULL) {
142     delete fScanObj;
143     fScanObj=NULL;
144   }
145 }
146 //_________________________________________________________________________
147 AliITSOnlineSPDscanAnalyzer& AliITSOnlineSPDscanAnalyzer::operator=(const AliITSOnlineSPDscanAnalyzer& handle) {
148   // assignment operator, only copies the filename and params (not the processed data)
149   if (this!=&handle) {
150     for (UInt_t hs=0; hs<6; hs++) {
151       for (UInt_t chipNr=0; chipNr<11; chipNr++) {
152         if (fMeanMultiplicity[hs][chipNr]!=NULL) {
153           delete fMeanMultiplicity[hs][chipNr];
154         }
155         if (fHitEventEfficiency[hs][chipNr]!=NULL) {
156           delete fHitEventEfficiency[hs][chipNr];
157         }
158       }
159     }
160     if (fTriggers!=NULL) {
161       delete fTriggers;
162       fTriggers=NULL;
163     }
164
165     DeleteUniformityHistograms();
166
167     if (fScanObj!=NULL) {
168       delete fScanObj;
169       fScanObj=NULL;
170     }
171    
172     fFileName=handle.fFileName;
173     fOverWrite=handle.fOverWrite;
174     fNoiseThreshold=handle.fNoiseThreshold;
175     fNoiseMinimumEvents=handle.fNoiseMinimumEvents;
176     fMinNrStepsBeforeIncrease=handle.fMinNrStepsBeforeIncrease;
177     fMinIncreaseFromBaseLine=handle.fMinIncreaseFromBaseLine;
178     fStepDownDacSafe=handle.fStepDownDacSafe;
179     fMaxBaseLineLevel=handle.fMaxBaseLineLevel;
180
181     for (UInt_t chipNr=0; chipNr<11; chipNr++) {
182       for (UInt_t hs=0; hs<6; hs++) {
183         fMeanMultiplicity[hs][chipNr]=NULL;
184         fHitEventEfficiency[hs][chipNr]=NULL;
185       }
186     }
187     for (UInt_t mod=0; mod<240; mod++) {
188       fbModuleScanned[mod]=kFALSE;
189     }
190
191     fHandler=NULL;
192     
193     fType=99;
194     fDacId=99;
195
196     Init();    
197   }
198   return *this;
199 }
200 //_________________________________________________________________________
201 void AliITSOnlineSPDscanAnalyzer::Init(Bool_t readFromGridFile) {
202   // first checks type of container and then initializes container obj
203   if (!readFromGridFile) {
204     FILE* fp0 = fopen(fFileName.Data(), "r");
205     if (fp0 == NULL) {
206       return;
207     }
208     else {
209       fclose(fp0);
210     }
211   }
212
213   fScanObj = new AliITSOnlineSPDscan(fFileName.Data(),readFromGridFile);
214   fType = fScanObj->GetType();
215   delete fScanObj;
216
217   // init container
218   switch(fType) {
219   case kUNIMA:
220   case kNOISE:
221     fScanObj = new AliITSOnlineSPDscanSingle(fFileName.Data(),readFromGridFile);
222     break;
223   case kMINTH:
224   case kDAC:
225   case kDELAY:
226     fScanObj = new AliITSOnlineSPDscanMultiple(fFileName.Data(),readFromGridFile);
227     fDacId = ((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacId();
228     break;
229   case kMEANTH:
230     fScanObj = new AliITSOnlineSPDscanMeanTh(fFileName.Data(),readFromGridFile);
231     fDacId = ((AliITSOnlineSPDscanMeanTh*)fScanObj)->GetDacId();
232     break;
233   default:
234     Error("AliITSOnlineSPDscanAnalyzer::Init","Type %d not defined!",fType);
235     fScanObj=NULL;
236     return;
237     break;
238   }
239
240 }
241 //_________________________________________________________________________
242 void AliITSOnlineSPDscanAnalyzer::SetParam(const Char_t *pname, const Char_t *pval) {
243   // set a parameter
244   TString name = pname;
245   TString val = pval;
246   if (name.CompareTo("fOverWrite")==0) {
247     if (val.CompareTo("YES")==0 || val.CompareTo("1")==0) {
248       fOverWrite = kTRUE;
249     }
250     else fOverWrite = kFALSE;
251   }
252   else if (name.CompareTo("fNoiseThreshold")==0) {
253     fNoiseThreshold = val.Atof();
254   }
255   else if (name.CompareTo("fNoiseMinimumEvents")==0) {
256     fNoiseMinimumEvents = val.Atoi();
257   }
258   else if (name.CompareTo("fMinNrStepsBeforeIncrease")==0) {
259     fMinNrStepsBeforeIncrease = val.Atoi();
260   }
261   else if (name.CompareTo("fMinIncreaseFromBaseLine")==0) {
262     fMinIncreaseFromBaseLine = val.Atof();
263   }
264   else if (name.CompareTo("fStepDownDacSafe")==0) {
265     fStepDownDacSafe = val.Atoi();
266   }
267   else if (name.CompareTo("fMaxBaseLineLevel")==0) {
268     fMaxBaseLineLevel = val.Atof();
269   }
270   else {
271     Error("AliITSOnlineSPDscanAnalyzer::SetParam","Parameter %s in configuration file unknown.",name.Data());
272   }
273 }
274 //_________________________________________________________________________
275 void AliITSOnlineSPDscanAnalyzer::ReadParamsFromLocation(const Char_t *dirName) {
276   // opens file (default name) in dir dirName and reads parameters from it
277   TString paramsFileName = Form("%s/standal_params.txt",dirName);
278   ifstream paramsFile;
279   paramsFile.open(paramsFileName, ifstream::in);
280   if (paramsFile.fail()) {
281     printf("No config file (%s) present. Using default tuning parameters.\n",paramsFileName.Data());
282   }
283   else {
284     while(1) {
285       Char_t paramN[50];
286       Char_t paramV[50];
287       paramsFile >> paramN;
288       if (paramsFile.eof()) break;
289       paramsFile >> paramV;
290       SetParam(paramN,paramV);
291       if (paramsFile.eof()) break;
292     }
293     paramsFile.close();
294   }
295 }
296 //_________________________________________________________________________
297 Bool_t AliITSOnlineSPDscanAnalyzer::IsChipPresent(UInt_t hs, UInt_t chipNr) {
298   // is the chip present?
299   if (fScanObj==NULL) {
300     Warning("AliITSOnlineSPDscanAnalyzer::IsChipPresent","No data!");
301     return kFALSE;
302   }
303   return fScanObj->GetChipPresent(hs,chipNr);
304 }
305 //_________________________________________________________________________
306 Bool_t AliITSOnlineSPDscanAnalyzer::ProcessDeadPixels(/*Char_t *oldcalibDir*/) {
307   // process dead pixel data, for uniformity scan, 
308   // NB: This will not be the general way of finding dead pixels.
309   if (fScanObj==NULL) {
310     Warning("AliITSOnlineSPDscanAnalyzer::ProcessDeadPixels","No data!");
311     return kFALSE;
312   }
313   // should be type kUNIMA
314   if (fType!=kUNIMA) {
315     Warning("AliITSOnlineSPDscanAnalyzer::ProcessDeadPixels","Dead pixels only for scan type %d.",kUNIMA);
316     return kFALSE;
317   }
318   // handler should be initialized
319   if (fHandler==NULL) {
320     Error("AliITSOnlineSPDscanAnalyzer::ProcessDeadPixels","Calibration handler is not initialized!");
321     return kFALSE;
322   }
323
324   UInt_t routerNr = fScanObj->GetRouterNr();
325   UInt_t rowStart = fScanObj->GetRowStart();
326   UInt_t rowEnd   = fScanObj->GetRowEnd();
327   for (UInt_t hs=0; hs<6; hs++) {
328     for (UInt_t chipNr=0; chipNr<10; chipNr++) {
329       if (fScanObj->GetChipPresent(hs,chipNr) && fScanObj->GetAverageMultiplicity(0,hs,chipNr)>0) { // check the status of the chippresent parameter in the mood header!!!!!!!!!!!!!!!!!!!!!!!!!!!!
330         if (fOverWrite) {fHandler->ResetDeadForChip(routerNr,hs,chipNr);}
331         for (UInt_t col=0; col<32; col++) {
332           for (UInt_t row=rowStart; row<=rowEnd; row++) {
333             if (col!=1 && col!=9 && col!=17 && col!=25) { //exclude test columns!!!
334               if (fScanObj->GetHits(0,hs,chipNr,col,row)==0) {
335                 fHandler->SetDeadPixel(routerNr,hs,chipNr,col,row);
336               }
337             }
338           }
339         }
340       }
341     }
342   }
343   return kTRUE;
344 }
345 //_________________________________________________________________________
346 Bool_t AliITSOnlineSPDscanAnalyzer::ProcessUniformity() {
347   // process uniformity scan data (thanks to Roberta Ferretti for providing this method)
348   if (fScanObj==NULL) {
349     Warning("AliITSOnlineSPDscanAnalyzer::ProcessUniformity","No data!");
350     return kFALSE;
351   }
352   // should be type kUNIMA
353   if (fType!=kUNIMA) {
354     Warning("AliITSOnlineSPDscanAnalyzer::ProcessUniformity","Only for scan type %d.",kUNIMA);
355     return kFALSE;
356   }
357
358   CreateUniformityHistograms(); // create all histograms that will be filled here
359
360   //  UInt_t routerNr = fScanObj->GetRouterNr();
361   UInt_t rowStart = fScanObj->GetRowStart();
362   UInt_t rowEnd   = fScanObj->GetRowEnd();
363   UInt_t nrTriggers = fScanObj->GetTriggers(0)/(rowEnd-rowStart+1);
364
365   Float_t pixel100=0;
366   Float_t zeri=0;
367   Float_t pixelN=0;
368   UInt_t numChipsActive=0;
369
370   for (UInt_t hs=0; hs<6; hs++) {
371     Float_t pixel100hs=0;
372     Float_t zerihs=0;
373     Float_t pixelNhs=0;
374     UInt_t numChipsActiveHS=0;
375
376     for (UInt_t chipNr=0; chipNr<10; chipNr++) {
377       Float_t pixel100chip=0;
378       Float_t zerichip=0;
379       Float_t pixelNchip=0;
380
381       if (fScanObj->GetChipPresent(hs,chipNr)) { // check the status of the chippresent parameter in the mood header!!!!!!!!!!!!!!!!!!!!!!!!!!!!
382         numChipsActive++;
383         numChipsActiveHS++;
384
385         for (UInt_t col=0; col<32; col++) {
386           for (UInt_t row=rowStart; row<=rowEnd; row++) {
387             if (col!=1 && col!=9 && col!=17 && col!=25) { //exclude test columns!!!
388             
389               if (fScanObj->GetHits(0,hs,chipNr,col,row)==nrTriggers) {   
390                         pixel100++;
391                         pixel100hs++;
392                         pixel100chip++;
393               }
394               if (fScanObj->GetHits(0,hs,chipNr,col,row)==0) {
395                         zeri++;
396                         zerihs++;
397                         zerichip++;
398               }
399               if (fScanObj->GetHits(0,hs,chipNr,col,row)>nrTriggers) {    
400                         pixelN++;
401                         pixelNhs++;
402                         pixelNchip++;
403               }
404             }
405           }
406         }
407       
408         Float_t tPeffChip=(pixel100chip/(28*(rowEnd-rowStart+1)))*100;
409         fTPeffChip[hs]->Fill(chipNr,tPeffChip);
410         
411         Float_t deadPixelChip=(zerichip/(28*(rowEnd-rowStart+1)))*100;
412         fDeadPixelChip[hs]->Fill(chipNr,deadPixelChip);
413         
414         Float_t noisyPixelChip=(pixelNchip/(28*(rowEnd-rowStart+1)))*100;
415         fNoisyPixelChip[hs]->Fill(chipNr,noisyPixelChip);
416       }
417     }
418     
419     Float_t tPeffHS=(pixel100hs/(28*numChipsActiveHS*(rowEnd-rowStart+1)))*100;
420     fTPeffHS->Fill(hs,tPeffHS);
421     
422     Float_t deadPixelHS=(zerihs/(28*numChipsActiveHS*(rowEnd-rowStart+1)))*100;
423     fDeadPixelHS->Fill(hs,deadPixelHS);
424     
425     Float_t noisyPixelHS=(pixelNhs/(28*numChipsActiveHS*(rowEnd-rowStart+1)))*100;
426     fNoisyPixelHS->Fill(hs,noisyPixelHS);
427   }
428   
429   fTPeff=(pixel100/(28*numChipsActive*(rowEnd-rowStart+1)))*100;
430   fDeadPixel=(zeri/(28*numChipsActive*(rowEnd-rowStart+1)))*100;
431   fNoisyPixel=(pixelN/(28*numChipsActive*(rowEnd-rowStart+1)))*100;
432
433   return kTRUE;
434 }
435 //_________________________________________________________________________
436 void AliITSOnlineSPDscanAnalyzer::CreateUniformityHistograms() {
437   // create uniformity histograms to be filled by "ProcessUniformity" method
438   DeleteUniformityHistograms(); // make sure no old histograms are lying around...
439   UInt_t eq = GetRouterNr();
440   TString label;
441
442   label = Form("Ratio of 'Good' Pixels Per HS (eq %d)",eq);
443   fTPeffHS = new TH1F(label.Data(),label.Data(),6,-0.5,5.5);
444   fTPeffHS->SetXTitle("hs");
445   fTPeffHS->SetYTitle("ratio [%]");
446   fTPeffHS->SetFillColor(kBlue);
447   fTPeffHS->SetStats(0);
448
449   label = Form("Ratio of 'Dead' Pixels Per HS (eq %d)",eq);
450   fDeadPixelHS = new TH1F(label.Data(),label.Data(),6,-0.5,5.5);
451   fDeadPixelHS->SetXTitle("hs");
452   fDeadPixelHS->SetYTitle("ratio [%]");
453   fDeadPixelHS->SetFillColor(kBlue);
454   fDeadPixelHS->SetStats(0);
455
456   label = Form("Ratio of 'Noisy' Pixels Per HS (eq %d)",eq);
457   fNoisyPixelHS = new TH1F(label.Data(),label.Data(),6,-0.5,5.5);
458   fNoisyPixelHS->SetXTitle("hs");
459   fNoisyPixelHS->SetYTitle("ratio [%]");
460   fNoisyPixelHS->SetFillColor(kBlue);
461   fNoisyPixelHS->SetStats(0);
462
463   for (UInt_t hs=0; hs<6; hs++) {
464     label = Form("Ratio of 'Good' Pixels Per Chip (eq %d, hs %d)",eq,hs);
465     fTPeffChip[hs] = new TH1F(label.Data(),label.Data(),10,-0.5,9.5);
466     fTPeffChip[hs]->SetXTitle("chip");
467     fTPeffChip[hs]->SetYTitle("ratio [%]");
468     fTPeffChip[hs]->SetFillColor(kBlue);
469     fTPeffChip[hs]->SetStats(0);
470
471     label = Form("Ratio of 'Dead' Pixels Per Chip (eq %d, hs %d)",eq,hs);
472     fDeadPixelChip[hs] = new TH1F(label.Data(),label.Data(),10,-0.5,9.5);
473     fDeadPixelChip[hs]->SetXTitle("chip");
474     fDeadPixelChip[hs]->SetYTitle("ratio [%]");
475     fDeadPixelChip[hs]->SetFillColor(kBlue);
476     fDeadPixelChip[hs]->SetStats(0);
477
478     label = Form("Ratio of 'Noisy' Pixels Per Chip (eq %d, hs %d)",eq,hs);
479     fNoisyPixelChip[hs] = new TH1F(label.Data(),label.Data(),10,-0.5,9.5);
480     fNoisyPixelChip[hs]->SetXTitle("chip");
481     fNoisyPixelChip[hs]->SetYTitle("ratio [%]");
482     fNoisyPixelChip[hs]->SetFillColor(kBlue);
483     fNoisyPixelChip[hs]->SetStats(0);
484   }
485
486 }
487 //_________________________________________________________________________
488 void AliITSOnlineSPDscanAnalyzer::DeleteUniformityHistograms() {
489   // remove uniformity histograms if they are created
490   if (fTPeffHS!=NULL) {
491     delete fTPeffHS;
492     fTPeffHS=NULL;
493   }
494   if (fDeadPixelHS!=NULL) {
495     delete fDeadPixelHS;
496     fDeadPixelHS=NULL;
497   }
498   if (fNoisyPixelHS!=NULL) {
499     delete fNoisyPixelHS;
500     fNoisyPixelHS=NULL;
501   }
502   for (UInt_t hs=0; hs<6; hs++) {
503     if (fTPeffChip[hs]!=NULL) {
504       delete fTPeffChip[hs];
505       fTPeffChip[hs]=NULL;
506     }
507     if (fDeadPixelChip[hs]!=NULL) {
508       delete fDeadPixelChip[hs];
509       fDeadPixelChip[hs]=NULL;
510     }
511     if (fNoisyPixelChip[hs]!=NULL) {
512       delete fNoisyPixelChip[hs];
513       fNoisyPixelChip[hs]=NULL;
514     }
515   }
516 }
517 //_________________________________________________________________________
518 Bool_t AliITSOnlineSPDscanAnalyzer::ProcessNoisyPixels(/*Char_t *oldcalibDir*/) {
519   // process noisy pixel data
520   if (fScanObj==NULL) {
521     Warning("AliITSOnlineSPDscanAnalyzer::ProcessNoisyPixels","No data!");
522     return kFALSE;
523   }
524   // should be type kNOISE
525   if (fType != kNOISE) {
526     Warning("AliITSOnlineSPDscanAnalyzer::ProcessNoisyPixels","Noisy pixels only for scan type %d.",kNOISE);
527     return kFALSE;
528   }
529   // handler should be initialized
530   if (fHandler==NULL) {
531     Error("AliITSOnlineSPDscanAnalyzer::ProcessNoisyPixels","Calibration handler is not initialized!");
532     return kFALSE;
533   }
534   // check if enough statistics
535   if (fScanObj->GetTriggers(0)<fNoiseMinimumEvents) {
536     Warning("AliITSOnlineSPDscanAnalyzer::ProcessNoisyPixels","Process noisy: Too few events.");
537     return kFALSE;
538   }
539
540   UInt_t routerNr = fScanObj->GetRouterNr();
541   for (UInt_t hs=0; hs<6; hs++) {
542     for (UInt_t chipNr=0; chipNr<10; chipNr++) {
543       if (fScanObj->GetChipPresent(hs,chipNr)) { // check the status of the chippresent parameter in the mood header!!!!!!!!!!!!!!!!!!!!!!!!!!!!
544         if (fOverWrite) {fHandler->ResetNoisyForChip(routerNr,hs,chipNr);}
545         for (UInt_t col=0; col<32; col++) {
546           for (UInt_t row=0; row<256; row++) {
547             if (fScanObj->GetHitsEfficiency(0,hs,chipNr,col,row)>fNoiseThreshold) {
548               fHandler->SetNoisyPixel(routerNr,hs,chipNr,col,row);
549             }
550           }
551         }
552       }
553     }
554   }
555   return kTRUE;
556 }
557 //_________________________________________________________________________
558 Int_t AliITSOnlineSPDscanAnalyzer::GetDelay(UInt_t hs, UInt_t chipNr) {
559   // get delay
560   if (hs>=6 || chipNr>10) return -1;
561   if (fScanObj==NULL) {
562     Warning("AliITSOnlineSPDscanAnalyzer::GetDelay","No data!");
563     return -1;
564   }
565   // should be type kDELAY or kDAC with id 42 (delay_ctrl)
566   if (fType!=kDELAY && (fType!=kDAC || fDacId!=42)) {
567     Warning("AliITSOnlineSPDscanAnalyzer::GetDelay","Delay only for scan type %d or %d and dac_id 42.",kDELAY,kDAC);
568     return -1;
569   }
570   if (fMeanMultiplicity[hs][chipNr]==NULL) {
571     if (!ProcessMeanMultiplicity()) {
572       return -1;
573     }
574   }
575
576   UInt_t maxStep=0;
577   Float_t maxVal=0;
578   for (UInt_t step=0; step<fScanObj->GetNSteps(); step++) {
579     Double_t thisDac;
580     Double_t thisMult;
581     fMeanMultiplicity[hs][chipNr]->GetPoint(step,thisDac,thisMult);
582     if (thisMult > maxVal) {
583       maxVal = thisMult;
584       maxStep = step;
585     }
586   }
587
588   if (maxVal>0) {
589     return ((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(maxStep);
590   }
591   else {
592     return -1;
593   }
594
595 }
596 //_________________________________________________________________________
597 Int_t AliITSOnlineSPDscanAnalyzer::GetNrNoisyUnima(UInt_t hs, UInt_t chipNr) {
598   // in case of a uniformity scan, returns the nr of noisy pixels, (here > 200 hits)
599   if (hs>=6 || chipNr>10) return -1;
600   if (fScanObj==NULL) {
601     Error("AliITSOnlineSPDscanAnalyzer::GetNrNoisyUnima","No data!");
602     return kFALSE;
603   }
604   // should be type kUNIMA
605   if (fType != kUNIMA) {
606     Error("AliITSOnlineSPDscanAnalyzer::GetNrNoisyUnima","Noisy pixels Unima only for scan type %d.",kUNIMA);
607     return kFALSE;
608   }
609   if (fScanObj->GetTriggers(0)!=25600) {
610     Error("AliITSOnlineSPDscanAnalyzer::GetNrNoisyUnima","Process noisy unima: Incorrect number of events (!=25600.");
611     return kFALSE;
612   }
613
614   Int_t nrNoisy=0;
615   if (fScanObj->GetChipPresent(hs,chipNr)) { // check the status of the chippresent parameter in the mood header!!!!!!!!!!!!!!!!!!!!!!!!!!!!
616     for (UInt_t col=0; col<32; col++) {
617       for (UInt_t row=0; row<256; row++) {
618         if (fScanObj->GetHits(0,hs,chipNr,col,row)>200) {
619           nrNoisy++;
620         }
621       }
622     }
623   }
624   else {
625     return -1;
626   }
627   return nrNoisy;
628 }
629 //_________________________________________________________________________
630 Int_t AliITSOnlineSPDscanAnalyzer::FindLastMinThDac(UInt_t hs, UInt_t chipNr) {
631   // returns dac value where fMinIncreaseFromBaseLine reached
632   if (hs>=6 || chipNr>10) return -1;
633   if (fMeanMultiplicity[hs][chipNr]==NULL) {
634     if (!ProcessMeanMultiplicity()) {
635       return -1;
636     }
637   }
638   Double_t firstVal, dummy1;
639   fMeanMultiplicity[hs][chipNr]->GetPoint(0,dummy1,firstVal);
640   UInt_t step=0;
641   while (step<fScanObj->GetNSteps()-1) {
642     Double_t graphVal, dummy2;
643     fMeanMultiplicity[hs][chipNr]->GetPoint(step+1,dummy2,graphVal);
644     if (graphVal>firstVal+fMinIncreaseFromBaseLine) break;
645     step++;
646   }
647   if (step==fScanObj->GetNSteps()-1) return -1;
648   return ((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(step);
649 }
650
651 Int_t AliITSOnlineSPDscanAnalyzer::FindClosestLowerStep(Float_t dacValueInput) {
652   // returns step closest (lower) to a dacvalue 
653   UInt_t step=0;
654   while (step<fScanObj->GetNSteps()-1) {
655     Int_t dacVal = ((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(step+1);
656     if (dacVal>=dacValueInput) break;
657     step++;
658   }
659   return step;
660 }
661 //_________________________________________________________________________
662 Float_t AliITSOnlineSPDscanAnalyzer::GetCompareLine(UInt_t step, UInt_t hs, UInt_t chipNr, Float_t basePar2) {
663   // returns value to compare mean mult with (when finding min th)
664   if (hs>=6 || chipNr>10) return -1;
665   if (step<fMinNrStepsBeforeIncrease) return -1;
666   Float_t baseLine = basePar2;
667   if (baseLine<0) baseLine=0;
668   Float_t baseAdd;
669   Double_t baseM=0;
670   Double_t baseS=0;
671   Double_t d,m;
672   for (UInt_t st=1;st<2*step/3;st++) { // skip first point...
673     fMeanMultiplicity[hs][chipNr]->GetPoint(st,d,m);
674     baseM+=m-baseLine;
675     baseS+=(m-baseLine)*(m-baseLine);
676   }
677   baseAdd=2*sqrt( baseS/(2*step/3-1) - (baseM/(2*step/3-1))*(baseM/(2*step/3-1)) );
678   baseAdd+=0.03; // magic number
679   if (baseAdd>fMinIncreaseFromBaseLine) baseAdd=fMinIncreaseFromBaseLine;
680   return baseLine + baseAdd;
681 }
682
683 Int_t AliITSOnlineSPDscanAnalyzer::GetMinTh(UInt_t hs, UInt_t chipNr) {
684   // calculates and returns the minimum threshold
685   if (hs>=6 || chipNr>10) return -1;
686   if (fScanObj==NULL) {
687     Error("AliITSOnlineSPDscanAnalyzer::GetMinTh","No data!");
688     return -1;
689   }
690   // should be type  kMINTH  or  kDAC with id 39 (pre_vth)
691   if (fType!=kMINTH && (fType!=kDAC || fDacId!=39)) {
692     Error("AliITSOnlineSPDscanAnalyzer::GetMinTh","MinTh only for scan type %d OR %d with dac_id 39.",kMINTH,kDAC);
693     return -1;
694   }
695   if (fMeanMultiplicity[hs][chipNr]==NULL) {
696     if (!ProcessMeanMultiplicity()) {
697       return -1;
698     }
699   }
700
701   Int_t lastDac = FindLastMinThDac(hs,chipNr);
702   if (lastDac==-1) {
703     Warning("AliITSOnlineSPDscanAnalyzer::GetMinTh","HS%d, Chip%d: Increase of Mean Multiplicity by %1.2f never reached.",hs,chipNr,fMinIncreaseFromBaseLine);
704     return -1;
705   }
706
707   Int_t minDac = ((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(0);
708   TString funcName = Form("Fit minth func HS%d CHIP%d",hs,chipNr);
709   TF1 *minThFunc = new TF1(funcName.Data(),itsSpdErrorf,100,500,3);
710   minThFunc->SetParameter(0,lastDac+10);
711   minThFunc->SetParameter(1,2);
712   minThFunc->SetParameter(2,0);
713   minThFunc->SetParName(0,"Mean");
714   minThFunc->SetParName(1,"Sigma");
715   minThFunc->SetParName(2,"BaseLine");
716   minThFunc->SetLineWidth(1);
717   if (fMeanMultiplicity[hs][chipNr]==NULL) {
718     if (!ProcessMeanMultiplicity()) {
719       return -1;
720     }
721   }
722   fMeanMultiplicity[hs][chipNr]->Fit(funcName,"Q0","",minDac,lastDac);
723
724   //  Double_t mean = fMinThFunc[hs][chipNr]->GetParameter(0);
725   //  Double_t sigma = fMinThFunc[hs][chipNr]->GetParameter(1);
726   Double_t baseLine = minThFunc->GetParameter(2);
727   delete minThFunc;
728
729   if (baseLine>fMaxBaseLineLevel) {
730     Warning("AliITSOnlineSPDscanAnalyzer::GetMinTh","HS%d, Chip%d: BaseLine too large (%1.2f>%1.2f).",hs,chipNr,baseLine,fMaxBaseLineLevel);
731     return -1;
732   }
733   UInt_t step=FindClosestLowerStep(lastDac);
734   Float_t compareLine = GetCompareLine(step,hs,chipNr,baseLine);
735   if (compareLine==-1) {
736     Warning("AliITSOnlineSPDscanAnalyzer::GetMinTh","HS%d, Chip%d: Not enough steps (%d<%d) before increase to get a compare line.",hs,chipNr,step,fMinNrStepsBeforeIncrease);
737     return -1;
738   }
739
740   Double_t mult, dummy;
741   mult=1000;
742   while (mult > compareLine && step>0) {
743     fMeanMultiplicity[hs][chipNr]->GetPoint(step,dummy,mult);
744     step--;
745   }
746   Int_t minth = ((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(step+1)-fStepDownDacSafe;
747
748   if (step>0) {
749     return minth;
750   }
751   else {
752     Warning("AliITSOnlineSPDscanAnalyzer::GetMinTh","HS%d, Chip%d: Did not find a point below the compare line (%f).",hs,chipNr,compareLine);
753     return -1;
754   }
755 }
756 //_________________________________________________________________________
757 Int_t AliITSOnlineSPDscanAnalyzer::GetMeanTh(UInt_t hs, UInt_t chipNr) {
758   // calculates and returns the mean threshold
759   if (hs>=6 || chipNr>10) return -1;
760   if (fScanObj==NULL) {
761     Error("AliITSOnlineSPDscanAnalyzer::GetMeanTh","No data!");
762     return -1;
763   }
764   // should be type  kMEANTH  or  kDAC with id 39
765   if (fType!=kMEANTH && (fType!=kDAC || fDacId!=105)) {
766     Error("AliITSOnlineSPDscanAnalyzer::GetMeanTh","MeanTh only for scan type %d OR %d with dac_id 105.",kMEANTH,kDAC);
767     return -1;
768   }
769   if (fHitEventEfficiency[hs][chipNr]==NULL) {
770    printf("processing hit efficiency \n");
771     if (!ProcessHitEventEfficiency()) {
772       printf("...not processed!!\n");
773       return -1;
774     }
775   }
776   Double_t x,y;
777   fHitEventEfficiency[hs][chipNr]->GetPoint(fHitEventEfficiency[hs][chipNr]->GetN()-1,x,y);
778   UInt_t min = (UInt_t)x;
779   fHitEventEfficiency[hs][chipNr]->GetPoint(0,x,y);
780   UInt_t max = (UInt_t)x;
781
782   TString funcName = Form("Fit meanth func HS%d CHIP%d",hs,chipNr);
783   TF1 *minThFunc = new TF1(funcName.Data(),itsSpdScurveForMeanTh,min,max,3);
784   minThFunc->SetParameter(0,3000);
785   minThFunc->SetParameter(1,100);
786   minThFunc->SetParameter(2,0.4);
787   minThFunc->SetParName(0,"Mean");
788   minThFunc->SetParName(1,"Sigma");
789   minThFunc->SetParName(2,"Half");
790   minThFunc->SetLineWidth(1);
791  fHitEventEfficiency[hs][chipNr]->Fit(funcName,"Q","",min,max);
792
793   Double_t value = minThFunc->GetParameter(0);
794   TLine *ly = new TLine(min,0.5,value,0.5); ly->SetLineStyle(6);
795   ly->Draw("same");
796   TLine *lx = new TLine(value,0.,value,0.5);
797   lx->SetLineStyle(6);
798   lx->Draw("same");
799   delete minThFunc;
800
801   return (UInt_t)value;
802 }
803
804 //_________________________________________________________________________
805 Bool_t AliITSOnlineSPDscanAnalyzer::ProcessMeanMultiplicity() {
806   // process mean multiplicity data
807   if (fScanObj==NULL) {
808     Error("AliITSOnlineSPDscanAnalyzer::ProcessMeanMultiplicity","No data!");
809     return kFALSE;
810   }
811   for (UInt_t step=0; step<fScanObj->GetNSteps(); step++) {
812     for (UInt_t hs=0; hs<6; hs++) {
813       for (UInt_t chipNr=0; chipNr<11; chipNr++) {
814         //      if (fScanObj->GetChipPresent(hs,chipNr)) { // check the status of the chippresent parameter in the mood header!!!!!!!!!!!!!!!!!!!!!!!!!!!!
815         if (step==0) {
816           if (fMeanMultiplicity[hs][chipNr]!=NULL) {
817             delete fMeanMultiplicity[hs][chipNr];
818           }
819           fMeanMultiplicity[hs][chipNr] = new TGraph();
820         }
821         Float_t multiplMean=fScanObj->GetAverageMultiplicity(step,hs,chipNr);
822         if (fType==kMINTH || fType==kMEANTH || fType==kDAC || fType==kDELAY) {
823           fMeanMultiplicity[hs][chipNr]->SetPoint(step,((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(step),multiplMean);
824         }
825         else {
826           fMeanMultiplicity[hs][chipNr]->SetPoint(step,0,multiplMean);
827         }
828       }
829       //      }
830     }
831   }
832   return kTRUE;
833 }
834 //_________________________________________________________________________
835 TGraph* AliITSOnlineSPDscanAnalyzer::GetMeanMultiplicityG(UInt_t hs, UInt_t chipNr) {
836   // returns mean multiplicity graph
837   if (hs>=6 || chipNr>10) return NULL;
838   if (fMeanMultiplicity[hs][chipNr]==NULL) {
839     if (!ProcessMeanMultiplicity()) {
840       return NULL;
841     }
842   }
843   return fMeanMultiplicity[hs][chipNr];
844 }
845 //_________________________________________________________________________
846 Bool_t AliITSOnlineSPDscanAnalyzer::ProcessHitEventEfficiency() {
847   // process hit event efficiency
848   if (fScanObj==NULL) {
849     Error("AliITSOnlineSPDscanAnalyzer::ProcessHitEventEfficiency","No data!");
850     return kFALSE;
851   }
852   for (UInt_t step=0; step<fScanObj->GetNSteps(); step++) {
853     for (UInt_t hs=0; hs<6; hs++) {
854       for (UInt_t chipNr=0; chipNr<11; chipNr++) {
855         //      if (fScanObj->GetChipPresent(hs,chipNr)) { // check the status of the chippresent parameter in the mood header!!!!!!!!!!!!!!!!!!!!!!!!!!!!
856         if (step==0) {
857           if (fHitEventEfficiency[hs][chipNr]!=NULL) {
858             delete fHitEventEfficiency[hs][chipNr];
859           }
860           fHitEventEfficiency[hs][chipNr] = new TGraph();
861         }
862         Float_t efficiency=fScanObj->GetHitEventsEfficiency(step,hs,chipNr);
863         if (fType==kMINTH || fType==kMEANTH || fType==kDAC || fType==kDELAY) {
864           fHitEventEfficiency[hs][chipNr]->SetPoint(step,((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(step),efficiency);
865         }
866         else {
867           fHitEventEfficiency[hs][chipNr]->SetPoint(step,0,efficiency);
868         }
869       }
870       //      }
871     }
872   }
873   return kTRUE;
874 }
875 //_________________________________________________________________________
876 TGraph* AliITSOnlineSPDscanAnalyzer::GetHitEventEfficiencyG(UInt_t hs, UInt_t chipNr) {
877   // returns hit event efficiency graph
878   if (hs>=6 || chipNr>10) return NULL;
879   if (fHitEventEfficiency[hs][chipNr]==NULL) {
880     if (!ProcessHitEventEfficiency()) {
881       return NULL;
882     }
883   }
884   return fHitEventEfficiency[hs][chipNr];
885 }
886 //_________________________________________________________________________
887 Bool_t AliITSOnlineSPDscanAnalyzer::ProcessNrTriggers() {
888   // process nr of triggers data
889   if (fScanObj==NULL) {
890     Error("AliITSOnlineSPDscanAnalyzer::ProcessNrTriggers","No data!");
891     return kFALSE;
892   }
893   for (UInt_t step=0; step<fScanObj->GetNSteps(); step++) {
894     if (step==0) {
895       if (fTriggers!=NULL) {
896         delete fTriggers;
897       }
898       fTriggers = new TGraph();
899     }
900     if (fType==kMINTH || fType==kMEANTH || fType==kDAC || fType==kDELAY) {
901       fTriggers->SetPoint(step,((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(step),fScanObj->GetTriggers(step));
902     }
903     else {
904       fTriggers->SetPoint(step,0,fScanObj->GetTriggers(step));
905     }
906   }
907   return kTRUE;
908 }
909 //_________________________________________________________________________
910 TGraph* AliITSOnlineSPDscanAnalyzer::GetNrTriggersG() {
911   // returns nr of triggers graph
912   if (fTriggers==NULL) {
913     if (!ProcessNrTriggers()) {
914       return NULL;
915     }
916   }
917   return fTriggers;
918 }
919 //_________________________________________________________________________
920 Bool_t AliITSOnlineSPDscanAnalyzer::GetHalfStavePresent(UInt_t hs) {
921   // returns half stave present info
922   if (hs<6 && fScanObj!=NULL) {
923     Int_t chipstatus=0;
924     for (Int_t chip=0; chip<10; chip++) {
925       chipstatus+=fScanObj->GetChipPresent(hs,chip);
926     }
927     if (chipstatus>0) return kTRUE;
928   }
929   return kFALSE;
930 }
931 //_________________________________________________________________________
932 UInt_t AliITSOnlineSPDscanAnalyzer::GetRouterNr() {
933   // returns the router nr of scan obj
934   if (fScanObj!=NULL) return fScanObj->GetRouterNr(); 
935   else return 99;
936 }
937 //_________________________________________________________________________
938 TH2F* AliITSOnlineSPDscanAnalyzer::GetHitMapTot(UInt_t step) {
939   // creates and returns a pointer to a hitmap histo (half sector style a la spdmood)
940   if (fScanObj==NULL) {
941     Error("AliITSOnlineSPDscanAnalyzer::GetHitMapTot","No data!");
942     return NULL;
943   }
944   TString histoname;
945   if (fType==kMINTH || fType==kMEANTH || fType==kDAC) {
946     histoname = Form("Router %d , DAC %d",GetRouterNr(),((AliITSOnlineSPDscanMultiple*)fScanObj)->GetDacValue(step));
947   }
948   else {
949     histoname = Form("Router %d ",GetRouterNr());
950   }
951   TH2F* fHitMapTot = new TH2F(histoname.Data(),histoname.Data(),32*10,-0.5,32*10-0.5,256*6,-0.5,256*6-0.5);
952   fHitMapTot->SetNdivisions(-10,"X");
953   fHitMapTot->SetNdivisions(-006,"Y");
954   fHitMapTot->SetTickLength(0,"X");
955   fHitMapTot->SetTickLength(0,"Y");
956   fHitMapTot->GetXaxis()->SetLabelColor(gStyle->GetCanvasColor());
957   fHitMapTot->GetYaxis()->SetLabelColor(gStyle->GetCanvasColor());
958   for (UInt_t hs=0; hs<6; hs++) {
959     for (UInt_t chipNr=0; chipNr<10; chipNr++) {
960       for (UInt_t col=0; col<32; col++) {
961         for (UInt_t row=0; row<256; row++) {
962           if(hs<2) fHitMapTot->Fill(chipNr*32+(31-col),(5-hs)*256+row,fScanObj->GetHits(step,hs,chipNr,col,row));
963           else fHitMapTot->Fill(chipNr*32+col,(5-hs)*256+row,fScanObj->GetHits(step,hs,chipNr,col,row));
964         }
965       }
966     }
967   }
968   return fHitMapTot;
969 }
970 //_________________________________________________________________________
971 TH2F* AliITSOnlineSPDscanAnalyzer::GetHitMapChip(UInt_t step, UInt_t hs, UInt_t chip) {
972   // creates and returns a pointer to a hitmap histo (chip style a la spdmood)
973   if (fScanObj==NULL) {
974     Error("AliITSOnlineSPDscanAnalyzer::GetHitMapChip","No data!");
975     return NULL;
976   }
977
978   TString histoName;
979   TString histoTitle;
980   histoName = Form("fChipHisto_%d_%d_%d", GetRouterNr(), hs, chip);
981   histoTitle = Form("Eq ID %d, Half Stave %d, Chip %d", GetRouterNr(), hs, chip);
982
983   TH2F *returnHisto = new TH2F(histoName.Data(), histoTitle.Data(), 32, -0.5, 31.5, 256, -0.5, 255.5);
984   returnHisto->SetMinimum(0);
985   for (UInt_t col=0; col<32; col++) {
986     for (UInt_t row=0; row<256; row++) {
987       if(hs<2) returnHisto->Fill(31-col,row,fScanObj->GetHits(step,hs,chip,col,row));
988       else returnHisto->Fill(col,row,fScanObj->GetHits(step,hs,chip,col,row));
989     }
990   }
991
992   return returnHisto;
993 }
994