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