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 fIndexes = AliTPCROC::Instance()->GetRowIndexes(fSector);
76 fData = new Float_t[fNChannels];
77 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata] = 0.;
80 //_____________________________________________________________________________
81 AliTPCCalROC::AliTPCCalROC(const AliTPCCalROC &c)
90 // AliTPCCalROC copy constructor
93 fNChannels = AliTPCROC::Instance()->GetNChannels(fSector);
94 fNRows = AliTPCROC::Instance()->GetNRows(fSector);
95 fIndexes = AliTPCROC::Instance()->GetRowIndexes(fSector);
97 fData = new Float_t[fNChannels];
98 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata] = c.fData[idata];
100 //____________________________________________________________________________
101 AliTPCCalROC & AliTPCCalROC::operator =(const AliTPCCalROC & param)
104 // assignment operator - dummy
111 //_____________________________________________________________________________
112 AliTPCCalROC::~AliTPCCalROC()
115 // AliTPCCalROC destructor
125 void AliTPCCalROC::Streamer(TBuffer &R__b)
128 // Stream an object of class AliTPCCalROC.
130 if (R__b.IsReading()) {
131 AliTPCCalROC::Class()->ReadBuffer(R__b, this);
132 fIndexes = AliTPCROC::Instance()->GetRowIndexes(fSector);
134 AliTPCCalROC::Class()->WriteBuffer(R__b,this);
141 void AliTPCCalROC::Add(Float_t c1){
143 // add c1 to each channel of the ROC
145 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata]+=c1;
149 void AliTPCCalROC::Multiply(Float_t c1){
151 // multiply each channel of the ROC with c1
153 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata]*=c1;
157 void AliTPCCalROC::Add(const AliTPCCalROC * roc, Double_t c1){
159 // multiply AliTPCCalROC roc by c1 and add each channel to the coresponing channel in the ROC
162 for (UInt_t idata = 0; idata< fNChannels; idata++){
163 fData[idata]+=roc->fData[idata]*c1;
168 void AliTPCCalROC::Multiply(const AliTPCCalROC* roc) {
170 // multiply each channel of the ROC with the coresponding channel of 'roc'
173 for (UInt_t idata = 0; idata< fNChannels; idata++){
174 fData[idata]*=roc->fData[idata];
179 void AliTPCCalROC::Divide(const AliTPCCalROC* roc) {
181 // divide each channel of the ROC by the coresponding channel of 'roc'
184 Float_t kEpsilon=0.00000000000000001;
185 for (UInt_t idata = 0; idata< fNChannels; idata++){
186 if (TMath::Abs(roc->fData[idata])>kEpsilon)
187 fData[idata]/=roc->fData[idata];
191 Double_t AliTPCCalROC::GetMean(AliTPCCalROC* outlierROC) {
193 // returns the mean value of the ROC
194 // pads with value != 0 in outlierROC are not used for the calculation
195 // return 0 if no data is accepted by the outlier cuts
197 if (!outlierROC) return TMath::Mean(fNChannels, fData);
198 Double_t *ddata = new Double_t[fNChannels];
200 for (UInt_t i=0;i<fNChannels;i++) {
201 if (!(outlierROC->GetValue(i))) {
202 ddata[NPoints]= fData[i];
208 mean = TMath::Mean(NPoints, ddata);
213 Double_t AliTPCCalROC::GetMedian(AliTPCCalROC* outlierROC) {
215 // returns the median value of the ROC
216 // pads with value != 0 in outlierROC are not used for the calculation
217 // return 0 if no data is accepted by the outlier cuts
219 if (!outlierROC) return TMath::Median(fNChannels, fData);
220 Double_t *ddata = new Double_t[fNChannels];
222 for (UInt_t i=0;i<fNChannels;i++) {
223 if (!(outlierROC->GetValue(i))) {
224 ddata[NPoints]= fData[i];
230 median = TMath::Median(NPoints, ddata);
235 Double_t AliTPCCalROC::GetRMS(AliTPCCalROC* outlierROC) {
237 // returns the RMS value of the ROC
238 // pads with value != 0 in outlierROC are not used for the calculation
239 // return 0 if no data is accepted by the outlier cuts
241 if (!outlierROC) return TMath::RMS(fNChannels, fData);
242 Double_t *ddata = new Double_t[fNChannels];
244 for (UInt_t i=0;i<fNChannels;i++) {
245 if (!(outlierROC->GetValue(i))) {
246 ddata[NPoints]= fData[i];
252 rms = TMath::RMS(NPoints, ddata);
257 Double_t AliTPCCalROC::GetLTM(Double_t *sigma, Double_t fraction, AliTPCCalROC* outlierROC){
259 // returns the LTM and sigma
260 // pads with value != 0 in outlierROC are not used for the calculation
261 // return 0 if no data is accepted by the outlier cuts
263 Double_t *ddata = new Double_t[fNChannels];
265 for (UInt_t i=0;i<fNChannels;i++) {
266 if (!outlierROC || !(outlierROC->GetValue(i))) {
267 ddata[NPoints]= fData[i];
272 Double_t ltm =0, lsigma=0;
274 Int_t hh = TMath::Min(TMath::Nint(fraction *NPoints), Int_t(NPoints));
275 AliMathBase::EvaluateUni(NPoints,ddata, ltm, lsigma, hh);
276 if (sigma) *sigma=lsigma;
283 TH1F * AliTPCCalROC::MakeHisto1D(Float_t min, Float_t max,Int_t type){
286 // type -1 = user defined range
287 // 0 = nsigma cut nsigma=min
288 // 1 = delta cut around median delta=min
293 Float_t mean = GetMean();
294 Float_t sigma = GetRMS();
295 Float_t nsigma = TMath::Abs(min);
296 min = mean-nsigma*sigma;
297 max = mean+nsigma*sigma;
301 Float_t mean = GetMedian();
308 // LTM mean +- nsigma
311 Float_t mean = GetLTM(&sigma,max);
318 sprintf(name,"%s ROC 1D%d",GetTitle(),fSector);
319 TH1F * his = new TH1F(name,name,100, min,max);
320 for (UInt_t irow=0; irow<fNRows; irow++){
321 UInt_t npads = (Int_t)GetNPads(irow);
322 for (UInt_t ipad=0; ipad<npads; ipad++){
323 his->Fill(GetValue(irow,ipad));
331 TH2F * AliTPCCalROC::MakeHisto2D(Float_t min, Float_t max,Int_t type){
334 // type -1 = user defined range
335 // 0 = nsigma cut nsigma=min
336 // 1 = delta cut around median delta=min
341 Float_t mean = GetMean();
342 Float_t sigma = GetRMS();
343 Float_t nsigma = TMath::Abs(min);
344 min = mean-nsigma*sigma;
345 max = mean+nsigma*sigma;
349 Float_t mean = GetMedian();
356 Float_t mean = GetLTM(&sigma,max);
364 for (UInt_t irow=0; irow<fNRows; irow++){
365 if (GetNPads(irow)>maxPad) maxPad = GetNPads(irow);
368 sprintf(name,"%s ROC%d",GetTitle(),fSector);
369 TH2F * his = new TH2F(name,name,fNRows+10,-5, fNRows+5, maxPad+10, -(Int_t(maxPad/2))-5, maxPad/2+5);
370 for (UInt_t irow=0; irow<fNRows; irow++){
371 UInt_t npads = (Int_t)GetNPads(irow);
372 for (UInt_t ipad=0; ipad<npads; ipad++){
373 his->Fill(irow+0.5,Int_t(ipad)-Int_t(npads/2)+0.5,GetValue(irow,ipad));
376 his->SetMaximum(max);
377 his->SetMinimum(min);
381 TH2F * AliTPCCalROC::MakeHistoOutliers(Float_t delta, Float_t fraction, Int_t type){
383 // Make Histogram with outliers
384 // mode = 0 - sigma cut used
385 // mode = 1 - absolute cut used
386 // fraction - fraction of values used to define sigma
387 // delta - in mode 0 - nsigma cut
388 // mode 1 - delta cut
391 Float_t mean = GetLTM(&sigma,fraction);
392 if (type==0) delta*=sigma;
394 for (UInt_t irow=0; irow<fNRows; irow++){
395 if (GetNPads(irow)>maxPad) maxPad = GetNPads(irow);
399 sprintf(name,"%s ROC Outliers%d",GetTitle(),fSector);
400 TH2F * his = new TH2F(name,name,fNRows+10,-5, fNRows+5, maxPad+10, -(Int_t(maxPad/2))-5, maxPad/2+5);
401 for (UInt_t irow=0; irow<fNRows; irow++){
402 UInt_t npads = (Int_t)GetNPads(irow);
403 for (UInt_t ipad=0; ipad<npads; ipad++){
404 if (TMath::Abs(GetValue(irow,ipad)-mean)>delta)
405 his->Fill(irow+0.5,Int_t(ipad)-Int_t(npads/2)+0.5,1);
413 void AliTPCCalROC::Draw(Option_t* opt){
415 // create histogram with values and draw it
420 if (option.Contains("1D")){
433 void AliTPCCalROC::Test() {
435 // example function to show functionality and test AliTPCCalROC
438 Float_t kEpsilon=0.00001;
440 // create CalROC with defined values
441 AliTPCCalROC roc0(0);
442 for (UInt_t irow = 0; irow <roc0.GetNrows(); irow++){
443 for (UInt_t ipad = 0; ipad <roc0.GetNPads(irow); ipad++){
444 Float_t value = irow+ipad/1000.;
445 roc0.SetValue(irow,ipad,value);
449 // copy CalROC, readout values and compare to original
450 AliTPCCalROC roc1(roc0);
451 for (UInt_t irow = 0; irow <roc1.GetNrows(); irow++){
452 for (UInt_t ipad = 0; ipad <roc1.GetNPads(irow); ipad++){
453 Float_t value = irow+ipad/1000.;
454 if (roc1.GetValue(irow,ipad)!=value){
455 printf("Read/Write error\trow=%d\tpad=%d\n",irow,ipad);
460 // write original CalROC to file
461 TFile f("calcTest.root","recreate");
463 AliTPCCalROC * roc2 = (AliTPCCalROC*)f.Get("Roc0");
466 // read CalROC from file and compare to original CalROC
467 for (UInt_t irow = 0; irow <roc0.GetNrows(); irow++){
468 if (roc0.GetNPads(irow)!=roc2->GetNPads(irow))
469 printf("NPads - Read/Write error\trow=%d\n",irow);
470 for (UInt_t ipad = 0; ipad <roc1.GetNPads(irow); ipad++){
471 Float_t value = irow+ipad/1000.;
472 if (roc2->GetValue(irow,ipad)!=value){
473 printf("Read/Write error\trow=%d\tpad=%d\n",irow,ipad);
483 AliTPCCalROC roc3(roc0);
485 for (UInt_t irow = 0; irow <roc3.GetNrows(); irow++){
486 for (UInt_t ipad = 0; ipad <roc3.GetNPads(irow); ipad++){
487 Float_t value = irow+ipad/1000. + 1.5;
488 if (TMath::Abs(roc3.GetValue(irow,ipad)-value) > kEpsilon){
489 printf("Add constant - error\trow=%d\tpad=%d\n",irow,ipad);
494 // Add another CalROC
495 AliTPCCalROC roc4(roc0);
496 roc4.Add(&roc0, -1.5);
497 for (UInt_t irow = 0; irow <roc4.GetNrows(); irow++){
498 for (UInt_t ipad = 0; ipad <roc4.GetNPads(irow); ipad++){
499 Float_t value = irow+ipad/1000. - 1.5 * (irow+ipad/1000.);
500 if (TMath::Abs(roc4.GetValue(irow,ipad)-value) > kEpsilon){
501 printf("Add CalROC - error\trow=%d\tpad=%d\n",irow,ipad);
506 // Multiply with constant
507 AliTPCCalROC roc5(roc0);
509 for (UInt_t irow = 0; irow <roc5.GetNrows(); irow++){
510 for (UInt_t ipad = 0; ipad <roc5.GetNPads(irow); ipad++){
511 Float_t value = (irow+ipad/1000.) * (-1.4);
512 if (TMath::Abs(roc5.GetValue(irow,ipad)-value) > kEpsilon){
513 printf("Multiply with constant - error\trow=%d\tpad=%d\n",irow,ipad);
518 // Multiply another CalROC
519 AliTPCCalROC roc6(roc0);
520 roc6.Multiply(&roc0);
521 for (UInt_t irow = 0; irow <roc6.GetNrows(); irow++){
522 for (UInt_t ipad = 0; ipad <roc6.GetNPads(irow); ipad++){
523 Float_t value = (irow+ipad/1000.) * (irow+ipad/1000.);
524 if (TMath::Abs(roc6.GetValue(irow,ipad)-value) > kEpsilon*100){
525 printf("Multiply with CalROC - error\trow=%d\tpad=%d\n",irow,ipad);
532 AliTPCCalROC roc7(roc0);
534 for (UInt_t irow = 0; irow <roc7.GetNrows(); irow++){
535 for (UInt_t ipad = 0; ipad <roc7.GetNPads(irow); ipad++){
537 if (irow+ipad == 0) value = roc0.GetValue(irow,ipad);
538 if (TMath::Abs(roc7.GetValue(irow,ipad)-value) > kEpsilon){
539 printf("Multiply with CalROC - error\trow=%d\tpad=%d\n",irow,ipad);
548 // create CalROC with defined values
550 AliTPCCalROC sroc0(0);
551 for (UInt_t ichannel = 0; ichannel < sroc0.GetNchannels(); ichannel++){
552 Float_t value = rnd.Gaus(10., 2.);
553 sroc0.SetValue(ichannel,value);
556 printf("Mean (should be close to 10): %f\n", sroc0.GetMean());
557 printf("RMS (should be close to 2): %f\n", sroc0.GetRMS());
558 printf("Median (should be close to 10): %f\n", sroc0.GetMedian());
559 printf("LTM (should be close to 10): %f\n", sroc0.GetLTM());
561 //AliTPCCalROC* sroc1 = sroc0.LocalFit(4, 8);
565 // std::cout << TMath::Abs(roc5.GetValue(irow,ipad)-value) << std::endl;
569 AliTPCCalROC * AliTPCCalROC::LocalFit(Int_t rowRadius, Int_t padRadius, AliTPCCalROC* ROCoutliers, Bool_t robust, Double_t chi2Threshold, Double_t robustFraction) {
571 // MakeLocalFit - smoothing
572 // returns a AliTPCCalROC with smoothed data
573 // rowRadius and padRadius specifies a window around a given pad.
574 // The data of this window are fitted with a parabolic function.
575 // This function is evaluated at the pad's position.
576 // At the edges the window is shifted, so that the specified pad is not anymore in the center of the window.
577 // rowRadius - radius - rows to be used for smoothing
578 // padradius - radius - pads to be used for smoothing
579 // ROCoutlier - map of outliers - pads not to be used for local smoothing
580 // robust - robust method of fitting - (much slower)
581 // chi2Threshold: Threshold for chi2 when EvalRobust is called
582 // robustFraction: Fraction of data that will be used in EvalRobust
584 AliTPCCalROC * ROCfitted = new AliTPCCalROC(fSector);
585 TLinearFitter fitterQ(6,"hyp5");
586 // TLinearFitter fitterQ(6,"x0++x1++x2++x3++x4++x5");
587 fitterQ.StoreData(kTRUE);
588 for (UInt_t row=0; row < GetNrows(); row++) {
589 //std::cout << "Entering row " << row << " of " << GetNrows() << " @ sector "<< fSector << " for local fitting... "<< std::endl;
590 for (UInt_t pad=0; pad < GetNPads(row); pad++)
591 ROCfitted->SetValue(row, pad, GetNeighbourhoodValue(&fitterQ, row, pad, rowRadius, padRadius, ROCoutliers, robust, chi2Threshold, robustFraction));
597 Double_t AliTPCCalROC::GetNeighbourhoodValue(TLinearFitter* fitterQ, Int_t row, Int_t pad, Int_t rRadius, Int_t pRadius, AliTPCCalROC* ROCoutliers, Bool_t robust, Double_t chi2Threshold, Double_t robustFraction) {
599 // AliTPCCalROC::GetNeighbourhoodValue - smoothing - PRIVATE
600 // in this function the fit for LocalFit is done
603 fitterQ->ClearPoints();
604 TVectorD fitParam(6);
608 Float_t lPad[3] = {0};
609 Float_t localXY[3] = {0};
611 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
612 tpcROCinstance->GetPositionLocal(fSector, row, pad, lPad); // calculate position lPad by pad and row number
614 TArrayI *neighbourhoodRows = 0;
615 TArrayI *neighbourhoodPads = 0;
617 //std::cerr << "Trying to get neighbourhood for row " << row << ", pad " << pad << std::endl;
618 GetNeighbourhood(neighbourhoodRows, neighbourhoodPads, row, pad, rRadius, pRadius);
619 //std::cerr << "Got neighbourhood for row " << row << ", pad " << pad << std::endl;
622 for (Int_t i=0; i < (2*rRadius+1)*(2*pRadius+1); i++) {
623 r = neighbourhoodRows->At(i);
624 p = neighbourhoodPads->At(i);
625 if (r == -1 || p == -1) continue; // window is bigger than ROC
626 tpcROCinstance->GetPositionLocal(fSector, r, p, localXY); // calculate position localXY by pad and row number
627 dlx = lPad[0] - localXY[0];
628 dly = lPad[1] - localXY[1];
635 if (!ROCoutliers || ROCoutliers->GetValue(r,p) != 1) {
636 fitterQ->AddPoint(&xx[1], GetValue(r, p), 1);
641 delete neighbourhoodRows;
642 delete neighbourhoodPads;
644 if (npoints < 0.5 * ((2*rRadius+1)*(2*pRadius+1)) ) {
645 // std::cerr << "Too few data points for fitting @ row " << row << ", pad " << pad << " in sector " << fSector << std::endl;
646 return 0.; // for diagnostic
649 fitterQ->GetParameters(fitParam);
651 if (robust) chi2Q = fitterQ->GetChisquare()/(npoints-6.);
652 //if (robust) chi2Q = fitterQ->GetChisquare()/(npoints-6.);
653 if (robust && chi2Q > chi2Threshold) {
654 //std::cout << "robust fitter called... " << std::endl;
655 fitterQ->EvalRobust(robustFraction);
656 fitterQ->GetParameters(fitParam);
658 Double_t value = fitParam[0];
660 //if (value < 0) std::cerr << "negative fit-value " << value << " in sector "<< this->fSector << " @ row: " << row << " and pad: " << pad << ", with fitter Chi2 = " << chi2Q << std::endl;
667 void AliTPCCalROC::GetNeighbourhood(TArrayI* &rowArray, TArrayI* &padArray, Int_t row, Int_t pad, Int_t rRadius, Int_t pRadius) {
669 // AliTPCCalROC::GetNeighbourhood - PRIVATE
670 // in this function the window for LocalFit is determined
672 rowArray = new TArrayI((2*rRadius+1)*(2*pRadius+1));
673 padArray = new TArrayI((2*rRadius+1)*(2*pRadius+1));
675 Int_t rmin = row - rRadius;
676 UInt_t rmax = row + rRadius;
678 // if window goes out of ROC
683 if (rmax >= GetNrows()) {
684 rmin = rmin - (rmax - GetNrows()+1);
685 rmax = GetNrows() - 1;
686 if (rmin < 0 ) rmin = 0; // if the window is bigger than the ROC
692 for (UInt_t r = rmin; r <= rmax; r++) {
693 pmin = pad - pRadius;
694 pmax = pad + pRadius;
699 if (pmax >= (Int_t)GetNPads(r)) {
700 pmin = pmin - (pmax - GetNPads(r)+1);
701 pmax = GetNPads(r) - 1;
702 if (pmin < 0 ) pmin = 0; // if the window is bigger than the ROC
704 for (Int_t p = pmin; p <= pmax; p++) {
710 for (Int_t j = i; j < rowArray->GetSize(); j++){ // unused padArray-entries, in the case that the window is bigger than the ROC
711 //std::cout << "trying to write -1" << std::endl;
714 //std::cout << "writing -1" << std::endl;
720 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){
722 // Makes a GlobalFit for the given secotr and return fit-parameters, covariance and chi2
723 // The origin of the fit function is the center of the ROC!
724 // fitType == 0: fit plane function
725 // fitType == 1: fit parabolic function
726 // ROCoutliers - pads with value !=0 are not used in fitting procedure
727 // chi2Threshold: Threshold for chi2 when EvalRobust is called
728 // robustFraction: Fraction of data that will be used in EvalRobust
729 // err: error of the data points
731 TLinearFitter* fitterG = 0;
735 fitterG = new TLinearFitter (6,"x0++x1++x2++x3++x4++x5");
736 fitParam.ResizeTo(6);
737 covMatrix.ResizeTo(6,6);
739 fitterG = new TLinearFitter(3,"x0++x1++x2");
740 fitParam.ResizeTo(3);
741 covMatrix.ResizeTo(3,3);
743 fitterG->StoreData(kTRUE);
744 fitterG->ClearPoints();
748 Float_t centerPad[3] = {0};
749 Float_t localXY[3] = {0};
751 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
752 tpcROCinstance->GetPositionLocal(fSector, GetNrows()/2, GetNPads(GetNrows()/2)/2, centerPad); // calculate center of ROC
754 // loop over all channels and read data into fitterG
755 for (UInt_t irow = 0; irow < GetNrows(); irow++) {
756 for (UInt_t ipad = 0; ipad < GetNPads(irow); ipad++) {
758 if (ROCoutliers && ROCoutliers->GetValue(irow, ipad) != 0) continue;
759 tpcROCinstance->GetPositionLocal(fSector, irow, ipad, localXY); // calculate position localXY by pad and row number
760 dlx = localXY[0] - centerPad[0];
761 dly = localXY[1] - centerPad[1];
769 fitterG->AddPoint(xx, GetValue(irow, ipad), err);
772 if(npoints>10) { // make sure there is something to fit
774 fitterG->GetParameters(fitParam);
775 fitterG->GetCovarianceMatrix(covMatrix);
777 chi2 = fitterG->GetChisquare()/(npoints-6.);
778 else chi2 = fitterG->GetChisquare()/(npoints-3.);
779 if (robust && chi2 > chi2Threshold) {
780 // std::cout << "robust fitter called... " << std::endl;
781 fitterG->EvalRobust(robustFraction);
782 fitterG->GetParameters(fitParam);
785 // set parameteres to 0
786 Int_t nParameters = 3;
790 for(Int_t i = 0; i < nParameters; i++)
798 AliTPCCalROC* AliTPCCalROC::CreateGlobalFitCalROC(TVectorD &fitParam, Int_t sector){
800 // Create ROC with global fit parameters
801 // The origin of the fit function is the center of the ROC!
802 // loop over all channels, write fit values into new ROC and return it
805 Float_t centerPad[3] = {0};
806 Float_t localXY[3] = {0};
807 AliTPCCalROC * ROCfitted = new AliTPCCalROC(sector);
808 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
809 tpcROCinstance->GetPositionLocal(sector, ROCfitted->GetNrows()/2, ROCfitted->GetNPads(ROCfitted->GetNrows()/2)/2, centerPad); // calculate center of ROC
811 if (fitParam.GetNoElements() == 6) fitType = 1;
814 if (fitType == 1) { // parabolic fit
815 for (UInt_t irow = 0; irow < ROCfitted->GetNrows(); irow++) {
816 for (UInt_t ipad = 0; ipad < ROCfitted->GetNPads(irow); ipad++) {
817 tpcROCinstance->GetPositionLocal(sector, irow, ipad, localXY); // calculate position localXY by pad and row number
818 dlx = localXY[0] - centerPad[0];
819 dly = localXY[1] - centerPad[1];
820 value = fitParam[0] + fitParam[1]*dlx + fitParam[2]*dly + fitParam[3]*dlx*dlx + fitParam[4]*dly*dly + fitParam[5]*dlx*dly;
821 ROCfitted->SetValue(irow, ipad, value);
826 for (UInt_t irow = 0; irow < ROCfitted->GetNrows(); irow++) {
827 for (UInt_t ipad = 0; ipad < ROCfitted->GetNPads(irow); ipad++) {
828 tpcROCinstance->GetPositionLocal(sector, irow, ipad, localXY); // calculate position localXY by pad and row number
829 dlx = localXY[0] - centerPad[0];
830 dly = localXY[1] - centerPad[1];
831 value = fitParam[0] + fitParam[1]*dlx + fitParam[2]*dly;
832 ROCfitted->SetValue(irow, ipad, value);