1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 ///////////////////////////////////////////////////////////////////////////////
19 // Calibration base class for a single ROC //
20 // Contains one float value per pad //
21 // mapping of the pads taken form AliTPCROC //
23 ///////////////////////////////////////////////////////////////////////////////
37 #include "AliTPCCalROC.h"
38 #include "AliMathBase.h"
40 #include "TRandom3.h" // only needed by the AliTPCCalROCTest() method
42 ClassImp(AliTPCCalROC)
45 //_____________________________________________________________________________
46 AliTPCCalROC::AliTPCCalROC()
55 // Default constructor
60 //_____________________________________________________________________________
61 AliTPCCalROC::AliTPCCalROC(UInt_t sector)
70 // Constructor that initializes a given sector
73 fNChannels = AliTPCROC::Instance()->GetNChannels(fSector);
74 fNRows = AliTPCROC::Instance()->GetNRows(fSector);
75 fkIndexes = AliTPCROC::Instance()->GetRowIndexes(fSector);
76 fData = new Float_t[fNChannels];
78 // for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata] = 0.;
81 //_____________________________________________________________________________
82 AliTPCCalROC::AliTPCCalROC(const AliTPCCalROC &c)
91 // AliTPCCalROC copy constructor
94 fNChannels = AliTPCROC::Instance()->GetNChannels(fSector);
95 fNRows = AliTPCROC::Instance()->GetNRows(fSector);
96 fkIndexes = AliTPCROC::Instance()->GetRowIndexes(fSector);
98 fData = new Float_t[fNChannels];
99 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata] = c.fData[idata];
101 //____________________________________________________________________________
102 AliTPCCalROC & AliTPCCalROC::operator =(const AliTPCCalROC & param)
105 // assignment operator - dummy
107 if (this == ¶m) return (*this);
108 fSector = param.fSector;
109 fNChannels = AliTPCROC::Instance()->GetNChannels(fSector);
110 fNRows = AliTPCROC::Instance()->GetNRows(fSector);
111 fkIndexes = AliTPCROC::Instance()->GetRowIndexes(fSector);
113 if (fData) delete [] fData;
114 fData = new Float_t[fNChannels];
115 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata] = param.fData[idata];
120 //_____________________________________________________________________________
121 AliTPCCalROC::~AliTPCCalROC()
124 // AliTPCCalROC destructor
134 void AliTPCCalROC::Streamer(TBuffer &R__b)
137 // Stream an object of class AliTPCCalROC.
139 if (R__b.IsReading()) {
140 AliTPCCalROC::Class()->ReadBuffer(R__b, this);
141 fkIndexes = AliTPCROC::Instance()->GetRowIndexes(fSector);
143 AliTPCCalROC::Class()->WriteBuffer(R__b,this);
150 Bool_t AliTPCCalROC::MedianFilter(Int_t deltaRow, Int_t deltaPad, AliTPCCalROC* outlierROC, Bool_t doEdge){
153 // Modify content of the object - raplace value by median in neighorhood
155 Float_t *newBuffer=new Float_t[fNChannels] ;
156 Double_t *cacheBuffer=new Double_t[fNChannels];
158 for (Int_t iRow=0; iRow< Int_t(fNRows); iRow++){
159 Int_t nPads=GetNPads(iRow); // number of rows in current row
160 for (Int_t iPad=0; iPad<nPads; iPad++){
161 Double_t value=GetValue(iRow,iPad);
164 for (Int_t dRow=-deltaRow; dRow<=deltaRow; dRow++){
165 Int_t jRow=iRow+dRow; //take the row - mirror if ouside of range
167 if (jRow<0) sign0=-1.;
168 if (UInt_t(jRow)>=fNRows) sign0=-1.;
169 Int_t jPads= GetNPads(iRow+sign0*dRow);
170 Int_t offset=(nPads-jPads)/2.;
172 for (Int_t dPad=-deltaPad; dPad<=deltaPad; dPad++){
174 Int_t jPad=iPad+dPad+offset; //take the pad - mirror if ouside of range
177 if (jPad>=jPads) sign=-1;
180 jPad=iPad-dPad+offset;
181 if (!doEdge) continue;
183 if (IsInRange(UInt_t(kRow),jPad)){
184 Bool_t isOutlier=(outlierROC==NULL)?kFALSE:outlierROC->GetValue(kRow,jPad)>0;
186 cacheBuffer[counter]=sign*(GetValue(kRow,jPad)-value);
192 newBuffer[fkIndexes[iRow]+iPad]=0.;
193 if (counter>1) newBuffer[fkIndexes[iRow]+iPad] = TMath::Median(counter, cacheBuffer)+value;
196 memcpy(fData, newBuffer,GetNchannels()*sizeof(Float_t));
198 delete []cacheBuffer;
204 Bool_t AliTPCCalROC::LTMFilter(Int_t deltaRow, Int_t deltaPad, Float_t fraction, Int_t type, AliTPCCalROC* outlierROC, Bool_t doEdge){
209 // Modify content of the class
210 // write LTM mean or median
211 if (fraction<0 || fraction>1) return kFALSE;
212 Float_t *newBuffer=new Float_t[fNChannels] ;
213 Double_t *cacheBuffer=new Double_t[fNChannels];
215 for (Int_t iRow=0; iRow< Int_t(fNRows); iRow++){
216 Int_t nPads=GetNPads(iRow); // number of rows in current row
217 for (Int_t iPad=0; iPad<nPads; iPad++){
218 Double_t value=GetValue(iRow,iPad);
221 for (Int_t dRow=-deltaRow; dRow<=deltaRow; dRow++){
222 Int_t jRow=iRow+dRow; //take the row - mirror if ouside of range
224 if (jRow<0) sign0=-1.;
225 if (UInt_t(jRow)>=fNRows) sign0=-1.;
226 Int_t jPads= GetNPads(iRow+sign0*dRow);
227 Int_t offset=(nPads-jPads)/2.;
229 for (Int_t dPad=-deltaPad; dPad<=deltaPad; dPad++){
231 Int_t jPad=iPad+dPad+offset; //take the pad - mirror if ouside of range
234 if (jPad>=jPads) sign=-1;
236 if (!doEdge) continue;
238 jPad=iPad-dPad+offset;
240 if (IsInRange(UInt_t(kRow),jPad)){
241 Bool_t isOutlier=(outlierROC==NULL)?kFALSE:outlierROC->GetValue(kRow,jPad)>0;
243 cacheBuffer[counter]=sign*(GetValue(kRow,jPad)-value);
249 Double_t mean=0,rms=0;
250 if (TMath::Nint(fraction*Double_t(counter))>1 ) AliMathBase::EvaluateUni(counter,cacheBuffer,mean,rms,TMath::Nint(fraction*Double_t(counter)));
252 newBuffer[fkIndexes[iRow]+iPad] = (type==0)? mean:rms;
255 memcpy(fData, newBuffer,GetNchannels()*sizeof(Float_t));
257 delete []cacheBuffer;
261 Bool_t AliTPCCalROC::Convolute(Double_t sigmaPad, Double_t sigmaRow, AliTPCCalROC*outlierROC, TF1 */*fpad*/, TF1 */*frow*/){
263 // convolute the calibration with function fpad,frow
264 // in range +-4 sigma
266 Float_t *newBuffer=new Float_t[fNChannels] ;
268 for (Int_t iRow=0; iRow< Int_t(fNRows); iRow++){
269 Int_t nPads=GetNPads(iRow); // number of rows in current row
270 for (Int_t iPad=0; iPad<nPads; iPad++){
271 Int_t jRow0=TMath::Max(TMath::Nint(iRow-sigmaRow*4.),0);
272 Int_t jRow1=TMath::Min(TMath::Nint(iRow+sigmaRow*4.),Int_t(fNRows));
273 Int_t jPad0=TMath::Max(TMath::Nint(iPad-sigmaPad*4.),0);
274 Int_t jPad1=TMath::Min(TMath::Nint(iRow+sigmaPad*4.),Int_t(nPads));
278 for (Int_t jRow=jRow0; jRow<=jRow1; jRow++){
279 for (Int_t jPad=jPad0; jPad<=jPad1; jPad++){
280 if (!IsInRange(jRow,jPad)) continue;
281 Bool_t isOutlier=(outlierROC==NULL)?kFALSE:outlierROC->GetValue(jRow,jPad)>0;
283 Double_t weight= TMath::Gaus(jPad,iPad,sigmaPad)*TMath::Gaus(jRow,iRow,sigmaRow);
284 sumCal+=weight*GetValue(jRow,jPad);
291 newBuffer[fkIndexes[iRow]+iPad]=sumCal;
295 memcpy(fData, newBuffer,GetNchannels()*sizeof(Float_t));
309 void AliTPCCalROC::Add(Float_t c1){
311 // add c1 to each channel of the ROC
313 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata]+=c1;
317 void AliTPCCalROC::Multiply(Float_t c1){
319 // multiply each channel of the ROC with c1
321 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata]*=c1;
325 void AliTPCCalROC::Add(const AliTPCCalROC * roc, Double_t c1){
327 // multiply AliTPCCalROC roc by c1 and add each channel to the coresponing channel in the ROC
331 for (UInt_t idata = 0; idata< fNChannels; idata++){
332 fData[idata]+=roc->fData[idata]*c1;
337 void AliTPCCalROC::Multiply(const AliTPCCalROC* roc) {
339 // multiply each channel of the ROC with the coresponding channel of 'roc'
343 for (UInt_t idata = 0; idata< fNChannels; idata++){
344 fData[idata]*=roc->fData[idata];
349 void AliTPCCalROC::Divide(const AliTPCCalROC* roc) {
351 // divide each channel of the ROC by the coresponding channel of 'roc'
355 Float_t kEpsilon=0.00000000000000001;
356 for (UInt_t idata = 0; idata< fNChannels; idata++){
357 if (TMath::Abs(roc->fData[idata])>kEpsilon)
358 fData[idata]/=roc->fData[idata];
362 void AliTPCCalROC::Reset()
367 memset(fData,0,sizeof(Float_t)*fNChannels); // set all values to 0
370 Double_t AliTPCCalROC::GetMean(AliTPCCalROC *const outlierROC) const {
372 // returns the mean value of the ROC
373 // pads with value != 0 in outlierROC are not used for the calculation
374 // return 0 if no data is accepted by the outlier cuts
376 if (!outlierROC) return TMath::Mean(fNChannels, fData);
377 Double_t *ddata = new Double_t[fNChannels];
379 for (UInt_t i=0;i<fNChannels;i++) {
380 if (!(outlierROC->GetValue(i))) {
381 ddata[nPoints]= fData[i];
387 mean = TMath::Mean(nPoints, ddata);
392 Double_t AliTPCCalROC::GetMedian(AliTPCCalROC *const outlierROC) const {
394 // returns the median value of the ROC
395 // pads with value != 0 in outlierROC are not used for the calculation
396 // return 0 if no data is accepted by the outlier cuts
398 if (!outlierROC) return TMath::Median(fNChannels, fData);
399 Double_t *ddata = new Double_t[fNChannels];
401 for (UInt_t i=0;i<fNChannels;i++) {
402 if (!(outlierROC->GetValue(i))) {
403 ddata[nPoints]= fData[i];
409 median = TMath::Median(nPoints, ddata);
414 Double_t AliTPCCalROC::GetRMS(AliTPCCalROC *const outlierROC) const {
416 // returns the RMS value of the ROC
417 // pads with value != 0 in outlierROC are not used for the calculation
418 // return 0 if no data is accepted by the outlier cuts
420 if (!outlierROC) return TMath::RMS(fNChannels, fData);
421 Double_t *ddata = new Double_t[fNChannels];
423 for (UInt_t i=0;i<fNChannels;i++) {
424 if (!(outlierROC->GetValue(i))) {
425 ddata[nPoints]= fData[i];
431 rms = TMath::RMS(nPoints, ddata);
436 Double_t AliTPCCalROC::GetLTM(Double_t *const sigma, Double_t fraction, AliTPCCalROC *const outlierROC){
438 // returns the LTM and sigma
439 // pads with value != 0 in outlierROC are not used for the calculation
440 // return 0 if no data is accepted by the outlier cuts
442 Double_t *ddata = new Double_t[fNChannels];
444 for (UInt_t i=0;i<fNChannels;i++) {
445 if (!outlierROC || !(outlierROC->GetValue(i))) {
446 ddata[nPoints]= fData[i];
451 Double_t ltm =0, lsigma=0;
453 Int_t hh = TMath::Min(TMath::Nint(fraction *nPoints), Int_t(nPoints));
454 AliMathBase::EvaluateUni(nPoints,ddata, ltm, lsigma, hh);
455 if (sigma) *sigma=lsigma;
462 TH1F * AliTPCCalROC::MakeHisto1D(Float_t min, Float_t max,Int_t type){
465 // type -1 = user defined range
466 // 0 = nsigma cut nsigma=min
467 // 1 = delta cut around median delta=min
472 Float_t mean = GetMean();
473 Float_t sigma = GetRMS();
474 Float_t nsigma = TMath::Abs(min);
475 min = mean-nsigma*sigma;
476 max = mean+nsigma*sigma;
480 Float_t mean = GetMedian();
487 // LTM mean +- nsigma
490 Float_t mean = GetLTM(&sigma,max);
496 TString name=Form("%s ROC 1D%d",GetTitle(),fSector);
497 TH1F * his = new TH1F(name.Data(),name.Data(),100, min,max);
498 for (UInt_t irow=0; irow<fNRows; irow++){
499 UInt_t npads = (Int_t)GetNPads(irow);
500 for (UInt_t ipad=0; ipad<npads; ipad++){
501 his->Fill(GetValue(irow,ipad));
509 TH2F * AliTPCCalROC::MakeHisto2D(Float_t min, Float_t max,Int_t type){
512 // type -1 = user defined range
513 // 0 = nsigma cut nsigma=min
514 // 1 = delta cut around median delta=min
519 Float_t mean = GetMean();
520 Float_t sigma = GetRMS();
521 Float_t nsigma = TMath::Abs(min);
522 min = mean-nsigma*sigma;
523 max = mean+nsigma*sigma;
527 Float_t mean = GetMedian();
534 Float_t mean = GetLTM(&sigma,max);
542 for (UInt_t irow=0; irow<fNRows; irow++){
543 if (GetNPads(irow)>maxPad) maxPad = GetNPads(irow);
546 TString name=Form("%s ROC%d",GetTitle(),fSector);
547 TH2F * his = new TH2F(name.Data(),name.Data(),fNRows+10,-5, fNRows+5, maxPad+10, -(Int_t(maxPad/2))-5, maxPad/2+5);
548 for (UInt_t irow=0; irow<fNRows; irow++){
549 UInt_t npads = (Int_t)GetNPads(irow);
550 for (UInt_t ipad=0; ipad<npads; ipad++){
551 his->Fill(irow+0.5,Int_t(ipad)-Int_t(npads/2)+0.5,GetValue(irow,ipad));
554 his->SetMaximum(max);
555 his->SetMinimum(min);
559 TH2F * AliTPCCalROC::MakeHistoOutliers(Float_t delta, Float_t fraction, Int_t type){
561 // Make Histogram with outliers
562 // mode = 0 - sigma cut used
563 // mode = 1 - absolute cut used
564 // fraction - fraction of values used to define sigma
565 // delta - in mode 0 - nsigma cut
566 // mode 1 - delta cut
569 Float_t mean = GetLTM(&sigma,fraction);
570 if (type==0) delta*=sigma;
572 for (UInt_t irow=0; irow<fNRows; irow++){
573 if (GetNPads(irow)>maxPad) maxPad = GetNPads(irow);
576 TString name=Form("%s ROC Outliers%d",GetTitle(),fSector);
577 TH2F * his = new TH2F(name.Data(),name.Data(),fNRows+10,-5, fNRows+5, maxPad+10, -(Int_t(maxPad/2))-5, maxPad/2+5);
578 for (UInt_t irow=0; irow<fNRows; irow++){
579 UInt_t npads = (Int_t)GetNPads(irow);
580 for (UInt_t ipad=0; ipad<npads; ipad++){
581 if (TMath::Abs(GetValue(irow,ipad)-mean)>delta)
582 his->Fill(irow+0.5,Int_t(ipad)-Int_t(npads/2)+0.5,1);
590 void AliTPCCalROC::Draw(Option_t* opt){
592 // create histogram with values and draw it
597 if (option.Contains("1D")){
610 void AliTPCCalROC::Test() {
612 // example function to show functionality and test AliTPCCalROC
615 Float_t kEpsilon=0.00001;
617 // create CalROC with defined values
618 AliTPCCalROC roc0(0);
619 for (UInt_t irow = 0; irow <roc0.GetNrows(); irow++){
620 for (UInt_t ipad = 0; ipad <roc0.GetNPads(irow); ipad++){
621 Float_t value = irow+ipad/1000.;
622 roc0.SetValue(irow,ipad,value);
626 // copy CalROC, readout values and compare to original
627 AliTPCCalROC roc1(roc0);
628 for (UInt_t irow = 0; irow <roc1.GetNrows(); irow++){
629 for (UInt_t ipad = 0; ipad <roc1.GetNPads(irow); ipad++){
630 Float_t value = irow+ipad/1000.;
631 if (roc1.GetValue(irow,ipad)!=value){
632 printf("Read/Write error\trow=%d\tpad=%d\n",irow,ipad);
637 // write original CalROC to file
638 TFile f("calcTest.root","recreate");
640 AliTPCCalROC * roc2 = (AliTPCCalROC*)f.Get("Roc0");
643 // read CalROC from file and compare to original CalROC
644 for (UInt_t irow = 0; irow <roc0.GetNrows(); irow++){
645 if (roc0.GetNPads(irow)!=roc2->GetNPads(irow))
646 printf("NPads - Read/Write error\trow=%d\n",irow);
647 for (UInt_t ipad = 0; ipad <roc1.GetNPads(irow); ipad++){
648 Float_t value = irow+ipad/1000.;
649 if (roc2->GetValue(irow,ipad)!=value){
650 printf("Read/Write error\trow=%d\tpad=%d\n",irow,ipad);
660 AliTPCCalROC roc3(roc0);
662 for (UInt_t irow = 0; irow <roc3.GetNrows(); irow++){
663 for (UInt_t ipad = 0; ipad <roc3.GetNPads(irow); ipad++){
664 Float_t value = irow+ipad/1000. + 1.5;
665 if (TMath::Abs(roc3.GetValue(irow,ipad)-value) > kEpsilon){
666 printf("Add constant - error\trow=%d\tpad=%d\n",irow,ipad);
671 // Add another CalROC
672 AliTPCCalROC roc4(roc0);
673 roc4.Add(&roc0, -1.5);
674 for (UInt_t irow = 0; irow <roc4.GetNrows(); irow++){
675 for (UInt_t ipad = 0; ipad <roc4.GetNPads(irow); ipad++){
676 Float_t value = irow+ipad/1000. - 1.5 * (irow+ipad/1000.);
677 if (TMath::Abs(roc4.GetValue(irow,ipad)-value) > kEpsilon){
678 printf("Add CalROC - error\trow=%d\tpad=%d\n",irow,ipad);
683 // Multiply with constant
684 AliTPCCalROC roc5(roc0);
686 for (UInt_t irow = 0; irow <roc5.GetNrows(); irow++){
687 for (UInt_t ipad = 0; ipad <roc5.GetNPads(irow); ipad++){
688 Float_t value = (irow+ipad/1000.) * (-1.4);
689 if (TMath::Abs(roc5.GetValue(irow,ipad)-value) > kEpsilon){
690 printf("Multiply with constant - error\trow=%d\tpad=%d\n",irow,ipad);
695 // Multiply another CalROC
696 AliTPCCalROC roc6(roc0);
697 roc6.Multiply(&roc0);
698 for (UInt_t irow = 0; irow <roc6.GetNrows(); irow++){
699 for (UInt_t ipad = 0; ipad <roc6.GetNPads(irow); ipad++){
700 Float_t value = (irow+ipad/1000.) * (irow+ipad/1000.);
701 if (TMath::Abs(roc6.GetValue(irow,ipad)-value) > kEpsilon*100){
702 printf("Multiply with CalROC - error\trow=%d\tpad=%d\n",irow,ipad);
709 AliTPCCalROC roc7(roc0);
711 for (UInt_t irow = 0; irow <roc7.GetNrows(); irow++){
712 for (UInt_t ipad = 0; ipad <roc7.GetNPads(irow); ipad++){
714 if (irow+ipad == 0) value = roc0.GetValue(irow,ipad);
715 if (TMath::Abs(roc7.GetValue(irow,ipad)-value) > kEpsilon){
716 printf("Multiply with CalROC - error\trow=%d\tpad=%d\n",irow,ipad);
725 // create CalROC with defined values
727 AliTPCCalROC sroc0(0);
728 for (UInt_t ichannel = 0; ichannel < sroc0.GetNchannels(); ichannel++){
729 Float_t value = rnd.Gaus(10., 2.);
730 sroc0.SetValue(ichannel,value);
733 printf("Mean (should be close to 10): %f\n", sroc0.GetMean());
734 printf("RMS (should be close to 2): %f\n", sroc0.GetRMS());
735 printf("Median (should be close to 10): %f\n", sroc0.GetMedian());
736 printf("LTM (should be close to 10): %f\n", sroc0.GetLTM());
738 //AliTPCCalROC* sroc1 = sroc0.LocalFit(4, 8);
742 // std::cout << TMath::Abs(roc5.GetValue(irow,ipad)-value) << std::endl;
746 AliTPCCalROC * AliTPCCalROC::LocalFit(Int_t rowRadius, Int_t padRadius, AliTPCCalROC* ROCoutliers, Bool_t robust, Double_t chi2Threshold, Double_t robustFraction) {
748 // MakeLocalFit - smoothing
749 // returns a AliTPCCalROC with smoothed data
750 // rowRadius and padRadius specifies a window around a given pad.
751 // The data of this window are fitted with a parabolic function.
752 // This function is evaluated at the pad's position.
753 // At the edges the window is shifted, so that the specified pad is not anymore in the center of the window.
754 // rowRadius - radius - rows to be used for smoothing
755 // padradius - radius - pads to be used for smoothing
756 // ROCoutlier - map of outliers - pads not to be used for local smoothing
757 // robust - robust method of fitting - (much slower)
758 // chi2Threshold: Threshold for chi2 when EvalRobust is called
759 // robustFraction: Fraction of data that will be used in EvalRobust
761 AliTPCCalROC * xROCfitted = new AliTPCCalROC(fSector);
762 TLinearFitter fitterQ(6,"hyp5");
763 // TLinearFitter fitterQ(6,"x0++x1++x2++x3++x4++x5");
764 fitterQ.StoreData(kTRUE);
765 for (UInt_t row=0; row < GetNrows(); row++) {
766 //std::cout << "Entering row " << row << " of " << GetNrows() << " @ sector "<< fSector << " for local fitting... "<< std::endl;
767 for (UInt_t pad=0; pad < GetNPads(row); pad++)
768 xROCfitted->SetValue(row, pad, GetNeighbourhoodValue(&fitterQ, row, pad, rowRadius, padRadius, ROCoutliers, robust, chi2Threshold, robustFraction));
774 Double_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) {
776 // AliTPCCalROC::GetNeighbourhoodValue - smoothing - PRIVATE
777 // in this function the fit for LocalFit is done
780 fitterQ->ClearPoints();
781 TVectorD fitParam(6);
785 Float_t lPad[3] = {0};
786 Float_t localXY[3] = {0};
788 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
789 tpcROCinstance->GetPositionLocal(fSector, row, pad, lPad); // calculate position lPad by pad and row number
791 TArrayI *neighbourhoodRows = 0;
792 TArrayI *neighbourhoodPads = 0;
794 //std::cerr << "Trying to get neighbourhood for row " << row << ", pad " << pad << std::endl;
795 GetNeighbourhood(neighbourhoodRows, neighbourhoodPads, row, pad, rRadius, pRadius);
796 //std::cerr << "Got neighbourhood for row " << row << ", pad " << pad << std::endl;
799 for (Int_t i=0; i < (2*rRadius+1)*(2*pRadius+1); i++) {
800 r = neighbourhoodRows->At(i);
801 p = neighbourhoodPads->At(i);
802 if (r == -1 || p == -1) continue; // window is bigger than ROC
803 tpcROCinstance->GetPositionLocal(fSector, r, p, localXY); // calculate position localXY by pad and row number
804 dlx = lPad[0] - localXY[0];
805 dly = lPad[1] - localXY[1];
812 if (!ROCoutliers || ROCoutliers->GetValue(r,p) != 1) {
813 fitterQ->AddPoint(&xx[1], GetValue(r, p), 1);
818 delete neighbourhoodRows;
819 delete neighbourhoodPads;
821 if (npoints < 0.5 * ((2*rRadius+1)*(2*pRadius+1)) ) {
822 // std::cerr << "Too few data points for fitting @ row " << row << ", pad " << pad << " in sector " << fSector << std::endl;
823 return 0.; // for diagnostic
826 fitterQ->GetParameters(fitParam);
828 if (robust) chi2Q = fitterQ->GetChisquare()/(npoints-6.);
829 //if (robust) chi2Q = fitterQ->GetChisquare()/(npoints-6.);
830 if (robust && chi2Q > chi2Threshold) {
831 //std::cout << "robust fitter called... " << std::endl;
832 fitterQ->EvalRobust(robustFraction);
833 fitterQ->GetParameters(fitParam);
835 Double_t value = fitParam[0];
837 //if (value < 0) std::cerr << "negative fit-value " << value << " in sector "<< this->fSector << " @ row: " << row << " and pad: " << pad << ", with fitter Chi2 = " << chi2Q << std::endl;
844 void AliTPCCalROC::GetNeighbourhood(TArrayI* &rowArray, TArrayI* &padArray, Int_t row, Int_t pad, Int_t rRadius, Int_t pRadius) {
846 // AliTPCCalROC::GetNeighbourhood - PRIVATE
847 // in this function the window for LocalFit is determined
849 rowArray = new TArrayI((2*rRadius+1)*(2*pRadius+1));
850 padArray = new TArrayI((2*rRadius+1)*(2*pRadius+1));
852 Int_t rmin = row - rRadius;
853 UInt_t rmax = row + rRadius;
855 // if window goes out of ROC
860 if (rmax >= GetNrows()) {
861 rmin = rmin - (rmax - GetNrows()+1);
862 rmax = GetNrows() - 1;
863 if (rmin < 0 ) rmin = 0; // if the window is bigger than the ROC
869 for (UInt_t r = rmin; r <= rmax; r++) {
870 pmin = pad - pRadius;
871 pmax = pad + pRadius;
876 if (pmax >= (Int_t)GetNPads(r)) {
877 pmin = pmin - (pmax - GetNPads(r)+1);
878 pmax = GetNPads(r) - 1;
879 if (pmin < 0 ) pmin = 0; // if the window is bigger than the ROC
881 for (Int_t p = pmin; p <= pmax; p++) {
887 for (Int_t j = i; j < rowArray->GetSize(); j++){ // unused padArray-entries, in the case that the window is bigger than the ROC
888 //std::cout << "trying to write -1" << std::endl;
891 //std::cout << "writing -1" << std::endl;
897 void 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){
899 // Makes a GlobalFit for the given secotr and return fit-parameters, covariance and chi2
900 // The origin of the fit function is the center of the ROC!
901 // fitType == 0: fit plane function
902 // fitType == 1: fit parabolic function
903 // ROCoutliers - pads with value !=0 are not used in fitting procedure
904 // chi2Threshold: Threshold for chi2 when EvalRobust is called
905 // robustFraction: Fraction of data that will be used in EvalRobust
906 // err: error of the data points
908 TLinearFitter* fitterG = 0;
912 fitterG = new TLinearFitter (6,"x0++x1++x2++x3++x4++x5");
913 fitParam.ResizeTo(6);
914 covMatrix.ResizeTo(6,6);
916 fitterG = new TLinearFitter(3,"x0++x1++x2");
917 fitParam.ResizeTo(3);
918 covMatrix.ResizeTo(3,3);
920 fitterG->StoreData(kTRUE);
921 fitterG->ClearPoints();
925 Float_t centerPad[3] = {0};
926 Float_t localXY[3] = {0};
928 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
929 tpcROCinstance->GetPositionLocal(fSector, GetNrows()/2, GetNPads(GetNrows()/2)/2, centerPad); // calculate center of ROC
931 // loop over all channels and read data into fitterG
932 for (UInt_t irow = 0; irow < GetNrows(); irow++) {
933 for (UInt_t ipad = 0; ipad < GetNPads(irow); ipad++) {
935 if (ROCoutliers && ROCoutliers->GetValue(irow, ipad) != 0) continue;
936 tpcROCinstance->GetPositionLocal(fSector, irow, ipad, localXY); // calculate position localXY by pad and row number
937 dlx = localXY[0] - centerPad[0];
938 dly = localXY[1] - centerPad[1];
946 fitterG->AddPoint(xx, GetValue(irow, ipad), err);
949 if(npoints>10) { // make sure there is something to fit
951 fitterG->GetParameters(fitParam);
952 fitterG->GetCovarianceMatrix(covMatrix);
954 chi2 = fitterG->GetChisquare()/(npoints-6.);
955 else chi2 = fitterG->GetChisquare()/(npoints-3.);
956 if (robust && chi2 > chi2Threshold) {
957 // std::cout << "robust fitter called... " << std::endl;
958 fitterG->EvalRobust(robustFraction);
959 fitterG->GetParameters(fitParam);
962 // set parameteres to 0
963 Int_t nParameters = 3;
967 for(Int_t i = 0; i < nParameters; i++)
975 AliTPCCalROC* AliTPCCalROC::CreateGlobalFitCalROC(TVectorD &fitParam, Int_t sector){
977 // Create ROC with global fit parameters
978 // The origin of the fit function is the center of the ROC!
979 // loop over all channels, write fit values into new ROC and return it
982 Float_t centerPad[3] = {0};
983 Float_t localXY[3] = {0};
984 AliTPCCalROC * xROCfitted = new AliTPCCalROC(sector);
985 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
986 tpcROCinstance->GetPositionLocal(sector, xROCfitted->GetNrows()/2, xROCfitted->GetNPads(xROCfitted->GetNrows()/2)/2, centerPad); // calculate center of ROC
988 if (fitParam.GetNoElements() == 6) fitType = 1;
991 if (fitType == 1) { // parabolic fit
992 for (UInt_t irow = 0; irow < xROCfitted->GetNrows(); irow++) {
993 for (UInt_t ipad = 0; ipad < xROCfitted->GetNPads(irow); ipad++) {
994 tpcROCinstance->GetPositionLocal(sector, irow, ipad, localXY); // calculate position localXY by pad and row number
995 dlx = localXY[0] - centerPad[0];
996 dly = localXY[1] - centerPad[1];
997 value = fitParam[0] + fitParam[1]*dlx + fitParam[2]*dly + fitParam[3]*dlx*dlx + fitParam[4]*dly*dly + fitParam[5]*dlx*dly;
998 xROCfitted->SetValue(irow, ipad, value);
1002 else { // linear fit
1003 for (UInt_t irow = 0; irow < xROCfitted->GetNrows(); irow++) {
1004 for (UInt_t ipad = 0; ipad < xROCfitted->GetNPads(irow); ipad++) {
1005 tpcROCinstance->GetPositionLocal(sector, irow, ipad, localXY); // calculate position localXY by pad and row number
1006 dlx = localXY[0] - centerPad[0];
1007 dly = localXY[1] - centerPad[1];
1008 value = fitParam[0] + fitParam[1]*dlx + fitParam[2]*dly;
1009 xROCfitted->SetValue(irow, ipad, value);