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 #include "TTreeStream.h"
33 ClassImp(AliTPCCalPad)
35 //_____________________________________________________________________________
36 AliTPCCalPad::AliTPCCalPad():TNamed()
39 // AliTPCCalPad default constructor
42 for (Int_t isec = 0; isec < kNsec; isec++) {
48 //_____________________________________________________________________________
49 AliTPCCalPad::AliTPCCalPad(const Text_t *name, const Text_t *title)
53 // AliTPCCalPad constructor
55 for (Int_t isec = 0; isec < kNsec; isec++) {
56 fROC[isec] = new AliTPCCalROC(isec);
61 //_____________________________________________________________________________
62 AliTPCCalPad::AliTPCCalPad(const AliTPCCalPad &c):TNamed(c)
65 // AliTPCCalPad copy constructor
68 ((AliTPCCalPad &) c).Copy(*this);
72 //_____________________________________________________________________________
73 AliTPCCalPad::AliTPCCalPad(TObjArray * array):TNamed()
76 // AliTPCCalPad default constructor
79 for (Int_t isec = 0; isec < kNsec; isec++) {
80 fROC[isec] = (AliTPCCalROC *)array->At(isec);
86 ///_____________________________________________________________________________
87 AliTPCCalPad::~AliTPCCalPad()
90 // AliTPCCalPad destructor
93 for (Int_t isec = 0; isec < kNsec; isec++) {
102 //_____________________________________________________________________________
103 AliTPCCalPad &AliTPCCalPad::operator=(const AliTPCCalPad &c)
106 // Assignment operator
109 if (this != &c) ((AliTPCCalPad &) c).Copy(*this);
114 //_____________________________________________________________________________
115 void AliTPCCalPad::Copy(TObject &c) const
121 for (Int_t isec = 0; isec < kNsec; isec++) {
123 fROC[isec]->Copy(*((AliTPCCalPad &) c).fROC[isec]);
129 //_____________________________________________________________________________
130 void AliTPCCalPad::Add(Float_t c1)
133 // add constant for all channels of all ROCs
136 for (Int_t isec = 0; isec < kNsec; isec++) {
143 //_____________________________________________________________________________
144 void AliTPCCalPad::Multiply(Float_t c1)
147 // multiply constant for all channels of all ROCs
149 for (Int_t isec = 0; isec < kNsec; isec++) {
151 fROC[isec]->Multiply(c1);
156 //_____________________________________________________________________________
157 void AliTPCCalPad::Add(const AliTPCCalPad * pad, Double_t c1)
160 // add calpad channel by channel multiplied by c1 - all ROCs
162 for (Int_t isec = 0; isec < kNsec; isec++) {
164 fROC[isec]->Add(pad->GetCalROC(isec),c1);
169 //_____________________________________________________________________________
170 void AliTPCCalPad::Multiply(const AliTPCCalPad * pad)
173 // multiply calpad channel by channel - all ROCs
175 for (Int_t isec = 0; isec < kNsec; isec++) {
177 fROC[isec]->Multiply(pad->GetCalROC(isec));
182 //_____________________________________________________________________________
183 void AliTPCCalPad::Divide(const AliTPCCalPad * pad)
186 // divide calpad channel by channel - all ROCs
188 for (Int_t isec = 0; isec < kNsec; isec++) {
190 fROC[isec]->Divide(pad->GetCalROC(isec));
195 //_____________________________________________________________________________
196 TGraph * AliTPCCalPad::MakeGraph(Int_t type, Float_t ratio){
202 for (Int_t i=0;i<72;i++) if (fROC[i]) npoints++;
203 TGraph * graph = new TGraph(npoints);
205 for (Int_t isec=0;isec<72;isec++){
206 if (!fROC[isec]) continue;
207 if (type==0) graph->SetPoint(npoints,isec,fROC[isec]->GetMean());
208 if (type==1) graph->SetPoint(npoints,isec,fROC[isec]->GetMedian());
209 if (type==2) graph->SetPoint(npoints,isec,fROC[isec]->GetLTM(0,ratio));
213 graph->GetXaxis()->SetTitle("Sector");
215 graph->GetYaxis()->SetTitle("Mean");
216 graph->SetMarkerStyle(22);
219 graph->GetYaxis()->SetTitle("Median");
220 graph->SetMarkerStyle(22);
223 graph->GetYaxis()->SetTitle(Form("Mean%f",ratio));
224 graph->SetMarkerStyle(24);
230 //_____________________________________________________________________________
231 Double_t AliTPCCalPad::GetMeanRMS(Double_t &rms)
234 // Calculate mean an RMS of all rocs
236 Double_t sum = 0, sum2 = 0, n=0, val=0;
237 for (Int_t isec = 0; isec < kNsec; isec++) {
238 AliTPCCalROC *calRoc = fROC[isec];
240 for (UInt_t irow=0; irow<calRoc->GetNrows(); irow++){
241 for (UInt_t ipad=0; ipad<calRoc->GetNPads(irow); ipad++){
242 val = calRoc->GetValue(irow,ipad);
252 Double_t mean = sum*n1;
253 rms = TMath::Sqrt(TMath::Abs(sum2*n1-mean*mean));
258 //_____________________________________________________________________________
259 Double_t AliTPCCalPad::GetMean()
262 // return mean of the mean of all ROCs
266 for (Int_t isec = 0; isec < kNsec; isec++) {
267 AliTPCCalROC *calRoc = fROC[isec];
269 arr[n] = calRoc->GetMean();
273 return TMath::Mean(n,arr);
276 //_____________________________________________________________________________
277 Double_t AliTPCCalPad::GetRMS()
280 // return mean of the RMS of all ROCs
284 for (Int_t isec = 0; isec < kNsec; isec++) {
285 AliTPCCalROC *calRoc = fROC[isec];
287 arr[n] = calRoc->GetRMS();
291 return TMath::Mean(n,arr);
294 //_____________________________________________________________________________
295 Double_t AliTPCCalPad::GetMedian()
298 // return mean of the median of all ROCs
302 for (Int_t isec = 0; isec < kNsec; isec++) {
303 AliTPCCalROC *calRoc = fROC[isec];
305 arr[n] = calRoc->GetMedian();
309 return TMath::Mean(n,arr);
312 //_____________________________________________________________________________
313 Double_t AliTPCCalPad::GetLTM(Double_t *sigma, Double_t fraction)
316 // return mean of the LTM and sigma of all ROCs
318 Double_t arrm[kNsec];
319 Double_t arrs[kNsec];
323 for (Int_t isec = 0; isec < kNsec; isec++) {
324 AliTPCCalROC *calRoc = fROC[isec];
326 if ( sigma ) sTemp=arrs+n;
327 arrm[n] = calRoc->GetLTM(sTemp,fraction);
331 if ( sigma ) *sigma = TMath::Mean(n,arrs);
332 return TMath::Mean(n,arrm);
335 //_____________________________________________________________________________
336 TH1F * AliTPCCalPad::MakeHisto1D(Float_t min, Float_t max,Int_t type){
339 // type -1 = user defined range
340 // 0 = nsigma cut nsigma=min
344 Float_t mean = GetMean();
345 Float_t sigma = GetRMS();
346 Float_t nsigma = TMath::Abs(min);
347 min = mean-nsigma*sigma;
348 max = mean+nsigma*sigma;
352 Float_t mean = GetMedian();
359 // LTM mean +- nsigma
362 Float_t mean = GetLTM(&sigma,max);
369 sprintf(name,"%s Pad 1D",GetTitle());
370 TH1F * his = new TH1F(name,name,100, min,max);
371 for (Int_t isec = 0; isec < kNsec; isec++) {
373 for (UInt_t irow=0; irow<fROC[isec]->GetNrows(); irow++){
374 UInt_t npads = (Int_t)fROC[isec]->GetNPads(irow);
375 for (UInt_t ipad=0; ipad<npads; ipad++){
376 his->Fill(fROC[isec]->GetValue(irow,ipad));
384 //_____________________________________________________________________________
385 TH2F *AliTPCCalPad::MakeHisto2D(Int_t side){
388 // side - specify the side A = 0 C = 1
389 // type - used types of determination of boundaries in z
390 Float_t kEpsilon = 0.000000000001;
391 TH2F * his = new TH2F(GetName(), GetName(), 250,-250,250,250,-250,250);
392 AliTPCROC * roc = AliTPCROC::Instance();
393 for (Int_t isec=0; isec<72; isec++){
394 if (side==0 && isec%36>=18) continue;
395 if (side>0 && isec%36<18) continue;
397 AliTPCCalROC * calRoc = fROC[isec];
398 for (UInt_t irow=0; irow<calRoc->GetNrows(); irow++)
399 for (UInt_t ipad=0; ipad<calRoc->GetNPads(irow); ipad++)
400 if (TMath::Abs(calRoc->GetValue(irow,ipad))>kEpsilon){
402 roc->GetPositionGlobal(isec,irow,ipad,xyz);
403 Int_t binx = 1+TMath::Nint((xyz[0]+250.)*0.5);
404 Int_t biny = 1+TMath::Nint((xyz[1]+250.)*0.5);
405 Float_t value = calRoc->GetValue(irow,ipad);
406 his->SetBinContent(binx,biny,value);
410 his->SetXTitle("x (cm)");
411 his->SetYTitle("y (cm)");
418 void AliTPCCalPad::MakeTree(const char * fileName, TObjArray * array) {
420 // Write tree with all available information
422 TTreeSRedirector cstream(fileName);
423 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
425 TString* names = new TString[array->GetEntries()];
426 for (Int_t ivalue = 0; ivalue < array->GetEntries(); ivalue++)
427 names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
429 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
431 // get statistic for given sector
433 TVectorF median(array->GetEntries());
434 TVectorF mean(array->GetEntries());
435 TVectorF rms(array->GetEntries());
436 TVectorF ltm(array->GetEntries());
438 TVectorF *vectorArray = new TVectorF[array->GetEntries()];
439 for (Int_t ivalue = 0; ivalue < array->GetEntries(); ivalue++)
440 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
442 for (Int_t ivalue = 0; ivalue < array->GetEntries(); ivalue++) {
443 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
444 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
445 median[ivalue] = calROC->GetMedian();
446 mean[ivalue] = calROC->GetMean();
447 rms[ivalue] = calROC->GetRMS();
448 ltm[ivalue] = calROC->GetLTM();
452 // fill vectors of variable per pad
454 TVectorF *posArray = new TVectorF[6];
455 for (Int_t ivalue = 0; ivalue < 6; ivalue++)
456 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
458 Float_t posG[3] = {0};
459 Float_t posL[3] = {0};
461 for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
462 for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
463 tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
464 tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
465 posArray[0][ichannel] = irow;
466 posArray[1][ichannel] = ipad;
467 posArray[2][ichannel] = posL[0];
468 posArray[3][ichannel] = posL[1];
469 posArray[4][ichannel] = posG[0];
470 posArray[5][ichannel] = posG[1];
472 // loop over array containing AliTPCCalPads
473 for (Int_t ivalue = 0; ivalue < array->GetEntries(); ivalue++) {
474 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
475 (vectorArray[ivalue])[ichannel] = calPad->GetCalROC(isector)->GetValue(irow, ipad);
481 cstream << "calPads" <<
482 "sector=" << isector;
484 for (Int_t ivalue = 0; ivalue < array->GetEntries(); ivalue++) {
485 cstream << "calPads" <<
486 (Char_t*)((names[ivalue] + "Median=").Data()) << median[ivalue] <<
487 (Char_t*)((names[ivalue] + "Mean=").Data()) << mean[ivalue] <<
488 (Char_t*)((names[ivalue] + "RMS=").Data()) << rms[ivalue] <<
489 (Char_t*)((names[ivalue] + "LTM=").Data()) << ltm[ivalue];
492 for (Int_t ivalue = 0; ivalue < array->GetEntries(); ivalue++) {
493 cstream << "calPads" <<
494 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
497 cstream << "calPads" <<
498 "row.=" << &posArray[0] <<
499 "pad.=" << &posArray[1] <<
500 "lx.=" << &posArray[2] <<
501 "ly.=" << &posArray[3] <<
502 "gx.=" << &posArray[4] <<
503 "gy.=" << &posArray[5];
505 cstream << "calPads" <<
509 delete[] vectorArray;