]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDCalibraVector.cxx
Bug fix for like sign V0s
[u/mrichter/AliRoot.git] / TRD / AliTRDCalibraVector.cxx
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 // AliTRDCalibraVector                                                    //       
21 //                                                                        //   
22 // This class is for the vector method of the TRD calibration.            //
23 //                                                                        //
24 // Author:                                                                //
25 //   R. Bailhache (R.Bailhache@gsi.de)                                    //
26 //                                                                        //
27 ////////////////////////////////////////////////////////////////////////////
28
29 #include <TGraphErrors.h>
30 #include <TH1F.h>
31 #include <TObjArray.h>
32 #include <TMath.h>
33 #include <TDirectory.h>
34 #include <TROOT.h>
35 #include <TFile.h>
36
37 #include "AliLog.h"
38
39 #include "AliTRDCalibraVector.h"
40 #include "AliTRDCommonParam.h"
41 #include "TArrayF.h"
42 #include "TArrayI.h"
43
44 ClassImp(AliTRDCalibraVector)
45
46 //______________________________________________________________________________________
47 AliTRDCalibraVector::AliTRDCalibraVector()
48   :TObject()
49   ,fNameCH("CH2d")
50   ,fNamePH("PH2d")
51   ,fNamePRF("PRF2d")
52   ,fDetectorPH(-1)
53   ,fDetectorCH(-1)
54   ,fDetectorPRF(-1)
55   ,fNumberBinCharge(0)
56   ,fNumberBinPRF(0)
57   ,fTimeMax(0)
58   ,fPRFRange(1.5)
59 {
60   //
61   // Default constructor
62   //
63
64   for (Int_t idet = 0; idet < 540; idet++){
65     
66     fPHEntries[idet]=new TArrayI();
67     fPHMean[idet]=new TArrayF();
68     fPHSquares[idet]=new TArrayF();
69
70     fPRFEntries[idet]=new TArrayI();
71     fPRFMean[idet]=new TArrayF();
72     fPRFSquares[idet]=new TArrayF();
73
74
75     fCHEntries[idet]=new TArrayI();
76     
77   }
78   
79   for(Int_t k = 0; k < 3; k++){
80     fDetCha0[k] = 0;
81     fDetCha2[k] = 0;
82   }
83  
84 }
85 //______________________________________________________________________________________
86 AliTRDCalibraVector::AliTRDCalibraVector(const AliTRDCalibraVector &c)
87   :TObject(c)
88   ,fNameCH("CH2d")
89   ,fNamePH("PH2d")
90   ,fNamePRF("PRF2d")
91   ,fDetectorPH(-1)
92   ,fDetectorCH(-1)
93   ,fDetectorPRF(-1)
94   ,fNumberBinCharge(c.fNumberBinCharge)
95   ,fNumberBinPRF(c.fNumberBinPRF)
96   ,fTimeMax(c.fTimeMax)
97   ,fPRFRange(c.fPRFRange)
98 {
99   //
100   // Copy constructor
101   //
102   for (Int_t idet = 0; idet < 540; idet++){
103     
104     const TArrayI *phEntries  = (TArrayI*)c.fPHEntries[idet];
105     const TArrayF *phMean     = (TArrayF*)c.fPHMean[idet];
106     const TArrayF *phSquares  = (TArrayF*)c.fPHSquares[idet];
107
108     const TArrayI *prfEntries  = (TArrayI*)c.fPRFEntries[idet];
109     const TArrayF *prfMean     = (TArrayF*)c.fPRFMean[idet];
110     const TArrayF *prfSquares  = (TArrayF*)c.fPRFSquares[idet];
111
112     const TArrayI *chEntries  = (TArrayI*)c.fCHEntries[idet];
113
114     if ( phEntries != 0x0 )  fPHEntries[idet]=new TArrayI(*phEntries);
115     if ( phMean != 0x0 )     fPHMean[idet]=new TArrayF(*phMean);
116     if ( phSquares != 0x0 )  fPHSquares[idet]=new TArrayF(*phSquares);
117
118     if ( prfEntries != 0x0 )  fPRFEntries[idet]=new TArrayI(*prfEntries);
119     if ( prfMean != 0x0 )     fPRFMean[idet]=new TArrayF(*prfMean);
120     if ( prfSquares != 0x0 )  fPRFSquares[idet]=new TArrayF(*prfSquares);
121
122
123     if ( chEntries != 0x0 )   fCHEntries[idet]=new TArrayI(*chEntries);
124
125   }
126
127   for(Int_t k = 0; k < 3; k++){
128     fDetCha0[k] = c.fDetCha0[k];
129     fDetCha2[k] = c.fDetCha2[k];
130   }
131
132    
133 }
134 //_____________________________________________________________________
135 AliTRDCalibraVector& AliTRDCalibraVector::operator = (const  AliTRDCalibraVector &source)
136 {
137   //
138   // assignment operator
139   //
140   if (&source == this) return *this;
141   new (this) AliTRDCalibraVector(source);
142
143   return *this;
144 }
145 //____________________________________________________________________________________
146 AliTRDCalibraVector::~AliTRDCalibraVector()
147 {
148   //
149   // AliTRDCalibraVector destructor
150   //
151
152   if(fPHEntries) delete fPHEntries;
153   if(fPHMean) delete fPHMean;
154   if(fPHSquares) delete fPHSquares;
155   if(fPRFEntries) delete fPRFEntries;
156   if(fPRFMean) delete fPRFMean;
157   if(fPRFSquares) delete fPRFSquares;
158   if(fCHEntries) delete fCHEntries;
159
160 }
161 //_____________________________________________________________________________
162 Int_t AliTRDCalibraVector::SearchBin(Float_t value, Int_t i) const
163 {
164   //
165   // Search the bin
166   //
167
168   Int_t reponse        = 0;
169   Float_t fbinmin      = 0;
170   Float_t fbinmax      = value;
171   Int_t fNumberOfBin   = -1;
172
173   switch(i)
174     {
175     case 0:
176       fbinmax      = 300.0;
177       fbinmin      = 0.0;
178       fNumberOfBin = fNumberBinCharge;
179       break;
180
181     case 2:
182       fbinmax      =   TMath::Abs(fPRFRange);
183       fbinmin      =  -TMath::Abs(fPRFRange);
184       fNumberOfBin = fNumberBinPRF;
185       break;
186       
187     default: 
188       return -1;
189     }
190   
191   // Return -1 if out
192   if ((value >= fbinmax) || 
193       (value <  fbinmin)) {
194     return -1;
195   }
196   else {
197     reponse = (Int_t) ((fNumberOfBin*(value-fbinmin)) / (fbinmax-fbinmin));
198   }
199
200   return reponse;
201
202 }
203 //_____________________________________________________________________________
204 Bool_t AliTRDCalibraVector::UpdateVectorCH(Int_t det, Int_t group, Float_t value)
205 {
206   //
207   // Fill the vector CH   
208   //
209
210   // Search bin
211   Int_t bin = SearchBin(value,0);
212   // Out
213   if (bin == -1) {
214     return kFALSE; 
215   }
216
217
218
219   if(fDetectorCH != det){
220     fCHEntries[det] = ((TArrayI *)GetCHEntries(det,kTRUE));
221   }
222
223   Int_t entries  = ((TArrayI *)fCHEntries[det])->At(group*fNumberBinCharge+bin);
224   
225   Int_t entriesn = entries+1;
226   ((TArrayI *)fCHEntries[det])->AddAt(entriesn,group*fNumberBinCharge+bin);
227   
228   fDetectorCH = det;
229
230  
231   return kTRUE;
232
233 }
234 //_____________________________________________________________________________
235 Bool_t AliTRDCalibraVector::UpdateVectorPRF(Int_t det, Int_t group, Float_t x, Float_t y)
236 {
237   //
238   // Fill the vector PRF
239   //
240
241   // Search bin
242   Int_t bin = SearchBin(x,2);
243   // Out
244   if (bin == -1) {
245     return kFALSE; 
246   }
247
248   
249   if(fDetectorPRF != det){
250     fPRFEntries[det] = ((TArrayI *)GetPRFEntries(det,kTRUE));
251     fPRFMean[det]    = ((TArrayF *)GetPRFMean(det,kTRUE));
252     fPRFSquares[det] = ((TArrayF *)GetPRFSquares(det,kTRUE));
253   }
254
255   Int_t entries  = ((TArrayI *)fPRFEntries[det])->At((Int_t)group*fNumberBinPRF+bin);
256   Float_t mean   = ((TArrayI *)fPRFMean[det])->At((Int_t)group*fNumberBinPRF+bin);
257   Float_t square = ((TArrayI *)fPRFSquares[det])->At((Int_t)group*fNumberBinPRF+bin);
258   
259   Int_t entriesn = entries+1;
260   ((TArrayI *)fPRFEntries[det])->AddAt(entriesn,(Int_t)group*fNumberBinPRF+bin);
261   Float_t meann = (mean*((Float_t)entries)+y)/((Float_t)entriesn);
262   ((TArrayF *)fPRFMean[det])->AddAt(meann,(Int_t)group*fNumberBinPRF+bin);
263   Float_t squaren = ((square*((Float_t)entries))+(y*y))/((Float_t)entriesn);
264   ((TArrayF *)fPRFSquares[det])->AddAt(squaren,(Int_t)group*fNumberBinPRF+bin);
265
266   
267   fDetectorPRF = det;
268
269   return kTRUE;
270   
271 }
272 //_____________________________________________________________________________
273 Bool_t AliTRDCalibraVector::UpdateVectorPH(Int_t det, Int_t group, Int_t time, Float_t value)
274 {
275   //
276   // Fill the vector PH  
277   //
278
279   // Search bin
280   Int_t bin = time;
281   // Out
282   if ((bin <         0) || 
283       (bin >= fTimeMax)) {
284     return kFALSE; 
285   }
286
287
288   if(fDetectorPH != det){
289     fPHEntries[det] = ((TArrayI *)GetPHEntries(det,kTRUE));
290     fPHMean[det]    = ((TArrayF *)GetPHMean(det,kTRUE));
291     fPHSquares[det] = ((TArrayF *)GetPHSquares(det,kTRUE));
292   }
293
294   Int_t entries  = ((TArrayI *)fPHEntries[det])->At(group*fTimeMax+bin);
295   Float_t mean   = ((TArrayI *)fPHMean[det])->At(group*fTimeMax+bin);
296   Float_t square = ((TArrayI *)fPHSquares[det])->At(group*fTimeMax+bin);
297   
298   Int_t entriesn = entries+1;
299   ((TArrayI *)fPHEntries[det])->AddAt(entriesn,group*fTimeMax+bin);
300   Float_t meann = (mean*((Float_t)entries)+value)/((Float_t)entriesn);
301   ((TArrayF *)fPHMean[det])->AddAt(meann,group*fTimeMax+bin);
302   Float_t squaren = ((square*((Float_t)entries))+(value*value))/((Float_t)entriesn);
303   ((TArrayF *)fPHSquares[det])->AddAt(squaren,group*fTimeMax+bin);
304
305   
306   fDetectorPH = det;
307
308   return kTRUE;
309   
310 }
311 //__________________________________________________________________________________
312 Bool_t AliTRDCalibraVector::Add(AliTRDCalibraVector *calvect)
313 {
314   //
315   // Add a other AliTRCalibraVector to this one
316   //
317
318   // Check compatibility
319   if(fNumberBinCharge != calvect->GetNumberBinCharge()) return kFALSE;
320   if(fNumberBinPRF    != calvect->GetNumberBinPRF()) return kFALSE;
321   if(fPRFRange        != calvect->GetPRFRange()) return kFALSE;
322   if(fTimeMax         != calvect->GetTimeMax()) return kFALSE;
323   for(Int_t k = 0; k < 3; k++){
324     if(fDetCha0[k] != calvect->GetDetCha0(k)) return kFALSE;
325     if(fDetCha2[k] != calvect->GetDetCha2(k)) return kFALSE;
326   }
327
328   //printf("pass0!\n");
329
330   // Add
331   for (Int_t idet = 0; idet < 540; idet++){
332     
333     const TArrayI *phEntriesvect   = calvect->GetPHEntries(idet);
334     const TArrayF *phMeanvect      = calvect->GetPHMean(idet);
335     const TArrayF *phSquaresvect   = calvect->GetPHSquares(idet);
336     
337     const TArrayI *prfEntriesvect  = calvect->GetPRFEntries(idet);
338     const TArrayF *prfMeanvect     = calvect->GetPRFMean(idet);
339     const TArrayF *prfSquaresvect  = calvect->GetPRFSquares(idet);
340     
341     const TArrayI *chEntriesvect   = calvect->GetCHEntries(idet);
342
343     //printf("idet %d!\n",idet);
344
345     //printf("phEntriesvect %d\n",(Bool_t) phEntriesvect);
346     //printf("phMeanvect %d\n",(Bool_t) phMeanvect);
347     //printf("phSquaresvect %d\n",(Bool_t) phSquaresvect);
348
349     //printf("prfEntriesvect %d\n",(Bool_t) prfEntriesvect);
350     //printf("prfMeanvect %d\n",(Bool_t) prfMeanvect);
351     //printf("prfSquaresvect %d\n",(Bool_t) prfSquaresvect);
352
353     //printf("chEntriesvect %d\n",(Bool_t) chEntriesvect);
354
355     if ( phEntriesvect != 0x0 ){
356       //Take the stuff
357       fPHEntries[idet] = ((TArrayI *)GetPHEntries(idet,kTRUE));
358       fPHMean[idet]    = ((TArrayF *)GetPHMean(idet,kTRUE));
359       fPHSquares[idet] = ((TArrayF *)GetPHSquares(idet,kTRUE));
360       Int_t total = ((TArrayI *)fPHEntries[idet])->GetSize();
361       // Add
362       for(Int_t k = 0; k < total; k++){
363         Int_t entries  = ((TArrayI *)fPHEntries[idet])->At(k);
364         Float_t mean   = ((TArrayF *)fPHMean[idet])->At(k);
365         Float_t square = ((TArrayF *)fPHSquares[idet])->At(k);
366   
367         Int_t entriesn = entries+((TArrayI *)phEntriesvect)->At(k);
368         if(entriesn <= 0) continue;
369         ((TArrayI *)fPHEntries[idet])->AddAt(entriesn,k);
370         Float_t meann = (mean*((Float_t)entries)+((TArrayF *)phMeanvect)->At(k)*((Float_t)((TArrayI *)phEntriesvect)->At(k)))/((Float_t)entriesn);
371         ((TArrayF *)fPHMean[idet])->AddAt(meann,k);
372         Float_t sq      = ((TArrayF *)phSquaresvect)->At(k)*((Float_t)((TArrayI *)phEntriesvect)->At(k));
373         Float_t squaren = ((square*((Float_t)entries))+sq)/((Float_t)entriesn);
374         ((TArrayF *)fPHSquares[idet])->AddAt(squaren,k);
375         //printf("test ph!\n");
376       }
377     }     
378
379     if ( prfEntriesvect != 0x0 ){
380       //Take the stuff
381       fPRFEntries[idet] = ((TArrayI *)GetPRFEntries(idet,kTRUE));
382       fPRFMean[idet]    = ((TArrayF *)GetPRFMean(idet,kTRUE));
383       fPRFSquares[idet] = ((TArrayF *)GetPRFSquares(idet,kTRUE));
384       Int_t total = fPRFEntries[idet]->GetSize();
385       // Add
386       for(Int_t k = 0; k < total; k++){
387         Int_t entries  = ((TArrayI *)fPRFEntries[idet])->At(k);
388         Float_t mean   = ((TArrayF *)fPRFMean[idet])->At(k);
389         Float_t square = ((TArrayF *)fPRFSquares[idet])->At(k);
390
391         //printf("entries0 %d\n",entries);
392         //printf("mean0 %f\n",mean);
393         //printf("square0 %f\n",square);
394
395         //printf("entries1 %d\n",prfEntriesvect->At(k));
396         //printf("mean1 %f\n",prfMeanvect->At(k));
397         //printf("square1 %f\n",prfSquaresvect->At(k));
398   
399
400         //printf("entries0 size %d\n",fPRFEntries->GetSize());
401         //printf("mean0 size %d\n",fPRFMean->GetSize());
402         //printf("squares0 size %d\n",fPRFSquares->GetSize());
403
404         //printf("entries1 size %d\n",prfEntriesvect->GetSize());
405         //printf("mean1 size %d\n",prfMeanvect->GetSize());
406         //printf("squares1 size %d\n",prfSquaresvect->GetSize());
407
408
409         Int_t entriesn = entries+((TArrayI *)prfEntriesvect)->At(k);
410         if(entriesn <= 0) continue;
411         ((TArrayI *)fPRFEntries[idet])->AddAt(entriesn,k);
412         Float_t meann = (mean*((Float_t)entries)+((TArrayF *)prfMeanvect)->At(k)*((Float_t)((TArrayI *)prfEntriesvect)->At(k)))/((Float_t)entriesn);
413         ((TArrayF *)fPRFMean[idet])->AddAt(meann,k);
414         Float_t sq      = ((TArrayF *)prfSquaresvect)->At(k)*((Float_t)((TArrayI *)prfEntriesvect)->At(k));
415         Float_t squaren = ((square*((Float_t)entries))+sq)/((Float_t)entriesn);
416         ((TArrayF *)fPRFSquares[idet])->AddAt(squaren,k);
417         //printf("test prf!\n");
418       }
419     }
420
421     if ( chEntriesvect != 0x0 ){
422       //Take the stuff
423       fCHEntries[idet] = ((TArrayI *)GetCHEntries(idet,kTRUE));
424       Int_t total = fCHEntries[idet]->GetSize();
425       //if(idet == 180) printf("total %d\n",total);
426       // Add
427       for(Int_t k = 0; k < total; k++){
428         Int_t entries  = ((TArrayI *)fCHEntries[idet])->At(k);
429         Int_t entriesn = entries+((TArrayI *)chEntriesvect)->At(k);
430         //if((idet == 180) && ((entries != 0) || (entriesn != 0))) printf("for k %d we have entries %d and entriesn %d\n",k,entries,entriesn);
431         if(entriesn <= 0) continue;
432         ((TArrayI *)fCHEntries[idet])->AddAt(entriesn,k);
433       }
434       //printf("test ch!\n");
435     }           
436   }
437   
438   return kTRUE;
439
440 //_____________________________________________________________________________
441 TGraphErrors *AliTRDCalibraVector::ConvertVectorPHTGraphErrors(Int_t det, Int_t group
442                                         , const Char_t * name)
443 {
444   //
445   // Convert the fVectorPHMean, fVectorPHSquares and fVectorPHEntries in TGraphErrors
446   //
447
448   // Take the info
449   fPHEntries[det] = ((TArrayI *)GetPHEntries(det,kTRUE));
450   fPHMean[det]    = ((TArrayF *)GetPHMean(det,kTRUE));
451   fPHSquares[det] = ((TArrayF *)GetPHSquares(det,kTRUE));
452   
453
454   // Init the stuff
455   TGraphErrors *histo;
456   Float_t sf = 10.0;
457   AliTRDCommonParam *parCom = AliTRDCommonParam::Instance();
458   if (!parCom) {
459     AliInfo("Could not get CommonParam, take the default 10MHz");
460   }
461   sf = parCom->GetSamplingFrequency();
462   // Axis
463   Double_t *x;
464   Double_t *y;
465   Double_t *ex;
466   Double_t *ey;
467   Double_t step = 0.0;
468   Double_t min  = 0.0;
469   x  = new Double_t[fTimeMax]; // Xaxis
470   y  = new Double_t[fTimeMax]; // Sum/entries
471   ex = new Double_t[fTimeMax]; // Nentries
472   ey = new Double_t[fTimeMax]; // Sum of square/nentries
473   step = 1.0 / sf;
474   min  = 0.0;
475   Int_t offset = group*fTimeMax;
476   
477   // Fill histo
478   for (Int_t k = 0; k < fTimeMax; k++) {
479     x[k]  = min + k*step;
480     y[k]  = 0.0;
481     ex[k] = 0.0;
482     ey[k] = 0.0;   
483     Int_t bin = offset+k;
484     // Fill only if there is more than 0 something
485     if (fPHEntries[det]->At(bin) > 0) {
486       ex[k] = ((TArrayI *)fPHEntries[det])->At(bin);
487       y[k]  = ((TArrayF *)fPHMean[det])->At(bin);
488       ey[k] = ((TArrayF *)fPHSquares[det])->At(bin);
489     }
490   }
491
492   // Define the TGraphErrors
493   histo = new TGraphErrors(fTimeMax,x,y,ex,ey);
494   histo->SetTitle(name); 
495   return histo;
496
497 }
498 //_____________________________________________________________________________
499 TGraphErrors *AliTRDCalibraVector::ConvertVectorPRFTGraphErrors(Int_t det, Int_t group
500                                         , const Char_t * name)
501 {
502   //
503   // Convert the fVectorPRFMean, fVectorPRFSquares and fVectorPRFEntries in TGraphErrors
504   //
505
506   // Take the info
507   fPRFEntries[det] = ((TArrayI *)GetPRFEntries(det,kTRUE));
508   fPRFMean[det]    = ((TArrayF *)GetPRFMean(det,kTRUE));
509   fPRFSquares[det] = ((TArrayF *)GetPRFSquares(det,kTRUE));
510   
511
512   // Init the stuff
513   TGraphErrors *histo;
514   // Axis
515   Double_t *x;
516   Double_t *y;
517   Double_t *ex;
518   Double_t *ey;
519   Double_t step = 0.0;
520   Double_t min  = 0.0;
521   x  = new Double_t[fNumberBinPRF]; // Xaxis
522   y  = new Double_t[fNumberBinPRF]; // Sum/entries
523   ex = new Double_t[fNumberBinPRF]; // Nentries
524   ey = new Double_t[fNumberBinPRF]; // Sum of square/nentries
525   step = (2*TMath::Abs(fPRFRange)) / fNumberBinPRF;
526   min  = -TMath::Abs(fPRFRange) + step / 2.0;
527   Int_t offset = group*fNumberBinPRF;
528   //printf("number of total: %d\n",fNumberBinPRF);
529   // Fill histo
530   for (Int_t k = 0; k < fNumberBinPRF; k++) {
531     x[k]  = min + k*step;
532     y[k]  = 0.0;
533     ex[k] = 0.0;
534     ey[k] = 0.0;
535     Int_t bin = offset+k;
536     // Fill only if there is more than 0 something
537     if (fPRFEntries[det]->At(bin) > 0) {
538       ex[k] = ((TArrayF *)fPRFEntries[det])->At(bin);
539       y[k]  = ((TArrayF *)fPRFMean[det])->At(bin);
540       ey[k] = ((TArrayF *)fPRFSquares[det])->At(bin);
541     }
542     //printf("Number of entries %f for %d\n",ex[k],k);
543   }
544
545   // Define the TGraphErrors
546   histo = new TGraphErrors(fNumberBinPRF,x,y,ex,ey);
547   histo->SetTitle(name); 
548   return histo;
549
550 }  
551 //_____________________________________________________________________________
552 TH1F *AliTRDCalibraVector::ConvertVectorCHHisto(Int_t det, Int_t group
553                                               , const Char_t * name)
554 {
555   //
556   // Convert the fVectorCHEntries in TH1F
557   //
558
559   // Take the info
560   fCHEntries[det] = ((TArrayI *)GetCHEntries(det,kTRUE));
561   
562   // Init the stuff
563   TH1F *histo = new TH1F(name,name,fNumberBinCharge,0,300);
564   histo->Sumw2();
565   Int_t offset = group*fNumberBinCharge;
566   // Fill histo
567   for (Int_t k = 0; k < fNumberBinCharge; k++) {
568     Int_t bin = offset+k;
569     histo->SetBinContent(k+1,((TArrayI *)fCHEntries[det])->At(bin));
570     histo->SetBinError(k+1,TMath::Sqrt(TMath::Abs(((TArrayI *)fCHEntries[det])->At(bin))));
571   }
572   
573   return histo;
574
575
576 //_____________________________________________________________________
577 TArrayI* AliTRDCalibraVector::GetPHEntries(Int_t det
578                                               , Bool_t force) /*FOLD00*/
579 {
580     //
581     // return pointer to Carge ROC Calibration
582     // if force is true create a new histogram if it doesn't exist allready
583     //
584     TArrayI**arr = &fPHEntries[0];
585     return GetEntriesPH(det, arr, force);
586 }
587 //_____________________________________________________________________
588 TArrayI* AliTRDCalibraVector::GetPRFEntries(Int_t det
589                                                , Bool_t force) /*FOLD00*/
590 {
591     //
592     // return pointer to Carge ROC Calibration
593     // if force is true create a new histogram if it doesn't exist allready
594     //
595     TArrayI**arr = &fPRFEntries[0];
596     return GetEntriesPRF(det, arr, force);
597 }
598 //_____________________________________________________________________
599 TArrayI* AliTRDCalibraVector::GetCHEntries(Int_t det
600                                               , Bool_t force) /*FOLD00*/
601 {
602     //
603     // return pointer to Carge ROC Calibration
604     // if force is true create a new histogram if it doesn't exist allready
605     //
606     TArrayI**arr = &fCHEntries[0];
607     return GetEntriesCH(det, arr, force);
608 }
609 //_____________________________________________________________________
610 TArrayF* AliTRDCalibraVector::GetPHMean(Int_t det
611                                            , Bool_t force) /*FOLD00*/
612 {
613     //
614     // return pointer to Carge ROC Calibration
615     // if force is true create a new histogram if it doesn't exist allready
616     //
617     TArrayF**arr = &fPHMean[0];
618     return GetMeanSquaresPH(det, arr, force);
619 }
620 //_____________________________________________________________________
621 TArrayF* AliTRDCalibraVector::GetPHSquares(Int_t det
622                                               , Bool_t force) /*FOLD00*/
623 {
624     //
625     // return pointer to Carge ROC Calibration
626     // if force is true create a new histogram if it doesn't exist allready
627     //
628     TArrayF**arr = &fPHSquares[0];
629     return GetMeanSquaresPH(det, arr, force);
630 }
631 //_____________________________________________________________________
632 TArrayF* AliTRDCalibraVector::GetPRFMean(Int_t det
633                                             , Bool_t force) /*FOLD00*/
634 {
635     //
636     // return pointer to Carge ROC Calibration
637     // if force is true create a new histogram if it doesn't exist allready
638     //
639     TArrayF**arr = &fPRFMean[0];
640     return GetMeanSquaresPRF(det, arr, force);
641 }
642 //_____________________________________________________________________
643 TArrayF* AliTRDCalibraVector::GetPRFSquares(Int_t det
644                                                , Bool_t force) /*FOLD00*/
645 {
646     //
647     // return pointer to Carge ROC Calibration
648     // if force is true create a new histogram if it doesn't exist allready
649     //
650     TArrayF**arr = &fPRFSquares[0];
651     return GetMeanSquaresPRF(det, arr, force);
652 }
653 //_____________________________________________________________________
654 TArrayI* AliTRDCalibraVector::GetEntriesCH(Int_t det
655                                               , TArrayI** arr
656                                               , Bool_t force) /*FOLD00*/
657 {
658     //
659     // return pointer to TArrayI Entries
660     // if force is true create a new TArrayI if it doesn't exist allready
661     //
662   if ( !force || (((TArrayI *)arr[det])->GetSize()>0))
663         return (TArrayI*)arr[det];
664
665     // if we are forced and TArrayI doesn't yes exist create it
666   Int_t stack  = GetStack(det);
667   Int_t ngroup = 0;
668   if(stack == 2) ngroup = fDetCha2[0]*fNumberBinCharge;
669   else ngroup = fDetCha0[0]*fNumberBinCharge;
670   // init
671   ((TArrayI *)arr[det])->Set(ngroup);
672   for(Int_t k = 0; k < ngroup; k++){
673     ((TArrayI *)arr[det])->AddAt(0,k);
674   }
675   return (TArrayI*)arr[det];
676 }
677 //_____________________________________________________________________
678 TArrayI* AliTRDCalibraVector::GetEntriesPRF(Int_t det
679                                                , TArrayI** arr
680                                                , Bool_t force) /*FOLD00*/
681 {
682     //
683     // return pointer to TArrayI Entries
684     // if force is true create a new TArrayI if it doesn't exist allready
685     //
686   if ( !force || (((TArrayI *)arr[det])->GetSize()>0))
687         return (TArrayI*)arr[det];
688
689   // if we are forced and TArrayI doesn't yes exist create it
690   Int_t stack  = GetStack(det);
691   Int_t ngroup = 0;
692   if(stack == 2) ngroup = fDetCha2[2]*fNumberBinPRF;
693   else ngroup = fDetCha0[2]*fNumberBinPRF;
694   // init
695   ((TArrayI *)arr[det])->Set(ngroup);
696   for(Int_t k = 0; k < ngroup; k++){
697     ((TArrayI *)arr[det])->AddAt(0,k);
698   }
699   return (TArrayI*)arr[det];
700
701 }
702 //_____________________________________________________________________
703 TArrayI* AliTRDCalibraVector::GetEntriesPH(Int_t det
704                                               , TArrayI** arr
705                                               , Bool_t force) /*FOLD00*/
706 {
707     //
708     // return pointer to TArrayI Entries
709     // if force is true create a new TArrayI if it doesn't exist allready
710     //
711     if ( !force || (((TArrayI *)arr[det])->GetSize()>0))
712         return (TArrayI*)arr[det];
713
714     // if we are forced and TArrayI doesn't yes exist create it
715     Int_t stack  = GetStack(det);
716     Int_t ngroup = 0;
717     if(stack == 2) ngroup = fDetCha2[1]*fTimeMax;
718     else ngroup = fDetCha0[1]*fTimeMax;
719     // init
720     ((TArrayI *)arr[det])->Set(ngroup);
721     for(Int_t k = 0; k < ngroup; k++){
722       ((TArrayI *)arr[det])->AddAt(0,k);
723     }
724     return (TArrayI*)arr[det];
725
726 }
727 //_____________________________________________________________________
728 TArrayF* AliTRDCalibraVector::GetMeanSquaresPH(Int_t det
729                                                   , TArrayF** arr
730                                                   , Bool_t force) /*FOLD00*/
731 {
732     //
733     // return pointer to TArrayF Mean or Squares
734     // if force is true create a new TArrayF if it doesn't exist allready
735     //
736     if ( !force || (((TArrayF *)arr[det])->GetSize()>0))
737         return (TArrayF*)arr[det];
738
739     // if we are forced and TArrayF doesn't yes exist create it
740     Int_t stack  = GetStack(det);
741     Int_t ngroup = 0;
742     if(stack == 2) ngroup = fDetCha2[1]*fTimeMax;
743     else ngroup = fDetCha0[1]*fTimeMax;
744     // init
745     ((TArrayF *)arr[det])->Set(ngroup);
746     for(Int_t k = 0; k < ngroup; k++){
747       ((TArrayF *)arr[det])->AddAt(0.0,k);
748     }
749     return ((TArrayF *)arr[det]);
750 }
751 //_____________________________________________________________________
752 TArrayF* AliTRDCalibraVector::GetMeanSquaresPRF(Int_t det
753                                                    , TArrayF** arr
754                                                    , Bool_t force) /*FOLD00*/
755 {
756     //
757     // return pointer to TArrayF Mean or Squares
758     // if force is true create a new TArrayF if it doesn't exist allready
759     //
760   if ( !force || (((TArrayF *)arr[det])->GetSize()>0))
761     return arr[det];
762   
763   // if we are forced and TArrayF doesn't yes exist create it
764   Int_t stack  = GetStack(det);
765   Int_t ngroup = 0;
766   if(stack == 2) ngroup = fDetCha2[2]*fNumberBinPRF;
767   else ngroup = fDetCha0[2]*fNumberBinPRF;
768   // init
769   ((TArrayF *)arr[det])->Set(ngroup);
770   for(Int_t k = 0; k < ngroup; k++){
771     ((TArrayF *)arr[det])->AddAt(0.0,k);
772   }
773   return ((TArrayF *)arr[det]);
774 }
775 //_____________________________________________________________________________
776 Int_t AliTRDCalibraVector::GetStack(Int_t d) const
777 {
778   //
779   // Reconstruct the stack number from the detector number
780   //
781
782   return ((Int_t) (d % 30) / 6);
783
784 }
785