ROOT logo
/**************************************************************************
 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
 *                                                                        *
 * Author: The ALICE Off-line Project.                                    *
 * Contributors are mentioned in the code where appropriate.              *
 *                                                                        *
 * Permission to use, copy, modify and distribute this software and its   *
 * documentation strictly for non-commercial purposes is hereby granted   *
 * without fee, provided that the above copyright notice appears in all   *
 * copies and that both the copyright notice and this permission notice   *
 * appear in the supporting documentation. The authors make no claims     *
 * about the suitability of this software for any purpose. It is          *
 * provided "as is" without express or implied warranty.                  *
 **************************************************************************/

/* $Id$ */

////////////////////////////////////////////////////////////////////////////
//                                                                        //   
// AliTRDCalibraVector                                                    //       
//                                                                        //   
// This class is for the vector method of the TRD calibration.            //
//                                                                        //
// Author:                                                                //
//   R. Bailhache (R.Bailhache@gsi.de)                                    //
//                                                                        //
////////////////////////////////////////////////////////////////////////////

#include <TGraphErrors.h>
#include <TH1F.h>
#include <TObjArray.h>
#include <TObject.h>
#include <TMath.h>
#include <TDirectory.h>
#include <TROOT.h>
#include <TFile.h>
#include <TString.h>

#include "AliLog.h"

#include "AliTRDCalibraVector.h"
#include "AliTRDCommonParam.h"
#include "AliTRDCalibraMode.h"
#include "AliTRDPhInfo.h"
#include "AliTRDEntriesInfo.h"
#include "AliTRDPrfInfo.h"
#include "AliTRDgeometry.h"

ClassImp(AliTRDCalibraVector)

