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];
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 fkIndexes = 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
106 if (this == ¶m) return (*this);
107 fSector = param.fSector;
108 fNChannels = AliTPCROC::Instance()->GetNChannels(fSector);
109 fNRows = AliTPCROC::Instance()->GetNRows(fSector);
110 fkIndexes = AliTPCROC::Instance()->GetRowIndexes(fSector);
112 if (fData) delete [] fData;
113 fData = new Float_t[fNChannels];
114 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata] = param.fData[idata];
119 //_____________________________________________________________________________
120 AliTPCCalROC::~AliTPCCalROC()
123 // AliTPCCalROC destructor
133 void AliTPCCalROC::Streamer(TBuffer &R__b)
136 // Stream an object of class AliTPCCalROC.
138 if (R__b.IsReading()) {
139 AliTPCCalROC::Class()->ReadBuffer(R__b, this);
140 fkIndexes = AliTPCROC::Instance()->GetRowIndexes(fSector);
142 AliTPCCalROC::Class()->WriteBuffer(R__b,this);
149 void AliTPCCalROC::Add(Float_t c1){
151 // add c1 to each channel of the ROC
153 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata]+=c1;
157 void AliTPCCalROC::Multiply(Float_t c1){
159 // multiply each channel of the ROC with c1
161 for (UInt_t idata = 0; idata< fNChannels; idata++) fData[idata]*=c1;
165 void AliTPCCalROC::Add(const AliTPCCalROC * roc, Double_t c1){
167 // multiply AliTPCCalROC roc by c1 and add each channel to the coresponing channel in the ROC
170 for (UInt_t idata = 0; idata< fNChannels; idata++){
171 fData[idata]+=roc->fData[idata]*c1;
176 void AliTPCCalROC::Multiply(const AliTPCCalROC* roc) {
178 // multiply each channel of the ROC with the coresponding channel of 'roc'
181 for (UInt_t idata = 0; idata< fNChannels; idata++){
182 fData[idata]*=roc->fData[idata];
187 void AliTPCCalROC::Divide(const AliTPCCalROC* roc) {
189 // divide each channel of the ROC by the coresponding channel of 'roc'
192 Float_t kEpsilon=0.00000000000000001;
193 for (UInt_t idata = 0; idata< fNChannels; idata++){
194 if (TMath::Abs(roc->fData[idata])>kEpsilon)
195 fData[idata]/=roc->fData[idata];
199 Double_t AliTPCCalROC::GetMean(AliTPCCalROC *const outlierROC) const {
201 // returns the mean value of the ROC
202 // pads with value != 0 in outlierROC are not used for the calculation
203 // return 0 if no data is accepted by the outlier cuts
205 if (!outlierROC) return TMath::Mean(fNChannels, fData);
206 Double_t *ddata = new Double_t[fNChannels];
208 for (UInt_t i=0;i<fNChannels;i++) {
209 if (!(outlierROC->GetValue(i))) {
210 ddata[nPoints]= fData[i];
216 mean = TMath::Mean(nPoints, ddata);
221 Double_t AliTPCCalROC::GetMedian(AliTPCCalROC *const outlierROC) const {
223 // returns the median value of the ROC
224 // pads with value != 0 in outlierROC are not used for the calculation
225 // return 0 if no data is accepted by the outlier cuts
227 if (!outlierROC) return TMath::Median(fNChannels, fData);
228 Double_t *ddata = new Double_t[fNChannels];
230 for (UInt_t i=0;i<fNChannels;i++) {
231 if (!(outlierROC->GetValue(i))) {
232 ddata[nPoints]= fData[i];
238 median = TMath::Median(nPoints, ddata);
243 Double_t AliTPCCalROC::GetRMS(AliTPCCalROC *const outlierROC) const {
245 // returns the RMS value of the ROC
246 // pads with value != 0 in outlierROC are not used for the calculation
247 // return 0 if no data is accepted by the outlier cuts
249 if (!outlierROC) return TMath::RMS(fNChannels, fData);
250 Double_t *ddata = new Double_t[fNChannels];
252 for (UInt_t i=0;i<fNChannels;i++) {
253 if (!(outlierROC->GetValue(i))) {
254 ddata[nPoints]= fData[i];
260 rms = TMath::RMS(nPoints, ddata);
265 Double_t AliTPCCalROC::GetLTM(Double_t *const sigma, Double_t fraction, AliTPCCalROC *const outlierROC){
267 // returns the LTM and sigma
268 // pads with value != 0 in outlierROC are not used for the calculation
269 // return 0 if no data is accepted by the outlier cuts
271 Double_t *ddata = new Double_t[fNChannels];
273 for (UInt_t i=0;i<fNChannels;i++) {
274 if (!outlierROC || !(outlierROC->GetValue(i))) {
275 ddata[nPoints]= fData[i];
280 Double_t ltm =0, lsigma=0;
282 Int_t hh = TMath::Min(TMath::Nint(fraction *nPoints), Int_t(nPoints));
283 AliMathBase::EvaluateUni(nPoints,ddata, ltm, lsigma, hh);
284 if (sigma) *sigma=lsigma;
291 TH1F * AliTPCCalROC::MakeHisto1D(Float_t min, Float_t max,Int_t type){
294 // type -1 = user defined range
295 // 0 = nsigma cut nsigma=min
296 // 1 = delta cut around median delta=min
301 Float_t mean = GetMean();
302 Float_t sigma = GetRMS();
303 Float_t nsigma = TMath::Abs(min);
304 min = mean-nsigma*sigma;
305 max = mean+nsigma*sigma;
309 Float_t mean = GetMedian();
316 // LTM mean +- nsigma
319 Float_t mean = GetLTM(&sigma,max);
325 TString name=Form("%s ROC 1D%d",GetTitle(),fSector);
326 TH1F * his = new TH1F(name.Data(),name.Data(),100, min,max);
327 for (UInt_t irow=0; irow<fNRows; irow++){
328 UInt_t npads = (Int_t)GetNPads(irow);
329 for (UInt_t ipad=0; ipad<npads; ipad++){
330 his->Fill(GetValue(irow,ipad));
338 TH2F * AliTPCCalROC::MakeHisto2D(Float_t min, Float_t max,Int_t type){
341 // type -1 = user defined range
342 // 0 = nsigma cut nsigma=min
343 // 1 = delta cut around median delta=min
348 Float_t mean = GetMean();
349 Float_t sigma = GetRMS();
350 Float_t nsigma = TMath::Abs(min);
351 min = mean-nsigma*sigma;
352 max = mean+nsigma*sigma;
356 Float_t mean = GetMedian();
363 Float_t mean = GetLTM(&sigma,max);
371 for (UInt_t irow=0; irow<fNRows; irow++){
372 if (GetNPads(irow)>maxPad) maxPad = GetNPads(irow);
375 TString name=Form("%s ROC%d",GetTitle(),fSector);
376 TH2F * his = new TH2F(name.Data(),name.Data(),fNRows+10,-5, fNRows+5, maxPad+10, -(Int_t(maxPad/2))-5, maxPad/2+5);
377 for (UInt_t irow=0; irow<fNRows; irow++){
378 UInt_t npads = (Int_t)GetNPads(irow);
379 for (UInt_t ipad=0; ipad<npads; ipad++){
380 his->Fill(irow+0.5,Int_t(ipad)-Int_t(npads/2)+0.5,GetValue(irow,ipad));
383 his->SetMaximum(max);
384 his->SetMinimum(min);
388 TH2F * AliTPCCalROC::MakeHistoOutliers(Float_t delta, Float_t fraction, Int_t type){
390 // Make Histogram with outliers
391 // mode = 0 - sigma cut used
392 // mode = 1 - absolute cut used
393 // fraction - fraction of values used to define sigma
394 // delta - in mode 0 - nsigma cut
395 // mode 1 - delta cut
398 Float_t mean = GetLTM(&sigma,fraction);
399 if (type==0) delta*=sigma;
401 for (UInt_t irow=0; irow<fNRows; irow++){
402 if (GetNPads(irow)>maxPad) maxPad = GetNPads(irow);
405 TString name=Form("%s ROC Outliers%d",GetTitle(),fSector);
406 TH2F * his = new TH2F(name.Data(),name.Data(),fNRows+10,-5, fNRows+5, maxPad+10, -(Int_t(maxPad/2))-5, maxPad/2+5);
407 for (UInt_t irow=0; irow<fNRows; irow++){
408 UInt_t npads = (Int_t)GetNPads(irow);
409 for (UInt_t ipad=0; ipad<npads; ipad++){
410 if (TMath::Abs(GetValue(irow,ipad)-mean)>delta)
411 his->Fill(irow+0.5,Int_t(ipad)-Int_t(npads/2)+0.5,1);
419 void AliTPCCalROC::Draw(Option_t* opt){
421 // create histogram with values and draw it
426 if (option.Contains("1D")){
439 void AliTPCCalROC::Test() {
441 // example function to show functionality and test AliTPCCalROC
444 Float_t kEpsilon=0.00001;
446 // create CalROC with defined values
447 AliTPCCalROC roc0(0);
448 for (UInt_t irow = 0; irow <roc0.GetNrows(); irow++){
449 for (UInt_t ipad = 0; ipad <roc0.GetNPads(irow); ipad++){
450 Float_t value = irow+ipad/1000.;
451 roc0.SetValue(irow,ipad,value);
455 // copy CalROC, readout values and compare to original
456 AliTPCCalROC roc1(roc0);
457 for (UInt_t irow = 0; irow <roc1.GetNrows(); irow++){
458 for (UInt_t ipad = 0; ipad <roc1.GetNPads(irow); ipad++){
459 Float_t value = irow+ipad/1000.;
460 if (roc1.GetValue(irow,ipad)!=value){
461 printf("Read/Write error\trow=%d\tpad=%d\n",irow,ipad);
466 // write original CalROC to file
467 TFile f("calcTest.root","recreate");
469 AliTPCCalROC * roc2 = (AliTPCCalROC*)f.Get("Roc0");
472 // read CalROC from file and compare to original CalROC
473 for (UInt_t irow = 0; irow <roc0.GetNrows(); irow++){
474 if (roc0.GetNPads(irow)!=roc2->GetNPads(irow))
475 printf("NPads - Read/Write error\trow=%d\n",irow);
476 for (UInt_t ipad = 0; ipad <roc1.GetNPads(irow); ipad++){
477 Float_t value = irow+ipad/1000.;
478 if (roc2->GetValue(irow,ipad)!=value){
479 printf("Read/Write error\trow=%d\tpad=%d\n",irow,ipad);
489 AliTPCCalROC roc3(roc0);
491 for (UInt_t irow = 0; irow <roc3.GetNrows(); irow++){
492 for (UInt_t ipad = 0; ipad <roc3.GetNPads(irow); ipad++){
493 Float_t value = irow+ipad/1000. + 1.5;
494 if (TMath::Abs(roc3.GetValue(irow,ipad)-value) > kEpsilon){
495 printf("Add constant - error\trow=%d\tpad=%d\n",irow,ipad);
500 // Add another CalROC
501 AliTPCCalROC roc4(roc0);
502 roc4.Add(&roc0, -1.5);
503 for (UInt_t irow = 0; irow <roc4.GetNrows(); irow++){
504 for (UInt_t ipad = 0; ipad <roc4.GetNPads(irow); ipad++){
505 Float_t value = irow+ipad/1000. - 1.5 * (irow+ipad/1000.);
506 if (TMath::Abs(roc4.GetValue(irow,ipad)-value) > kEpsilon){
507 printf("Add CalROC - error\trow=%d\tpad=%d\n",irow,ipad);
512 // Multiply with constant
513 AliTPCCalROC roc5(roc0);
515 for (UInt_t irow = 0; irow <roc5.GetNrows(); irow++){
516 for (UInt_t ipad = 0; ipad <roc5.GetNPads(irow); ipad++){
517 Float_t value = (irow+ipad/1000.) * (-1.4);
518 if (TMath::Abs(roc5.GetValue(irow,ipad)-value) > kEpsilon){
519 printf("Multiply with constant - error\trow=%d\tpad=%d\n",irow,ipad);
524 // Multiply another CalROC
525 AliTPCCalROC roc6(roc0);
526 roc6.Multiply(&roc0);
527 for (UInt_t irow = 0; irow <roc6.GetNrows(); irow++){
528 for (UInt_t ipad = 0; ipad <roc6.GetNPads(irow); ipad++){
529 Float_t value = (irow+ipad/1000.) * (irow+ipad/1000.);
530 if (TMath::Abs(roc6.GetValue(irow,ipad)-value) > kEpsilon*100){
531 printf("Multiply with CalROC - error\trow=%d\tpad=%d\n",irow,ipad);
538 AliTPCCalROC roc7(roc0);
540 for (UInt_t irow = 0; irow <roc7.GetNrows(); irow++){
541 for (UInt_t ipad = 0; ipad <roc7.GetNPads(irow); ipad++){
543 if (irow+ipad == 0) value = roc0.GetValue(irow,ipad);
544 if (TMath::Abs(roc7.GetValue(irow,ipad)-value) > kEpsilon){
545 printf("Multiply with CalROC - error\trow=%d\tpad=%d\n",irow,ipad);
554 // create CalROC with defined values
556 AliTPCCalROC sroc0(0);
557 for (UInt_t ichannel = 0; ichannel < sroc0.GetNchannels(); ichannel++){
558 Float_t value = rnd.Gaus(10., 2.);
559 sroc0.SetValue(ichannel,value);
562 printf("Mean (should be close to 10): %f\n", sroc0.GetMean());
563 printf("RMS (should be close to 2): %f\n", sroc0.GetRMS());
564 printf("Median (should be close to 10): %f\n", sroc0.GetMedian());
565 printf("LTM (should be close to 10): %f\n", sroc0.GetLTM());
567 //AliTPCCalROC* sroc1 = sroc0.LocalFit(4, 8);
571 // std::cout << TMath::Abs(roc5.GetValue(irow,ipad)-value) << std::endl;
575 AliTPCCalROC * AliTPCCalROC::LocalFit(Int_t rowRadius, Int_t padRadius, AliTPCCalROC* ROCoutliers, Bool_t robust, Double_t chi2Threshold, Double_t robustFraction) {
577 // MakeLocalFit - smoothing
578 // returns a AliTPCCalROC with smoothed data
579 // rowRadius and padRadius specifies a window around a given pad.
580 // The data of this window are fitted with a parabolic function.
581 // This function is evaluated at the pad's position.
582 // At the edges the window is shifted, so that the specified pad is not anymore in the center of the window.
583 // rowRadius - radius - rows to be used for smoothing
584 // padradius - radius - pads to be used for smoothing
585 // ROCoutlier - map of outliers - pads not to be used for local smoothing
586 // robust - robust method of fitting - (much slower)
587 // chi2Threshold: Threshold for chi2 when EvalRobust is called
588 // robustFraction: Fraction of data that will be used in EvalRobust
590 AliTPCCalROC * xROCfitted = new AliTPCCalROC(fSector);
591 TLinearFitter fitterQ(6,"hyp5");
592 // TLinearFitter fitterQ(6,"x0++x1++x2++x3++x4++x5");
593 fitterQ.StoreData(kTRUE);
594 for (UInt_t row=0; row < GetNrows(); row++) {
595 //std::cout << "Entering row " << row << " of " << GetNrows() << " @ sector "<< fSector << " for local fitting... "<< std::endl;
596 for (UInt_t pad=0; pad < GetNPads(row); pad++)
597 xROCfitted->SetValue(row, pad, GetNeighbourhoodValue(&fitterQ, row, pad, rowRadius, padRadius, ROCoutliers, robust, chi2Threshold, robustFraction));
603 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) {
605 // AliTPCCalROC::GetNeighbourhoodValue - smoothing - PRIVATE
606 // in this function the fit for LocalFit is done
609 fitterQ->ClearPoints();
610 TVectorD fitParam(6);
614 Float_t lPad[3] = {0};
615 Float_t localXY[3] = {0};
617 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
618 tpcROCinstance->GetPositionLocal(fSector, row, pad, lPad); // calculate position lPad by pad and row number
620 TArrayI *neighbourhoodRows = 0;
621 TArrayI *neighbourhoodPads = 0;
623 //std::cerr << "Trying to get neighbourhood for row " << row << ", pad " << pad << std::endl;
624 GetNeighbourhood(neighbourhoodRows, neighbourhoodPads, row, pad, rRadius, pRadius);
625 //std::cerr << "Got neighbourhood for row " << row << ", pad " << pad << std::endl;
628 for (Int_t i=0; i < (2*rRadius+1)*(2*pRadius+1); i++) {
629 r = neighbourhoodRows->At(i);
630 p = neighbourhoodPads->At(i);
631 if (r == -1 || p == -1) continue; // window is bigger than ROC
632 tpcROCinstance->GetPositionLocal(fSector, r, p, localXY); // calculate position localXY by pad and row number
633 dlx = lPad[0] - localXY[0];
634 dly = lPad[1] - localXY[1];
641 if (!ROCoutliers || ROCoutliers->GetValue(r,p) != 1) {
642 fitterQ->AddPoint(&xx[1], GetValue(r, p), 1);
647 delete neighbourhoodRows;
648 delete neighbourhoodPads;
650 if (npoints < 0.5 * ((2*rRadius+1)*(2*pRadius+1)) ) {
651 // std::cerr << "Too few data points for fitting @ row " << row << ", pad " << pad << " in sector " << fSector << std::endl;
652 return 0.; // for diagnostic
655 fitterQ->GetParameters(fitParam);
657 if (robust) chi2Q = fitterQ->GetChisquare()/(npoints-6.);
658 //if (robust) chi2Q = fitterQ->GetChisquare()/(npoints-6.);
659 if (robust && chi2Q > chi2Threshold) {
660 //std::cout << "robust fitter called... " << std::endl;
661 fitterQ->EvalRobust(robustFraction);
662 fitterQ->GetParameters(fitParam);
664 Double_t value = fitParam[0];
666 //if (value < 0) std::cerr << "negative fit-value " << value << " in sector "<< this->fSector << " @ row: " << row << " and pad: " << pad << ", with fitter Chi2 = " << chi2Q << std::endl;
673 void AliTPCCalROC::GetNeighbourhood(TArrayI* &rowArray, TArrayI* &padArray, Int_t row, Int_t pad, Int_t rRadius, Int_t pRadius) {
675 // AliTPCCalROC::GetNeighbourhood - PRIVATE
676 // in this function the window for LocalFit is determined
678 rowArray = new TArrayI((2*rRadius+1)*(2*pRadius+1));
679 padArray = new TArrayI((2*rRadius+1)*(2*pRadius+1));
681 Int_t rmin = row - rRadius;
682 UInt_t rmax = row + rRadius;
684 // if window goes out of ROC
689 if (rmax >= GetNrows()) {
690 rmin = rmin - (rmax - GetNrows()+1);
691 rmax = GetNrows() - 1;
692 if (rmin < 0 ) rmin = 0; // if the window is bigger than the ROC
698 for (UInt_t r = rmin; r <= rmax; r++) {
699 pmin = pad - pRadius;
700 pmax = pad + pRadius;
705 if (pmax >= (Int_t)GetNPads(r)) {
706 pmin = pmin - (pmax - GetNPads(r)+1);
707 pmax = GetNPads(r) - 1;
708 if (pmin < 0 ) pmin = 0; // if the window is bigger than the ROC
710 for (Int_t p = pmin; p <= pmax; p++) {
716 for (Int_t j = i; j < rowArray->GetSize(); j++){ // unused padArray-entries, in the case that the window is bigger than the ROC
717 //std::cout << "trying to write -1" << std::endl;
720 //std::cout << "writing -1" << std::endl;
726 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){
728 // Makes a GlobalFit for the given secotr and return fit-parameters, covariance and chi2
729 // The origin of the fit function is the center of the ROC!
730 // fitType == 0: fit plane function
731 // fitType == 1: fit parabolic function
732 // ROCoutliers - pads with value !=0 are not used in fitting procedure
733 // chi2Threshold: Threshold for chi2 when EvalRobust is called
734 // robustFraction: Fraction of data that will be used in EvalRobust
735 // err: error of the data points
737 TLinearFitter* fitterG = 0;
741 fitterG = new TLinearFitter (6,"x0++x1++x2++x3++x4++x5");
742 fitParam.ResizeTo(6);
743 covMatrix.ResizeTo(6,6);
745 fitterG = new TLinearFitter(3,"x0++x1++x2");
746 fitParam.ResizeTo(3);
747 covMatrix.ResizeTo(3,3);
749 fitterG->StoreData(kTRUE);
750 fitterG->ClearPoints();
754 Float_t centerPad[3] = {0};
755 Float_t localXY[3] = {0};
757 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
758 tpcROCinstance->GetPositionLocal(fSector, GetNrows()/2, GetNPads(GetNrows()/2)/2, centerPad); // calculate center of ROC
760 // loop over all channels and read data into fitterG
761 for (UInt_t irow = 0; irow < GetNrows(); irow++) {
762 for (UInt_t ipad = 0; ipad < GetNPads(irow); ipad++) {
764 if (ROCoutliers && ROCoutliers->GetValue(irow, ipad) != 0) continue;
765 tpcROCinstance->GetPositionLocal(fSector, irow, ipad, localXY); // calculate position localXY by pad and row number
766 dlx = localXY[0] - centerPad[0];
767 dly = localXY[1] - centerPad[1];
775 fitterG->AddPoint(xx, GetValue(irow, ipad), err);
778 if(npoints>10) { // make sure there is something to fit
780 fitterG->GetParameters(fitParam);
781 fitterG->GetCovarianceMatrix(covMatrix);
783 chi2 = fitterG->GetChisquare()/(npoints-6.);
784 else chi2 = fitterG->GetChisquare()/(npoints-3.);
785 if (robust && chi2 > chi2Threshold) {
786 // std::cout << "robust fitter called... " << std::endl;
787 fitterG->EvalRobust(robustFraction);
788 fitterG->GetParameters(fitParam);
791 // set parameteres to 0
792 Int_t nParameters = 3;
796 for(Int_t i = 0; i < nParameters; i++)
804 AliTPCCalROC* AliTPCCalROC::CreateGlobalFitCalROC(TVectorD &fitParam, Int_t sector){
806 // Create ROC with global fit parameters
807 // The origin of the fit function is the center of the ROC!
808 // loop over all channels, write fit values into new ROC and return it
811 Float_t centerPad[3] = {0};
812 Float_t localXY[3] = {0};
813 AliTPCCalROC * xROCfitted = new AliTPCCalROC(sector);
814 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
815 tpcROCinstance->GetPositionLocal(sector, xROCfitted->GetNrows()/2, xROCfitted->GetNPads(xROCfitted->GetNrows()/2)/2, centerPad); // calculate center of ROC
817 if (fitParam.GetNoElements() == 6) fitType = 1;
820 if (fitType == 1) { // parabolic fit
821 for (UInt_t irow = 0; irow < xROCfitted->GetNrows(); irow++) {
822 for (UInt_t ipad = 0; ipad < xROCfitted->GetNPads(irow); ipad++) {
823 tpcROCinstance->GetPositionLocal(sector, irow, ipad, localXY); // calculate position localXY by pad and row number
824 dlx = localXY[0] - centerPad[0];
825 dly = localXY[1] - centerPad[1];
826 value = fitParam[0] + fitParam[1]*dlx + fitParam[2]*dly + fitParam[3]*dlx*dlx + fitParam[4]*dly*dly + fitParam[5]*dlx*dly;
827 xROCfitted->SetValue(irow, ipad, value);
832 for (UInt_t irow = 0; irow < xROCfitted->GetNrows(); irow++) {
833 for (UInt_t ipad = 0; ipad < xROCfitted->GetNPads(irow); ipad++) {
834 tpcROCinstance->GetPositionLocal(sector, irow, ipad, localXY); // calculate position localXY by pad and row number
835 dlx = localXY[0] - centerPad[0];
836 dly = localXY[1] - centerPad[1];
837 value = fitParam[0] + fitParam[1]*dlx + fitParam[2]*dly;
838 xROCfitted->SetValue(irow, ipad, value);