]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/Base/AliTPCCalROC.cxx
ATO-119 - solving merge conflict
[u/mrichter/AliRoot.git] / TPC / Base / AliTPCCalROC.cxx
CommitLineData
07627591 1/**************************************************************************
2 * Copyright(c) 1998-1999, 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
17///////////////////////////////////////////////////////////////////////////////
18// //
a6d2bd0c 19// Calibration base class for a single ROC //
20// Contains one float value per pad //
c5bbaa2c 21// mapping of the pads taken form AliTPCROC //
07627591 22// //
23///////////////////////////////////////////////////////////////////////////////
24
72fbbc82 25//
26// ROOT includes
27//
07627591 28#include "TMath.h"
c5bbaa2c 29#include "TClass.h"
30#include "TFile.h"
184bcc16 31#include "TH1F.h"
2e9bedc9 32#include "TH2F.h"
72fbbc82 33#include "TArrayI.h"
9cc76bba 34
72fbbc82 35//
36//
37#include "AliTPCCalROC.h"
184bcc16 38#include "AliMathBase.h"
e0bc0200 39
ca5dca67 40#include "TRandom3.h" // only needed by the AliTPCCalROCTest() method
41
07627591 42ClassImp(AliTPCCalROC)
07627591 43
44
45//_____________________________________________________________________________
179c6296 46AliTPCCalROC::AliTPCCalROC()
f691a5e2 47 :TNamed(),
179c6296 48 fSector(0),
49 fNChannels(0),
50 fNRows(0),
a5ade541 51 fkIndexes(0),
179c6296 52 fData(0)
07627591 53{
54 //
55 // Default constructor
56 //
179c6296 57
07627591 58}
59
60//_____________________________________________________________________________
179c6296 61AliTPCCalROC::AliTPCCalROC(UInt_t sector)
f691a5e2 62 :TNamed(),
179c6296 63 fSector(0),
64 fNChannels(0),
65 fNRows(0),
a5ade541 66 fkIndexes(0),
179c6296 67 fData(0)
07627591 68{
69 //
70 // Constructor that initializes a given sector
71 //
07627591 72 fSector = sector;
c5bbaa2c 73 fNChannels = AliTPCROC::Instance()->GetNChannels(fSector);
74 fNRows = AliTPCROC::Instance()->GetNRows(fSector);
a5ade541 75 fkIndexes = AliTPCROC::Instance()->GetRowIndexes(fSector);
c5bbaa2c 76 fData = new Float_t[fNChannels];
77 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata] = 0.;
07627591 78}
79
80//_____________________________________________________________________________
179c6296 81AliTPCCalROC::AliTPCCalROC(const AliTPCCalROC &c)
f691a5e2 82 :TNamed(c),
179c6296 83 fSector(0),
84 fNChannels(0),
85 fNRows(0),
a5ade541 86 fkIndexes(0),
179c6296 87 fData(0)
07627591 88{
89 //
90 // AliTPCCalROC copy constructor
91 //
92 fSector = c.fSector;
c5bbaa2c 93 fNChannels = AliTPCROC::Instance()->GetNChannels(fSector);
94 fNRows = AliTPCROC::Instance()->GetNRows(fSector);
a5ade541 95 fkIndexes = AliTPCROC::Instance()->GetRowIndexes(fSector);
c5bbaa2c 96 //
97 fData = new Float_t[fNChannels];
98 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata] = c.fData[idata];
07627591 99}
179c6296 100//____________________________________________________________________________
101AliTPCCalROC & AliTPCCalROC::operator =(const AliTPCCalROC & param)
102{
103 //
104 // assignment operator - dummy
105 //
04420071 106 if (this == &param) return (*this);
77ba0580 107 fSector = param.fSector;
108 fNChannels = AliTPCROC::Instance()->GetNChannels(fSector);
109 fNRows = AliTPCROC::Instance()->GetNRows(fSector);
110 fkIndexes = AliTPCROC::Instance()->GetRowIndexes(fSector);
111 //
df0035f8 112 if (fData) delete [] fData;
77ba0580 113 fData = new Float_t[fNChannels];
114 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata] = param.fData[idata];
179c6296 115 return (*this);
116}
117
07627591 118
119//_____________________________________________________________________________
120AliTPCCalROC::~AliTPCCalROC()
121{
122 //
123 // AliTPCCalROC destructor
124 //
07627591 125 if (fData) {
126 delete [] fData;
127 fData = 0;
128 }
129}
130
c5bbaa2c 131
132
133void AliTPCCalROC::Streamer(TBuffer &R__b)
134{
a6d2bd0c 135 //
c5bbaa2c 136 // Stream an object of class AliTPCCalROC.
a6d2bd0c 137 //
c5bbaa2c 138 if (R__b.IsReading()) {
139 AliTPCCalROC::Class()->ReadBuffer(R__b, this);
a5ade541 140 fkIndexes = AliTPCROC::Instance()->GetRowIndexes(fSector);
c5bbaa2c 141 } else {
142 AliTPCCalROC::Class()->WriteBuffer(R__b,this);
143 }
144}
145
ff42e0c6 146
147//
148
b1ebb72e 149Bool_t AliTPCCalROC::MedianFilter(Int_t deltaRow, Int_t deltaPad, AliTPCCalROC* outlierROC, Bool_t doEdge){
150 //){
4c066a96 151 //
ff42e0c6 152 // Modify content of the object - raplace value by median in neighorhood
153 //
4c066a96 154 Float_t *newBuffer=new Float_t[fNChannels] ;
ff42e0c6 155 Double_t *cacheBuffer=new Double_t[fNChannels];
4c066a96 156 //
ff42e0c6 157 for (Int_t iRow=0; iRow< Int_t(fNRows); iRow++){
4c066a96 158 Int_t nPads=GetNPads(iRow); // number of rows in current row
159 for (Int_t iPad=0; iPad<nPads; iPad++){
ff42e0c6 160 Double_t value=GetValue(iRow,iPad);
4c066a96 161 Int_t counter=0;
ff42e0c6 162 //
4c066a96 163 for (Int_t dRow=-deltaRow; dRow<=deltaRow; dRow++){
ff42e0c6 164 Int_t jRow=iRow+dRow; //take the row - mirror if ouside of range
165 Float_t sign0=1.;
166 if (jRow<0) sign0=-1.;
167 if (UInt_t(jRow)>=fNRows) sign0=-1.;
168 Int_t jPads= GetNPads(iRow+sign0*dRow);
4c066a96 169 Int_t offset=(nPads-jPads)/2.;
ff42e0c6 170 //
4c066a96 171 for (Int_t dPad=-deltaPad; dPad<=deltaPad; dPad++){
ff42e0c6 172 Float_t sign=sign0;
173 Int_t jPad=iPad+dPad+offset; //take the pad - mirror if ouside of range
174 Int_t kRow=jRow;
175 if (jPad<0) sign=-1;
176 if (jPad>=jPads) sign=-1;
177 if (sign<0){
178 kRow=iRow-dRow;
179 jPad=iPad-dPad+offset;
b1ebb72e 180 if (!doEdge) continue;
ff42e0c6 181 }
182 if (IsInRange(UInt_t(kRow),jPad)){
b1ebb72e 183 Bool_t isOutlier=(outlierROC==NULL)?kFALSE:outlierROC->GetValue(kRow,jPad)>0;
184 if (!isOutlier){
185 cacheBuffer[counter]=sign*(GetValue(kRow,jPad)-value);
186 counter++;
187 }
ff42e0c6 188 }
4c066a96 189 }
190 }
b1ebb72e 191 newBuffer[fkIndexes[iRow]+iPad]=0.;
192 if (counter>1) newBuffer[fkIndexes[iRow]+iPad] = TMath::Median(counter, cacheBuffer)+value;
4c066a96 193 }
194 }
195 memcpy(fData, newBuffer,GetNchannels()*sizeof(Float_t));
196 delete []newBuffer;
197 delete []cacheBuffer;
198 return kTRUE;
199}
200
ff42e0c6 201
202
b1ebb72e 203Bool_t AliTPCCalROC::LTMFilter(Int_t deltaRow, Int_t deltaPad, Float_t fraction, Int_t type, AliTPCCalROC* outlierROC, Bool_t doEdge){
204 //){
4c066a96 205 //
206 //
207 // //
208 // Modify content of the class
ff42e0c6 209 // write LTM mean or median
4c066a96 210 if (fraction<0 || fraction>1) return kFALSE;
211 Float_t *newBuffer=new Float_t[fNChannels] ;
212 Double_t *cacheBuffer=new Double_t[fNChannels];
213 //
ff42e0c6 214 for (Int_t iRow=0; iRow< Int_t(fNRows); iRow++){
4c066a96 215 Int_t nPads=GetNPads(iRow); // number of rows in current row
216 for (Int_t iPad=0; iPad<nPads; iPad++){
ff42e0c6 217 Double_t value=GetValue(iRow,iPad);
4c066a96 218 Int_t counter=0;
ff42e0c6 219 //
4c066a96 220 for (Int_t dRow=-deltaRow; dRow<=deltaRow; dRow++){
ff42e0c6 221 Int_t jRow=iRow+dRow; //take the row - mirror if ouside of range
222 Float_t sign0=1.;
223 if (jRow<0) sign0=-1.;
224 if (UInt_t(jRow)>=fNRows) sign0=-1.;
225 Int_t jPads= GetNPads(iRow+sign0*dRow);
4c066a96 226 Int_t offset=(nPads-jPads)/2.;
ff42e0c6 227 //
4c066a96 228 for (Int_t dPad=-deltaPad; dPad<=deltaPad; dPad++){
ff42e0c6 229 Float_t sign=sign0;
230 Int_t jPad=iPad+dPad+offset; //take the pad - mirror if ouside of range
231 Int_t kRow=jRow;
232 if (jPad<0) sign=-1;
233 if (jPad>=jPads) sign=-1;
234 if (sign<0){
b1ebb72e 235 if (!doEdge) continue;
ff42e0c6 236 kRow=iRow-dRow;
237 jPad=iPad-dPad+offset;
238 }
239 if (IsInRange(UInt_t(kRow),jPad)){
b1ebb72e 240 Bool_t isOutlier=(outlierROC==NULL)?kFALSE:outlierROC->GetValue(kRow,jPad)>0;
241 if (!isOutlier){
242 cacheBuffer[counter]=sign*(GetValue(kRow,jPad)-value);
243 counter++;
244 }
ff42e0c6 245 }
4c066a96 246 }
247 }
b1ebb72e 248 Double_t mean=0,rms=0;
249 if (TMath::Nint(fraction*Double_t(counter))>1 ) AliMathBase::EvaluateUni(counter,cacheBuffer,mean,rms,TMath::Nint(fraction*Double_t(counter)));
ff42e0c6 250 mean+=value;
4c066a96 251 newBuffer[fkIndexes[iRow]+iPad] = (type==0)? mean:rms;
252 }
253 }
254 memcpy(fData, newBuffer,GetNchannels()*sizeof(Float_t));
255 delete []newBuffer;
256 delete []cacheBuffer;
257 return kTRUE;
258}
259
0a0c96eb 260Bool_t AliTPCCalROC::Convolute(Double_t sigmaPad, Double_t sigmaRow, AliTPCCalROC*outlierROC, TF1 */*fpad*/, TF1 */*frow*/){
261 //
262 // convolute the calibration with function fpad,frow
263 // in range +-4 sigma
264
265 Float_t *newBuffer=new Float_t[fNChannels] ;
266 //
267 for (Int_t iRow=0; iRow< Int_t(fNRows); iRow++){
268 Int_t nPads=GetNPads(iRow); // number of rows in current row
269 for (Int_t iPad=0; iPad<nPads; iPad++){
270 Int_t jRow0=TMath::Max(TMath::Nint(iRow-sigmaRow*4.),0);
271 Int_t jRow1=TMath::Min(TMath::Nint(iRow+sigmaRow*4.),Int_t(fNRows));
272 Int_t jPad0=TMath::Max(TMath::Nint(iPad-sigmaPad*4.),0);
273 Int_t jPad1=TMath::Min(TMath::Nint(iRow+sigmaPad*4.),Int_t(nPads));
274 //
275 Double_t sumW=0;
276 Double_t sumCal=0;
277 for (Int_t jRow=jRow0; jRow<=jRow1; jRow++){
278 for (Int_t jPad=jPad0; jPad<=jPad1; jPad++){
279 if (!IsInRange(jRow,jPad)) continue;
280 Bool_t isOutlier=(outlierROC==NULL)?kFALSE:outlierROC->GetValue(jRow,jPad)>0;
281 if (!isOutlier){
282 Double_t weight= TMath::Gaus(jPad,iPad,sigmaPad)*TMath::Gaus(jRow,iRow,sigmaRow);
283 sumCal+=weight*GetValue(jRow,jPad);
284 sumW+=weight;
285 }
286 }
287 }
288 if (sumW>0){
289 sumCal/=sumW;
290 newBuffer[fkIndexes[iRow]+iPad]=sumCal;
291 }
292 }
293 }
294 memcpy(fData, newBuffer,GetNchannels()*sizeof(Float_t));
295 delete []newBuffer;
7f4692a6 296 return kTRUE;
0a0c96eb 297}
4c066a96 298
92e31745 299
300//
301
7b83b139 302
92e31745 303
304
a6d2bd0c 305
306// algebra fuctions:
43b569c6 307
308void AliTPCCalROC::Add(Float_t c1){
309 //
a6d2bd0c 310 // add c1 to each channel of the ROC
43b569c6 311 //
312 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata]+=c1;
313}
a6d2bd0c 314
315
43b569c6 316void AliTPCCalROC::Multiply(Float_t c1){
317 //
a6d2bd0c 318 // multiply each channel of the ROC with c1
43b569c6 319 //
320 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata]*=c1;
321}
322
a6d2bd0c 323
43b569c6 324void AliTPCCalROC::Add(const AliTPCCalROC * roc, Double_t c1){
325 //
a6d2bd0c 326 // multiply AliTPCCalROC roc by c1 and add each channel to the coresponing channel in the ROC
327 // - pad by pad
43b569c6 328 //
329 for (UInt_t idata = 0; idata< fNChannels; idata++){
330 fData[idata]+=roc->fData[idata]*c1;
331 }
332}
333
334
335void AliTPCCalROC::Multiply(const AliTPCCalROC* roc) {
336 //
a6d2bd0c 337 // multiply each channel of the ROC with the coresponding channel of 'roc'
338 // - pad by pad -
43b569c6 339 //
340 for (UInt_t idata = 0; idata< fNChannels; idata++){
341 fData[idata]*=roc->fData[idata];
342 }
343}
344
345
346void AliTPCCalROC::Divide(const AliTPCCalROC* roc) {
347 //
a6d2bd0c 348 // divide each channel of the ROC by the coresponding channel of 'roc'
349 // - pad by pad -
43b569c6 350 //
351 Float_t kEpsilon=0.00000000000000001;
352 for (UInt_t idata = 0; idata< fNChannels; idata++){
353 if (TMath::Abs(roc->fData[idata])>kEpsilon)
354 fData[idata]/=roc->fData[idata];
355 }
356}
357
a5ade541 358Double_t AliTPCCalROC::GetMean(AliTPCCalROC *const outlierROC) const {
a6d2bd0c 359 //
360 // returns the mean value of the ROC
361 // pads with value != 0 in outlierROC are not used for the calculation
9cc76bba 362 // return 0 if no data is accepted by the outlier cuts
a6d2bd0c 363 //
ca5dca67 364 if (!outlierROC) return TMath::Mean(fNChannels, fData);
365 Double_t *ddata = new Double_t[fNChannels];
a5ade541 366 Int_t nPoints = 0;
ca5dca67 367 for (UInt_t i=0;i<fNChannels;i++) {
368 if (!(outlierROC->GetValue(i))) {
a5ade541 369 ddata[nPoints]= fData[i];
370 nPoints++;
ca5dca67 371 }
372 }
9cc76bba 373 Double_t mean = 0;
a5ade541 374 if(nPoints>0)
375 mean = TMath::Mean(nPoints, ddata);
ca5dca67 376 delete [] ddata;
377 return mean;
378}
43b569c6 379
a5ade541 380Double_t AliTPCCalROC::GetMedian(AliTPCCalROC *const outlierROC) const {
a6d2bd0c 381 //
382 // returns the median value of the ROC
383 // pads with value != 0 in outlierROC are not used for the calculation
9cc76bba 384 // return 0 if no data is accepted by the outlier cuts
a6d2bd0c 385 //
ca5dca67 386 if (!outlierROC) return TMath::Median(fNChannels, fData);
387 Double_t *ddata = new Double_t[fNChannels];
a5ade541 388 Int_t nPoints = 0;
ca5dca67 389 for (UInt_t i=0;i<fNChannels;i++) {
390 if (!(outlierROC->GetValue(i))) {
a5ade541 391 ddata[nPoints]= fData[i];
392 nPoints++;
ca5dca67 393 }
394 }
9cc76bba 395 Double_t median = 0;
a5ade541 396 if(nPoints>0)
397 median = TMath::Median(nPoints, ddata);
ca5dca67 398 delete [] ddata;
9cc76bba 399 return median;
ca5dca67 400}
43b569c6 401
a5ade541 402Double_t AliTPCCalROC::GetRMS(AliTPCCalROC *const outlierROC) const {
a6d2bd0c 403 //
404 // returns the RMS value of the ROC
405 // pads with value != 0 in outlierROC are not used for the calculation
9cc76bba 406 // return 0 if no data is accepted by the outlier cuts
a6d2bd0c 407 //
ca5dca67 408 if (!outlierROC) return TMath::RMS(fNChannels, fData);
409 Double_t *ddata = new Double_t[fNChannels];
a5ade541 410 Int_t nPoints = 0;
ca5dca67 411 for (UInt_t i=0;i<fNChannels;i++) {
412 if (!(outlierROC->GetValue(i))) {
a5ade541 413 ddata[nPoints]= fData[i];
414 nPoints++;
ca5dca67 415 }
416 }
9cc76bba 417 Double_t rms = 0;
a5ade541 418 if(nPoints>0)
419 rms = TMath::RMS(nPoints, ddata);
ca5dca67 420 delete [] ddata;
9cc76bba 421 return rms;
ca5dca67 422}
c5bbaa2c 423
a5ade541 424Double_t AliTPCCalROC::GetLTM(Double_t *const sigma, Double_t fraction, AliTPCCalROC *const outlierROC){
184bcc16 425 //
a6d2bd0c 426 // returns the LTM and sigma
427 // pads with value != 0 in outlierROC are not used for the calculation
9cc76bba 428 // return 0 if no data is accepted by the outlier cuts
184bcc16 429 //
430 Double_t *ddata = new Double_t[fNChannels];
a5ade541 431 UInt_t nPoints = 0;
ca5dca67 432 for (UInt_t i=0;i<fNChannels;i++) {
433 if (!outlierROC || !(outlierROC->GetValue(i))) {
a5ade541 434 ddata[nPoints]= fData[i];
435 nPoints++;
ca5dca67 436 }
437 }
9cc76bba 438
439 Double_t ltm =0, lsigma=0;
a5ade541 440 if(nPoints>0) {
441 Int_t hh = TMath::Min(TMath::Nint(fraction *nPoints), Int_t(nPoints));
442 AliMathBase::EvaluateUni(nPoints,ddata, ltm, lsigma, hh);
9cc76bba 443 if (sigma) *sigma=lsigma;
444 }
445
184bcc16 446 delete [] ddata;
9cc76bba 447 return ltm;
184bcc16 448}
c5bbaa2c 449
184bcc16 450TH1F * AliTPCCalROC::MakeHisto1D(Float_t min, Float_t max,Int_t type){
2e9bedc9 451 //
184bcc16 452 // make 1D histo
453 // type -1 = user defined range
454 // 0 = nsigma cut nsigma=min
455 // 1 = delta cut around median delta=min
a6d2bd0c 456 //
184bcc16 457 if (type>=0){
458 if (type==0){
459 // nsigma range
460 Float_t mean = GetMean();
461 Float_t sigma = GetRMS();
462 Float_t nsigma = TMath::Abs(min);
463 min = mean-nsigma*sigma;
464 max = mean+nsigma*sigma;
465 }
466 if (type==1){
467 // fixed range
468 Float_t mean = GetMedian();
469 Float_t delta = min;
470 min = mean-delta;
471 max = mean+delta;
472 }
473 if (type==2){
474 //
475 // LTM mean +- nsigma
476 //
477 Double_t sigma;
478 Float_t mean = GetLTM(&sigma,max);
479 sigma*=min;
480 min = mean-sigma;
481 max = mean+sigma;
482 }
483 }
5dbad769 484 TString name=Form("%s ROC 1D%d",GetTitle(),fSector);
485 TH1F * his = new TH1F(name.Data(),name.Data(),100, min,max);
184bcc16 486 for (UInt_t irow=0; irow<fNRows; irow++){
487 UInt_t npads = (Int_t)GetNPads(irow);
e0bc0200 488 for (UInt_t ipad=0; ipad<npads; ipad++){
184bcc16 489 his->Fill(GetValue(irow,ipad));
490 }
491 }
492 return his;
493}
494
495
496
497TH2F * AliTPCCalROC::MakeHisto2D(Float_t min, Float_t max,Int_t type){
2e9bedc9 498 //
184bcc16 499 // make 2D histo
500 // type -1 = user defined range
501 // 0 = nsigma cut nsigma=min
502 // 1 = delta cut around median delta=min
a6d2bd0c 503 //
184bcc16 504 if (type>=0){
505 if (type==0){
506 // nsigma range
507 Float_t mean = GetMean();
508 Float_t sigma = GetRMS();
509 Float_t nsigma = TMath::Abs(min);
510 min = mean-nsigma*sigma;
511 max = mean+nsigma*sigma;
512 }
513 if (type==1){
514 // fixed range
515 Float_t mean = GetMedian();
516 Float_t delta = min;
517 min = mean-delta;
518 max = mean+delta;
519 }
520 if (type==2){
521 Double_t sigma;
522 Float_t mean = GetLTM(&sigma,max);
523 sigma*=min;
524 min = mean-sigma;
525 max = mean+sigma;
526
527 }
528 }
2e9bedc9 529 UInt_t maxPad = 0;
530 for (UInt_t irow=0; irow<fNRows; irow++){
531 if (GetNPads(irow)>maxPad) maxPad = GetNPads(irow);
532 }
5dbad769 533
534 TString name=Form("%s ROC%d",GetTitle(),fSector);
535 TH2F * his = new TH2F(name.Data(),name.Data(),fNRows+10,-5, fNRows+5, maxPad+10, -(Int_t(maxPad/2))-5, maxPad/2+5);
2e9bedc9 536 for (UInt_t irow=0; irow<fNRows; irow++){
537 UInt_t npads = (Int_t)GetNPads(irow);
e0bc0200 538 for (UInt_t ipad=0; ipad<npads; ipad++){
2e9bedc9 539 his->Fill(irow+0.5,Int_t(ipad)-Int_t(npads/2)+0.5,GetValue(irow,ipad));
540 }
541 }
184bcc16 542 his->SetMaximum(max);
543 his->SetMinimum(min);
544 return his;
545}
546
547TH2F * AliTPCCalROC::MakeHistoOutliers(Float_t delta, Float_t fraction, Int_t type){
548 //
549 // Make Histogram with outliers
550 // mode = 0 - sigma cut used
551 // mode = 1 - absolute cut used
552 // fraction - fraction of values used to define sigma
553 // delta - in mode 0 - nsigma cut
554 // mode 1 - delta cut
a6d2bd0c 555 //
184bcc16 556 Double_t sigma;
557 Float_t mean = GetLTM(&sigma,fraction);
558 if (type==0) delta*=sigma;
559 UInt_t maxPad = 0;
560 for (UInt_t irow=0; irow<fNRows; irow++){
561 if (GetNPads(irow)>maxPad) maxPad = GetNPads(irow);
562 }
563
5dbad769 564 TString name=Form("%s ROC Outliers%d",GetTitle(),fSector);
565 TH2F * his = new TH2F(name.Data(),name.Data(),fNRows+10,-5, fNRows+5, maxPad+10, -(Int_t(maxPad/2))-5, maxPad/2+5);
184bcc16 566 for (UInt_t irow=0; irow<fNRows; irow++){
567 UInt_t npads = (Int_t)GetNPads(irow);
e0bc0200 568 for (UInt_t ipad=0; ipad<npads; ipad++){
184bcc16 569 if (TMath::Abs(GetValue(irow,ipad)-mean)>delta)
570 his->Fill(irow+0.5,Int_t(ipad)-Int_t(npads/2)+0.5,1);
571 }
572 }
573 return his;
574}
575
576
577
578void AliTPCCalROC::Draw(Option_t* opt){
579 //
580 // create histogram with values and draw it
581 //
582 TH1 * his=0;
583 TString option=opt;
584 option.ToUpper();
585 if (option.Contains("1D")){
586 his = MakeHisto1D();
587 }
588 else{
589 his = MakeHisto2D();
590 }
2e9bedc9 591 his->Draw(option);
592}
593
594
184bcc16 595
596
597
ca5dca67 598void AliTPCCalROC::Test() {
c5bbaa2c 599 //
ca5dca67 600 // example function to show functionality and test AliTPCCalROC
c5bbaa2c 601 //
ca5dca67 602
603 Float_t kEpsilon=0.00001;
604
605 // create CalROC with defined values
c5bbaa2c 606 AliTPCCalROC roc0(0);
607 for (UInt_t irow = 0; irow <roc0.GetNrows(); irow++){
608 for (UInt_t ipad = 0; ipad <roc0.GetNPads(irow); ipad++){
609 Float_t value = irow+ipad/1000.;
610 roc0.SetValue(irow,ipad,value);
611 }
612 }
ca5dca67 613
614 // copy CalROC, readout values and compare to original
c5bbaa2c 615 AliTPCCalROC roc1(roc0);
616 for (UInt_t irow = 0; irow <roc1.GetNrows(); irow++){
617 for (UInt_t ipad = 0; ipad <roc1.GetNPads(irow); ipad++){
618 Float_t value = irow+ipad/1000.;
619 if (roc1.GetValue(irow,ipad)!=value){
ca5dca67 620 printf("Read/Write error\trow=%d\tpad=%d\n",irow,ipad);
c5bbaa2c 621 }
622 }
ca5dca67 623 }
624
625 // write original CalROC to file
c5bbaa2c 626 TFile f("calcTest.root","recreate");
627 roc0.Write("Roc0");
628 AliTPCCalROC * roc2 = (AliTPCCalROC*)f.Get("Roc0");
629 f.Close();
ca5dca67 630
631 // read CalROC from file and compare to original CalROC
c5bbaa2c 632 for (UInt_t irow = 0; irow <roc0.GetNrows(); irow++){
633 if (roc0.GetNPads(irow)!=roc2->GetNPads(irow))
634 printf("NPads - Read/Write error\trow=%d\n",irow);
635 for (UInt_t ipad = 0; ipad <roc1.GetNPads(irow); ipad++){
636 Float_t value = irow+ipad/1000.;
637 if (roc2->GetValue(irow,ipad)!=value){
ca5dca67 638 printf("Read/Write error\trow=%d\tpad=%d\n",irow,ipad);
c5bbaa2c 639 }
640 }
ca5dca67 641 }
642
643 //
644 // Algebra Tests
645 //
646
647 // Add constant
648 AliTPCCalROC roc3(roc0);
649 roc3.Add(1.5);
650 for (UInt_t irow = 0; irow <roc3.GetNrows(); irow++){
651 for (UInt_t ipad = 0; ipad <roc3.GetNPads(irow); ipad++){
652 Float_t value = irow+ipad/1000. + 1.5;
653 if (TMath::Abs(roc3.GetValue(irow,ipad)-value) > kEpsilon){
654 printf("Add constant - error\trow=%d\tpad=%d\n",irow,ipad);
655 }
656 }
657 }
658
659 // Add another CalROC
660 AliTPCCalROC roc4(roc0);
661 roc4.Add(&roc0, -1.5);
662 for (UInt_t irow = 0; irow <roc4.GetNrows(); irow++){
663 for (UInt_t ipad = 0; ipad <roc4.GetNPads(irow); ipad++){
664 Float_t value = irow+ipad/1000. - 1.5 * (irow+ipad/1000.);
665 if (TMath::Abs(roc4.GetValue(irow,ipad)-value) > kEpsilon){
666 printf("Add CalROC - error\trow=%d\tpad=%d\n",irow,ipad);
667 }
668 }
669 }
670
671 // Multiply with constant
672 AliTPCCalROC roc5(roc0);
673 roc5.Multiply(-1.4);
674 for (UInt_t irow = 0; irow <roc5.GetNrows(); irow++){
675 for (UInt_t ipad = 0; ipad <roc5.GetNPads(irow); ipad++){
676 Float_t value = (irow+ipad/1000.) * (-1.4);
677 if (TMath::Abs(roc5.GetValue(irow,ipad)-value) > kEpsilon){
678 printf("Multiply with constant - error\trow=%d\tpad=%d\n",irow,ipad);
679 }
680 }
681 }
682
683 // Multiply another CalROC
684 AliTPCCalROC roc6(roc0);
685 roc6.Multiply(&roc0);
686 for (UInt_t irow = 0; irow <roc6.GetNrows(); irow++){
687 for (UInt_t ipad = 0; ipad <roc6.GetNPads(irow); ipad++){
688 Float_t value = (irow+ipad/1000.) * (irow+ipad/1000.);
689 if (TMath::Abs(roc6.GetValue(irow,ipad)-value) > kEpsilon*100){
690 printf("Multiply with CalROC - error\trow=%d\tpad=%d\n",irow,ipad);
691 }
692 }
693 }
694
695
696 // Divide by CalROC
697 AliTPCCalROC roc7(roc0);
698 roc7.Divide(&roc0);
699 for (UInt_t irow = 0; irow <roc7.GetNrows(); irow++){
700 for (UInt_t ipad = 0; ipad <roc7.GetNPads(irow); ipad++){
701 Float_t value = 1.;
702 if (irow+ipad == 0) value = roc0.GetValue(irow,ipad);
703 if (TMath::Abs(roc7.GetValue(irow,ipad)-value) > kEpsilon){
704 printf("Multiply with CalROC - error\trow=%d\tpad=%d\n",irow,ipad);
705 }
706 }
707 }
708
709 //
710 // Statistics Test
711 //
712
713 // create CalROC with defined values
714 TRandom3 rnd(0);
715 AliTPCCalROC sroc0(0);
716 for (UInt_t ichannel = 0; ichannel < sroc0.GetNchannels(); ichannel++){
717 Float_t value = rnd.Gaus(10., 2.);
718 sroc0.SetValue(ichannel,value);
719 }
720
721 printf("Mean (should be close to 10): %f\n", sroc0.GetMean());
722 printf("RMS (should be close to 2): %f\n", sroc0.GetRMS());
723 printf("Median (should be close to 10): %f\n", sroc0.GetMedian());
724 printf("LTM (should be close to 10): %f\n", sroc0.GetLTM());
725
726 //AliTPCCalROC* sroc1 = sroc0.LocalFit(4, 8);
727
728 //delete sroc1;
729
730// std::cout << TMath::Abs(roc5.GetValue(irow,ipad)-value) << std::endl;
c5bbaa2c 731}
732
72fbbc82 733
a6d2bd0c 734AliTPCCalROC * AliTPCCalROC::LocalFit(Int_t rowRadius, Int_t padRadius, AliTPCCalROC* ROCoutliers, Bool_t robust, Double_t chi2Threshold, Double_t robustFraction) {
72fbbc82 735 //
736 // MakeLocalFit - smoothing
a6d2bd0c 737 // returns a AliTPCCalROC with smoothed data
738 // rowRadius and padRadius specifies a window around a given pad.
739 // The data of this window are fitted with a parabolic function.
740 // This function is evaluated at the pad's position.
741 // At the edges the window is shifted, so that the specified pad is not anymore in the center of the window.
72fbbc82 742 // rowRadius - radius - rows to be used for smoothing
743 // padradius - radius - pads to be used for smoothing
744 // ROCoutlier - map of outliers - pads not to be used for local smoothing
745 // robust - robust method of fitting - (much slower)
a6d2bd0c 746 // chi2Threshold: Threshold for chi2 when EvalRobust is called
747 // robustFraction: Fraction of data that will be used in EvalRobust
748 //
a5ade541 749 AliTPCCalROC * xROCfitted = new AliTPCCalROC(fSector);
a6d2bd0c 750 TLinearFitter fitterQ(6,"hyp5");
751 // TLinearFitter fitterQ(6,"x0++x1++x2++x3++x4++x5");
72fbbc82 752 fitterQ.StoreData(kTRUE);
f116e22c 753 for (UInt_t row=0; row < GetNrows(); row++) {
72fbbc82 754 //std::cout << "Entering row " << row << " of " << GetNrows() << " @ sector "<< fSector << " for local fitting... "<< std::endl;
f116e22c 755 for (UInt_t pad=0; pad < GetNPads(row); pad++)
a5ade541 756 xROCfitted->SetValue(row, pad, GetNeighbourhoodValue(&fitterQ, row, pad, rowRadius, padRadius, ROCoutliers, robust, chi2Threshold, robustFraction));
72fbbc82 757 }
a5ade541 758 return xROCfitted;
72fbbc82 759}
760
761
a5ade541 762Double_t AliTPCCalROC::GetNeighbourhoodValue(TLinearFitter* fitterQ, Int_t row, Int_t pad, Int_t rRadius, Int_t pRadius, AliTPCCalROC *const ROCoutliers, Bool_t robust, Double_t chi2Threshold, Double_t robustFraction) {
a6d2bd0c 763 //
764 // AliTPCCalROC::GetNeighbourhoodValue - smoothing - PRIVATE
765 // in this function the fit for LocalFit is done
72fbbc82 766 //
72fbbc82 767
768 fitterQ->ClearPoints();
769 TVectorD fitParam(6);
770 Int_t npoints=0;
771 Double_t xx[6];
772 Float_t dlx, dly;
773 Float_t lPad[3] = {0};
774 Float_t localXY[3] = {0};
775
776 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
777 tpcROCinstance->GetPositionLocal(fSector, row, pad, lPad); // calculate position lPad by pad and row number
778
779 TArrayI *neighbourhoodRows = 0;
780 TArrayI *neighbourhoodPads = 0;
ca5dca67 781
782 //std::cerr << "Trying to get neighbourhood for row " << row << ", pad " << pad << std::endl;
72fbbc82 783 GetNeighbourhood(neighbourhoodRows, neighbourhoodPads, row, pad, rRadius, pRadius);
ca5dca67 784 //std::cerr << "Got neighbourhood for row " << row << ", pad " << pad << std::endl;
72fbbc82 785
786 Int_t r, p;
787 for (Int_t i=0; i < (2*rRadius+1)*(2*pRadius+1); i++) {
788 r = neighbourhoodRows->At(i);
789 p = neighbourhoodPads->At(i);
ca5dca67 790 if (r == -1 || p == -1) continue; // window is bigger than ROC
72fbbc82 791 tpcROCinstance->GetPositionLocal(fSector, r, p, localXY); // calculate position localXY by pad and row number
792 dlx = lPad[0] - localXY[0];
793 dly = lPad[1] - localXY[1];
a6d2bd0c 794 //xx[0] = 1;
72fbbc82 795 xx[1] = dlx;
796 xx[2] = dly;
797 xx[3] = dlx*dlx;
798 xx[4] = dly*dly;
799 xx[5] = dlx*dly;
ca5dca67 800 if (!ROCoutliers || ROCoutliers->GetValue(r,p) != 1) {
a6d2bd0c 801 fitterQ->AddPoint(&xx[1], GetValue(r, p), 1);
72fbbc82 802 npoints++;
803 }
804 }
ca5dca67 805
806 delete neighbourhoodRows;
807 delete neighbourhoodPads;
808
72fbbc82 809 if (npoints < 0.5 * ((2*rRadius+1)*(2*pRadius+1)) ) {
810 // std::cerr << "Too few data points for fitting @ row " << row << ", pad " << pad << " in sector " << fSector << std::endl;
811 return 0.; // for diagnostic
812 }
813 fitterQ->Eval();
814 fitterQ->GetParameters(fitParam);
815 Float_t chi2Q = 0;
a6d2bd0c 816 if (robust) chi2Q = fitterQ->GetChisquare()/(npoints-6.);
72fbbc82 817 //if (robust) chi2Q = fitterQ->GetChisquare()/(npoints-6.);
a6d2bd0c 818 if (robust && chi2Q > chi2Threshold) {
72fbbc82 819 //std::cout << "robust fitter called... " << std::endl;
a6d2bd0c 820 fitterQ->EvalRobust(robustFraction);
72fbbc82 821 fitterQ->GetParameters(fitParam);
822 }
823 Double_t value = fitParam[0];
824
72fbbc82 825 //if (value < 0) std::cerr << "negative fit-value " << value << " in sector "<< this->fSector << " @ row: " << row << " and pad: " << pad << ", with fitter Chi2 = " << chi2Q << std::endl;
72fbbc82 826 return value;
827}
828
829
830
831
832void AliTPCCalROC::GetNeighbourhood(TArrayI* &rowArray, TArrayI* &padArray, Int_t row, Int_t pad, Int_t rRadius, Int_t pRadius) {
833 //
a6d2bd0c 834 // AliTPCCalROC::GetNeighbourhood - PRIVATE
835 // in this function the window for LocalFit is determined
72fbbc82 836 //
837 rowArray = new TArrayI((2*rRadius+1)*(2*pRadius+1));
838 padArray = new TArrayI((2*rRadius+1)*(2*pRadius+1));
72fbbc82 839
840 Int_t rmin = row - rRadius;
f116e22c 841 UInt_t rmax = row + rRadius;
72fbbc82 842
843 // if window goes out of ROC
844 if (rmin < 0) {
845 rmax = rmax - rmin;
846 rmin = 0;
847 }
848 if (rmax >= GetNrows()) {
849 rmin = rmin - (rmax - GetNrows()+1);
850 rmax = GetNrows() - 1;
851 if (rmin < 0 ) rmin = 0; // if the window is bigger than the ROC
852 }
853
854 Int_t pmin, pmax;
855 Int_t i = 0;
856
f116e22c 857 for (UInt_t r = rmin; r <= rmax; r++) {
72fbbc82 858 pmin = pad - pRadius;
859 pmax = pad + pRadius;
860 if (pmin < 0) {
861 pmax = pmax - pmin;
862 pmin = 0;
863 }
9389f9a4 864 if (pmax >= (Int_t)GetNPads(r)) {
72fbbc82 865 pmin = pmin - (pmax - GetNPads(r)+1);
866 pmax = GetNPads(r) - 1;
867 if (pmin < 0 ) pmin = 0; // if the window is bigger than the ROC
868 }
869 for (Int_t p = pmin; p <= pmax; p++) {
ca5dca67 870 (*rowArray)[i] = r;
871 (*padArray)[i] = p;
72fbbc82 872 i++;
873 }
874 }
875 for (Int_t j = i; j < rowArray->GetSize(); j++){ // unused padArray-entries, in the case that the window is bigger than the ROC
876 //std::cout << "trying to write -1" << std::endl;
ca5dca67 877 (*rowArray)[j] = -1;
878 (*padArray)[j] = -1;
72fbbc82 879 //std::cout << "writing -1" << std::endl;
ca5dca67 880 }
72fbbc82 881}
882
883
f116e22c 884
b233e5e9 885void AliTPCCalROC::GlobalFit(const AliTPCCalROC* ROCoutliers, Bool_t robust, TVectorD &fitParam, TMatrixD &covMatrix, Float_t & chi2, Int_t fitType, Double_t chi2Threshold, Double_t robustFraction, Double_t err){
a6d2bd0c 886 //
887 // Makes a GlobalFit for the given secotr and return fit-parameters, covariance and chi2
888 // The origin of the fit function is the center of the ROC!
889 // fitType == 0: fit plane function
890 // fitType == 1: fit parabolic function
891 // ROCoutliers - pads with value !=0 are not used in fitting procedure
892 // chi2Threshold: Threshold for chi2 when EvalRobust is called
893 // robustFraction: Fraction of data that will be used in EvalRobust
b233e5e9 894 // err: error of the data points
f116e22c 895 //
f116e22c 896 TLinearFitter* fitterG = 0;
897 Double_t xx[6];
898
b233e5e9 899 if (fitType == 1) {
f116e22c 900 fitterG = new TLinearFitter (6,"x0++x1++x2++x3++x4++x5");
b233e5e9 901 fitParam.ResizeTo(6);
902 covMatrix.ResizeTo(6,6);
903 } else {
f116e22c 904 fitterG = new TLinearFitter(3,"x0++x1++x2");
b233e5e9 905 fitParam.ResizeTo(3);
906 covMatrix.ResizeTo(3,3);
907 }
92e56aeb 908 fitterG->StoreData(kTRUE);
f116e22c 909 fitterG->ClearPoints();
910 Int_t npoints=0;
911
912 Float_t dlx, dly;
913 Float_t centerPad[3] = {0};
914 Float_t localXY[3] = {0};
915
916 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
917 tpcROCinstance->GetPositionLocal(fSector, GetNrows()/2, GetNPads(GetNrows()/2)/2, centerPad); // calculate center of ROC
918
919 // loop over all channels and read data into fitterG
b233e5e9 920 for (UInt_t irow = 0; irow < GetNrows(); irow++) {
921 for (UInt_t ipad = 0; ipad < GetNPads(irow); ipad++) {
922 // fill fitterG
923 if (ROCoutliers && ROCoutliers->GetValue(irow, ipad) != 0) continue;
924 tpcROCinstance->GetPositionLocal(fSector, irow, ipad, localXY); // calculate position localXY by pad and row number
925 dlx = localXY[0] - centerPad[0];
926 dly = localXY[1] - centerPad[1];
927 xx[0] = 1;
928 xx[1] = dlx;
929 xx[2] = dly;
930 xx[3] = dlx*dlx;
931 xx[4] = dly*dly;
932 xx[5] = dlx*dly;
933 npoints++;
934 fitterG->AddPoint(xx, GetValue(irow, ipad), err);
f116e22c 935 }
936 }
9cc76bba 937 if(npoints>10) { // make sure there is something to fit
938 fitterG->Eval();
f116e22c 939 fitterG->GetParameters(fitParam);
9cc76bba 940 fitterG->GetCovarianceMatrix(covMatrix);
941 if (fitType == 1)
942 chi2 = fitterG->GetChisquare()/(npoints-6.);
943 else chi2 = fitterG->GetChisquare()/(npoints-3.);
944 if (robust && chi2 > chi2Threshold) {
945 // std::cout << "robust fitter called... " << std::endl;
946 fitterG->EvalRobust(robustFraction);
947 fitterG->GetParameters(fitParam);
948 }
949 } else {
950 // set parameteres to 0
951 Int_t nParameters = 3;
952 if (fitType == 1)
953 nParameters = 6;
954
955 for(Int_t i = 0; i < nParameters; i++)
956 fitParam[i] = 0;
f116e22c 957 }
9cc76bba 958
f116e22c 959 delete fitterG;
960}
961
962
92e56aeb 963AliTPCCalROC* AliTPCCalROC::CreateGlobalFitCalROC(TVectorD &fitParam, Int_t sector){
f116e22c 964 //
965 // Create ROC with global fit parameters
a6d2bd0c 966 // The origin of the fit function is the center of the ROC!
967 // loop over all channels, write fit values into new ROC and return it
f116e22c 968 //
969 Float_t dlx, dly;
970 Float_t centerPad[3] = {0};
971 Float_t localXY[3] = {0};
a5ade541 972 AliTPCCalROC * xROCfitted = new AliTPCCalROC(sector);
f116e22c 973 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
a5ade541 974 tpcROCinstance->GetPositionLocal(sector, xROCfitted->GetNrows()/2, xROCfitted->GetNPads(xROCfitted->GetNrows()/2)/2, centerPad); // calculate center of ROC
f116e22c 975 Int_t fitType = 1;
976 if (fitParam.GetNoElements() == 6) fitType = 1;
977 else fitType = 0;
978 Double_t value = 0;
979 if (fitType == 1) { // parabolic fit
a5ade541 980 for (UInt_t irow = 0; irow < xROCfitted->GetNrows(); irow++) {
981 for (UInt_t ipad = 0; ipad < xROCfitted->GetNPads(irow); ipad++) {
f116e22c 982 tpcROCinstance->GetPositionLocal(sector, irow, ipad, localXY); // calculate position localXY by pad and row number
b233e5e9 983 dlx = localXY[0] - centerPad[0];
984 dly = localXY[1] - centerPad[1];
f116e22c 985 value = fitParam[0] + fitParam[1]*dlx + fitParam[2]*dly + fitParam[3]*dlx*dlx + fitParam[4]*dly*dly + fitParam[5]*dlx*dly;
a5ade541 986 xROCfitted->SetValue(irow, ipad, value);
f116e22c 987 }
988 }
989 }
990 else { // linear fit
a5ade541 991 for (UInt_t irow = 0; irow < xROCfitted->GetNrows(); irow++) {
992 for (UInt_t ipad = 0; ipad < xROCfitted->GetNPads(irow); ipad++) {
f116e22c 993 tpcROCinstance->GetPositionLocal(sector, irow, ipad, localXY); // calculate position localXY by pad and row number
b233e5e9 994 dlx = localXY[0] - centerPad[0];
995 dly = localXY[1] - centerPad[1];
f116e22c 996 value = fitParam[0] + fitParam[1]*dlx + fitParam[2]*dly;
a5ade541 997 xROCfitted->SetValue(irow, ipad, value);
f116e22c 998 }
999 }
1000 }
a5ade541 1001 return xROCfitted;
f116e22c 1002}
1003