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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // TPC calibration class for parameters which saved per pad //
22 ///////////////////////////////////////////////////////////////////////////////
24 #include "AliTPCCalPad.h"
25 #include "AliTPCCalROC.h"
26 #include <TObjArray.h>
31 ClassImp(AliTPCCalPad)
33 //_____________________________________________________________________________
34 AliTPCCalPad::AliTPCCalPad():TNamed()
37 // AliTPCCalPad default constructor
40 for (Int_t isec = 0; isec < kNsec; isec++) {
46 //_____________________________________________________________________________
47 AliTPCCalPad::AliTPCCalPad(const Text_t *name, const Text_t *title)
51 // AliTPCCalPad constructor
53 for (Int_t isec = 0; isec < kNsec; isec++) {
54 fROC[isec] = new AliTPCCalROC(isec);
59 //_____________________________________________________________________________
60 AliTPCCalPad::AliTPCCalPad(const AliTPCCalPad &c):TNamed(c)
63 // AliTPCCalPad copy constructor
66 ((AliTPCCalPad &) c).Copy(*this);
70 //_____________________________________________________________________________
71 AliTPCCalPad::AliTPCCalPad(TObjArray * array):TNamed()
74 // AliTPCCalPad default constructor
77 for (Int_t isec = 0; isec < kNsec; isec++) {
78 fROC[isec] = (AliTPCCalROC *)array->At(isec);
84 ///_____________________________________________________________________________
85 AliTPCCalPad::~AliTPCCalPad()
88 // AliTPCCalPad destructor
91 for (Int_t isec = 0; isec < kNsec; isec++) {
100 //_____________________________________________________________________________
101 AliTPCCalPad &AliTPCCalPad::operator=(const AliTPCCalPad &c)
104 // Assignment operator
107 if (this != &c) ((AliTPCCalPad &) c).Copy(*this);
112 //_____________________________________________________________________________
113 void AliTPCCalPad::Copy(TObject &c) const
119 for (Int_t isec = 0; isec < kNsec; isec++) {
121 fROC[isec]->Copy(*((AliTPCCalPad &) c).fROC[isec]);
127 //_____________________________________________________________________________
128 void AliTPCCalPad::Add(Float_t c1)
131 // add constant for all channels of all ROCs
134 for (Int_t isec = 0; isec < kNsec; isec++) {
141 //_____________________________________________________________________________
142 void AliTPCCalPad::Multiply(Float_t c1)
145 // multiply constant for all channels of all ROCs
147 for (Int_t isec = 0; isec < kNsec; isec++) {
149 fROC[isec]->Multiply(c1);
154 //_____________________________________________________________________________
155 void AliTPCCalPad::Add(const AliTPCCalPad * pad, Double_t c1)
158 // add calpad channel by channel multiplied by c1 - all ROCs
160 for (Int_t isec = 0; isec < kNsec; isec++) {
162 fROC[isec]->Add(pad->GetCalROC(isec),c1);
167 //_____________________________________________________________________________
168 void AliTPCCalPad::Multiply(const AliTPCCalPad * pad)
171 // multiply calpad channel by channel - all ROCs
173 for (Int_t isec = 0; isec < kNsec; isec++) {
175 fROC[isec]->Multiply(pad->GetCalROC(isec));
180 //_____________________________________________________________________________
181 void AliTPCCalPad::Divide(const AliTPCCalPad * pad)
184 // divide calpad channel by channel - all ROCs
186 for (Int_t isec = 0; isec < kNsec; isec++) {
188 fROC[isec]->Divide(pad->GetCalROC(isec));
193 //_____________________________________________________________________________
194 TGraph * AliTPCCalPad::MakeGraph(Int_t type, Float_t ratio){
200 for (Int_t i=0;i<72;i++) if (fROC[i]) npoints++;
201 TGraph * graph = new TGraph(npoints);
203 for (Int_t isec=0;isec<72;isec++){
204 if (!fROC[isec]) continue;
205 if (type==0) graph->SetPoint(npoints,isec,fROC[isec]->GetMean());
206 if (type==1) graph->SetPoint(npoints,isec,fROC[isec]->GetMedian());
207 if (type==2) graph->SetPoint(npoints,isec,fROC[isec]->GetLTM(0,ratio));
211 graph->GetXaxis()->SetTitle("Sector");
213 graph->GetYaxis()->SetTitle("Mean");
214 graph->SetMarkerStyle(22);
217 graph->GetYaxis()->SetTitle("Median");
218 graph->SetMarkerStyle(22);
221 graph->GetYaxis()->SetTitle(Form("Mean%f",ratio));
222 graph->SetMarkerStyle(24);
228 //_____________________________________________________________________________
229 Double_t AliTPCCalPad::GetMeanRMS(Double_t &rms)
232 // Calculate mean an RMS of all rocs
234 Double_t sum = 0, sum2 = 0, n=0, val=0;
235 for (Int_t isec = 0; isec < kNsec; isec++) {
236 AliTPCCalROC *calRoc = fROC[isec];
238 for (UInt_t irow=0; irow<calRoc->GetNrows(); irow++){
239 for (UInt_t ipad=0; ipad<calRoc->GetNPads(irow); ipad++){
240 val = calRoc->GetValue(irow,ipad);
250 Double_t mean = sum*n1;
251 rms = TMath::Sqrt(TMath::Abs(sum2*n1-mean*mean));
256 //_____________________________________________________________________________
257 Double_t AliTPCCalPad::GetMean()
260 // return mean of the mean of all ROCs
264 for (Int_t isec = 0; isec < kNsec; isec++) {
265 AliTPCCalROC *calRoc = fROC[isec];
267 arr[n] = calRoc->GetMean();
271 return TMath::Mean(n,arr);
274 //_____________________________________________________________________________
275 Double_t AliTPCCalPad::GetRMS()
278 // return mean of the RMS of all ROCs
282 for (Int_t isec = 0; isec < kNsec; isec++) {
283 AliTPCCalROC *calRoc = fROC[isec];
285 arr[n] = calRoc->GetRMS();
289 return TMath::Mean(n,arr);
292 //_____________________________________________________________________________
293 Double_t AliTPCCalPad::GetMedian()
296 // return mean of the median of all ROCs
300 for (Int_t isec = 0; isec < kNsec; isec++) {
301 AliTPCCalROC *calRoc = fROC[isec];
303 arr[n] = calRoc->GetMedian();
307 return TMath::Mean(n,arr);
310 //_____________________________________________________________________________
311 Double_t AliTPCCalPad::GetLTM(Double_t *sigma, Double_t fraction)
314 // return mean of the LTM and sigma of all ROCs
316 Double_t arrm[kNsec];
317 Double_t arrs[kNsec];
321 for (Int_t isec = 0; isec < kNsec; isec++) {
322 AliTPCCalROC *calRoc = fROC[isec];
324 if ( sigma ) sTemp=arrs+n;
325 arrm[n] = calRoc->GetLTM(sTemp,fraction);
329 if ( sigma ) *sigma = TMath::Mean(n,arrs);
330 return TMath::Mean(n,arrm);
333 //_____________________________________________________________________________
334 TH1F * AliTPCCalPad::MakeHisto1D(Float_t min, Float_t max,Int_t type){
337 // type -1 = user defined range
338 // 0 = nsigma cut nsigma=min
342 Float_t mean = GetMean();
343 Float_t sigma = GetRMS();
344 Float_t nsigma = TMath::Abs(min);
345 min = mean-nsigma*sigma;
346 max = mean+nsigma*sigma;
350 Float_t mean = GetMedian();
357 // LTM mean +- nsigma
360 Float_t mean = GetLTM(&sigma,max);
367 sprintf(name,"%s Pad 1D",GetTitle());
368 TH1F * his = new TH1F(name,name,100, min,max);
369 for (Int_t isec = 0; isec < kNsec; isec++) {
371 for (UInt_t irow=0; irow<fROC[isec]->GetNrows(); irow++){
372 UInt_t npads = (Int_t)fROC[isec]->GetNPads(irow);
373 for (UInt_t ipad=0; ipad<npads; ipad++){
374 his->Fill(fROC[isec]->GetValue(irow,ipad));
382 //_____________________________________________________________________________
383 TH2F *AliTPCCalPad::MakeHisto2D(Int_t side){
386 // side - specify the side A = 0 C = 1
387 // type - used types of determination of boundaries in z
388 Float_t kEpsilon = 0.000000000001;
389 TH2F * his = new TH2F(GetName(), GetName(), 250,-250,250,250,-250,250);
390 AliTPCROC * roc = AliTPCROC::Instance();
391 for (Int_t isec=0; isec<72; isec++){
392 if (side==0 && isec%36>=18) continue;
393 if (side>0 && isec%36<18) continue;
395 AliTPCCalROC * calRoc = fROC[isec];
396 for (UInt_t irow=0; irow<calRoc->GetNrows(); irow++)
397 for (UInt_t ipad=0; ipad<calRoc->GetNPads(irow); ipad++)
398 if (TMath::Abs(calRoc->GetValue(irow,ipad))>kEpsilon){
400 roc->GetPositionGlobal(isec,irow,ipad,xyz);
401 Int_t binx = 1+TMath::Nint((xyz[0]+250.)*0.5);
402 Int_t biny = 1+TMath::Nint((xyz[1]+250.)*0.5);
403 Float_t value = calRoc->GetValue(irow,ipad);
404 his->SetBinContent(binx,biny,value);
408 his->SetXTitle("x (cm)");
409 his->SetYTitle("y (cm)");