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