1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 ////////////////////////////////////////////////////////////////////////////
20 // AliTRDCalibraVector //
22 // This class is for the vector method of the TRD calibration. //
25 // R. Bailhache (R.Bailhache@gsi.de) //
27 ////////////////////////////////////////////////////////////////////////////
29 #include <TGraphErrors.h>
31 #include <TObjArray.h>
33 #include <TDirectory.h>
39 #include "AliTRDCalibraVector.h"
40 #include "AliTRDCommonParam.h"
41 #include "AliTRDarrayF.h"
42 #include "AliTRDarrayI.h"
44 ClassImp(AliTRDCalibraVector)
46 //______________________________________________________________________________________
47 AliTRDCalibraVector::AliTRDCalibraVector()
60 ,fVectorPHSquares(540)
61 ,fVectorPHEntries(540)
62 ,fVectorCHEntries(540)
64 ,fVectorPRFSquares(540)
65 ,fVectorPRFEntries(540)
72 // Default constructor
75 for(Int_t k = 0; k < 3; k++){
81 //______________________________________________________________________________________
82 AliTRDCalibraVector::AliTRDCalibraVector(const AliTRDCalibraVector &c)
95 ,fVectorPHSquares(540)
96 ,fVectorPHEntries(540)
97 ,fVectorCHEntries(540)
99 ,fVectorPRFSquares(540)
100 ,fVectorPRFEntries(540)
101 ,fNumberBinCharge(c.fNumberBinCharge)
102 ,fNumberBinPRF(c.fNumberBinPRF)
103 ,fTimeMax(c.fTimeMax)
104 ,fPRFRange(c.fPRFRange)
109 for (Int_t idet = 0; idet < 540; idet++){
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);
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);
119 const AliTRDarrayI *chEntries = (AliTRDarrayI*)c.fVectorCHEntries.UncheckedAt(idet);
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);
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);
130 if ( chEntries != 0x0 ) fVectorCHEntries.AddAt(new AliTRDarrayI(*chEntries), idet);
134 for(Int_t k = 0; k < 3; k++){
135 fDetCha0[k] = c.fDetCha0[k];
136 fDetCha2[k] = c.fDetCha2[k];
139 fVectorPHEntries.SetName(c.fVectorPHEntries.GetName());
140 fVectorCHEntries.SetName(c.fVectorCHEntries.GetName());
141 fVectorPRFEntries.SetName(c.fVectorPRFEntries.GetName());
144 //_____________________________________________________________________
145 AliTRDCalibraVector& AliTRDCalibraVector::operator = (const AliTRDCalibraVector &source)
148 // assignment operator
150 if (&source == this) return *this;
151 new (this) AliTRDCalibraVector(source);
155 //____________________________________________________________________________________
156 AliTRDCalibraVector::~AliTRDCalibraVector()
159 // AliTRDCalibraVector destructor
161 fVectorPHMean.Delete();
162 fVectorPHSquares.Delete();
163 fVectorPHEntries.Delete();
164 fVectorCHEntries.Delete();
165 fVectorPRFMean.Delete();
166 fVectorPRFSquares.Delete();
167 fVectorPRFEntries.Delete();
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;
178 //_____________________________________________________________________________
179 Int_t AliTRDCalibraVector::SearchBin(Float_t value, Int_t i) const
187 Float_t fbinmax = value;
188 Int_t fNumberOfBin = -1;
195 fNumberOfBin = fNumberBinCharge;
199 fbinmax = TMath::Abs(fPRFRange);
200 fbinmin = -TMath::Abs(fPRFRange);
201 fNumberOfBin = fNumberBinPRF;
209 if ((value >= fbinmax) ||
214 reponse = (Int_t) ((fNumberOfBin*(value-fbinmin)) / (fbinmax-fbinmin));
220 //_____________________________________________________________________________
221 Bool_t AliTRDCalibraVector::UpdateVectorCH(Int_t det, Int_t group, Float_t value)
224 // Fill the vector CH
228 Int_t bin = SearchBin(value,0);
236 if(fDetectorCH != det){
237 fCHEntries = ((AliTRDarrayI *)GetCHEntries(det,kTRUE));
240 Int_t entries = fCHEntries->At(group*fNumberBinCharge+bin);
242 Int_t entriesn = entries+1;
243 fCHEntries->AddAt(entriesn,group*fNumberBinCharge+bin);
251 //_____________________________________________________________________________
252 Bool_t AliTRDCalibraVector::UpdateVectorPRF(Int_t det, Int_t group, Float_t x, Float_t y)
255 // Fill the vector PRF
259 Int_t bin = SearchBin(x,2);
266 if(fDetectorPRF != det){
267 fPRFEntries = ((AliTRDarrayI *)GetPRFEntries(det,kTRUE));
268 fPRFMean = ((AliTRDarrayF *)GetPRFMean(det,kTRUE));
269 fPRFSquares = ((AliTRDarrayF *)GetPRFSquares(det,kTRUE));
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);
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);
289 //_____________________________________________________________________________
290 Bool_t AliTRDCalibraVector::UpdateVectorPH(Int_t det, Int_t group, Int_t time, Float_t value)
293 // Fill the vector PH
305 if(fDetectorPH != det){
306 fPHEntries = ((AliTRDarrayI *)GetPHEntries(det,kTRUE));
307 fPHMean = ((AliTRDarrayF *)GetPHMean(det,kTRUE));
308 fPHSquares = ((AliTRDarrayF *)GetPHSquares(det,kTRUE));
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);
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);
328 //__________________________________________________________________________________
329 Bool_t AliTRDCalibraVector::Add(AliTRDCalibraVector *calvect)
332 // Add a other AliTRCalibraVector to this one
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;
345 //printf("pass0!\n");
348 for (Int_t idet = 0; idet < 540; idet++){
350 const AliTRDarrayI *phEntriesvect = calvect->GetPHEntries(idet);
351 const AliTRDarrayF *phMeanvect = calvect->GetPHMean(idet);
352 const AliTRDarrayF *phSquaresvect = calvect->GetPHSquares(idet);
354 const AliTRDarrayI *prfEntriesvect = calvect->GetPRFEntries(idet);
355 const AliTRDarrayF *prfMeanvect = calvect->GetPRFMean(idet);
356 const AliTRDarrayF *prfSquaresvect = calvect->GetPRFSquares(idet);
358 const AliTRDarrayI *chEntriesvect = calvect->GetCHEntries(idet);
360 //printf("idet %d!\n",idet);
362 //printf("phEntriesvect %d\n",(Bool_t) phEntriesvect);
363 //printf("phMeanvect %d\n",(Bool_t) phMeanvect);
364 //printf("phSquaresvect %d\n",(Bool_t) phSquaresvect);
366 //printf("prfEntriesvect %d\n",(Bool_t) prfEntriesvect);
367 //printf("prfMeanvect %d\n",(Bool_t) prfMeanvect);
368 //printf("prfSquaresvect %d\n",(Bool_t) prfSquaresvect);
370 //printf("chEntriesvect %d\n",(Bool_t) chEntriesvect);
372 if ( phEntriesvect != 0x0 ){
374 fPHEntries = ((AliTRDarrayI *)GetPHEntries(idet,kTRUE));
375 fPHMean = ((AliTRDarrayF *)GetPHMean(idet,kTRUE));
376 fPHSquares = ((AliTRDarrayF *)GetPHSquares(idet,kTRUE));
377 Int_t total = fPHEntries->GetSize();
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);
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");
396 if ( prfEntriesvect != 0x0 ){
398 fPRFEntries = ((AliTRDarrayI *)GetPRFEntries(idet,kTRUE));
399 fPRFMean = ((AliTRDarrayF *)GetPRFMean(idet,kTRUE));
400 fPRFSquares = ((AliTRDarrayF *)GetPRFSquares(idet,kTRUE));
401 Int_t total = fPRFEntries->GetSize();
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);
408 //printf("entries0 %d\n",entries);
409 //printf("mean0 %f\n",mean);
410 //printf("square0 %f\n",square);
412 //printf("entries1 %d\n",prfEntriesvect->At(k));
413 //printf("mean1 %f\n",prfMeanvect->At(k));
414 //printf("square1 %f\n",prfSquaresvect->At(k));
417 //printf("entries0 size %d\n",fPRFEntries->GetSize());
418 //printf("mean0 size %d\n",fPRFMean->GetSize());
419 //printf("squares0 size %d\n",fPRFSquares->GetSize());
421 //printf("entries1 size %d\n",prfEntriesvect->GetSize());
422 //printf("mean1 size %d\n",prfMeanvect->GetSize());
423 //printf("squares1 size %d\n",prfSquaresvect->GetSize());
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");
438 if ( chEntriesvect != 0x0 ){
440 fCHEntries = ((AliTRDarrayI *)GetCHEntries(idet,kTRUE));
441 Int_t total = fCHEntries->GetSize();
442 //if(idet == 180) printf("total %d\n",total);
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);
451 //printf("test ch!\n");
457 //_____________________________________________________________________________
458 TGraphErrors *AliTRDCalibraVector::ConvertVectorPHTGraphErrors(Int_t det, Int_t group
459 , const Char_t * name)
462 // Convert the fVectorPHMean, fVectorPHSquares and fVectorPHEntries in TGraphErrors
466 fPHEntries = ((AliTRDarrayI *)GetPHEntries(det,kTRUE));
467 fPHMean = ((AliTRDarrayF *)GetPHMean(det,kTRUE));
468 fPHSquares = ((AliTRDarrayF *)GetPHSquares(det,kTRUE));
474 AliTRDCommonParam *parCom = AliTRDCommonParam::Instance();
476 AliInfo("Could not get CommonParam, take the default 10MHz");
478 sf = parCom->GetSamplingFrequency();
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
492 Int_t offset = group*fTimeMax;
495 for (Int_t k = 0; k < fTimeMax; k++) {
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);
509 // Define the TGraphErrors
510 histo = new TGraphErrors(fTimeMax,x,y,ex,ey);
511 histo->SetTitle(name);
515 //_____________________________________________________________________________
516 TGraphErrors *AliTRDCalibraVector::ConvertVectorPRFTGraphErrors(Int_t det, Int_t group
517 , const Char_t * name)
520 // Convert the fVectorPRFMean, fVectorPRFSquares and fVectorPRFEntries in TGraphErrors
524 fPRFEntries = ((AliTRDarrayI *)GetPRFEntries(det,kTRUE));
525 fPRFMean = ((AliTRDarrayF *)GetPRFMean(det,kTRUE));
526 fPRFSquares = ((AliTRDarrayF *)GetPRFSquares(det,kTRUE));
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);
547 for (Int_t k = 0; k < fNumberBinPRF; k++) {
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);
559 //printf("Number of entries %f for %d\n",ex[k],k);
562 // Define the TGraphErrors
563 histo = new TGraphErrors(fNumberBinPRF,x,y,ex,ey);
564 histo->SetTitle(name);
568 //_____________________________________________________________________________
569 TH1F *AliTRDCalibraVector::ConvertVectorCHHisto(Int_t det, Int_t group
570 , const Char_t * name)
573 // Convert the fVectorCHEntries in TH1F
577 fCHEntries = ((AliTRDarrayI *)GetCHEntries(det,kTRUE));
580 TH1F *histo = new TH1F(name,name,fNumberBinCharge,0,300);
582 Int_t offset = group*fNumberBinCharge;
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))));
593 //_____________________________________________________________________
594 AliTRDarrayI* AliTRDCalibraVector::GetPHEntries(Int_t det
595 , Bool_t force) /*FOLD00*/
598 // return pointer to Carge ROC Calibration
599 // if force is true create a new histogram if it doesn't exist allready
601 TObjArray *arr = &fVectorPHEntries;
602 return GetEntriesPH(det, arr, force);
604 //_____________________________________________________________________
605 AliTRDarrayI* AliTRDCalibraVector::GetPRFEntries(Int_t det
606 , Bool_t force) /*FOLD00*/
609 // return pointer to Carge ROC Calibration
610 // if force is true create a new histogram if it doesn't exist allready
612 TObjArray *arr = &fVectorPRFEntries;
613 return GetEntriesPRF(det, arr, force);
615 //_____________________________________________________________________
616 AliTRDarrayI* AliTRDCalibraVector::GetCHEntries(Int_t det
617 , Bool_t force) /*FOLD00*/
620 // return pointer to Carge ROC Calibration
621 // if force is true create a new histogram if it doesn't exist allready
623 TObjArray *arr = &fVectorCHEntries;
624 return GetEntriesCH(det, arr, force);
626 //_____________________________________________________________________
627 AliTRDarrayF* AliTRDCalibraVector::GetPHMean(Int_t det
628 , Bool_t force) /*FOLD00*/
631 // return pointer to Carge ROC Calibration
632 // if force is true create a new histogram if it doesn't exist allready
634 TObjArray *arr = &fVectorPHMean;
635 return GetMeanSquaresPH(det, arr, force);
637 //_____________________________________________________________________
638 AliTRDarrayF* AliTRDCalibraVector::GetPHSquares(Int_t det
639 , Bool_t force) /*FOLD00*/
642 // return pointer to Carge ROC Calibration
643 // if force is true create a new histogram if it doesn't exist allready
645 TObjArray *arr = &fVectorPHSquares;
646 return GetMeanSquaresPH(det, arr, force);
648 //_____________________________________________________________________
649 AliTRDarrayF* AliTRDCalibraVector::GetPRFMean(Int_t det
650 , Bool_t force) /*FOLD00*/
653 // return pointer to Carge ROC Calibration
654 // if force is true create a new histogram if it doesn't exist allready
656 TObjArray *arr = &fVectorPRFMean;
657 return GetMeanSquaresPRF(det, arr, force);
659 //_____________________________________________________________________
660 AliTRDarrayF* AliTRDCalibraVector::GetPRFSquares(Int_t det
661 , Bool_t force) /*FOLD00*/
664 // return pointer to Carge ROC Calibration
665 // if force is true create a new histogram if it doesn't exist allready
667 TObjArray *arr = &fVectorPRFSquares;
668 return GetMeanSquaresPRF(det, arr, force);
670 //_____________________________________________________________________
671 AliTRDarrayI* AliTRDCalibraVector::GetEntriesCH(Int_t det
673 , Bool_t force) /*FOLD00*/
676 // return pointer to AliTRDarrayI Entries
677 // if force is true create a new AliTRDarrayI if it doesn't exist allready
679 if ( !force || arr->UncheckedAt(det) )
680 return (AliTRDarrayI*)arr->UncheckedAt(det);
682 // if we are forced and AliTRDarrayI doesn't yes exist create it
683 AliTRDarrayI *croc = new AliTRDarrayI();
684 Int_t stack = GetStack(det);
686 if(stack == 2) ngroup = fDetCha2[0]*fNumberBinCharge;
687 else ngroup = fDetCha0[0]*fNumberBinCharge;
689 croc->Expand(ngroup);
690 for(Int_t k = 0; k < ngroup; k++){
693 arr->AddAt(croc,det);
696 //_____________________________________________________________________
697 AliTRDarrayI* AliTRDCalibraVector::GetEntriesPRF(Int_t det
699 , Bool_t force) /*FOLD00*/
702 // return pointer to AliTRDarrayI Entries
703 // if force is true create a new AliTRDarrayI if it doesn't exist allready
705 if ( !force || arr->UncheckedAt(det) )
706 return (AliTRDarrayI*)arr->UncheckedAt(det);
708 // if we are forced and AliTRDarrayI doesn't yes exist create it
709 AliTRDarrayI *croc = new AliTRDarrayI();
710 Int_t stack = GetStack(det);
712 if(stack == 2) ngroup = fDetCha2[2]*fNumberBinPRF;
713 else ngroup = fDetCha0[2]*fNumberBinPRF;
715 croc->Expand(ngroup);
716 for(Int_t k = 0; k < ngroup; k++){
719 arr->AddAt(croc,det);
723 //_____________________________________________________________________
724 AliTRDarrayI* AliTRDCalibraVector::GetEntriesPH(Int_t det
726 , Bool_t force) /*FOLD00*/
729 // return pointer to AliTRDarrayI Entries
730 // if force is true create a new AliTRDarrayI if it doesn't exist allready
732 if ( !force || arr->UncheckedAt(det) )
733 return (AliTRDarrayI*)arr->UncheckedAt(det);
735 // if we are forced and AliTRDarrayI doesn't yes exist create it
736 AliTRDarrayI *croc = new AliTRDarrayI();
737 Int_t stack = GetStack(det);
739 if(stack == 2) ngroup = fDetCha2[1]*fTimeMax;
740 else ngroup = fDetCha0[1]*fTimeMax;
742 croc->Expand(ngroup);
743 for(Int_t k = 0; k < ngroup; k++){
746 arr->AddAt(croc,det);
750 //_____________________________________________________________________
751 AliTRDarrayF* AliTRDCalibraVector::GetMeanSquaresPH(Int_t det
753 , Bool_t force) /*FOLD00*/
756 // return pointer to AliTRDarrayF Mean or Squares
757 // if force is true create a new AliTRDarrayF if it doesn't exist allready
759 if ( !force || arr->UncheckedAt(det) )
760 return (AliTRDarrayF*)arr->UncheckedAt(det);
762 // if we are forced and AliTRDarrayF doesn't yes exist create it
763 AliTRDarrayF *croc = new AliTRDarrayF();
764 Int_t stack = GetStack(det);
766 if(stack == 2) ngroup = fDetCha2[1]*fTimeMax;
767 else ngroup = fDetCha0[1]*fTimeMax;
769 croc->Expand(ngroup);
770 for(Int_t k = 0; k < ngroup; k++){
773 arr->AddAt(croc,det);
776 //_____________________________________________________________________
777 AliTRDarrayF* AliTRDCalibraVector::GetMeanSquaresPRF(Int_t det
779 , Bool_t force) /*FOLD00*/
782 // return pointer to AliTRDarrayF Mean or Squares
783 // if force is true create a new AliTRDarrayF if it doesn't exist allready
785 if ( !force || arr->UncheckedAt(det) )
786 return (AliTRDarrayF*)arr->UncheckedAt(det);
788 // if we are forced and AliTRDarrayF doesn't yes exist create it
789 AliTRDarrayF *croc = new AliTRDarrayF();
790 Int_t stack = GetStack(det);
792 if(stack == 2) ngroup = fDetCha2[2]*fNumberBinPRF;
793 else ngroup = fDetCha0[2]*fNumberBinPRF;
795 croc->Expand(ngroup);
796 for(Int_t k = 0; k < ngroup; k++){
799 arr->AddAt(croc,det);
802 //_____________________________________________________________________________
803 Int_t AliTRDCalibraVector::GetStack(Int_t d) const
806 // Reconstruct the stack number from the detector number
809 return ((Int_t) (d % 30) / 6);