//______________________________________________________________________________________
AliTRDCalibraVector::AliTRDCalibraVector()
  :TObject()
  ,fModeCH(0)
  ,fModePH(0)
  ,fModePRF(0)
  ,fNbGroupPRF(0)
  ,fDetectorPH(-1)
  ,fDetectorCH(-1)
  ,fDetectorPRF(-1)
  ,fStopFillCH(kFALSE)
  ,fHisto(0x0)
  ,fGraph(0x0)
  ,fCalVector(0x0)
  ,fNumberBinCharge(0)
  ,fNumberBinPRF(0)
  ,fTimeMax(0)
  ,fPRFRange(1.5)
{
  //
  // Default constructor
  //

  for (Int_t idet = 0; idet < 540; idet++){
    
    fPHEntries[idet]= 0x0;
    fPHMean[idet]= 0x0;
    fPHSquares[idet]= 0x0;

    fPRFEntries[idet]= 0x0;
    fPRFMean[idet]= 0x0;
    fPRFSquares[idet]= 0x0;


    fCHEntries[idet]= 0x0;
    
  }
  
  for(Int_t k = 0; k < 3; k++){
    fDetCha0[k] = 0;
    fDetCha2[k] = 0;
  }
 
}
//______________________________________________________________________________________
AliTRDCalibraVector::AliTRDCalibraVector(const AliTRDCalibraVector &c)
  :TObject(c)
  ,fModeCH(c.fModeCH)
  ,fModePH(c.fModePH)
  ,fModePRF(c.fModePRF)
  ,fNbGroupPRF(c.fNbGroupPRF)
  ,fDetectorPH(-1)
  ,fDetectorCH(-1)
  ,fDetectorPRF(-1)
  ,fStopFillCH(kFALSE)
  ,fHisto(0x0)
  ,fGraph(0x0)
  ,fCalVector(0x0)
  ,fNumberBinCharge(c.fNumberBinCharge)
  ,fNumberBinPRF(c.fNumberBinPRF)
  ,fTimeMax(c.fTimeMax)
  ,fPRFRange(c.fPRFRange)
{
  //
  // Copy constructor
  //
  
  for(Int_t k = 0; k < 3; k++){
    fDetCha0[k] = c.fDetCha0[k];
    fDetCha2[k] = c.fDetCha2[k];
  }

  for (Int_t idet = 0; idet < 540; idet++){
    
    const AliTRDEntriesInfo *phEntries  = (AliTRDEntriesInfo*)c.fPHEntries[idet];
    const AliTRDPhInfo      *phMean     = (AliTRDPhInfo *)c.fPHMean[idet];
    const AliTRDPhInfo      *phSquares  = (AliTRDPhInfo *)c.fPHSquares[idet];

    const AliTRDEntriesInfo *prfEntries = (AliTRDEntriesInfo*)c.fPRFEntries[idet];
    const AliTRDPrfInfo     *prfMean     = (AliTRDPrfInfo *)c.fPRFMean[idet];
    const AliTRDPrfInfo     *prfSquares  = (AliTRDPrfInfo *)c.fPRFSquares[idet];

    const AliTRDEntriesInfo *chEntries  = (AliTRDEntriesInfo*)c.fCHEntries[idet];

    if ( chEntries != 0x0 ) fCHEntries[idet] = new AliTRDEntriesInfo(*chEntries);
    
    if ( phEntries != 0x0 ) {
      fPHMean[idet]    = new AliTRDPhInfo(*phMean);
      fPHSquares[idet] = new AliTRDPhInfo(*phSquares);
      fPHEntries[idet] = new AliTRDEntriesInfo(*phEntries);
    }

    if ( prfEntries != 0x0 ) {
      fPRFEntries[idet] = new AliTRDEntriesInfo(*prfEntries);
      fPRFMean[idet] = new AliTRDPrfInfo(*prfMean);
      fPRFSquares[idet] = new AliTRDPrfInfo(*prfSquares);
    }
    
  }
   
}
//_____________________________________________________________________
AliTRDCalibraVector& AliTRDCalibraVector::operator = (const  AliTRDCalibraVector &source)
{
  //
  // assignment operator
  //
  if (&source == this) return *this;
  new (this) AliTRDCalibraVector(source);

  return *this;
}
//____________________________________________________________________________________
AliTRDCalibraVector::~AliTRDCalibraVector()
{
  //
  // AliTRDCalibraVector destructor
  //

  for (Int_t i=0; i<540; i++) {
    delete fPHEntries[i];
    delete fPHMean[i];
    delete fPHSquares[i];
    delete fPRFEntries[i];
    delete fPRFMean[i];
    delete fPRFSquares[i];
    delete fCHEntries[i];
  }

  if(fHisto) delete fHisto;
  if(fGraph) delete fGraph;
  if(fCalVector) delete fCalVector;

}
//_____________________________________________________________________________
Long64_t AliTRDCalibraVector::Merge(const TCollection* list) 
{
  // Merge list of objects (needed by PROOF)

  if (!list)
    return 0;
  
  if (list->IsEmpty())
    return 1;
  
  TIterator* iter = list->MakeIterator();
  TObject* obj = 0;
  
  // collection of generated histograms
  Int_t count=0;
  while((obj = iter->Next()) != 0) 
    {
      AliTRDCalibraVector* entry = dynamic_cast<AliTRDCalibraVector*>(obj);
      if (entry == 0) continue; 
      
      if(this->Add(entry)) count++;
    
    }
  
  return count;
}
//_____________________________________________________________________________
void AliTRDCalibraVector::TestInit(Int_t i, Int_t detmax)
{
  //
  // Init to see the size
  //

  for(Int_t det = 0; det < detmax; det++){

    if(i==2) {
      
      fPRFEntries[det] = ((AliTRDEntriesInfo *)GetPRFEntries(det,kTRUE));
      fPRFMean[det]    = ((AliTRDPrfInfo  *)GetPRFMean(det,kTRUE));
      fPRFSquares[det] = ((AliTRDPrfInfo  *)GetPRFSquares(det,kTRUE));

    }

    if(i==1) {

      fPHEntries[det] = ((AliTRDEntriesInfo *)GetPHEntries(det,kTRUE));
      fPHMean[det]    = ((AliTRDPhInfo  *)GetPHMean(det,kTRUE));
      fPHSquares[det] = ((AliTRDPhInfo  *)GetPHSquares(det,kTRUE));
      
    }
    
    if(i==0) fCHEntries[det] = ((AliTRDEntriesInfo *)GetCHEntries(det,kTRUE));

  }

}
//_____________________________________________________________________________
Int_t AliTRDCalibraVector::SearchBin(Float_t value, Int_t i) const
{
  //
  // Search the bin
  //

  Int_t reponse        = 0;
  Float_t fbinmin      = 0;
  Float_t fbinmax      = value;
  Int_t fNumberOfBin   = -1;

  switch(i)
    {
    case 0:
      fbinmax      = 300.0;
      fbinmin      = 0.0;
      fNumberOfBin = fNumberBinCharge;
      break;

    case 2:
      fbinmax      =   TMath::Abs(fPRFRange);
      fbinmin      =  -TMath::Abs(fPRFRange);
      fNumberOfBin = fNumberBinPRF;
      break;
      
    default: 
      return -1;
    }
  
  // Return -1 if out
  if ((value >= fbinmax) || 
      (value <  fbinmin)) {
    return -1;
  }
  else {
    reponse = (Int_t) ((fNumberOfBin*(value-fbinmin)) / (fbinmax-fbinmin));
  }

  return reponse;

}
//_____________________________________________________________________________
Bool_t AliTRDCalibraVector::UpdateVectorCH(Int_t det, Int_t group, Float_t value)
{
  //
  // Fill the vector CH   
  //

  // Search bin
  Int_t bin = SearchBin(value,0);
  // Out
  if (bin == -1) {
    return kFALSE; 
  }



  if(fDetectorCH != det){
    fCHEntries[det] = ((AliTRDEntriesInfo *)GetCHEntries(det,kTRUE));
  }

  Int_t entries  = ((AliTRDEntriesInfo *)fCHEntries[det])->At(group*fNumberBinCharge+bin);
  
  Int_t entriesn = entries+1;

  if(entriesn > 65535) {
    fStopFillCH = kTRUE;
    return kTRUE;
  }

  ((AliTRDEntriesInfo *)fCHEntries[det])->AddAt(entriesn,group*fNumberBinCharge+bin);
  
  fDetectorCH = det;

 
  return kTRUE;

}
//_____________________________________________________________________________
Bool_t AliTRDCalibraVector::UpdateVectorPRF(Int_t det, Int_t group, Float_t x, Float_t y)
{
  //
  // Fill the vector PRF
  //

  // Search bin
  Int_t bin = SearchBin(x,2);
  // Out
  if (bin == -1) {
    return kFALSE; 
  }

  
  if(fDetectorPRF != det){
    fPRFEntries[det] = ((AliTRDEntriesInfo *)GetPRFEntries(det,kTRUE));
    fPRFMean[det]    = ((AliTRDPrfInfo  *)GetPRFMean(det,kTRUE));
    fPRFSquares[det] = ((AliTRDPrfInfo  *)GetPRFSquares(det,kTRUE));
  }

  Int_t entries  = ((AliTRDEntriesInfo *)fPRFEntries[det])->At((Int_t)group*fNumberBinPRF+bin);
  Float_t mean   = ((AliTRDPrfInfo  *)fPRFMean[det])->At((Int_t)group*fNumberBinPRF+bin);
  Float_t square = ((AliTRDPrfInfo  *)fPRFSquares[det])->At((Int_t)group*fNumberBinPRF+bin);
  
  Int_t entriesn = entries+1;

  if(entriesn > 65535) return kTRUE;
  
  ((AliTRDEntriesInfo *)fPRFEntries[det])->AddAt(entriesn,(Int_t)group*fNumberBinPRF+bin);
  Float_t meann = (mean*((Float_t)entries)+y)/((Float_t)entriesn);
  ((AliTRDPrfInfo *)fPRFMean[det])->AddAt(meann,(Int_t)group*fNumberBinPRF+bin);
  Float_t squaren = ((square*((Float_t)entries))+(y*y))/((Float_t)entriesn);
  ((AliTRDPrfInfo *)fPRFSquares[det])->AddAt(squaren,(Int_t)group*fNumberBinPRF+bin);

  
  fDetectorPRF = det;

  return kTRUE;
  
}
//_____________________________________________________________________________
Bool_t AliTRDCalibraVector::UpdateVectorPH(Int_t det, Int_t group, Int_t time, Float_t value)
{
  //
  // Fill the vector PH  
  //

  // Search bin
  Int_t bin = time;
  // Out
  if ((bin <         0) || 
      (bin >= fTimeMax)) {
    return kFALSE; 
  }


  if(fDetectorPH != det){
    fPHEntries[det] = ((AliTRDEntriesInfo *)GetPHEntries(det,kTRUE));
    fPHMean[det]    = ((AliTRDPhInfo  *)GetPHMean(det,kTRUE));
    fPHSquares[det] = ((AliTRDPhInfo  *)GetPHSquares(det,kTRUE));
  }

  Int_t entries  = ((AliTRDEntriesInfo *)fPHEntries[det])->At(group*fTimeMax+bin);
  Float_t mean   = ((AliTRDPhInfo  *)fPHMean[det])->At(group*fTimeMax+bin);
  Float_t square = ((AliTRDPhInfo  *)fPHSquares[det])->AtS(group*fTimeMax+bin);
  
  Int_t entriesn = entries+1;
  Float_t meann = (mean*((Float_t)entries)+value)/((Float_t)entriesn);
  Float_t squaren = ((square*((Float_t)entries))+(value*value))/((Float_t)entriesn);
  
  if(entriesn > 65535) return kTRUE;
  //printf("meann %f, squaren %f\n",meann,squaren);
  if((meann > 3000.0) || (meann < 0.0) || (squaren > (3000.0*3000.0)) || (squaren < 0.0)) return kFALSE;

  ((AliTRDEntriesInfo *)fPHEntries[det])->AddAt(entriesn,group*fTimeMax+bin);
  ((AliTRDPhInfo *)fPHMean[det])->AddAt(meann,group*fTimeMax+bin);
  ((AliTRDPhInfo *)fPHSquares[det])->AddAtS(squaren,group*fTimeMax+bin);
  
  fDetectorPH = det;

  return kTRUE;
  
}
//_____________________________________________________________________________
Bool_t AliTRDCalibraVector::FillVectorCH(Int_t det, Int_t group, Int_t bin, Int_t entries)
{
  //
  // Fill the vector CH   
  //

  if(entries > 65535) return kFALSE;

  if(fDetectorCH != det){
    fCHEntries[det] = ((AliTRDEntriesInfo *)GetCHEntries(det,kTRUE));
  }

  ((AliTRDEntriesInfo *)fCHEntries[det])->AddAt(entries,group*fNumberBinCharge+bin);
  
  fDetectorCH = det;
  
  
  return kTRUE;
  
}
//_____________________________________________________________________________
Bool_t AliTRDCalibraVector::FillVectorPRF(Int_t det, Int_t group, Int_t bin, Int_t entries, Float_t mean, Float_t square)
{
  //
  // Fill the vector PRF
  //

  if((entries > 65535) || (mean > 1.0) || (mean < 0.0) || (square > 1.0) || (square < 0.0)) return kFALSE;

  if(fDetectorPRF != det){
    fPRFEntries[det] = ((AliTRDEntriesInfo *)GetPRFEntries(det,kTRUE));
    fPRFMean[det]    = ((AliTRDPrfInfo  *)GetPRFMean(det,kTRUE));
    fPRFSquares[det] = ((AliTRDPrfInfo  *)GetPRFSquares(det,kTRUE));
  }

  ((AliTRDEntriesInfo *)fPRFEntries[det])->AddAt(entries,(Int_t)group*fNumberBinPRF+bin);
  ((AliTRDPrfInfo  *)fPRFMean[det])->AddAt(mean,(Int_t)group*fNumberBinPRF+bin);
  ((AliTRDPrfInfo  *)fPRFSquares[det])->AddAt(square,(Int_t)group*fNumberBinPRF+bin);

  
  fDetectorPRF = det;

  return kTRUE;
  
}
//_____________________________________________________________________________
Bool_t AliTRDCalibraVector::FillVectorPH(Int_t det, Int_t group, Int_t bin, Int_t entries, Float_t mean, Float_t square)
{
  //
  // Fill the vector PH  
  //

  if((entries > 65535) || (mean > 3000.0) || (mean < 0.0) || (square > (3000.0*3000.0)) || (square < 0.0)) return kFALSE;


  if(fDetectorPH != det){
    fPHEntries[det] = ((AliTRDEntriesInfo *)GetPHEntries(det,kTRUE));
    fPHMean[det]    = ((AliTRDPhInfo  *)GetPHMean(det,kTRUE));
    fPHSquares[det] = ((AliTRDPhInfo  *)GetPHSquares(det,kTRUE));
  }

  ((AliTRDEntriesInfo *)fPHEntries[det])->AddAt(entries,group*fTimeMax+bin);
  ((AliTRDPhInfo  *)fPHMean[det])->AddAt(mean,group*fTimeMax+bin);
  ((AliTRDPhInfo  *)fPHSquares[det])->AddAtS(square,group*fTimeMax+bin);
  
  fDetectorPH = det;

  return kTRUE;
  
}
//__________________________________________________________________________________
Bool_t AliTRDCalibraVector::Add(AliTRDCalibraVector *calvect)
{
  //
  // Add a other AliTRCalibraVector to this one
  //

  Bool_t result = kTRUE;

  // Check compatibility
  if(fNumberBinCharge != calvect->GetNumberBinCharge()) return kFALSE;
  if(fNumberBinPRF    != calvect->GetNumberBinPRF()) return kFALSE;
  if(fPRFRange        != calvect->GetPRFRange()) return kFALSE;
  if(fTimeMax         != calvect->GetTimeMax()) return kFALSE;
  for(Int_t k = 0; k < 3; k++){
    if(fDetCha0[k] != calvect->GetDetCha0(k)) return kFALSE;
    if(fDetCha2[k] != calvect->GetDetCha2(k)) return kFALSE;
  }

  //printf("All ok for variables before adding!\n"); 

  // Add
  for (Int_t idet = 0; idet < 540; idet++){

    //printf("Detector %d\n",idet);
    
    const AliTRDEntriesInfo *phEntriesvect    = (AliTRDEntriesInfo *) calvect->GetPHEntries(idet);
    const AliTRDPhInfo      *phMeanvect       = (AliTRDPhInfo *) calvect->GetPHMean(idet);
    const AliTRDPhInfo      *phSquaresvect    = (AliTRDPhInfo *) calvect->GetPHSquares(idet);
    
    const AliTRDEntriesInfo *prfEntriesvect   = (AliTRDEntriesInfo *) calvect->GetPRFEntries(idet);
    const AliTRDPrfInfo    *prfMeanvect       = (AliTRDPrfInfo *) calvect->GetPRFMean(idet);
    const AliTRDPrfInfo    *prfSquaresvect    = (AliTRDPrfInfo *) calvect->GetPRFSquares(idet);
    
    const AliTRDEntriesInfo *chEntriesvect    = (AliTRDEntriesInfo *) calvect->GetCHEntries(idet);

    if ( phEntriesvect != 0x0 ){
      //Take the stuff
      fPHEntries[idet] = ((AliTRDEntriesInfo *)GetPHEntries(idet,kTRUE));
      fPHMean[idet]    = ((AliTRDPhInfo *)GetPHMean(idet,kTRUE));
      fPHSquares[idet] = ((AliTRDPhInfo *)GetPHSquares(idet,kTRUE));

      Int_t total = fPHEntries[idet]->GetSize();
      //printf("Total size PH %d\n",total);
      // Add
      for(Int_t k = 0; k < total; k++){
	Int_t entriesv  = ((AliTRDEntriesInfo *)phEntriesvect)->At(k);
	Float_t meanv   = ((AliTRDPhInfo  *)phMeanvect)->At(k);
	Float_t squarev = ((AliTRDPhInfo  *)phSquaresvect)->AtS(k);
	
	Int_t entries  = ((AliTRDEntriesInfo *)fPHEntries[idet])->At(k);
	Float_t mean   = ((AliTRDPhInfo  *)fPHMean[idet])->At(k);
	Float_t square = ((AliTRDPhInfo  *)fPHSquares[idet])->AtS(k);
  
	Int_t entriesn = entries+entriesv;
	Float_t meann = (mean*((Float_t)entries)+meanv*((Float_t)entriesv))/((Float_t)entriesn);
	Float_t squaren = ((square*((Float_t)entries))+(squarev*((Float_t)entriesv)))/((Float_t)entriesn);
	
	if((entriesn > 0) && (entriesn <= 65535) && (meann >= 0) && (meann < 3000.0) && (squaren >= 0.0) && (squaren < (3000.0*3000.0))) {
	
	  ((AliTRDEntriesInfo *)fPHEntries[idet])->AddAt(entriesn,k);
	  ((AliTRDPhInfo *)fPHMean[idet])->AddAt(meann,k);
	  ((AliTRDPhInfo *)fPHSquares[idet])->AddAtS(squaren,k);
      
	}
      }
    }     

    if ( prfEntriesvect != 0x0 ){
      //Take the stuff
      fPRFEntries[idet] = ((AliTRDEntriesInfo *)GetPRFEntries(idet,kTRUE));
      fPRFMean[idet]    = ((AliTRDPrfInfo  *)GetPRFMean(idet,kTRUE));
      fPRFSquares[idet] = ((AliTRDPrfInfo  *)GetPRFSquares(idet,kTRUE));
      
      Int_t total = fPRFEntries[idet]->GetSize(); 
      //Int_t total0 = fPRFMean[idet]->GetSize(); 
      //Int_t total1 = fPRFSquares[idet]->GetSize(); 
      //printf("Total size PRF %d\n",total);     
      //printf("Total0 size PRF %d\n",total0);     
      //printf("Total1 size PRF %d\n",total1);     
      // Add
      for(Int_t k = 0; k < total; k++){


	Int_t entries  = ((AliTRDEntriesInfo *)fPRFEntries[idet])->At(k);
	Float_t mean   = ((AliTRDPrfInfo  *)fPRFMean[idet])->At(k);
	Float_t square = ((AliTRDPrfInfo  *)fPRFSquares[idet])->At(k);

	Int_t entriesv  = ((AliTRDEntriesInfo *)prfEntriesvect)->At(k);
	Float_t meanv   = ((AliTRDPrfInfo  *)prfMeanvect)->At(k);
	Float_t squarev = ((AliTRDPrfInfo  *)prfSquaresvect)->At(k);

	Int_t entriesn = entries + entriesv;
	
	if((entriesn > 0) && (entriesn <= 65535)) {
	  
	  ((AliTRDEntriesInfo *)fPRFEntries[idet])->AddAt(entriesn,k);
	  
	  Float_t meann = (mean*((Float_t)entries)+meanv*((Float_t)entriesv))/((Float_t)entriesn);
	  //printf("test0\n");
	  ((AliTRDPrfInfo *)fPRFMean[idet])->AddAt(meann,k);
	  //printf("test1\n");
	  
	  Float_t squaren = ((square*((Float_t)entries))+(squarev*((Float_t)entriesv)))/((Float_t)entriesn);
	  //printf("test2\n");
	  ((AliTRDPrfInfo *)fPRFSquares[idet])->AddAt(squaren,k);
	  //printf("test3\n");
		  
	}
      }
    }

    if ( chEntriesvect != 0x0 ){
      //Take the stuff
      fCHEntries[idet] = ((AliTRDEntriesInfo *)GetCHEntries(idet,kTRUE));
      //printf("TestAdd\n");
      fStopFillCH = ((AliTRDEntriesInfo *)fCHEntries[idet])->TestAdd((AliTRDEntriesInfo *)chEntriesvect);
      //
      if(!fStopFillCH) {
	fStopFillCH = kTRUE;
	result = kFALSE;
      }
      else {
	
	((AliTRDEntriesInfo *)fCHEntries[idet])->Add(chEntriesvect);
	//printf("Add\n");
      }
    }           
  }
  
  return result;
}
//_____________________________________________________________________________________________________________________
AliTRDCalibraVector *AliTRDCalibraVector::AddStatsPerDetectorCH()
{
  //
  // Create a AliTRDCalibraVector detector wise
  //

  // Use a AliTRDCalibraMode to navigate in the calibration groups
  AliTRDCalibraMode calibMode = AliTRDCalibraMode();
  calibMode.SetNz(0,GetNz(0));
  calibMode.SetNrphi(0,GetNrphi(0));
  if(((calibMode.GetNz(0) == 100) && (calibMode.GetNrphi(0) == 100)) || ((calibMode.GetNz(0) == 10) && (calibMode.GetNrphi(0) == 10))) return 0x0;
  
  Int_t nybins = 6*4*18*fDetCha0[0]+6*18*fDetCha2[0];
  Int_t nxbins = fNumberBinCharge;
  
  // Check 
  Int_t perChamber2 = 0;
  Int_t perChamber0 = 0;
  calibMode.ModePadCalibration(2,0);
  calibMode.ModePadFragmentation(0,2,0,0);
  calibMode.SetDetChamb2(0);
  perChamber2 = (Int_t) calibMode.GetDetChamb2(0);
  calibMode.ModePadCalibration(0,0);
  calibMode.ModePadFragmentation(0,0,0,0);
  calibMode.SetDetChamb0(0);
  perChamber0 = (Int_t) calibMode.GetDetChamb0(0);
  if(nybins != (6*18*perChamber2+6*4*18*perChamber0)) return 0x0;
    
  // Create calvector 
  if(!fCalVector) fCalVector = new AliTRDCalibraVector();
  else{ 
    fCalVector->~AliTRDCalibraVector();
    new(fCalVector) AliTRDCalibraVector();
  }
  fCalVector->SetNumberBinCharge(nxbins);
  fCalVector->SetDetCha0(0,1);
  fCalVector->SetDetCha2(0,1);
  fCalVector->SetNzNrphi(0,0,0);
 
  
  for(Int_t det = 0; det < 540; det++){
    
    // Take
    AliTRDEntriesInfo *entriesch = (AliTRDEntriesInfo *)GetCHEntries(det,kFALSE);
    if(!entriesch) continue;  
  
    // Number of groups
    Int_t numberofgroup = 0;
    if(AliTRDgeometry::GetStack(det) == 2) numberofgroup = perChamber2;
    else numberofgroup = perChamber0;
  
    // Check if one can merge calibration groups for this chamber
    // entries is the number of entries in each bin after adding the different the calibration group in the detector
    fStopFillCH = kFALSE;
    Int_t firstnumberofgroup = -1;
    Int_t entries[500];
    for(Int_t k = 0; k < nxbins; k++){
      entries[k] = 0;
    }
    // Loop over group in the detector
    for(Int_t k = 0; k < numberofgroup; k++){
      // Loop over bins
      for(Int_t nx = 0; nx < nxbins; nx++) {
	Int_t binnumber = k*nxbins+nx;
	entries[nx] += entriesch->At(binnumber);
	// as soon as one bin is over threshold stop 
	if(!fStopFillCH) {
	  if(entries[nx] > 65535) {
	    firstnumberofgroup = k;
	    fStopFillCH = kTRUE;
	  }
	}
	else continue;
      }
    }
    if(fStopFillCH && (firstnumberofgroup == 0)) return 0x0;
    if(!fStopFillCH) firstnumberofgroup = numberofgroup;
    
    // Now add up to possible 
    for(Int_t k = 0; k < nxbins; k++){
      entries[k] = 0;
    }
    for(Int_t k = 0; k < firstnumberofgroup; k++){
      for(Int_t nx = 0; nx < nxbins; nx++) {
	Int_t binnumber = k*nxbins+nx;
	entries[nx] += entriesch->At(binnumber);
      }
    }

    // Finally fill
    for(Int_t nx = 0; nx < nxbins; nx++){
      fCalVector->FillVectorCH(det,0,nx,(Int_t)entries[nx]);  
    }
  }
  
  return fCalVector;
} 
//_____________________________________________________________________________________________________________________
AliTRDCalibraVector *AliTRDCalibraVector::AddStatsPerDetectorPH()
{
  //
  // Create a AliTRDCalibraVector detector wise
  //

  AliTRDCalibraMode calibMode = AliTRDCalibraMode();
  calibMode.SetNz(1,GetNz(1));
  calibMode.SetNrphi(1,GetNrphi(1));
  if(((calibMode.GetNz(1) == 100) && (calibMode.GetNrphi(1) == 100)) || ((calibMode.GetNz(1) == 10) && (calibMode.GetNrphi(1) == 10))) return 0x0;
  

  // Check
  Int_t nybins = 6*4*18*fDetCha0[1]+6*18*fDetCha2[1];
  Int_t nxbins = fTimeMax;
 
  Int_t perChamber2 = 0;
  Int_t perChamber0 = 0;
  calibMode.ModePadCalibration(2,1);
  calibMode.ModePadFragmentation(0,2,0,1);
  calibMode.SetDetChamb2(1);
  perChamber2 = (Int_t) calibMode.GetDetChamb2(1);
  calibMode.ModePadCalibration(0,1);
  calibMode.ModePadFragmentation(0,0,0,1);
  calibMode.SetDetChamb0(1);
  perChamber0 = (Int_t) calibMode.GetDetChamb0(1);
  
  if(nybins != (6*18*perChamber2+6*4*18*perChamber0)) return 0x0;
  
  // Create calvector 
  if(!fCalVector) fCalVector = new AliTRDCalibraVector();
  else{ 
    fCalVector->~AliTRDCalibraVector();
    new(fCalVector) AliTRDCalibraVector();
  }
  fCalVector->SetTimeMax(nxbins);
  fCalVector->SetDetCha0(1,1);
  fCalVector->SetDetCha2(1,1);
  fCalVector->SetNzNrphi(1,0,0);
  
 
  for(Int_t det = 0; det < 540; det++){
    
    // Take
    AliTRDEntriesInfo *entriesph = (AliTRDEntriesInfo *)GetPHEntries(det,kFALSE);
    if(!entriesph) continue;
    AliTRDPhInfo      *meanph    = (AliTRDPhInfo *)GetPHMean(det,kFALSE);
    AliTRDPhInfo      *squaresph = (AliTRDPhInfo *)GetPHSquares(det,kFALSE);

    // Number of groups
    Int_t numberofgroup = 0;
    if(AliTRDgeometry::GetStack(det) == 2) numberofgroup = perChamber2;
    else numberofgroup = perChamber0;
    
    // PH
    for(Int_t nx = 0; nx < nxbins; nx++) {
      
      Double_t entries = 0.0;
      Double_t sumw2   = 0.0;
      Double_t sumw    = 0.0;

      // Sum the contributions of the different calibration group in the detector      
      for(Int_t k = 0; k < numberofgroup; k++){
	  
	Int_t binnumber = k*nxbins+nx;	  
	
	Int_t entriesv  = ((AliTRDEntriesInfo *)entriesph)->At(binnumber);
       	Float_t sumw2v  = ((AliTRDPhInfo *)squaresph)->AtS(binnumber)*entriesv;
	Float_t sumwv   = ((AliTRDPhInfo *)meanph)->At(binnumber)*entriesv;
	
	
	if(((entries+entriesv) > 65535) || ((entries+entriesv) <= 0)) continue;

	entries += entriesv;
	sumw2   += sumw2v;
	sumw    += sumwv;
      
      }

      if(entries > 0) {
	sumw2 = sumw2/((Float_t)entries);
	sumw  = sumw/((Float_t)entries);
      }
      
      fCalVector->FillVectorPH(det,0,nx,(Int_t)entries,(Float_t)sumw,(Float_t)sumw2);
    }
  }

  return fCalVector;
  
} 
//_____________________________________________________________________________________________________________________
AliTRDCalibraVector *AliTRDCalibraVector::AddStatsPerDetectorPRF()
{
  //
  // Create a AliTRDCalibraVector detector wise
  //

  AliTRDCalibraMode calibMode = AliTRDCalibraMode();
  calibMode.SetNz(2,GetNz(2));
  calibMode.SetNrphi(2,GetNrphi(2));
  if(((calibMode.GetNz(2) == 100) && (calibMode.GetNrphi(2) == 100)) || ((calibMode.GetNz(2) == 10) && (calibMode.GetNrphi(2) == 10))) return 0x0;

  // Check  
  Int_t nybins =  6*4*18*fDetCha0[2]+ 6*18*fDetCha2[2];
  Int_t nxbins = fNumberBinPRF;
  
  Int_t perChamber2 = 0;
  Int_t perChamber0 = 0;
  calibMode.ModePadCalibration(2,2);
  calibMode.ModePadFragmentation(0,2,0,2);
  calibMode.SetDetChamb2(2);
  perChamber2 = (Int_t) calibMode.GetDetChamb2(2);
  calibMode.ModePadCalibration(0,2);
  calibMode.ModePadFragmentation(0,0,0,2);
  calibMode.SetDetChamb0(2);
  perChamber0 = (Int_t) calibMode.GetDetChamb0(2);
  
  if(nybins != (6*18*perChamber2+6*4*18*perChamber0)) return 0x0;
    
  // Create calvector 
  if(!fCalVector) fCalVector = new AliTRDCalibraVector();
  else{ 
    fCalVector->~AliTRDCalibraVector();
    new(fCalVector) AliTRDCalibraVector();
  }
  fCalVector->SetNumberBinPRF(nxbins);
  fCalVector->SetDetCha0(2,1);
  fCalVector->SetDetCha2(2,1);
  fCalVector->SetNzNrphi(2,0,0);
  fCalVector->SetNbGroupPRF(fNbGroupPRF);

  
  for(Int_t det = 0; det < 540; det++){
    
    // Take
    AliTRDEntriesInfo *entriesprf = (AliTRDEntriesInfo *) GetPRFEntries(det,kFALSE);
    if(!entriesprf) continue;
    AliTRDPrfInfo     *meanprf    = (AliTRDPrfInfo *) GetPRFMean(det,kFALSE);
    AliTRDPrfInfo     *squaresprf = (AliTRDPrfInfo *) GetPRFSquares(det,kFALSE);

    // Number of groups
    Int_t numberofgroup = 0;
    if(AliTRDgeometry::GetStack(det) == 2) numberofgroup = perChamber2;
    else numberofgroup = perChamber0;
    
    for(Int_t nx = 0; nx < nxbins; nx++) {
      
      Double_t entries = 0.0;
      Double_t sumw2   = 0.0;
      Double_t sumw    = 0.0;
      
      // Sum the contributions of the different groups in the detector for one bin
      for(Int_t k = 0; k < numberofgroup; k++){
	  
	Int_t binnumber = k*nxbins+nx;	  

	Int_t entriesv  = ((AliTRDEntriesInfo *)entriesprf)->At(binnumber);
       	Float_t sumw2v  = ((AliTRDPrfInfo *)squaresprf)->At(binnumber)*entriesv;
	Float_t sumwv   = ((AliTRDPrfInfo *)meanprf)->At(binnumber)*entriesv;
	
	if(((entries+entriesv) > 65535) || ((entries+entriesv) <= 0)) continue;

	entries += entriesv;
	sumw2   += sumw2v;
	sumw    += sumwv;
      
      }

      if(entries > 0) {
	sumw2 = sumw2/((Float_t)entries);
	sumw  = sumw/((Float_t)entries);
      }
      
      fCalVector->FillVectorPRF(det,0,nx,(Int_t)entries,(Float_t)sumw,(Float_t)sumw2);
      
    }
  }

  return fCalVector;
}
//_______________________________________________________________________________
Bool_t AliTRDCalibraVector::FindTheMaxEntries(Int_t i, Int_t &detectormax, Int_t &groupmax)
{
  //
  // Find detectormax and groupmax with the biggest number of entries
  //

  // Coverity
  //Int_t numberofTB = 0;
  //if(i==0) numberofTB = (Int_t) GetNumberBinCharge();
  //if(i==1) numberofTB = GetTimeMax();
  //if(i==2) numberofTB = GetNumberBinPRF();
  if((i!=0) && (i!=1) && (i!=2)) AliInfo("Didn't understand i");


  // Init
  Double_t entries [540];
  for(Int_t idet = 0; idet < 540; idet++){
    entries[idet] = 0.0;
  }

  AliTRDEntriesInfo *entriesd = 0x0;
  // Take the number of entries per detector
  for(Int_t idet = 0; idet < 540; idet++){
 
    if(i==0) entriesd = (AliTRDEntriesInfo *) GetCHEntries(idet,kFALSE);
    if(i==1) entriesd = (AliTRDEntriesInfo *) GetPHEntries(idet,kFALSE);
    if(i==2) entriesd = (AliTRDEntriesInfo *) GetPRFEntries(idet,kFALSE);
    if(!entriesd) continue;

    entries[idet] = entriesd->GetSum();
    
  }

  // Search detector max
  Double_t max = -10;
  detectormax = -1;
  for(Int_t idet = 0; idet < 540; idet++){
    if(entries[idet] > max) {
      max = entries[idet];
      detectormax = idet;
    }
  }
  if((max == 0.0) || (detectormax <0.0)) return kFALSE;

  // Search group max
  if(i==0) entriesd = (AliTRDEntriesInfo *) GetCHEntries(detectormax,kFALSE);
  if(i==1) entriesd = (AliTRDEntriesInfo *) GetPHEntries(detectormax,kFALSE);
  if(i==2) entriesd = (AliTRDEntriesInfo *) GetPRFEntries(detectormax,kFALSE);  
  if(!entriesd) return kFALSE;
  // Number of groups
  Int_t numberofgroup = 0;
  if(AliTRDgeometry::GetStack(detectormax) == 2) numberofgroup = fDetCha2[i];
  else numberofgroup = fDetCha0[i];
  // Init
  Double_t nbgroup [2304];
  for(Int_t k = 0; k < 2304; k++){
    nbgroup[k] = 0.0;
  }
  Int_t nxbins = 0;
  if(i==0) nxbins = fNumberBinCharge;
  if(i==1) nxbins = fTimeMax;
  if(i==2) nxbins = fNumberBinPRF;
  // Compute the number of entries per group
  for(Int_t k = 0; k < numberofgroup; k++){
    for(Int_t nx = 0; nx < nxbins; nx++) {
      Int_t binnumber = k*nxbins+nx;	  
      nbgroup[k] += ((AliTRDEntriesInfo  *)entriesd)->At(binnumber);
    }
  }
  max = -10.0;
  groupmax = -1;
  for(Int_t io = 0; io < numberofgroup; io++){
    if(nbgroup[io] > max){
      max = nbgroup[io];
      groupmax = io;
    }
  }
  if((max == 0.0) || (groupmax < 0.0) || (groupmax >= numberofgroup)) return kFALSE;

  return kTRUE;

}
//_____________________________________________________________________________
TGraphErrors *AliTRDCalibraVector::ConvertVectorPHTGraphErrors(Int_t det, Int_t group , const Char_t * name)
{
  //
  // Convert the fVectorPHMean, fVectorPHSquares and fVectorPHEntries in TGraphErrors
  //

  // Take the info
  fPHEntries[det] = ((AliTRDEntriesInfo  *)GetPHEntries(det,kTRUE));
  fPHMean[det]    = ((AliTRDPhInfo *)GetPHMean(det,kTRUE));
  fPHSquares[det] = ((AliTRDPhInfo *)GetPHSquares(det,kTRUE));
  

  // Axis
  Float_t sf = 10.0;
  AliTRDCommonParam *parCom = AliTRDCommonParam::Instance();
  if (parCom) {
    sf = parCom->GetSamplingFrequency();
  }
  else {
    AliInfo("Could not get CommonParam, take the default 10MHz");
  }
  // Axis
  Double_t x[35];  // Xaxis
  Double_t y[35];  // Sum/entries
  Double_t ex[35]; // Nentries
  Double_t ey[35]; // Sum of square/nentries
  Double_t step = 0.0;
  Double_t min  = 0.0;
  if(sf > 0.0) step = 1.0 / sf;
  min  = 0.0;
  Int_t offset = group*fTimeMax;
  
  // Fill histo
  for (Int_t k = 0; k < fTimeMax; k++) {
    x[k]  = min + k*step;
    y[k]  = 0.0;
    ex[k] = 0.0;
    ey[k] = 0.0;   
    Int_t bin = offset+k;
    // Fill only if there is more than 0 something
    if (((AliTRDEntriesInfo *)fPHEntries[det])->At(bin) > 0) {
      ex[k] = ((AliTRDEntriesInfo *)fPHEntries[det])->At(bin);
      y[k]  = ((AliTRDPhInfo *)fPHMean[det])->At(bin);
      ey[k] = ((AliTRDPhInfo *)fPHSquares[det])->AtS(bin);
    }
  }

  // Define the TGraphErrors
  if(!fGraph) fGraph = new TGraphErrors(fTimeMax,&x[0],&y[0],&ex[0],&ey[0]);
  else{ 
    fGraph->~TGraphErrors();
    new(fGraph) TGraphErrors(fTimeMax,&x[0],&y[0],&ex[0],&ey[0]);
  } 
  fGraph->SetTitle(name); 

  return fGraph;

}
//_____________________________________________________________________________
TGraphErrors *AliTRDCalibraVector::ConvertVectorPRFTGraphErrors(Int_t det, Int_t group , const Char_t * name)
{
  //
  // Convert the fVectorPRFMean, fVectorPRFSquares and fVectorPRFEntries in TGraphErrors
  //

  // Take the info
  fPRFEntries[det] = ((AliTRDEntriesInfo *)GetPRFEntries(det,kTRUE));
  fPRFMean[det]    = ((AliTRDPrfInfo     *)GetPRFMean(det,kTRUE));
  fPRFSquares[det] = ((AliTRDPrfInfo     *)GetPRFSquares(det,kTRUE));
  

  // Axis
  Double_t x[200];  // Xaxis
  Double_t y[200];  // Sum/entries
  Double_t ex[200]; //Nentries
  Double_t ey[200]; // Sum of square/nentries
  Double_t step = 0.0;
  Double_t min  = 0.0;
  if(fNumberBinPRF) step = (2*TMath::Abs(fPRFRange)) / fNumberBinPRF;
  min  = -TMath::Abs(fPRFRange) + step / 2.0;
  Int_t offset = group*fNumberBinPRF;
  //printf("number of total: %d\n",fNumberBinPRF);
  // Fill histo
  for (Int_t k = 0; k < fNumberBinPRF; k++) {
    x[k]  = min + k*step;
    y[k]  = 0.0;
    ex[k] = 0.0;
    ey[k] = 0.0;
    Int_t bin = offset+k;
    // Fill only if there is more than 0 something
    if (((AliTRDEntriesInfo *)fPRFEntries[det])->At(bin) > 0) {
      ex[k] = ((AliTRDEntriesInfo *)fPRFEntries[det])->At(bin);
      y[k]  = ((AliTRDPrfInfo *)fPRFMean[det])->At(bin);
      ey[k] = ((AliTRDPrfInfo *)fPRFSquares[det])->At(bin);
    }
    //printf("Number of entries %f for %d\n",ex[k],k);
  }

  // Define the TGraphErrors
  if(!fGraph) fGraph = new TGraphErrors(fNumberBinPRF,&x[0],&y[0],&ex[0],&ey[0]);
  else{ 
    fGraph->~TGraphErrors();
    new(fGraph) TGraphErrors(fNumberBinPRF,&x[0],&y[0],&ex[0],&ey[0]);
  }
  fGraph->SetTitle(name); 

  return fGraph;



}
//_____________________________________________________________________________
TH1F *AliTRDCalibraVector::CorrectTheError(const TGraphErrors *hist, Int_t &nbEntries)
{
  //
  // In the case of the vectors method the trees contains TGraphErrors for PH and PRF
  // to be able to add them after
  // We convert it to a TH1F to be able to applied the same fit function method
  // After having called this function you can not add the statistics anymore
  //

  Int_t nbins       = hist->GetN();
  Double_t *x       = hist->GetX();
  Double_t *entries = hist->GetEX();
  Double_t *mean    = hist->GetY();
  Double_t *square  = hist->GetEY();
  nbEntries   = 0;

  if (nbins < 2) {
    return 0x0; 
  }

  Double_t step     = x[1] - x[0]; 
  Double_t minvalue = x[0] - step/2;
  Double_t maxvalue = x[(nbins-1)] + step/2;

  if(!fHisto) fHisto = new TH1F("projcorrecterror","",nbins,minvalue,maxvalue);
  else{ 
    fHisto->~TH1F();
    new(fHisto) TH1F("projcorrecterror","",nbins,minvalue,maxvalue);
  }

  for (Int_t k = 0; k < nbins; k++) {
    fHisto->SetBinContent(k+1,mean[k]);
    if (entries[k] > 0.0) {
      nbEntries += (Int_t) entries[k];
      Double_t d = TMath::Abs(square[k] - (mean[k]*mean[k]));
      fHisto->SetBinError(k+1,TMath::Sqrt(d/entries[k]));
    }
    else {
      fHisto->SetBinError(k+1,0.0);
    }
  }

  return fHisto;
 
}  
//_____________________________________________________________________________
TH1F *AliTRDCalibraVector::ConvertVectorCHHisto(Int_t det, Int_t group, const Char_t * name)
{
  //
  // Convert the fVectorCHEntries in TH1F
  //

  // Take the info
  fCHEntries[det] = ((AliTRDEntriesInfo *)GetCHEntries(det,kTRUE));
  
  // Init the stuff
  if(!fHisto) fHisto = new TH1F(name,name,fNumberBinCharge,0,300);
  else{ 
    fHisto->~TH1F();
    new(fHisto) TH1F(name,name,fNumberBinCharge,0,300);
  }
  fHisto->Sumw2();
  Int_t offset = group*fNumberBinCharge;
  // Fill histo
  for (Int_t k = 0; k < fNumberBinCharge; k++) {
    Int_t bin = offset+k;
    fHisto->SetBinContent(k+1,((AliTRDEntriesInfo *)fCHEntries[det])->At(bin));
    fHisto->SetBinError(k+1,TMath::Sqrt(TMath::Abs(((AliTRDEntriesInfo *)fCHEntries[det])->At(bin))));
  }
  
  return fHisto;

} 
//_____________________________________________________________________
TObject* AliTRDCalibraVector::GetPHEntries(Int_t det
                                              , Bool_t force) /*FOLD00*/
{
    //
    // return pointer to Carge ROC Calibration
    // if force is true create a new histogram if it doesn't exist allready
    //
    AliTRDEntriesInfo**arr = &fPHEntries[0];
    return (TObject *) GetEntriesPH(det, arr, force);
}
//_____________________________________________________________________
TObject* AliTRDCalibraVector::GetPRFEntries(Int_t det
                                               , Bool_t force) /*FOLD00*/
{
    //
    // return pointer to Carge ROC Calibration
    // if force is true create a new histogram if it doesn't exist allready
    //
    AliTRDEntriesInfo**arr = &fPRFEntries[0];
    return (TObject *) GetEntriesPRF(det, arr, force);
}
//_____________________________________________________________________
TObject* AliTRDCalibraVector::GetCHEntries(Int_t det
                                              , Bool_t force) /*FOLD00*/
{
    //
    // return pointer to Carge ROC Calibration
    // if force is true create a new histogram if it doesn't exist allready
    //
    AliTRDEntriesInfo**arr = &fCHEntries[0];
    return (TObject *) GetEntriesCH(det, arr, force);
}
//_____________________________________________________________________
TObject* AliTRDCalibraVector::GetPHMean(Int_t det
                                           , Bool_t force) /*FOLD00*/
{
    //
    // return pointer to ROC Calibration
    // if force is true create a new array
    //
    AliTRDPhInfo**arr = &fPHMean[0];
    return (TObject *) GetMeanSquaresPH(det, arr, force);
}
//_____________________________________________________________________
TObject* AliTRDCalibraVector::GetPHSquares(Int_t det
                                              , Bool_t force) /*FOLD00*/
{
    //
    // return pointer to ROC Calibration
    // if force is true create a new array
    //
    AliTRDPhInfo**arr = &fPHSquares[0];
    return (TObject *)  GetMeanSquaresPH(det, arr, force);
}
//_____________________________________________________________________
TObject* AliTRDCalibraVector::GetPRFMean(Int_t det
                                            , Bool_t force) /*FOLD00*/
{
    //
    // return pointer to ROC Calibration
    // if force is true create a new array
    //
    AliTRDPrfInfo**arr = &fPRFMean[0];
    return (TObject *) GetMeanSquaresPRF(det, arr, force);
}
//_____________________________________________________________________
TObject* AliTRDCalibraVector::GetPRFSquares(Int_t det
                                               , Bool_t force) /*FOLD00*/
{
    //
    // return pointer to ROC Calibration
    // if force is true create a new array
    //
    AliTRDPrfInfo**arr = &fPRFSquares[0];
    return (TObject *) GetMeanSquaresPRF(det, arr, force);
}
//_____________________________________________________________________
AliTRDEntriesInfo* AliTRDCalibraVector::GetEntriesCH(Int_t det
                                               , AliTRDEntriesInfo** arr
                                               , Bool_t force) /*FOLD00*/
{
    //
    // return pointer to UShort_t array Entries
    // if force is true create a new UShort_t array if it doesn't exist allready
    //
  if ( (!force) || (arr[det]))
    return (AliTRDEntriesInfo*)arr[det];

  // if we are forced and TArrayI doesn't yes exist create it
  Int_t ngroup = GetTotalNumberOfBinsInDetector(det,0,fNumberBinCharge); 
  // init
  arr[det] = new AliTRDEntriesInfo(ngroup);
  
  return (AliTRDEntriesInfo*)arr[det];

}
//_____________________________________________________________________
AliTRDEntriesInfo* AliTRDCalibraVector::GetEntriesPRF(Int_t det
                                               , AliTRDEntriesInfo** arr
                                               , Bool_t force) /*FOLD00*/
{
    //
    // return pointer to UShort_t array Entries
    // if force is true create a new UShort_t array if it doesn't exist allready
    //
  if ( (!force) || (arr[det]))
    return (AliTRDEntriesInfo*)arr[det];

  // if we are forced and TArrayI doesn't yes exist create it
  Int_t ngroup = GetTotalNumberOfBinsInDetector(det,2,fNumberBinPRF); 
  // init
  arr[det] = new AliTRDEntriesInfo(ngroup);
  
  return (AliTRDEntriesInfo*)arr[det];

}
//_____________________________________________________________________
AliTRDEntriesInfo *AliTRDCalibraVector::GetEntriesPH(Int_t det
                                              , AliTRDEntriesInfo ** arr
                                              , Bool_t force) /*FOLD00*/
{
    //
    // return pointer to UShort_t array Entries
    // if force is true create a new UShort_t array if it doesn't exist allready
    //
    if ( (!force) || (arr[det]))
	return (AliTRDEntriesInfo *)arr[det];

    // if we are forced and UShort_t doesn't yet exist create it
    Int_t ngroup = GetTotalNumberOfBinsInDetector(det,1,fTimeMax); 
    // init
    arr[det] = new AliTRDEntriesInfo(ngroup);
    
    return (AliTRDEntriesInfo*)arr[det];
   
}
//_____________________________________________________________________
AliTRDPhInfo* AliTRDCalibraVector::GetMeanSquaresPH(Int_t det
                                                  , AliTRDPhInfo** arr
                                                  , Bool_t force) /*FOLD00*/
{
    //
    // return pointer to Float_t array Mean or Squares
    // if force is true create a new Float_t array if it doesn't exist allready
    //
    if ( (!force) || (arr[det]))
	return (AliTRDPhInfo*)arr[det];

    // if we are forced and Float_t array doesn't yes exist create it
    Int_t ngroup = GetTotalNumberOfBinsInDetector(det,1,fTimeMax); 
    // init
    arr[det] = new AliTRDPhInfo(ngroup);
    
    return ((AliTRDPhInfo *)arr[det]);
}
//_____________________________________________________________________
AliTRDPrfInfo* AliTRDCalibraVector::GetMeanSquaresPRF(Int_t det
                                                   , AliTRDPrfInfo** arr
                                                   , Bool_t force) /*FOLD00*/
{
    //
    // return pointer to Float_t array Mean or Squares
    // if force is true create a new array if it doesn't exist allready
    //
  if ( (!force) || (arr[det]))
    return arr[det];
  
  // if we are forced and the array doesn't yet exist create it
  Int_t ngroup = GetTotalNumberOfBinsInDetector(det,2,fNumberBinPRF); 
  // init
  arr[det] = new AliTRDPrfInfo(ngroup);
  
  return (AliTRDPrfInfo*)arr[det];
   
}
//_____________________________________________________________________________
Int_t AliTRDCalibraVector::GetTotalNumberOfBinsInDetector(Int_t det, Int_t i, Int_t nbBin) const
{
  //
  // Get the total number of bins (Nb of bins*Nb of groups) in the detector det for the group i
  //

  Int_t ngroup = 0;
  Int_t stack  = AliTRDgeometry::GetStack(det);
  if(stack == 2) ngroup = fDetCha2[i]*nbBin;
  else ngroup = fDetCha0[i]*nbBin;

  return ngroup;
 
}
//____________________________________________________________________________
Int_t AliTRDCalibraVector::GetNz(Int_t i) const
{
  //
  // Get Nz the granularity in row
  //

  Int_t nz = 0;
  if(i==0) nz = (Int_t)(fModeCH>>4);
  if(i==1) nz = (Int_t)(fModePH>>4);
  if(i==2) nz = (Int_t)(fModePRF>>4);
  
  return nz;

}
//____________________________________________________________________________
Int_t AliTRDCalibraVector::GetNrphi(Int_t i) const
{
  //
  // Get Nrphi the granularity in col
  //

  Int_t nrphi = 0;
  if(i==0) nrphi = (Int_t)(fModeCH&15);
  if(i==1) nrphi = (Int_t)(fModePH&15);
  if(i==2) nrphi = (Int_t)(fModePRF&15);
  
  return nrphi;

}
//_________________________________________________________________________________
TString AliTRDCalibraVector::GetNamePH() const
{
  //
  // Get the name of PH to know the granularity
  //
  
  Int_t nz = GetNz(1);
  Int_t nrphi = GetNrphi(1);

  TString name("Nz");
  name += nz;
  name += "Nrphi";
  name += nrphi;
  
  return name;

}   
//_________________________________________________________________________________
TString AliTRDCalibraVector::GetNameCH() const
{
  //
  // Get the name of CH to know the granularity
  //
  
  Int_t nz = GetNz(0);
  Int_t nrphi = GetNrphi(0);

  TString name("Nz");
  name += nz;
  name += "Nrphi";
  name += nrphi;
  
  return name;

}   
//_________________________________________________________________________________
TString AliTRDCalibraVector::GetNamePRF() const
{
  //
  // Get the name of PRF to know the granularity
  //

  Int_t nz = GetNz(2);
  Int_t nrphi = GetNrphi(2);
  
  TString name("Nz");
  name += nz;
  name += "Nrphi";
  name += nrphi;
  name += "Ngp";
  name += fNbGroupPRF;
  
  return name;

}
//____________________________________________________________________________
void AliTRDCalibraVector::SetNzNrphi(Int_t i, Int_t nz, Int_t nrphi) 
{
  //
  // Set NzNrphi for the granularity
  //
  
  if(i==0) {
    fModeCH = nz;
    fModeCH = fModeCH << 4;
    fModeCH |= nrphi;
  }
  if(i==1) {
    fModePH = nz;
    fModePH = fModePH << 4;
    fModePH |= nrphi;
  }
  if(i==2) {
    fModePRF = nz;
    fModePRF = fModePRF << 4;
    fModePRF |= nrphi;
  }
  
}

 AliTRDCalibraVector.cxx:1
 AliTRDCalibraVector.cxx:2
 AliTRDCalibraVector.cxx:3
 AliTRDCalibraVector.cxx:4
 AliTRDCalibraVector.cxx:5
 AliTRDCalibraVector.cxx:6
 AliTRDCalibraVector.cxx:7
 AliTRDCalibraVector.cxx:8
 AliTRDCalibraVector.cxx:9
 AliTRDCalibraVector.cxx:10
 AliTRDCalibraVector.cxx:11
 AliTRDCalibraVector.cxx:12
 AliTRDCalibraVector.cxx:13
 AliTRDCalibraVector.cxx:14
 AliTRDCalibraVector.cxx:15
 AliTRDCalibraVector.cxx:16
 AliTRDCalibraVector.cxx:17
 AliTRDCalibraVector.cxx:18
 AliTRDCalibraVector.cxx:19
 AliTRDCalibraVector.cxx:20
 AliTRDCalibraVector.cxx:21
 AliTRDCalibraVector.cxx:22
 AliTRDCalibraVector.cxx:23
 AliTRDCalibraVector.cxx:24
 AliTRDCalibraVector.cxx:25
 AliTRDCalibraVector.cxx:26
 AliTRDCalibraVector.cxx:27
 AliTRDCalibraVector.cxx:28
 AliTRDCalibraVector.cxx:29
 AliTRDCalibraVector.cxx:30
 AliTRDCalibraVector.cxx:31
 AliTRDCalibraVector.cxx:32
 AliTRDCalibraVector.cxx:33
 AliTRDCalibraVector.cxx:34
 AliTRDCalibraVector.cxx:35
 AliTRDCalibraVector.cxx:36
 AliTRDCalibraVector.cxx:37
 AliTRDCalibraVector.cxx:38
 AliTRDCalibraVector.cxx:39
 AliTRDCalibraVector.cxx:40
 AliTRDCalibraVector.cxx:41
 AliTRDCalibraVector.cxx:42
 AliTRDCalibraVector.cxx:43
 AliTRDCalibraVector.cxx:44
 AliTRDCalibraVector.cxx:45
 AliTRDCalibraVector.cxx:46
 AliTRDCalibraVector.cxx:47
 AliTRDCalibraVector.cxx:48
 AliTRDCalibraVector.cxx:49
 AliTRDCalibraVector.cxx:50
 AliTRDCalibraVector.cxx:51
 AliTRDCalibraVector.cxx:52
 AliTRDCalibraVector.cxx:53
 AliTRDCalibraVector.cxx:54
 AliTRDCalibraVector.cxx:55
 AliTRDCalibraVector.cxx:56
 AliTRDCalibraVector.cxx:57
 AliTRDCalibraVector.cxx:58
 AliTRDCalibraVector.cxx:59
 AliTRDCalibraVector.cxx:60
 AliTRDCalibraVector.cxx:61
 AliTRDCalibraVector.cxx:62
 AliTRDCalibraVector.cxx:63
 AliTRDCalibraVector.cxx:64
 AliTRDCalibraVector.cxx:65
 AliTRDCalibraVector.cxx:66
 AliTRDCalibraVector.cxx:67
 AliTRDCalibraVector.cxx:68
 AliTRDCalibraVector.cxx:69
 AliTRDCalibraVector.cxx:70
 AliTRDCalibraVector.cxx:71
 AliTRDCalibraVector.cxx:72
 AliTRDCalibraVector.cxx:73
 AliTRDCalibraVector.cxx:74
 AliTRDCalibraVector.cxx:75
 AliTRDCalibraVector.cxx:76
 AliTRDCalibraVector.cxx:77
 AliTRDCalibraVector.cxx:78
 AliTRDCalibraVector.cxx:79
 AliTRDCalibraVector.cxx:80
 AliTRDCalibraVector.cxx:81
 AliTRDCalibraVector.cxx:82
 AliTRDCalibraVector.cxx:83
 AliTRDCalibraVector.cxx:84
 AliTRDCalibraVector.cxx:85
 AliTRDCalibraVector.cxx:86
 AliTRDCalibraVector.cxx:87
 AliTRDCalibraVector.cxx:88
 AliTRDCalibraVector.cxx:89
 AliTRDCalibraVector.cxx:90
 AliTRDCalibraVector.cxx:91
 AliTRDCalibraVector.cxx:92
 AliTRDCalibraVector.cxx:93
 AliTRDCalibraVector.cxx:94
 AliTRDCalibraVector.cxx:95
 AliTRDCalibraVector.cxx:96
 AliTRDCalibraVector.cxx:97
 AliTRDCalibraVector.cxx:98
 AliTRDCalibraVector.cxx:99
 AliTRDCalibraVector.cxx:100
 AliTRDCalibraVector.cxx:101
 AliTRDCalibraVector.cxx:102
 AliTRDCalibraVector.cxx:103
 AliTRDCalibraVector.cxx:104
 AliTRDCalibraVector.cxx:105
 AliTRDCalibraVector.cxx:106
 AliTRDCalibraVector.cxx:107
 AliTRDCalibraVector.cxx:108
 AliTRDCalibraVector.cxx:109
 AliTRDCalibraVector.cxx:110
 AliTRDCalibraVector.cxx:111
 AliTRDCalibraVector.cxx:112
 AliTRDCalibraVector.cxx:113
 AliTRDCalibraVector.cxx:114
 AliTRDCalibraVector.cxx:115
 AliTRDCalibraVector.cxx:116
 AliTRDCalibraVector.cxx:117
 AliTRDCalibraVector.cxx:118
 AliTRDCalibraVector.cxx:119
 AliTRDCalibraVector.cxx:120
 AliTRDCalibraVector.cxx:121
 AliTRDCalibraVector.cxx:122
 AliTRDCalibraVector.cxx:123
 AliTRDCalibraVector.cxx:124
 AliTRDCalibraVector.cxx:125
 AliTRDCalibraVector.cxx:126
 AliTRDCalibraVector.cxx:127
 AliTRDCalibraVector.cxx:128
 AliTRDCalibraVector.cxx:129
 AliTRDCalibraVector.cxx:130
 AliTRDCalibraVector.cxx:131
 AliTRDCalibraVector.cxx:132
 AliTRDCalibraVector.cxx:133
 AliTRDCalibraVector.cxx:134
 AliTRDCalibraVector.cxx:135
 AliTRDCalibraVector.cxx:136
 AliTRDCalibraVector.cxx:137
 AliTRDCalibraVector.cxx:138
 AliTRDCalibraVector.cxx:139
 AliTRDCalibraVector.cxx:140
 AliTRDCalibraVector.cxx:141
 AliTRDCalibraVector.cxx:142
 AliTRDCalibraVector.cxx:143
 AliTRDCalibraVector.cxx:144
 AliTRDCalibraVector.cxx:145
 AliTRDCalibraVector.cxx:146
 AliTRDCalibraVector.cxx:147
 AliTRDCalibraVector.cxx:148
 AliTRDCalibraVector.cxx:149
 AliTRDCalibraVector.cxx:150
 AliTRDCalibraVector.cxx:151
 AliTRDCalibraVector.cxx:152
 AliTRDCalibraVector.cxx:153
 AliTRDCalibraVector.cxx:154
 AliTRDCalibraVector.cxx:155
 AliTRDCalibraVector.cxx:156
 AliTRDCalibraVector.cxx:157
 AliTRDCalibraVector.cxx:158
 AliTRDCalibraVector.cxx:159
 AliTRDCalibraVector.cxx:160
 AliTRDCalibraVector.cxx:161
 AliTRDCalibraVector.cxx:162
 AliTRDCalibraVector.cxx:163
 AliTRDCalibraVector.cxx:164
 AliTRDCalibraVector.cxx:165
 AliTRDCalibraVector.cxx:166
 AliTRDCalibraVector.cxx:167
 AliTRDCalibraVector.cxx:168
 AliTRDCalibraVector.cxx:169
 AliTRDCalibraVector.cxx:170
 AliTRDCalibraVector.cxx:171
 AliTRDCalibraVector.cxx:172
 AliTRDCalibraVector.cxx:173
 AliTRDCalibraVector.cxx:174
 AliTRDCalibraVector.cxx:175
 AliTRDCalibraVector.cxx:176
 AliTRDCalibraVector.cxx:177
 AliTRDCalibraVector.cxx:178
 AliTRDCalibraVector.cxx:179
 AliTRDCalibraVector.cxx:180
 AliTRDCalibraVector.cxx:181
 AliTRDCalibraVector.cxx:182
 AliTRDCalibraVector.cxx:183
 AliTRDCalibraVector.cxx:184
 AliTRDCalibraVector.cxx:185
 AliTRDCalibraVector.cxx:186
 AliTRDCalibraVector.cxx:187
 AliTRDCalibraVector.cxx:188
 AliTRDCalibraVector.cxx:189
 AliTRDCalibraVector.cxx:190
 AliTRDCalibraVector.cxx:191
 AliTRDCalibraVector.cxx:192
 AliTRDCalibraVector.cxx:193
 AliTRDCalibraVector.cxx:194
 AliTRDCalibraVector.cxx:195
 AliTRDCalibraVector.cxx:196
 AliTRDCalibraVector.cxx:197
 AliTRDCalibraVector.cxx:198
 AliTRDCalibraVector.cxx:199
 AliTRDCalibraVector.cxx:200
 AliTRDCalibraVector.cxx:201
 AliTRDCalibraVector.cxx:202
 AliTRDCalibraVector.cxx:203
 AliTRDCalibraVector.cxx:204
 AliTRDCalibraVector.cxx:205
 AliTRDCalibraVector.cxx:206
 AliTRDCalibraVector.cxx:207
 AliTRDCalibraVector.cxx:208
 AliTRDCalibraVector.cxx:209
 AliTRDCalibraVector.cxx:210
 AliTRDCalibraVector.cxx:211
 AliTRDCalibraVector.cxx:212
 AliTRDCalibraVector.cxx:213
 AliTRDCalibraVector.cxx:214
 AliTRDCalibraVector.cxx:215
 AliTRDCalibraVector.cxx:216
 AliTRDCalibraVector.cxx:217
 AliTRDCalibraVector.cxx:218
 AliTRDCalibraVector.cxx:219
 AliTRDCalibraVector.cxx:220
 AliTRDCalibraVector.cxx:221
 AliTRDCalibraVector.cxx:222
 AliTRDCalibraVector.cxx:223
 AliTRDCalibraVector.cxx:224
 AliTRDCalibraVector.cxx:225
 AliTRDCalibraVector.cxx:226
 AliTRDCalibraVector.cxx:227
 AliTRDCalibraVector.cxx:228
 AliTRDCalibraVector.cxx:229
 AliTRDCalibraVector.cxx:230
 AliTRDCalibraVector.cxx:231
 AliTRDCalibraVector.cxx:232
 AliTRDCalibraVector.cxx:233
 AliTRDCalibraVector.cxx:234
 AliTRDCalibraVector.cxx:235
 AliTRDCalibraVector.cxx:236
 AliTRDCalibraVector.cxx:237
 AliTRDCalibraVector.cxx:238
 AliTRDCalibraVector.cxx:239
 AliTRDCalibraVector.cxx:240
 AliTRDCalibraVector.cxx:241
 AliTRDCalibraVector.cxx:242
 AliTRDCalibraVector.cxx:243
 AliTRDCalibraVector.cxx:244
 AliTRDCalibraVector.cxx:245
 AliTRDCalibraVector.cxx:246
 AliTRDCalibraVector.cxx:247
 AliTRDCalibraVector.cxx:248
 AliTRDCalibraVector.cxx:249
 AliTRDCalibraVector.cxx:250
 AliTRDCalibraVector.cxx:251
 AliTRDCalibraVector.cxx:252
 AliTRDCalibraVector.cxx:253
 AliTRDCalibraVector.cxx:254
 AliTRDCalibraVector.cxx:255
 AliTRDCalibraVector.cxx:256
 AliTRDCalibraVector.cxx:257
 AliTRDCalibraVector.cxx:258
 AliTRDCalibraVector.cxx:259
 AliTRDCalibraVector.cxx:260
 AliTRDCalibraVector.cxx:261
 AliTRDCalibraVector.cxx:262
 AliTRDCalibraVector.cxx:263
 AliTRDCalibraVector.cxx:264
 AliTRDCalibraVector.cxx:265
 AliTRDCalibraVector.cxx:266
 AliTRDCalibraVector.cxx:267
 AliTRDCalibraVector.cxx:268
 AliTRDCalibraVector.cxx:269
 AliTRDCalibraVector.cxx:270
 AliTRDCalibraVector.cxx:271
 AliTRDCalibraVector.cxx:272
 AliTRDCalibraVector.cxx:273
 AliTRDCalibraVector.cxx:274
 AliTRDCalibraVector.cxx:275
 AliTRDCalibraVector.cxx:276
 AliTRDCalibraVector.cxx:277
 AliTRDCalibraVector.cxx:278
 AliTRDCalibraVector.cxx:279
 AliTRDCalibraVector.cxx:280
 AliTRDCalibraVector.cxx:281
 AliTRDCalibraVector.cxx:282
 AliTRDCalibraVector.cxx:283
 AliTRDCalibraVector.cxx:284
 AliTRDCalibraVector.cxx:285
 AliTRDCalibraVector.cxx:286
 AliTRDCalibraVector.cxx:287
 AliTRDCalibraVector.cxx:288
 AliTRDCalibraVector.cxx:289
 AliTRDCalibraVector.cxx:290
 AliTRDCalibraVector.cxx:291
 AliTRDCalibraVector.cxx:292
 AliTRDCalibraVector.cxx:293
 AliTRDCalibraVector.cxx:294
 AliTRDCalibraVector.cxx:295
 AliTRDCalibraVector.cxx:296
 AliTRDCalibraVector.cxx:297
 AliTRDCalibraVector.cxx:298
 AliTRDCalibraVector.cxx:299
 AliTRDCalibraVector.cxx:300
 AliTRDCalibraVector.cxx:301
 AliTRDCalibraVector.cxx:302
 AliTRDCalibraVector.cxx:303
 AliTRDCalibraVector.cxx:304
 AliTRDCalibraVector.cxx:305
 AliTRDCalibraVector.cxx:306
 AliTRDCalibraVector.cxx:307
 AliTRDCalibraVector.cxx:308
 AliTRDCalibraVector.cxx:309
 AliTRDCalibraVector.cxx:310
 AliTRDCalibraVector.cxx:311
 AliTRDCalibraVector.cxx:312
 AliTRDCalibraVector.cxx:313
 AliTRDCalibraVector.cxx:314
 AliTRDCalibraVector.cxx:315
 AliTRDCalibraVector.cxx:316
 AliTRDCalibraVector.cxx:317
 AliTRDCalibraVector.cxx:318
 AliTRDCalibraVector.cxx:319
 AliTRDCalibraVector.cxx:320
 AliTRDCalibraVector.cxx:321
 AliTRDCalibraVector.cxx:322
 AliTRDCalibraVector.cxx:323
 AliTRDCalibraVector.cxx:324
 AliTRDCalibraVector.cxx:325
 AliTRDCalibraVector.cxx:326
 AliTRDCalibraVector.cxx:327
 AliTRDCalibraVector.cxx:328
 AliTRDCalibraVector.cxx:329
 AliTRDCalibraVector.cxx:330
 AliTRDCalibraVector.cxx:331
 AliTRDCalibraVector.cxx:332
 AliTRDCalibraVector.cxx:333
 AliTRDCalibraVector.cxx:334
 AliTRDCalibraVector.cxx:335
 AliTRDCalibraVector.cxx:336
 AliTRDCalibraVector.cxx:337
 AliTRDCalibraVector.cxx:338
 AliTRDCalibraVector.cxx:339
 AliTRDCalibraVector.cxx:340
 AliTRDCalibraVector.cxx:341
 AliTRDCalibraVector.cxx:342
 AliTRDCalibraVector.cxx:343
 AliTRDCalibraVector.cxx:344
 AliTRDCalibraVector.cxx:345
 AliTRDCalibraVector.cxx:346
 AliTRDCalibraVector.cxx:347
 AliTRDCalibraVector.cxx:348
 AliTRDCalibraVector.cxx:349
 AliTRDCalibraVector.cxx:350
 AliTRDCalibraVector.cxx:351
 AliTRDCalibraVector.cxx:352
 AliTRDCalibraVector.cxx:353
 AliTRDCalibraVector.cxx:354
 AliTRDCalibraVector.cxx:355
 AliTRDCalibraVector.cxx:356
 AliTRDCalibraVector.cxx:357
 AliTRDCalibraVector.cxx:358
 AliTRDCalibraVector.cxx:359
 AliTRDCalibraVector.cxx:360
 AliTRDCalibraVector.cxx:361
 AliTRDCalibraVector.cxx:362
 AliTRDCalibraVector.cxx:363
 AliTRDCalibraVector.cxx:364
 AliTRDCalibraVector.cxx:365
 AliTRDCalibraVector.cxx:366
 AliTRDCalibraVector.cxx:367
 AliTRDCalibraVector.cxx:368
 AliTRDCalibraVector.cxx:369
 AliTRDCalibraVector.cxx:370
 AliTRDCalibraVector.cxx:371
 AliTRDCalibraVector.cxx:372
 AliTRDCalibraVector.cxx:373
 AliTRDCalibraVector.cxx:374
 AliTRDCalibraVector.cxx:375
 AliTRDCalibraVector.cxx:376
 AliTRDCalibraVector.cxx:377
 AliTRDCalibraVector.cxx:378
 AliTRDCalibraVector.cxx:379
 AliTRDCalibraVector.cxx:380
 AliTRDCalibraVector.cxx:381
 AliTRDCalibraVector.cxx:382
 AliTRDCalibraVector.cxx:383
 AliTRDCalibraVector.cxx:384
 AliTRDCalibraVector.cxx:385
 AliTRDCalibraVector.cxx:386
 AliTRDCalibraVector.cxx:387
 AliTRDCalibraVector.cxx:388
 AliTRDCalibraVector.cxx:389
 AliTRDCalibraVector.cxx:390
 AliTRDCalibraVector.cxx:391
 AliTRDCalibraVector.cxx:392
 AliTRDCalibraVector.cxx:393
 AliTRDCalibraVector.cxx:394
 AliTRDCalibraVector.cxx:395
 AliTRDCalibraVector.cxx:396
 AliTRDCalibraVector.cxx:397
 AliTRDCalibraVector.cxx:398
 AliTRDCalibraVector.cxx:399
 AliTRDCalibraVector.cxx:400
 AliTRDCalibraVector.cxx:401
 AliTRDCalibraVector.cxx:402
 AliTRDCalibraVector.cxx:403
 AliTRDCalibraVector.cxx:404
 AliTRDCalibraVector.cxx:405
 AliTRDCalibraVector.cxx:406
 AliTRDCalibraVector.cxx:407
 AliTRDCalibraVector.cxx:408
 AliTRDCalibraVector.cxx:409
 AliTRDCalibraVector.cxx:410
 AliTRDCalibraVector.cxx:411
 AliTRDCalibraVector.cxx:412
 AliTRDCalibraVector.cxx:413
 AliTRDCalibraVector.cxx:414
 AliTRDCalibraVector.cxx:415
 AliTRDCalibraVector.cxx:416
 AliTRDCalibraVector.cxx:417
 AliTRDCalibraVector.cxx:418
 AliTRDCalibraVector.cxx:419
 AliTRDCalibraVector.cxx:420
 AliTRDCalibraVector.cxx:421
 AliTRDCalibraVector.cxx:422
 AliTRDCalibraVector.cxx:423
 AliTRDCalibraVector.cxx:424
 AliTRDCalibraVector.cxx:425
 AliTRDCalibraVector.cxx:426
 AliTRDCalibraVector.cxx:427
 AliTRDCalibraVector.cxx:428
 AliTRDCalibraVector.cxx:429
 AliTRDCalibraVector.cxx:430
 AliTRDCalibraVector.cxx:431
 AliTRDCalibraVector.cxx:432
 AliTRDCalibraVector.cxx:433
 AliTRDCalibraVector.cxx:434
 AliTRDCalibraVector.cxx:435
 AliTRDCalibraVector.cxx:436
 AliTRDCalibraVector.cxx:437
 AliTRDCalibraVector.cxx:438
 AliTRDCalibraVector.cxx:439
 AliTRDCalibraVector.cxx:440
 AliTRDCalibraVector.cxx:441
 AliTRDCalibraVector.cxx:442
 AliTRDCalibraVector.cxx:443
 AliTRDCalibraVector.cxx:444
 AliTRDCalibraVector.cxx:445
 AliTRDCalibraVector.cxx:446
 AliTRDCalibraVector.cxx:447
 AliTRDCalibraVector.cxx:448
 AliTRDCalibraVector.cxx:449
 AliTRDCalibraVector.cxx:450
 AliTRDCalibraVector.cxx:451
 AliTRDCalibraVector.cxx:452
 AliTRDCalibraVector.cxx:453
 AliTRDCalibraVector.cxx:454
 AliTRDCalibraVector.cxx:455
 AliTRDCalibraVector.cxx:456
 AliTRDCalibraVector.cxx:457
 AliTRDCalibraVector.cxx:458
 AliTRDCalibraVector.cxx:459
 AliTRDCalibraVector.cxx:460
 AliTRDCalibraVector.cxx:461
 AliTRDCalibraVector.cxx:462
 AliTRDCalibraVector.cxx:463
 AliTRDCalibraVector.cxx:464
 AliTRDCalibraVector.cxx:465
 AliTRDCalibraVector.cxx:466
 AliTRDCalibraVector.cxx:467
 AliTRDCalibraVector.cxx:468
 AliTRDCalibraVector.cxx:469
 AliTRDCalibraVector.cxx:470
 AliTRDCalibraVector.cxx:471
 AliTRDCalibraVector.cxx:472
 AliTRDCalibraVector.cxx:473
 AliTRDCalibraVector.cxx:474
 AliTRDCalibraVector.cxx:475
 AliTRDCalibraVector.cxx:476
 AliTRDCalibraVector.cxx:477
 AliTRDCalibraVector.cxx:478
 AliTRDCalibraVector.cxx:479
 AliTRDCalibraVector.cxx:480
 AliTRDCalibraVector.cxx:481
 AliTRDCalibraVector.cxx:482
 AliTRDCalibraVector.cxx:483
 AliTRDCalibraVector.cxx:484
 AliTRDCalibraVector.cxx:485
 AliTRDCalibraVector.cxx:486
 AliTRDCalibraVector.cxx:487
 AliTRDCalibraVector.cxx:488
 AliTRDCalibraVector.cxx:489
 AliTRDCalibraVector.cxx:490
 AliTRDCalibraVector.cxx:491
 AliTRDCalibraVector.cxx:492
 AliTRDCalibraVector.cxx:493
 AliTRDCalibraVector.cxx:494
 AliTRDCalibraVector.cxx:495
 AliTRDCalibraVector.cxx:496
 AliTRDCalibraVector.cxx:497
 AliTRDCalibraVector.cxx:498
 AliTRDCalibraVector.cxx:499
 AliTRDCalibraVector.cxx:500
 AliTRDCalibraVector.cxx:501
 AliTRDCalibraVector.cxx:502
 AliTRDCalibraVector.cxx:503
 AliTRDCalibraVector.cxx:504
 AliTRDCalibraVector.cxx:505
 AliTRDCalibraVector.cxx:506
 AliTRDCalibraVector.cxx:507
 AliTRDCalibraVector.cxx:508
 AliTRDCalibraVector.cxx:509
 AliTRDCalibraVector.cxx:510
 AliTRDCalibraVector.cxx:511
 AliTRDCalibraVector.cxx:512
 AliTRDCalibraVector.cxx:513
 AliTRDCalibraVector.cxx:514
 AliTRDCalibraVector.cxx:515
 AliTRDCalibraVector.cxx:516
 AliTRDCalibraVector.cxx:517
 AliTRDCalibraVector.cxx:518
 AliTRDCalibraVector.cxx:519
 AliTRDCalibraVector.cxx:520
 AliTRDCalibraVector.cxx:521
 AliTRDCalibraVector.cxx:522
 AliTRDCalibraVector.cxx:523
 AliTRDCalibraVector.cxx:524
 AliTRDCalibraVector.cxx:525
 AliTRDCalibraVector.cxx:526
 AliTRDCalibraVector.cxx:527
 AliTRDCalibraVector.cxx:528
 AliTRDCalibraVector.cxx:529
 AliTRDCalibraVector.cxx:530
 AliTRDCalibraVector.cxx:531
 AliTRDCalibraVector.cxx:532
 AliTRDCalibraVector.cxx:533
 AliTRDCalibraVector.cxx:534
 AliTRDCalibraVector.cxx:535
 AliTRDCalibraVector.cxx:536
 AliTRDCalibraVector.cxx:537
 AliTRDCalibraVector.cxx:538
 AliTRDCalibraVector.cxx:539
 AliTRDCalibraVector.cxx:540
 AliTRDCalibraVector.cxx:541
 AliTRDCalibraVector.cxx:542
 AliTRDCalibraVector.cxx:543
 AliTRDCalibraVector.cxx:544
 AliTRDCalibraVector.cxx:545
 AliTRDCalibraVector.cxx:546
 AliTRDCalibraVector.cxx:547
 AliTRDCalibraVector.cxx:548
 AliTRDCalibraVector.cxx:549
 AliTRDCalibraVector.cxx:550
 AliTRDCalibraVector.cxx:551
 AliTRDCalibraVector.cxx:552
 AliTRDCalibraVector.cxx:553
 AliTRDCalibraVector.cxx:554
 AliTRDCalibraVector.cxx:555
 AliTRDCalibraVector.cxx:556
 AliTRDCalibraVector.cxx:557
 AliTRDCalibraVector.cxx:558
 AliTRDCalibraVector.cxx:559
 AliTRDCalibraVector.cxx:560
 AliTRDCalibraVector.cxx:561
 AliTRDCalibraVector.cxx:562
 AliTRDCalibraVector.cxx:563
 AliTRDCalibraVector.cxx:564
 AliTRDCalibraVector.cxx:565
 AliTRDCalibraVector.cxx:566
 AliTRDCalibraVector.cxx:567
 AliTRDCalibraVector.cxx:568
 AliTRDCalibraVector.cxx:569
 AliTRDCalibraVector.cxx:570
 AliTRDCalibraVector.cxx:571
 AliTRDCalibraVector.cxx:572
 AliTRDCalibraVector.cxx:573
 AliTRDCalibraVector.cxx:574
 AliTRDCalibraVector.cxx:575
 AliTRDCalibraVector.cxx:576
 AliTRDCalibraVector.cxx:577
 AliTRDCalibraVector.cxx:578
 AliTRDCalibraVector.cxx:579
 AliTRDCalibraVector.cxx:580
 AliTRDCalibraVector.cxx:581
 AliTRDCalibraVector.cxx:582
 AliTRDCalibraVector.cxx:583
 AliTRDCalibraVector.cxx:584
 AliTRDCalibraVector.cxx:585
 AliTRDCalibraVector.cxx:586
 AliTRDCalibraVector.cxx:587
 AliTRDCalibraVector.cxx:588
 AliTRDCalibraVector.cxx:589
 AliTRDCalibraVector.cxx:590
 AliTRDCalibraVector.cxx:591
 AliTRDCalibraVector.cxx:592
 AliTRDCalibraVector.cxx:593
 AliTRDCalibraVector.cxx:594
 AliTRDCalibraVector.cxx:595
 AliTRDCalibraVector.cxx:596
 AliTRDCalibraVector.cxx:597
 AliTRDCalibraVector.cxx:598
 AliTRDCalibraVector.cxx:599
 AliTRDCalibraVector.cxx:600
 AliTRDCalibraVector.cxx:601
 AliTRDCalibraVector.cxx:602
 AliTRDCalibraVector.cxx:603
 AliTRDCalibraVector.cxx:604
 AliTRDCalibraVector.cxx:605
 AliTRDCalibraVector.cxx:606
 AliTRDCalibraVector.cxx:607
 AliTRDCalibraVector.cxx:608
 AliTRDCalibraVector.cxx:609
 AliTRDCalibraVector.cxx:610
 AliTRDCalibraVector.cxx:611
 AliTRDCalibraVector.cxx:612
 AliTRDCalibraVector.cxx:613
 AliTRDCalibraVector.cxx:614
 AliTRDCalibraVector.cxx:615
 AliTRDCalibraVector.cxx:616
 AliTRDCalibraVector.cxx:617
 AliTRDCalibraVector.cxx:618
 AliTRDCalibraVector.cxx:619
 AliTRDCalibraVector.cxx:620
 AliTRDCalibraVector.cxx:621
 AliTRDCalibraVector.cxx:622
 AliTRDCalibraVector.cxx:623
 AliTRDCalibraVector.cxx:624
 AliTRDCalibraVector.cxx:625
 AliTRDCalibraVector.cxx:626
 AliTRDCalibraVector.cxx:627
 AliTRDCalibraVector.cxx:628
 AliTRDCalibraVector.cxx:629
 AliTRDCalibraVector.cxx:630
 AliTRDCalibraVector.cxx:631
 AliTRDCalibraVector.cxx:632
 AliTRDCalibraVector.cxx:633
 AliTRDCalibraVector.cxx:634
 AliTRDCalibraVector.cxx:635
 AliTRDCalibraVector.cxx:636
 AliTRDCalibraVector.cxx:637
 AliTRDCalibraVector.cxx:638
 AliTRDCalibraVector.cxx:639
 AliTRDCalibraVector.cxx:640
 AliTRDCalibraVector.cxx:641
 AliTRDCalibraVector.cxx:642
 AliTRDCalibraVector.cxx:643
 AliTRDCalibraVector.cxx:644
 AliTRDCalibraVector.cxx:645
 AliTRDCalibraVector.cxx:646
 AliTRDCalibraVector.cxx:647
 AliTRDCalibraVector.cxx:648
 AliTRDCalibraVector.cxx:649
 AliTRDCalibraVector.cxx:650
 AliTRDCalibraVector.cxx:651
 AliTRDCalibraVector.cxx:652
 AliTRDCalibraVector.cxx:653
 AliTRDCalibraVector.cxx:654
 AliTRDCalibraVector.cxx:655
 AliTRDCalibraVector.cxx:656
 AliTRDCalibraVector.cxx:657
 AliTRDCalibraVector.cxx:658
 AliTRDCalibraVector.cxx:659
 AliTRDCalibraVector.cxx:660
 AliTRDCalibraVector.cxx:661
 AliTRDCalibraVector.cxx:662
 AliTRDCalibraVector.cxx:663
 AliTRDCalibraVector.cxx:664
 AliTRDCalibraVector.cxx:665
 AliTRDCalibraVector.cxx:666
 AliTRDCalibraVector.cxx:667
 AliTRDCalibraVector.cxx:668
 AliTRDCalibraVector.cxx:669
 AliTRDCalibraVector.cxx:670
 AliTRDCalibraVector.cxx:671
 AliTRDCalibraVector.cxx:672
 AliTRDCalibraVector.cxx:673
 AliTRDCalibraVector.cxx:674
 AliTRDCalibraVector.cxx:675
 AliTRDCalibraVector.cxx:676
 AliTRDCalibraVector.cxx:677
 AliTRDCalibraVector.cxx:678
 AliTRDCalibraVector.cxx:679
 AliTRDCalibraVector.cxx:680
 AliTRDCalibraVector.cxx:681
 AliTRDCalibraVector.cxx:682
 AliTRDCalibraVector.cxx:683
 AliTRDCalibraVector.cxx:684
 AliTRDCalibraVector.cxx:685
 AliTRDCalibraVector.cxx:686
 AliTRDCalibraVector.cxx:687
 AliTRDCalibraVector.cxx:688
 AliTRDCalibraVector.cxx:689
 AliTRDCalibraVector.cxx:690
 AliTRDCalibraVector.cxx:691
 AliTRDCalibraVector.cxx:692
 AliTRDCalibraVector.cxx:693
 AliTRDCalibraVector.cxx:694
 AliTRDCalibraVector.cxx:695
 AliTRDCalibraVector.cxx:696
 AliTRDCalibraVector.cxx:697
 AliTRDCalibraVector.cxx:698
 AliTRDCalibraVector.cxx:699
 AliTRDCalibraVector.cxx:700
 AliTRDCalibraVector.cxx:701
 AliTRDCalibraVector.cxx:702
 AliTRDCalibraVector.cxx:703
 AliTRDCalibraVector.cxx:704
 AliTRDCalibraVector.cxx:705
 AliTRDCalibraVector.cxx:706
 AliTRDCalibraVector.cxx:707
 AliTRDCalibraVector.cxx:708
 AliTRDCalibraVector.cxx:709
 AliTRDCalibraVector.cxx:710
 AliTRDCalibraVector.cxx:711
 AliTRDCalibraVector.cxx:712
 AliTRDCalibraVector.cxx:713
 AliTRDCalibraVector.cxx:714
 AliTRDCalibraVector.cxx:715
 AliTRDCalibraVector.cxx:716
 AliTRDCalibraVector.cxx:717
 AliTRDCalibraVector.cxx:718
 AliTRDCalibraVector.cxx:719
 AliTRDCalibraVector.cxx:720
 AliTRDCalibraVector.cxx:721
 AliTRDCalibraVector.cxx:722
 AliTRDCalibraVector.cxx:723
 AliTRDCalibraVector.cxx:724
 AliTRDCalibraVector.cxx:725
 AliTRDCalibraVector.cxx:726
 AliTRDCalibraVector.cxx:727
 AliTRDCalibraVector.cxx:728
 AliTRDCalibraVector.cxx:729
 AliTRDCalibraVector.cxx:730
 AliTRDCalibraVector.cxx:731
 AliTRDCalibraVector.cxx:732
 AliTRDCalibraVector.cxx:733
 AliTRDCalibraVector.cxx:734
 AliTRDCalibraVector.cxx:735
 AliTRDCalibraVector.cxx:736
 AliTRDCalibraVector.cxx:737
 AliTRDCalibraVector.cxx:738
 AliTRDCalibraVector.cxx:739
 AliTRDCalibraVector.cxx:740
 AliTRDCalibraVector.cxx:741
 AliTRDCalibraVector.cxx:742
 AliTRDCalibraVector.cxx:743
 AliTRDCalibraVector.cxx:744
 AliTRDCalibraVector.cxx:745
 AliTRDCalibraVector.cxx:746
 AliTRDCalibraVector.cxx:747
 AliTRDCalibraVector.cxx:748
 AliTRDCalibraVector.cxx:749
 AliTRDCalibraVector.cxx:750
 AliTRDCalibraVector.cxx:751
 AliTRDCalibraVector.cxx:752
 AliTRDCalibraVector.cxx:753
 AliTRDCalibraVector.cxx:754
 AliTRDCalibraVector.cxx:755
 AliTRDCalibraVector.cxx:756
 AliTRDCalibraVector.cxx:757
 AliTRDCalibraVector.cxx:758
 AliTRDCalibraVector.cxx:759
 AliTRDCalibraVector.cxx:760
 AliTRDCalibraVector.cxx:761
 AliTRDCalibraVector.cxx:762
 AliTRDCalibraVector.cxx:763
 AliTRDCalibraVector.cxx:764
 AliTRDCalibraVector.cxx:765
 AliTRDCalibraVector.cxx:766
 AliTRDCalibraVector.cxx:767
 AliTRDCalibraVector.cxx:768
 AliTRDCalibraVector.cxx:769
 AliTRDCalibraVector.cxx:770
 AliTRDCalibraVector.cxx:771
 AliTRDCalibraVector.cxx:772
 AliTRDCalibraVector.cxx:773
 AliTRDCalibraVector.cxx:774
 AliTRDCalibraVector.cxx:775
 AliTRDCalibraVector.cxx:776
 AliTRDCalibraVector.cxx:777
 AliTRDCalibraVector.cxx:778
 AliTRDCalibraVector.cxx:779
 AliTRDCalibraVector.cxx:780
 AliTRDCalibraVector.cxx:781
 AliTRDCalibraVector.cxx:782
 AliTRDCalibraVector.cxx:783
 AliTRDCalibraVector.cxx:784
 AliTRDCalibraVector.cxx:785
 AliTRDCalibraVector.cxx:786
 AliTRDCalibraVector.cxx:787
 AliTRDCalibraVector.cxx:788
 AliTRDCalibraVector.cxx:789
 AliTRDCalibraVector.cxx:790
 AliTRDCalibraVector.cxx:791
 AliTRDCalibraVector.cxx:792
 AliTRDCalibraVector.cxx:793
 AliTRDCalibraVector.cxx:794
 AliTRDCalibraVector.cxx:795
 AliTRDCalibraVector.cxx:796
 AliTRDCalibraVector.cxx:797
 AliTRDCalibraVector.cxx:798
 AliTRDCalibraVector.cxx:799
 AliTRDCalibraVector.cxx:800
 AliTRDCalibraVector.cxx:801
 AliTRDCalibraVector.cxx:802
 AliTRDCalibraVector.cxx:803
 AliTRDCalibraVector.cxx:804
 AliTRDCalibraVector.cxx:805
 AliTRDCalibraVector.cxx:806
 AliTRDCalibraVector.cxx:807
 AliTRDCalibraVector.cxx:808
 AliTRDCalibraVector.cxx:809
 AliTRDCalibraVector.cxx:810
 AliTRDCalibraVector.cxx:811
 AliTRDCalibraVector.cxx:812
 AliTRDCalibraVector.cxx:813
 AliTRDCalibraVector.cxx:814
 AliTRDCalibraVector.cxx:815
 AliTRDCalibraVector.cxx:816
 AliTRDCalibraVector.cxx:817
 AliTRDCalibraVector.cxx:818
 AliTRDCalibraVector.cxx:819
 AliTRDCalibraVector.cxx:820
 AliTRDCalibraVector.cxx:821
 AliTRDCalibraVector.cxx:822
 AliTRDCalibraVector.cxx:823
 AliTRDCalibraVector.cxx:824
 AliTRDCalibraVector.cxx:825
 AliTRDCalibraVector.cxx:826
 AliTRDCalibraVector.cxx:827
 AliTRDCalibraVector.cxx:828
 AliTRDCalibraVector.cxx:829
 AliTRDCalibraVector.cxx:830
 AliTRDCalibraVector.cxx:831
 AliTRDCalibraVector.cxx:832
 AliTRDCalibraVector.cxx:833
 AliTRDCalibraVector.cxx:834
 AliTRDCalibraVector.cxx:835
 AliTRDCalibraVector.cxx:836
 AliTRDCalibraVector.cxx:837
 AliTRDCalibraVector.cxx:838
 AliTRDCalibraVector.cxx:839
 AliTRDCalibraVector.cxx:840
 AliTRDCalibraVector.cxx:841
 AliTRDCalibraVector.cxx:842
 AliTRDCalibraVector.cxx:843
 AliTRDCalibraVector.cxx:844
 AliTRDCalibraVector.cxx:845
 AliTRDCalibraVector.cxx:846
 AliTRDCalibraVector.cxx:847
 AliTRDCalibraVector.cxx:848
 AliTRDCalibraVector.cxx:849
 AliTRDCalibraVector.cxx:850
 AliTRDCalibraVector.cxx:851
 AliTRDCalibraVector.cxx:852
 AliTRDCalibraVector.cxx:853
 AliTRDCalibraVector.cxx:854
 AliTRDCalibraVector.cxx:855
 AliTRDCalibraVector.cxx:856
 AliTRDCalibraVector.cxx:857
 AliTRDCalibraVector.cxx:858
 AliTRDCalibraVector.cxx:859
 AliTRDCalibraVector.cxx:860
 AliTRDCalibraVector.cxx:861
 AliTRDCalibraVector.cxx:862
 AliTRDCalibraVector.cxx:863
 AliTRDCalibraVector.cxx:864
 AliTRDCalibraVector.cxx:865
 AliTRDCalibraVector.cxx:866
 AliTRDCalibraVector.cxx:867
 AliTRDCalibraVector.cxx:868
 AliTRDCalibraVector.cxx:869
 AliTRDCalibraVector.cxx:870
 AliTRDCalibraVector.cxx:871
 AliTRDCalibraVector.cxx:872
 AliTRDCalibraVector.cxx:873
 AliTRDCalibraVector.cxx:874
 AliTRDCalibraVector.cxx:875
 AliTRDCalibraVector.cxx:876
 AliTRDCalibraVector.cxx:877
 AliTRDCalibraVector.cxx:878
 AliTRDCalibraVector.cxx:879
 AliTRDCalibraVector.cxx:880
 AliTRDCalibraVector.cxx:881
 AliTRDCalibraVector.cxx:882
 AliTRDCalibraVector.cxx:883
 AliTRDCalibraVector.cxx:884
 AliTRDCalibraVector.cxx:885
 AliTRDCalibraVector.cxx:886
 AliTRDCalibraVector.cxx:887
 AliTRDCalibraVector.cxx:888
 AliTRDCalibraVector.cxx:889
 AliTRDCalibraVector.cxx:890
 AliTRDCalibraVector.cxx:891
 AliTRDCalibraVector.cxx:892
 AliTRDCalibraVector.cxx:893
 AliTRDCalibraVector.cxx:894
 AliTRDCalibraVector.cxx:895
 AliTRDCalibraVector.cxx:896
 AliTRDCalibraVector.cxx:897
 AliTRDCalibraVector.cxx:898
 AliTRDCalibraVector.cxx:899
 AliTRDCalibraVector.cxx:900
 AliTRDCalibraVector.cxx:901
 AliTRDCalibraVector.cxx:902
 AliTRDCalibraVector.cxx:903
 AliTRDCalibraVector.cxx:904
 AliTRDCalibraVector.cxx:905
 AliTRDCalibraVector.cxx:906
 AliTRDCalibraVector.cxx:907
 AliTRDCalibraVector.cxx:908
 AliTRDCalibraVector.cxx:909
 AliTRDCalibraVector.cxx:910
 AliTRDCalibraVector.cxx:911
 AliTRDCalibraVector.cxx:912
 AliTRDCalibraVector.cxx:913
 AliTRDCalibraVector.cxx:914
 AliTRDCalibraVector.cxx:915
 AliTRDCalibraVector.cxx:916
 AliTRDCalibraVector.cxx:917
 AliTRDCalibraVector.cxx:918
 AliTRDCalibraVector.cxx:919
 AliTRDCalibraVector.cxx:920
 AliTRDCalibraVector.cxx:921
 AliTRDCalibraVector.cxx:922
 AliTRDCalibraVector.cxx:923
 AliTRDCalibraVector.cxx:924
 AliTRDCalibraVector.cxx:925
 AliTRDCalibraVector.cxx:926
 AliTRDCalibraVector.cxx:927
 AliTRDCalibraVector.cxx:928
 AliTRDCalibraVector.cxx:929
 AliTRDCalibraVector.cxx:930
 AliTRDCalibraVector.cxx:931
 AliTRDCalibraVector.cxx:932
 AliTRDCalibraVector.cxx:933
 AliTRDCalibraVector.cxx:934
 AliTRDCalibraVector.cxx:935
 AliTRDCalibraVector.cxx:936
 AliTRDCalibraVector.cxx:937
 AliTRDCalibraVector.cxx:938
 AliTRDCalibraVector.cxx:939
 AliTRDCalibraVector.cxx:940
 AliTRDCalibraVector.cxx:941
 AliTRDCalibraVector.cxx:942
 AliTRDCalibraVector.cxx:943
 AliTRDCalibraVector.cxx:944
 AliTRDCalibraVector.cxx:945
 AliTRDCalibraVector.cxx:946
 AliTRDCalibraVector.cxx:947
 AliTRDCalibraVector.cxx:948
 AliTRDCalibraVector.cxx:949
 AliTRDCalibraVector.cxx:950
 AliTRDCalibraVector.cxx:951
 AliTRDCalibraVector.cxx:952
 AliTRDCalibraVector.cxx:953
 AliTRDCalibraVector.cxx:954
 AliTRDCalibraVector.cxx:955
 AliTRDCalibraVector.cxx:956
 AliTRDCalibraVector.cxx:957
 AliTRDCalibraVector.cxx:958
 AliTRDCalibraVector.cxx:959
 AliTRDCalibraVector.cxx:960
 AliTRDCalibraVector.cxx:961
 AliTRDCalibraVector.cxx:962
 AliTRDCalibraVector.cxx:963
 AliTRDCalibraVector.cxx:964
 AliTRDCalibraVector.cxx:965
 AliTRDCalibraVector.cxx:966
 AliTRDCalibraVector.cxx:967
 AliTRDCalibraVector.cxx:968
 AliTRDCalibraVector.cxx:969
 AliTRDCalibraVector.cxx:970
 AliTRDCalibraVector.cxx:971
 AliTRDCalibraVector.cxx:972
 AliTRDCalibraVector.cxx:973
 AliTRDCalibraVector.cxx:974
 AliTRDCalibraVector.cxx:975
 AliTRDCalibraVector.cxx:976
 AliTRDCalibraVector.cxx:977
 AliTRDCalibraVector.cxx:978
 AliTRDCalibraVector.cxx:979
 AliTRDCalibraVector.cxx:980
 AliTRDCalibraVector.cxx:981
 AliTRDCalibraVector.cxx:982
 AliTRDCalibraVector.cxx:983
 AliTRDCalibraVector.cxx:984
 AliTRDCalibraVector.cxx:985
 AliTRDCalibraVector.cxx:986
 AliTRDCalibraVector.cxx:987
 AliTRDCalibraVector.cxx:988
 AliTRDCalibraVector.cxx:989
 AliTRDCalibraVector.cxx:990
 AliTRDCalibraVector.cxx:991
 AliTRDCalibraVector.cxx:992
 AliTRDCalibraVector.cxx:993
 AliTRDCalibraVector.cxx:994
 AliTRDCalibraVector.cxx:995
 AliTRDCalibraVector.cxx:996
 AliTRDCalibraVector.cxx:997
 AliTRDCalibraVector.cxx:998
 AliTRDCalibraVector.cxx:999
 AliTRDCalibraVector.cxx:1000
 AliTRDCalibraVector.cxx:1001
 AliTRDCalibraVector.cxx:1002
 AliTRDCalibraVector.cxx:1003
 AliTRDCalibraVector.cxx:1004
 AliTRDCalibraVector.cxx:1005
 AliTRDCalibraVector.cxx:1006
 AliTRDCalibraVector.cxx:1007
 AliTRDCalibraVector.cxx:1008
 AliTRDCalibraVector.cxx:1009
 AliTRDCalibraVector.cxx:1010
 AliTRDCalibraVector.cxx:1011
 AliTRDCalibraVector.cxx:1012
 AliTRDCalibraVector.cxx:1013
 AliTRDCalibraVector.cxx:1014
 AliTRDCalibraVector.cxx:1015
 AliTRDCalibraVector.cxx:1016
 AliTRDCalibraVector.cxx:1017
 AliTRDCalibraVector.cxx:1018
 AliTRDCalibraVector.cxx:1019
 AliTRDCalibraVector.cxx:1020
 AliTRDCalibraVector.cxx:1021
 AliTRDCalibraVector.cxx:1022
 AliTRDCalibraVector.cxx:1023
 AliTRDCalibraVector.cxx:1024
 AliTRDCalibraVector.cxx:1025
 AliTRDCalibraVector.cxx:1026
 AliTRDCalibraVector.cxx:1027
 AliTRDCalibraVector.cxx:1028
 AliTRDCalibraVector.cxx:1029
 AliTRDCalibraVector.cxx:1030
 AliTRDCalibraVector.cxx:1031
 AliTRDCalibraVector.cxx:1032
 AliTRDCalibraVector.cxx:1033
 AliTRDCalibraVector.cxx:1034
 AliTRDCalibraVector.cxx:1035
 AliTRDCalibraVector.cxx:1036
 AliTRDCalibraVector.cxx:1037
 AliTRDCalibraVector.cxx:1038
 AliTRDCalibraVector.cxx:1039
 AliTRDCalibraVector.cxx:1040
 AliTRDCalibraVector.cxx:1041
 AliTRDCalibraVector.cxx:1042
 AliTRDCalibraVector.cxx:1043
 AliTRDCalibraVector.cxx:1044
 AliTRDCalibraVector.cxx:1045
 AliTRDCalibraVector.cxx:1046
 AliTRDCalibraVector.cxx:1047
 AliTRDCalibraVector.cxx:1048
 AliTRDCalibraVector.cxx:1049
 AliTRDCalibraVector.cxx:1050
 AliTRDCalibraVector.cxx:1051
 AliTRDCalibraVector.cxx:1052
 AliTRDCalibraVector.cxx:1053
 AliTRDCalibraVector.cxx:1054
 AliTRDCalibraVector.cxx:1055
 AliTRDCalibraVector.cxx:1056
 AliTRDCalibraVector.cxx:1057
 AliTRDCalibraVector.cxx:1058
 AliTRDCalibraVector.cxx:1059
 AliTRDCalibraVector.cxx:1060
 AliTRDCalibraVector.cxx:1061
 AliTRDCalibraVector.cxx:1062
 AliTRDCalibraVector.cxx:1063
 AliTRDCalibraVector.cxx:1064
 AliTRDCalibraVector.cxx:1065
 AliTRDCalibraVector.cxx:1066
 AliTRDCalibraVector.cxx:1067
 AliTRDCalibraVector.cxx:1068
 AliTRDCalibraVector.cxx:1069
 AliTRDCalibraVector.cxx:1070
 AliTRDCalibraVector.cxx:1071
 AliTRDCalibraVector.cxx:1072
 AliTRDCalibraVector.cxx:1073
 AliTRDCalibraVector.cxx:1074
 AliTRDCalibraVector.cxx:1075
 AliTRDCalibraVector.cxx:1076
 AliTRDCalibraVector.cxx:1077
 AliTRDCalibraVector.cxx:1078
 AliTRDCalibraVector.cxx:1079
 AliTRDCalibraVector.cxx:1080
 AliTRDCalibraVector.cxx:1081
 AliTRDCalibraVector.cxx:1082
 AliTRDCalibraVector.cxx:1083
 AliTRDCalibraVector.cxx:1084
 AliTRDCalibraVector.cxx:1085
 AliTRDCalibraVector.cxx:1086
 AliTRDCalibraVector.cxx:1087
 AliTRDCalibraVector.cxx:1088
 AliTRDCalibraVector.cxx:1089
 AliTRDCalibraVector.cxx:1090
 AliTRDCalibraVector.cxx:1091
 AliTRDCalibraVector.cxx:1092
 AliTRDCalibraVector.cxx:1093
 AliTRDCalibraVector.cxx:1094
 AliTRDCalibraVector.cxx:1095
 AliTRDCalibraVector.cxx:1096
 AliTRDCalibraVector.cxx:1097
 AliTRDCalibraVector.cxx:1098
 AliTRDCalibraVector.cxx:1099
 AliTRDCalibraVector.cxx:1100
 AliTRDCalibraVector.cxx:1101
 AliTRDCalibraVector.cxx:1102
 AliTRDCalibraVector.cxx:1103
 AliTRDCalibraVector.cxx:1104
 AliTRDCalibraVector.cxx:1105
 AliTRDCalibraVector.cxx:1106
 AliTRDCalibraVector.cxx:1107
 AliTRDCalibraVector.cxx:1108
 AliTRDCalibraVector.cxx:1109
 AliTRDCalibraVector.cxx:1110
 AliTRDCalibraVector.cxx:1111
 AliTRDCalibraVector.cxx:1112
 AliTRDCalibraVector.cxx:1113
 AliTRDCalibraVector.cxx:1114
 AliTRDCalibraVector.cxx:1115
 AliTRDCalibraVector.cxx:1116
 AliTRDCalibraVector.cxx:1117
 AliTRDCalibraVector.cxx:1118
 AliTRDCalibraVector.cxx:1119
 AliTRDCalibraVector.cxx:1120
 AliTRDCalibraVector.cxx:1121
 AliTRDCalibraVector.cxx:1122
 AliTRDCalibraVector.cxx:1123
 AliTRDCalibraVector.cxx:1124
 AliTRDCalibraVector.cxx:1125
 AliTRDCalibraVector.cxx:1126
 AliTRDCalibraVector.cxx:1127
 AliTRDCalibraVector.cxx:1128
 AliTRDCalibraVector.cxx:1129
 AliTRDCalibraVector.cxx:1130
 AliTRDCalibraVector.cxx:1131
 AliTRDCalibraVector.cxx:1132
 AliTRDCalibraVector.cxx:1133
 AliTRDCalibraVector.cxx:1134
 AliTRDCalibraVector.cxx:1135
 AliTRDCalibraVector.cxx:1136
 AliTRDCalibraVector.cxx:1137
 AliTRDCalibraVector.cxx:1138
 AliTRDCalibraVector.cxx:1139
 AliTRDCalibraVector.cxx:1140
 AliTRDCalibraVector.cxx:1141
 AliTRDCalibraVector.cxx:1142
 AliTRDCalibraVector.cxx:1143
 AliTRDCalibraVector.cxx:1144
 AliTRDCalibraVector.cxx:1145
 AliTRDCalibraVector.cxx:1146
 AliTRDCalibraVector.cxx:1147
 AliTRDCalibraVector.cxx:1148
 AliTRDCalibraVector.cxx:1149
 AliTRDCalibraVector.cxx:1150
 AliTRDCalibraVector.cxx:1151
 AliTRDCalibraVector.cxx:1152
 AliTRDCalibraVector.cxx:1153
 AliTRDCalibraVector.cxx:1154
 AliTRDCalibraVector.cxx:1155
 AliTRDCalibraVector.cxx:1156
 AliTRDCalibraVector.cxx:1157
 AliTRDCalibraVector.cxx:1158
 AliTRDCalibraVector.cxx:1159
 AliTRDCalibraVector.cxx:1160
 AliTRDCalibraVector.cxx:1161
 AliTRDCalibraVector.cxx:1162
 AliTRDCalibraVector.cxx:1163
 AliTRDCalibraVector.cxx:1164
 AliTRDCalibraVector.cxx:1165
 AliTRDCalibraVector.cxx:1166
 AliTRDCalibraVector.cxx:1167
 AliTRDCalibraVector.cxx:1168
 AliTRDCalibraVector.cxx:1169
 AliTRDCalibraVector.cxx:1170
 AliTRDCalibraVector.cxx:1171
 AliTRDCalibraVector.cxx:1172
 AliTRDCalibraVector.cxx:1173
 AliTRDCalibraVector.cxx:1174
 AliTRDCalibraVector.cxx:1175
 AliTRDCalibraVector.cxx:1176
 AliTRDCalibraVector.cxx:1177
 AliTRDCalibraVector.cxx:1178
 AliTRDCalibraVector.cxx:1179
 AliTRDCalibraVector.cxx:1180
 AliTRDCalibraVector.cxx:1181
 AliTRDCalibraVector.cxx:1182
 AliTRDCalibraVector.cxx:1183
 AliTRDCalibraVector.cxx:1184
 AliTRDCalibraVector.cxx:1185
 AliTRDCalibraVector.cxx:1186
 AliTRDCalibraVector.cxx:1187
 AliTRDCalibraVector.cxx:1188
 AliTRDCalibraVector.cxx:1189
 AliTRDCalibraVector.cxx:1190
 AliTRDCalibraVector.cxx:1191
 AliTRDCalibraVector.cxx:1192
 AliTRDCalibraVector.cxx:1193
 AliTRDCalibraVector.cxx:1194
 AliTRDCalibraVector.cxx:1195
 AliTRDCalibraVector.cxx:1196
 AliTRDCalibraVector.cxx:1197
 AliTRDCalibraVector.cxx:1198
 AliTRDCalibraVector.cxx:1199
 AliTRDCalibraVector.cxx:1200
 AliTRDCalibraVector.cxx:1201
 AliTRDCalibraVector.cxx:1202
 AliTRDCalibraVector.cxx:1203
 AliTRDCalibraVector.cxx:1204
 AliTRDCalibraVector.cxx:1205
 AliTRDCalibraVector.cxx:1206
 AliTRDCalibraVector.cxx:1207
 AliTRDCalibraVector.cxx:1208
 AliTRDCalibraVector.cxx:1209
 AliTRDCalibraVector.cxx:1210
 AliTRDCalibraVector.cxx:1211
 AliTRDCalibraVector.cxx:1212
 AliTRDCalibraVector.cxx:1213
 AliTRDCalibraVector.cxx:1214
 AliTRDCalibraVector.cxx:1215
 AliTRDCalibraVector.cxx:1216
 AliTRDCalibraVector.cxx:1217
 AliTRDCalibraVector.cxx:1218
 AliTRDCalibraVector.cxx:1219
 AliTRDCalibraVector.cxx:1220
 AliTRDCalibraVector.cxx:1221
 AliTRDCalibraVector.cxx:1222
 AliTRDCalibraVector.cxx:1223
 AliTRDCalibraVector.cxx:1224
 AliTRDCalibraVector.cxx:1225
 AliTRDCalibraVector.cxx:1226
 AliTRDCalibraVector.cxx:1227
 AliTRDCalibraVector.cxx:1228
 AliTRDCalibraVector.cxx:1229
 AliTRDCalibraVector.cxx:1230
 AliTRDCalibraVector.cxx:1231
 AliTRDCalibraVector.cxx:1232
 AliTRDCalibraVector.cxx:1233
 AliTRDCalibraVector.cxx:1234
 AliTRDCalibraVector.cxx:1235
 AliTRDCalibraVector.cxx:1236
 AliTRDCalibraVector.cxx:1237
 AliTRDCalibraVector.cxx:1238
 AliTRDCalibraVector.cxx:1239
 AliTRDCalibraVector.cxx:1240
 AliTRDCalibraVector.cxx:1241
 AliTRDCalibraVector.cxx:1242
 AliTRDCalibraVector.cxx:1243
 AliTRDCalibraVector.cxx:1244
 AliTRDCalibraVector.cxx:1245
 AliTRDCalibraVector.cxx:1246
 AliTRDCalibraVector.cxx:1247
 AliTRDCalibraVector.cxx:1248
 AliTRDCalibraVector.cxx:1249
 AliTRDCalibraVector.cxx:1250
 AliTRDCalibraVector.cxx:1251
 AliTRDCalibraVector.cxx:1252
 AliTRDCalibraVector.cxx:1253
 AliTRDCalibraVector.cxx:1254
 AliTRDCalibraVector.cxx:1255
 AliTRDCalibraVector.cxx:1256
 AliTRDCalibraVector.cxx:1257
 AliTRDCalibraVector.cxx:1258
 AliTRDCalibraVector.cxx:1259
 AliTRDCalibraVector.cxx:1260
 AliTRDCalibraVector.cxx:1261
 AliTRDCalibraVector.cxx:1262
 AliTRDCalibraVector.cxx:1263
 AliTRDCalibraVector.cxx:1264
 AliTRDCalibraVector.cxx:1265
 AliTRDCalibraVector.cxx:1266
 AliTRDCalibraVector.cxx:1267
 AliTRDCalibraVector.cxx:1268
 AliTRDCalibraVector.cxx:1269
 AliTRDCalibraVector.cxx:1270
 AliTRDCalibraVector.cxx:1271
 AliTRDCalibraVector.cxx:1272
 AliTRDCalibraVector.cxx:1273
 AliTRDCalibraVector.cxx:1274
 AliTRDCalibraVector.cxx:1275
 AliTRDCalibraVector.cxx:1276
 AliTRDCalibraVector.cxx:1277
 AliTRDCalibraVector.cxx:1278
 AliTRDCalibraVector.cxx:1279
 AliTRDCalibraVector.cxx:1280
 AliTRDCalibraVector.cxx:1281
 AliTRDCalibraVector.cxx:1282
 AliTRDCalibraVector.cxx:1283
 AliTRDCalibraVector.cxx:1284
 AliTRDCalibraVector.cxx:1285
 AliTRDCalibraVector.cxx:1286
 AliTRDCalibraVector.cxx:1287
 AliTRDCalibraVector.cxx:1288
 AliTRDCalibraVector.cxx:1289
 AliTRDCalibraVector.cxx:1290
 AliTRDCalibraVector.cxx:1291
 AliTRDCalibraVector.cxx:1292
 AliTRDCalibraVector.cxx:1293
 AliTRDCalibraVector.cxx:1294
 AliTRDCalibraVector.cxx:1295
 AliTRDCalibraVector.cxx:1296
 AliTRDCalibraVector.cxx:1297
 AliTRDCalibraVector.cxx:1298
 AliTRDCalibraVector.cxx:1299
 AliTRDCalibraVector.cxx:1300
 AliTRDCalibraVector.cxx:1301
 AliTRDCalibraVector.cxx:1302
 AliTRDCalibraVector.cxx:1303
 AliTRDCalibraVector.cxx:1304
 AliTRDCalibraVector.cxx:1305
 AliTRDCalibraVector.cxx:1306
 AliTRDCalibraVector.cxx:1307
 AliTRDCalibraVector.cxx:1308
 AliTRDCalibraVector.cxx:1309
 AliTRDCalibraVector.cxx:1310
 AliTRDCalibraVector.cxx:1311
 AliTRDCalibraVector.cxx:1312
 AliTRDCalibraVector.cxx:1313
 AliTRDCalibraVector.cxx:1314
 AliTRDCalibraVector.cxx:1315
 AliTRDCalibraVector.cxx:1316
 AliTRDCalibraVector.cxx:1317
 AliTRDCalibraVector.cxx:1318
 AliTRDCalibraVector.cxx:1319
 AliTRDCalibraVector.cxx:1320
 AliTRDCalibraVector.cxx:1321
 AliTRDCalibraVector.cxx:1322
 AliTRDCalibraVector.cxx:1323
 AliTRDCalibraVector.cxx:1324
 AliTRDCalibraVector.cxx:1325
 AliTRDCalibraVector.cxx:1326
 AliTRDCalibraVector.cxx:1327
 AliTRDCalibraVector.cxx:1328
 AliTRDCalibraVector.cxx:1329
 AliTRDCalibraVector.cxx:1330
 AliTRDCalibraVector.cxx:1331
 AliTRDCalibraVector.cxx:1332
 AliTRDCalibraVector.cxx:1333
 AliTRDCalibraVector.cxx:1334
 AliTRDCalibraVector.cxx:1335
 AliTRDCalibraVector.cxx:1336
 AliTRDCalibraVector.cxx:1337
 AliTRDCalibraVector.cxx:1338
 AliTRDCalibraVector.cxx:1339
 AliTRDCalibraVector.cxx:1340
 AliTRDCalibraVector.cxx:1341
 AliTRDCalibraVector.cxx:1342
 AliTRDCalibraVector.cxx:1343
 AliTRDCalibraVector.cxx:1344
 AliTRDCalibraVector.cxx:1345
 AliTRDCalibraVector.cxx:1346
 AliTRDCalibraVector.cxx:1347
 AliTRDCalibraVector.cxx:1348
 AliTRDCalibraVector.cxx:1349
 AliTRDCalibraVector.cxx:1350
 AliTRDCalibraVector.cxx:1351
 AliTRDCalibraVector.cxx:1352
 AliTRDCalibraVector.cxx:1353
 AliTRDCalibraVector.cxx:1354
 AliTRDCalibraVector.cxx:1355
 AliTRDCalibraVector.cxx:1356
 AliTRDCalibraVector.cxx:1357
 AliTRDCalibraVector.cxx:1358
 AliTRDCalibraVector.cxx:1359
 AliTRDCalibraVector.cxx:1360
 AliTRDCalibraVector.cxx:1361
 AliTRDCalibraVector.cxx:1362
 AliTRDCalibraVector.cxx:1363
 AliTRDCalibraVector.cxx:1364
 AliTRDCalibraVector.cxx:1365
 AliTRDCalibraVector.cxx:1366
 AliTRDCalibraVector.cxx:1367
 AliTRDCalibraVector.cxx:1368
 AliTRDCalibraVector.cxx:1369
 AliTRDCalibraVector.cxx:1370
 AliTRDCalibraVector.cxx:1371
 AliTRDCalibraVector.cxx:1372
 AliTRDCalibraVector.cxx:1373
 AliTRDCalibraVector.cxx:1374
 AliTRDCalibraVector.cxx:1375
 AliTRDCalibraVector.cxx:1376
 AliTRDCalibraVector.cxx:1377
 AliTRDCalibraVector.cxx:1378
 AliTRDCalibraVector.cxx:1379
 AliTRDCalibraVector.cxx:1380
 AliTRDCalibraVector.cxx:1381
 AliTRDCalibraVector.cxx:1382
 AliTRDCalibraVector.cxx:1383
 AliTRDCalibraVector.cxx:1384
 AliTRDCalibraVector.cxx:1385
 AliTRDCalibraVector.cxx:1386
 AliTRDCalibraVector.cxx:1387
 AliTRDCalibraVector.cxx:1388
 AliTRDCalibraVector.cxx:1389
 AliTRDCalibraVector.cxx:1390
 AliTRDCalibraVector.cxx:1391
 AliTRDCalibraVector.cxx:1392
 AliTRDCalibraVector.cxx:1393
 AliTRDCalibraVector.cxx:1394
 AliTRDCalibraVector.cxx:1395
 AliTRDCalibraVector.cxx:1396
 AliTRDCalibraVector.cxx:1397
 AliTRDCalibraVector.cxx:1398
 AliTRDCalibraVector.cxx:1399
 AliTRDCalibraVector.cxx:1400
 AliTRDCalibraVector.cxx:1401
 AliTRDCalibraVector.cxx:1402
 AliTRDCalibraVector.cxx:1403
 AliTRDCalibraVector.cxx:1404
 AliTRDCalibraVector.cxx:1405
 AliTRDCalibraVector.cxx:1406
 AliTRDCalibraVector.cxx:1407
 AliTRDCalibraVector.cxx:1408
 AliTRDCalibraVector.cxx:1409
 AliTRDCalibraVector.cxx:1410
 AliTRDCalibraVector.cxx:1411
 AliTRDCalibraVector.cxx:1412
 AliTRDCalibraVector.cxx:1413
 AliTRDCalibraVector.cxx:1414
 AliTRDCalibraVector.cxx:1415
 AliTRDCalibraVector.cxx:1416
 AliTRDCalibraVector.cxx:1417
 AliTRDCalibraVector.cxx:1418
 AliTRDCalibraVector.cxx:1419
 AliTRDCalibraVector.cxx:1420
 AliTRDCalibraVector.cxx:1421
 AliTRDCalibraVector.cxx:1422
 AliTRDCalibraVector.cxx:1423
 AliTRDCalibraVector.cxx:1424
 AliTRDCalibraVector.cxx:1425
 AliTRDCalibraVector.cxx:1426
 AliTRDCalibraVector.cxx:1427
 AliTRDCalibraVector.cxx:1428
 AliTRDCalibraVector.cxx:1429
 AliTRDCalibraVector.cxx:1430
 AliTRDCalibraVector.cxx:1431
 AliTRDCalibraVector.cxx:1432
 AliTRDCalibraVector.cxx:1433
 AliTRDCalibraVector.cxx:1434
 AliTRDCalibraVector.cxx:1435
 AliTRDCalibraVector.cxx:1436
 AliTRDCalibraVector.cxx:1437
 AliTRDCalibraVector.cxx:1438
 AliTRDCalibraVector.cxx:1439
 AliTRDCalibraVector.cxx:1440
 AliTRDCalibraVector.cxx:1441
 AliTRDCalibraVector.cxx:1442
 AliTRDCalibraVector.cxx:1443
 AliTRDCalibraVector.cxx:1444
 AliTRDCalibraVector.cxx:1445
 AliTRDCalibraVector.cxx:1446
 AliTRDCalibraVector.cxx:1447
 AliTRDCalibraVector.cxx:1448
 AliTRDCalibraVector.cxx:1449
 AliTRDCalibraVector.cxx:1450
 AliTRDCalibraVector.cxx:1451
 AliTRDCalibraVector.cxx:1452
 AliTRDCalibraVector.cxx:1453
 AliTRDCalibraVector.cxx:1454
 AliTRDCalibraVector.cxx:1455
 AliTRDCalibraVector.cxx:1456
 AliTRDCalibraVector.cxx:1457
 AliTRDCalibraVector.cxx:1458