]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCCalPad.cxx
Ignoring *.rootmap
[u/mrichter/AliRoot.git] / TPC / AliTPCCalPad.cxx
CommitLineData
07627591 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
16/* $Id$ */
17
18///////////////////////////////////////////////////////////////////////////////
19// //
20// TPC calibration class for parameters which saved per pad //
21// //
22///////////////////////////////////////////////////////////////////////////////
23
24#include "AliTPCCalPad.h"
25#include "AliTPCCalROC.h"
184bcc16 26#include <TObjArray.h>
27#include <TAxis.h>
28#include <TGraph.h>
200be8a6 29#include <TGraph2D.h>
30#include <TH2F.h>
586007f3 31#include "TTreeStream.h"
32
07627591 33ClassImp(AliTPCCalPad)
34
35//_____________________________________________________________________________
36AliTPCCalPad::AliTPCCalPad():TNamed()
37{
38 //
39 // AliTPCCalPad default constructor
40 //
41
42 for (Int_t isec = 0; isec < kNsec; isec++) {
43 fROC[isec] = 0;
44 }
45
46}
47
48//_____________________________________________________________________________
49AliTPCCalPad::AliTPCCalPad(const Text_t *name, const Text_t *title)
50 :TNamed(name,title)
51{
52 //
53 // AliTPCCalPad constructor
54 //
55 for (Int_t isec = 0; isec < kNsec; isec++) {
56 fROC[isec] = new AliTPCCalROC(isec);
57 }
58}
59
60
61//_____________________________________________________________________________
62AliTPCCalPad::AliTPCCalPad(const AliTPCCalPad &c):TNamed(c)
63{
64 //
65 // AliTPCCalPad copy constructor
66 //
67
68 ((AliTPCCalPad &) c).Copy(*this);
69
70}
71
184bcc16 72//_____________________________________________________________________________
73AliTPCCalPad::AliTPCCalPad(TObjArray * array):TNamed()
74{
75 //
76 // AliTPCCalPad default constructor
77 //
78
79 for (Int_t isec = 0; isec < kNsec; isec++) {
80 fROC[isec] = (AliTPCCalROC *)array->At(isec);
81 }
82
83}
84
85
07627591 86///_____________________________________________________________________________
87AliTPCCalPad::~AliTPCCalPad()
88{
89 //
90 // AliTPCCalPad destructor
91 //
92
93 for (Int_t isec = 0; isec < kNsec; isec++) {
94 if (fROC[isec]) {
95 delete fROC[isec];
96 fROC[isec] = 0;
97 }
98 }
99
100}
101
102//_____________________________________________________________________________
103AliTPCCalPad &AliTPCCalPad::operator=(const AliTPCCalPad &c)
104{
105 //
106 // Assignment operator
107 //
108
109 if (this != &c) ((AliTPCCalPad &) c).Copy(*this);
110 return *this;
111
112}
113
114//_____________________________________________________________________________
115void AliTPCCalPad::Copy(TObject &c) const
116{
117 //
118 // Copy function
119 //
120
121 for (Int_t isec = 0; isec < kNsec; isec++) {
122 if (fROC[isec]) {
123 fROC[isec]->Copy(*((AliTPCCalPad &) c).fROC[isec]);
124 }
125 }
126 TObject::Copy(c);
127}
184bcc16 128
90127643 129//_____________________________________________________________________________
130void AliTPCCalPad::Add(Float_t c1)
131{
132 //
133 // add constant for all channels of all ROCs
134 //
135
136 for (Int_t isec = 0; isec < kNsec; isec++) {
137 if (fROC[isec]){
138 fROC[isec]->Add(c1);
139 }
140 }
141}
142
143//_____________________________________________________________________________
144void AliTPCCalPad::Multiply(Float_t c1)
145{
146 //
147 // multiply constant for all channels of all ROCs
148 //
149 for (Int_t isec = 0; isec < kNsec; isec++) {
150 if (fROC[isec]){
151 fROC[isec]->Multiply(c1);
152 }
153 }
154}
155
156//_____________________________________________________________________________
157void AliTPCCalPad::Add(const AliTPCCalPad * pad, Double_t c1)
158{
159 //
160 // add calpad channel by channel multiplied by c1 - all ROCs
161 //
162 for (Int_t isec = 0; isec < kNsec; isec++) {
163 if (fROC[isec]){
164 fROC[isec]->Add(pad->GetCalROC(isec),c1);
165 }
166 }
167}
168
169//_____________________________________________________________________________
170void AliTPCCalPad::Multiply(const AliTPCCalPad * pad)
171{
172 //
173 // multiply calpad channel by channel - all ROCs
174 //
175 for (Int_t isec = 0; isec < kNsec; isec++) {
176 if (fROC[isec]){
177 fROC[isec]->Multiply(pad->GetCalROC(isec));
178 }
179 }
180}
181
182//_____________________________________________________________________________
183void AliTPCCalPad::Divide(const AliTPCCalPad * pad)
184{
185 //
186 // divide calpad channel by channel - all ROCs
187 //
188 for (Int_t isec = 0; isec < kNsec; isec++) {
189 if (fROC[isec]){
190 fROC[isec]->Divide(pad->GetCalROC(isec));
191 }
192 }
193}
194
195//_____________________________________________________________________________
184bcc16 196TGraph * AliTPCCalPad::MakeGraph(Int_t type, Float_t ratio){
197 //
198 // type=1 - mean
199 // 2 - median
200 // 3 - LTM
201 Int_t npoints = 0;
202 for (Int_t i=0;i<72;i++) if (fROC[i]) npoints++;
203 TGraph * graph = new TGraph(npoints);
204 npoints=0;
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));
210 npoints++;
211 }
212
213 graph->GetXaxis()->SetTitle("Sector");
214 if (type==0) {
215 graph->GetYaxis()->SetTitle("Mean");
216 graph->SetMarkerStyle(22);
217 }
218 if (type==1) {
219 graph->GetYaxis()->SetTitle("Median");
220 graph->SetMarkerStyle(22);
221 }
222 if (type==2) {
223 graph->GetYaxis()->SetTitle(Form("Mean%f",ratio));
224 graph->SetMarkerStyle(24);
225 }
226
227 return graph;
228}
200be8a6 229
90127643 230//_____________________________________________________________________________
231Double_t AliTPCCalPad::GetMeanRMS(Double_t &rms)
232{
233 //
234 // Calculate mean an RMS of all rocs
235 //
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];
239 if ( calRoc ){
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);
243 sum+=val;
244 sum2+=val*val;
245 n++;
246 }
247 }
248
249 }
250 }
251 Double_t n1 = 1./n;
252 Double_t mean = sum*n1;
253 rms = TMath::Sqrt(TMath::Abs(sum2*n1-mean*mean));
254 return mean;
255}
256
257
258//_____________________________________________________________________________
259Double_t AliTPCCalPad::GetMean()
260{
261 //
262 // return mean of the mean of all ROCs
263 //
264 Double_t arr[kNsec];
265 Int_t n=0;
266 for (Int_t isec = 0; isec < kNsec; isec++) {
267 AliTPCCalROC *calRoc = fROC[isec];
268 if ( calRoc ){
269 arr[n] = calRoc->GetMean();
270 n++;
271 }
272 }
273 return TMath::Mean(n,arr);
274}
275
276//_____________________________________________________________________________
277Double_t AliTPCCalPad::GetRMS()
278{
279 //
280 // return mean of the RMS of all ROCs
281 //
282 Double_t arr[kNsec];
283 Int_t n=0;
284 for (Int_t isec = 0; isec < kNsec; isec++) {
285 AliTPCCalROC *calRoc = fROC[isec];
286 if ( calRoc ){
287 arr[n] = calRoc->GetRMS();
288 n++;
289 }
290 }
291 return TMath::Mean(n,arr);
292}
293
294//_____________________________________________________________________________
295Double_t AliTPCCalPad::GetMedian()
296{
297 //
298 // return mean of the median of all ROCs
299 //
300 Double_t arr[kNsec];
301 Int_t n=0;
302 for (Int_t isec = 0; isec < kNsec; isec++) {
303 AliTPCCalROC *calRoc = fROC[isec];
304 if ( calRoc ){
305 arr[n] = calRoc->GetMedian();
306 n++;
307 }
308 }
309 return TMath::Mean(n,arr);
310}
311
312//_____________________________________________________________________________
313Double_t AliTPCCalPad::GetLTM(Double_t *sigma, Double_t fraction)
314{
315 //
316 // return mean of the LTM and sigma of all ROCs
317 //
318 Double_t arrm[kNsec];
319 Double_t arrs[kNsec];
320 Double_t *sTemp=0x0;
321 Int_t n=0;
322
323 for (Int_t isec = 0; isec < kNsec; isec++) {
324 AliTPCCalROC *calRoc = fROC[isec];
325 if ( calRoc ){
326 if ( sigma ) sTemp=arrs+n;
327 arrm[n] = calRoc->GetLTM(sTemp,fraction);
328 n++;
329 }
330 }
331 if ( sigma ) *sigma = TMath::Mean(n,arrs);
332 return TMath::Mean(n,arrm);
333}
334
335//_____________________________________________________________________________
336TH1F * AliTPCCalPad::MakeHisto1D(Float_t min, Float_t max,Int_t type){
337 //
338 // make 1D histo
339 // type -1 = user defined range
340 // 0 = nsigma cut nsigma=min
341 if (type>=0){
342 if (type==0){
343 // nsigma range
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;
349 }
350 if (type==1){
351 // fixed range
352 Float_t mean = GetMedian();
353 Float_t delta = min;
354 min = mean-delta;
355 max = mean+delta;
356 }
357 if (type==2){
358 //
359 // LTM mean +- nsigma
360 //
361 Double_t sigma;
362 Float_t mean = GetLTM(&sigma,max);
363 sigma*=min;
364 min = mean-sigma;
365 max = mean+sigma;
366 }
367 }
368 char name[1000];
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++) {
372 if (fROC[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));
377 }
378 }
379 }
380 }
381 return his;
382}
383
384//_____________________________________________________________________________
200be8a6 385TH2F *AliTPCCalPad::MakeHisto2D(Int_t side){
386 //
387 // Make 2D graph
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;
396 if (fROC[isec]){
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){
401 Float_t xyz[3];
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);
407 }
408 }
409 }
410 his->SetXTitle("x (cm)");
411 his->SetYTitle("y (cm)");
412 return his;
413}
414
415
416
586007f3 417
418void AliTPCCalPad::MakeTree(const char * fileName, TObjArray * array) {
419 //
420 // Write tree with all available information
421 //
422 TTreeSRedirector cstream(fileName);
423 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
424
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());
428
429 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
430 //
431 // get statistic for given sector
432 //
433 TVectorF median(array->GetEntries());
434 TVectorF mean(array->GetEntries());
435 TVectorF rms(array->GetEntries());
436 TVectorF ltm(array->GetEntries());
437
438 TVectorF *vectorArray = new TVectorF[array->GetEntries()];
439 for (Int_t ivalue = 0; ivalue < array->GetEntries(); ivalue++)
440 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
441
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();
449 }
450
451 //
452 // fill vectors of variable per pad
453 //
454 TVectorF *posArray = new TVectorF[6];
455 for (Int_t ivalue = 0; ivalue < 6; ivalue++)
456 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
457
458 Float_t posG[3] = {0};
459 Float_t posL[3] = {0};
460 Int_t ichannel = 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];
471
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);
476 }
477 ichannel++;
478 }
479 }
480
481 cstream << "calPads" <<
482 "sector=" << isector;
483
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];
490 }
491
492 for (Int_t ivalue = 0; ivalue < array->GetEntries(); ivalue++) {
493 cstream << "calPads" <<
494 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
495 }
496
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];
504
505 cstream << "calPads" <<
506 "\n";
507
508 delete[] posArray;
509 delete[] vectorArray;
510 }
511 delete[] names;
512}
513
514