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