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: AliTriggerRunScalers.cxx 22322 2007-11-22 11:43:14Z cvetan $ */

///////////////////////////////////////////////////////////////////////////////
//
// Class to define the ALICE Trigger Scalers Record per Run.
// ReadScalers(): read the txt file (rXXXX.cnt) provided by CTP and creates array of records.
// ConsistencyCheck(): checks if the last scaler added to record is consistent, 
// i.e. time increases and valued are not decreasing in time.
// 
//
//////////////////////////////////////////////////////////////////////////////

#include <stdlib.h>
#include <Riostream.h>

//#include <TObject.h>
//#include <TArrayC.h>
//#include <TFile.h>
#include <TClass.h>
#include <TSystem.h>
#include <TObjString.h>
#include <TObjArray.h>
#include <TMath.h>
#include <TGraphErrors.h>

#include "AliLog.h"  
#include "AliTimeStamp.h"
#include "AliTriggerScalers.h"
#include "AliTriggerScalersESD.h"
#include "AliTriggerScalersRecord.h"
#include "AliTriggerScalersRecordESD.h"
#include "AliTriggerRunScalers.h"
#include "AliTriggerConfiguration.h"
#include "AliTriggerClass.h"
#include "AliTriggerBCMask.h"

using std::endl;
using std::cout;
using std::ifstream;
ClassImp( AliTriggerRunScalers )

//_____________________________________________________________________________
AliTriggerRunScalers::AliTriggerRunScalers():
  TObject(),
  fVersion(0),
  fRunNumber(0),
  fnClasses(0),
  fClassIndex(),                    
  fScalersRecord(),
  fScalersRecordESD()
{
  // Default constructor
}
//______________________________________________________________________________
AliTriggerRunScalers::~AliTriggerRunScalers() 
{
 // Destructor
 fScalersRecord.SetOwner(); 
 fScalersRecord.Delete(); 
 fScalersRecordESD.SetOwner(); 
 fScalersRecordESD.Delete(); 
}
//_____________________________________________________________________________
void AliTriggerRunScalers::AddTriggerScalers( AliTriggerScalersRecord* scaler ) 
{ 
  // Add scaler and check consistency
  fScalersRecord.AddLast( scaler );
  UInt_t* overflow[6];
  for(Int_t i=0; i<6; i++) {
     overflow[i] = new UInt_t[fnClasses];
     for(Int_t j=0; j<fnClasses; j++) overflow[i][j] = 0;
  }
  if (AliTriggerRunScalers::ConsistencyCheck(fScalersRecord.GetEntriesFast()-1,kFALSE,overflow)){
   AliErrorClass("Trigger counters not in the right order or decreasing!");
  //  scaler->Print();
  //  fScalersRecord.Sort(); 
 }
}
//_____________________________________________________________________________

AliTriggerRunScalers::AliTriggerRunScalers(const AliTriggerRunScalers &run) :
 TObject(),
 fVersion(run.fVersion),
 fRunNumber(run.fRunNumber),
 fnClasses(run.fnClasses),
 fClassIndex(),                    
 fScalersRecord(),
 fScalersRecordESD()
{
// copy constructor
  fClassIndex.Set(run.fClassIndex.GetSize());
for (Int_t i = 0; i < run.fClassIndex.GetSize(); i++) {
    if (run.fClassIndex[i]) fClassIndex.AddAt(run.fClassIndex[i], i);
  }
for (Int_t i = 0; i < run.fScalersRecord.GetEntriesFast(); i++) {
    if (run.fScalersRecord[i]) fScalersRecord.Add(run.fScalersRecord[i]->Clone());
  }
for (Int_t i = 0; i < run.fScalersRecordESD.GetEntriesFast(); i++) {
    if (run.fScalersRecordESD[i]) fScalersRecordESD.Add(run.fScalersRecordESD[i]->Clone());
  }

}
//_____________________________________________________________________________
AliTriggerRunScalers &AliTriggerRunScalers::operator=(const AliTriggerRunScalers& run)
{
// assignment operator
if(&run == this) return *this;
((TObject *)this)->operator=(run);

fVersion = run.fVersion;
fRunNumber = run.fRunNumber;
fnClasses = run.fnClasses;

for (Int_t i = 0; i < run.fClassIndex.GetSize(); i++) {
    if (run.fClassIndex[i]) fClassIndex.AddAt(run.fClassIndex[i], i);
  }
for (Int_t i = 0; i < run.fScalersRecord.GetEntriesFast(); i++) {
    if (run.fScalersRecord[i]) fScalersRecord.Add(run.fScalersRecord[i]->Clone());
  }
for (Int_t i = 0; i < run.fScalersRecordESD.GetEntriesFast(); i++) {
    if (run.fScalersRecordESD[i]) fScalersRecordESD.Add(run.fScalersRecordESD[i]->Clone());
  }
return *this;
} 
//_____________________________________________________________________________
AliTriggerRunScalers* AliTriggerRunScalers::ReadScalers( TString & filename )
{
  // Read scalers from text file(.cnt) provided by CTP 
  // for given run and convert it to root format 
  if( gSystem->AccessPathName( filename.Data() ) ) {
     AliErrorClass( Form( "file (%s) not found", filename.Data() ) );
     return NULL;
  }

  ifstream *file = new ifstream ( filename.Data() );
  if (!*file) {
    AliErrorClass(Form("Error opening file (%s) !\n",filename.Data()));
    file->close();
    delete file;
    return NULL;
  }
  
  AliTriggerRunScalers* rScaler = new AliTriggerRunScalers();
  
  TString strLine;
  Bool_t verflag = kFALSE;
  Bool_t classflag = kFALSE;
  UChar_t nclass = 0;
  while (strLine.ReadLine(*file)) {
    if (strLine.BeginsWith("#")) continue;
    
    TObjArray *tokens = strLine.Tokenize(" \t");
    Int_t ntokens = tokens->GetEntriesFast();
    // 1st line, version, it is one number, 
    if (!verflag) {
      if (ntokens != 1) { 
        AliErrorClass( Form( "Error reading version number from (%s), line :%s\n", 
                              filename.Data() , strLine.Data() ) );  
	delete tokens;
        return NULL;
      }
  //    cout << "Version "<< ((TObjString*)tokens->At(0))->String().Atoi() << endl;
      rScaler->SetVersion( ((TObjString*)tokens->At(0))->String().Atoi() );
      verflag = kTRUE;
      delete tokens;
      continue;
    }
   
    // 2nd line, run number , number of classes, list of classes used in this partition

    if (!classflag) {
      if ( !((TObjString*)tokens->At(1))->String().IsDigit() ) {
        AliErrorClass( Form( "Error reading Run number from (%s)\n", filename.Data() )); 
      }
  //    cout << "Run Number " << ((TObjString*)tokens->At(0))->String().Atoi() << endl;
      rScaler->SetRunNumber( ((TObjString*)tokens->At(0))->String().Atoi() );
      nclass = (UChar_t)((TObjString*)tokens->At(1))->String().Atoi();
  //    cout << "Number of classes " << nclass << endl;
      rScaler->SetNumClasses( nclass );
      if ( nclass != ntokens - 2 ) {
        AliErrorClass( Form( "Error reading number of classes from (%s)\n", filename.Data() )); 
      }
      for (UChar_t i=0; i<nclass; ++i) {
        rScaler->SetClass( i, (Char_t)((TObjString*)tokens->At(2+i))->String().Atoi() );
      }
      classflag = kTRUE;
      delete tokens;
      continue;
    }
    
    // Records
    // Each record consists of 1+(number of classes in partition) lines
    //  1st line of record = time stamp (4 words)
    //                        1st word = ORBIT(24 bits)
    //                        2nd word = Period Counter (28 bit)
    //                        3rd word = seconds (32 bits) from epoch
    //                        4th word = microsecs (32 bits)
    //  other lines = 6 words of counters (L0 before,L0 after, ....)
    if (ntokens != 4 && ntokens !=5) { 
      AliErrorClass( Form( "Error reading timestamp from (%s): line (%s)", 
                            filename.Data(), strLine.Data() )); 
      return NULL;
    }

    UInt_t orbit     = strtoul(((TObjString*)tokens->At(0))->String(), NULL, 10);
    UInt_t period    = strtoul(((TObjString*)tokens->At(1))->String(), NULL, 10);
    UInt_t seconds   = strtoul(((TObjString*)tokens->At(2))->String(), NULL, 10);
    UInt_t microSecs = strtoul(((TObjString*)tokens->At(3))->String(), NULL, 10); 

    AliTriggerScalersRecord * rec = new AliTriggerScalersRecord();
    if(ntokens==5){
      UInt_t tgroup=strtoul(((TObjString*)tokens->At(4))->String(), NULL, 10);
      rec->SetTimeGroup(tgroup);
    }
    rec->SetTimeStamp( orbit, period, seconds, microSecs );
    TString strLine1;
    for (Int_t i=0; i<nclass; ++i) {
      strLine1.ReadLine(*file);
      if (strLine1.BeginsWith("#")) continue;
      TObjArray *tokens1 = strLine1.Tokenize(" \t");
      Int_t ntokens1 = tokens1->GetEntriesFast();
      if( ntokens1 != 6 ) {
        AliErrorClass( Form( "Error reading scalers from (%s): line (%s)", 
			     filename.Data(), strLine1.Data() ));
	delete rec;
	delete tokens1;
	return rScaler;
      }

      UInt_t lOCB = strtoul(((TObjString*)tokens1->At(0))->String(), NULL, 10);
      UInt_t lOCA = strtoul(((TObjString*)tokens1->At(1))->String(), NULL, 10);
      UInt_t l1CB = strtoul(((TObjString*)tokens1->At(2))->String(), NULL, 10);
      UInt_t l1CA = strtoul(((TObjString*)tokens1->At(3))->String(), NULL, 10);
      UInt_t l2CB = strtoul(((TObjString*)tokens1->At(4))->String(), NULL, 10);
      UInt_t l2CA = strtoul(((TObjString*)tokens1->At(5))->String(), NULL, 10);

      rScaler->GetClass(i);
      rec->AddTriggerScalers( rScaler->GetClass(i),
                              lOCB, lOCA, l1CB,
                              l1CA, l2CB, l2CA );

      delete tokens1;
    } 
    rScaler->AddTriggerScalers( rec );
    
    delete tokens;     
  }
  file->close();
  delete file;

  return  rScaler; 
}
  
//_____________________________________________________________________________
Int_t  AliTriggerRunScalers::FindNearestScalersRecord( const AliTimeStamp *stamp ) const
{
   // Find Trigger scaler record with the closest timestamp <= "stamp"
   // using a binary search. 
   // return the index in the array of records, if the timestamp 
   // is out of range return -1

   Int_t   base, position=-1, last, result = 0;
   Int_t op2 = 0;
   
   //fScalersRecord.Sort();

   base = 0;
   last = fScalersRecord.GetEntriesFast();

   while (last >= base) {
      position = (base+last) / 2;
      AliDebug(1, Form("position= %d   base= %d    last= %d  ",position,base,last));
      AliTriggerScalersRecord* rec = (AliTriggerScalersRecord*)fScalersRecord.At(position);
      if( rec && rec->GetTimeStamp()) op2 = 1;
      if( op2 && (result = stamp->Compare(rec->GetTimeStamp())) == 0  )
         return position;  // exact match 
      if (!op2 || result < 0)
         last = position-1;
      else
         base = position+1;
      op2 = 0;  
   }
   if( (position == 0 && result < 0) || position >= fScalersRecord.GetEntriesFast() ) 
    return -1;  // out of range
   else 
    return (result < 0 ) ? position-1 : position; // nearst < stamp   
}
//_____________________________________________________________________________
Int_t AliTriggerRunScalers::ConsistencyCheck(Int_t position,Bool_t correctOverflow, UInt_t** overflow)
{
   //Check if counters are consistent(increase). Example: lOCB(n) < lOCB(n+1) and lOCB > lOCA
   // scalers coding 0,1,2,3,4,5=0b,0a,1b,1a,2b,2a
   // returns: 
   //         1 = decresing time 
   //         2 = too big jump in scalers, looks like some readings are missing
   //         3 = (level+1) > (level)
   if (position == 0){
      if(correctOverflow){
        AliErrorClass("position=0\n");
        return 1;
      }else return 0; // to work correctly in AddScalers
   };
   UInt_t c2[6], c1[6];
   ULong64_t c64[6]; 
   Bool_t increase[6];  
   for(Int_t i=0;i<6;i++){increase[i]=0;}
   ULong64_t const max1 = 4294967295ul;  //32bit counters overflow after 4294967295
   ULong64_t const max2 = 1000000000ul;  //when counters overflow they seem to be decreasing. Assume decrease cannot be smaller than max2.

   AliTriggerScalersRecord* scalers2 = (AliTriggerScalersRecord*)fScalersRecord.At(position);
   AliTriggerScalersRecord* scalers1 = (AliTriggerScalersRecord*)fScalersRecord.At(position-1);
   if (scalers2->Compare((AliTriggerScalersRecord*)fScalersRecord.At(position-1)) == -1) return 1;

   AliTriggerScalersRecordESD* recESD = 0;
   if(correctOverflow){
     recESD = new AliTriggerScalersRecordESD();
     recESD->SetTimeStamp(scalers2->GetTimeStamp());
     recESD->SetTimeGroup(scalers2->GetTimeGroup());
   }
   for( Int_t ic=0; ic<fnClasses; ++ic ){
      TObjArray* scalersArray2 = (TObjArray*)scalers2->GetTriggerScalers();
      AliTriggerScalers* counters2 = (AliTriggerScalers*)scalersArray2->At(ic);
      UChar_t iclass = counters2->GetClassIndex();
      counters2->GetAllScalers(c2);
      TObjArray* scalersArray1 = (TObjArray*)scalers1->GetTriggerScalers();
      AliTriggerScalers* counters1 = (AliTriggerScalers*)scalersArray1->At(ic);
      counters1->GetAllScalers(c1);
      for(Int_t i=0;i<5;i++){
         if ( c2[i] >= c1[i] ) increase[i]=1;
         else if ( c2[i] < c1[i] && (c1[i] - c2[i]) > max2) overflow[i][ic]++;
         else return 2;
      }
      for(Int_t i=0;i<5;i++){
         if ((c2[i] - c1[i]) < (c2[i+1] - c1[i+1]) && increase[i] && increase[i+1] ) {
                 if ( ((c2[i+1] - c1[i+1]) - (c2[i] - c1[i])) < 16 ) {AliWarningClass("Trigger scaler Level[i+1] > Level[i]. Diff < 16!");}
                 else {
  		    cout << " 1 pos= " << position << "i= "<< i << endl;
		    return 3; 
		 }
	 }
         else if ( (max1 - c1[i]+c2[i]) < (c2[i+1] - c1[i+1]) && overflow[i][ic] && increase[i+1] ) {
                 if ( ((c2[i+1] - c1[i+1]) - (max1 - c1[i]+c2[i])) < 16 ) {AliWarningClass("Trigger scaler Level[i+1] > Level[i]. Diff < 16!");}
                 else{
  		    cout << " 2 pos= " << position << "i= "<< i << endl;
		    return 3; 
		 }
	 }
         else if ( (c2[i] - c1[i]) < (max1 - c1[i+1] + c2[i+1]) && increase[i] && overflow[i+1][ic] ) {
                 if ( ((max1 - c1[i+1] + c2[i+1]) - (c2[i] - c1[i])) < 16 ) {AliWarningClass("Trigger scaler Level[i+1] > Level[i]. Diff < 16!");}
                 else{
  		    cout << " 3 pos= " << position << "i= "<< i << endl;
		    return 3; 
		    }
	 }
         else if ( (max1 - c1[i] + c2[i] ) < (max1 - c1[i+1] + c2[i+1]) && overflow[i][ic] && overflow[i+1][ic] ) {
                 if ( ((max1 - c1[i+1] + c2[i+1]) - (max1 - c1[i] + c2[i] )) < 16 ) {AliWarningClass("Trigger scaler Level[i+1] > Level[i]. Diff < 16!");}
                 else{
  		    cout << " 4 pos= " << position << "i= "<< i << endl;
		    return 3;
		    }
       }
      }
      if(correctOverflow){ 
        for(Int_t i=0;i<6;i++){ c64[i]=c2[i]+max1*overflow[i][ic]; }
        AliTriggerScalersESD* s= new AliTriggerScalersESD(iclass,c64);
        recESD->AddTriggerScalers(s);
         }

 }
 if(correctOverflow)fScalersRecordESD.AddLast(recESD);
 return 0;
}
//____________________________________________________________________________
Int_t AliTriggerRunScalers::CorrectScalersOverflow()
{
 // Run over fScalersRecord, check overflow using CheckConsistency methos
 // and save corrected result in fScalersRecordESD.

 // Temporary fix for the OCDB entries written with v4-16-Release
 // where the wrong sorting was used
 fScalersRecord.Sort();
 UInt_t c1[6];
 ULong64_t c64[6];
 AliTriggerScalersRecordESD* recESD = new AliTriggerScalersRecordESD();
 // add 0
 AliTriggerScalersRecord* scalers = (AliTriggerScalersRecord*)fScalersRecord.At(0);
 recESD->SetTimeStamp(scalers->GetTimeStamp());
 recESD->SetTimeGroup(scalers->GetTimeGroup());
 for( Int_t ic=0; ic<fnClasses; ++ic ){
    TObjArray* scalersArray = (TObjArray*)scalers->GetTriggerScalers();
    AliTriggerScalers* counters = (AliTriggerScalers*)scalersArray->At(ic);
    counters->GetAllScalers(c1);
    UChar_t iclass = counters->GetClassIndex();
    for(Int_t i=0; i<6; i++)c64[i]=c1[i];
    AliTriggerScalersESD* s= new AliTriggerScalersESD(iclass,c64);
    recESD->AddTriggerScalers(s);
 }
 fScalersRecordESD.AddLast(recESD);

 UInt_t* overflow[6];
 for(Int_t i=0; i<6; i++) {
    overflow[i] = new UInt_t[fnClasses];
    for(Int_t j=0; j<fnClasses; j++) overflow[i][j] = 0;
 }

 for(Int_t i=1;i<fScalersRecord.GetEntriesFast(); i++){
  if(ConsistencyCheck(i,kTRUE,overflow)){
    fScalersRecord.At(i)->Print();
    fScalersRecord.At(i-1)->Print();
    fScalersRecordESD.SetOwner(); 
    fScalersRecordESD.Delete(); 
    AliErrorClass("Inconsistent scalers, they will not be provided.\n");
    return 1;
  }
 }
 if(fScalersRecordESD.GetEntriesFast() != fScalersRecord.GetEntriesFast()){
    AliErrorClass("Internal error: #scalers ESD != #scalers \n");
    return 1;
 }
 return 0;
}
//_____________________________________________________________________________
AliTriggerScalersESD* AliTriggerRunScalers::GetScalersForEventClass(const AliTimeStamp* stamp,const Int_t classIndex) const
{
 // Find scalers for event for class in fScalersRecordESD
 // Assumes that fScalerRecord = fScalerRecordESD
 Int_t position = FindNearestScalersRecord(stamp);
 if ( position == -1 ) { 
  AliErrorClass("Event AliTimeStamp out of range!");
  return 0; 
 }
 // check also position=max
 AliTriggerScalersRecordESD* scalrec1 = (AliTriggerScalersRecordESD*)fScalersRecordESD.At(position);
 AliTriggerScalersRecordESD* scalrec2 = (AliTriggerScalersRecordESD*)fScalersRecordESD.At(position+1);
 TObjArray* scalers1 = (TObjArray*)scalrec1->GetTriggerScalers();
 TObjArray* scalers2 = (TObjArray*)scalrec2->GetTriggerScalers();

 if(scalers1->GetEntriesFast() != fnClasses){
  AliErrorClass("Internal error: #classes in RecordESD != fnClasses\n");
  return 0; 
 }

 AliTriggerScalersESD *s1,*s2;
 for ( Int_t ic=0; ic < (Int_t)fnClasses; ++ic ){

      s1 = (AliTriggerScalersESD*)scalers1->At(ic);
      s2 = (AliTriggerScalersESD*)scalers2->At(ic);

      Bool_t classfound = (s1->GetClassIndex() == classIndex) && (s2->GetClassIndex() == classIndex);
      if(classfound){
        ULong64_t max = 16777216ul;
        AliTriggerScalersRecordESD* scalrec0 = (AliTriggerScalersRecordESD*)fScalersRecordESD.At(0);
        TObjArray* scalers0 = (TObjArray*)scalrec0->GetTriggerScalers();
        AliTriggerScalersESD *s0 = (AliTriggerScalersESD*)scalers0->At(ic);
	ULong64_t base[6],c1[6],c2[6],cint[6];
        ULong64_t orbit = max*(stamp->GetPeriod()) + stamp->GetOrbit();
	s0->GetAllScalers(base);
	s1->GetAllScalers(c1);
	s2->GetAllScalers(c2);
	ULong64_t orbit1 = max*(scalrec1->GetTimeStamp()->GetPeriod())+scalrec1->GetTimeStamp()->GetOrbit();
	ULong64_t orbit2 = max*(scalrec2->GetTimeStamp()->GetPeriod())+scalrec2->GetTimeStamp()->GetOrbit();
        for(Int_t i=0;i<6;i++){
	   Double_t slope=Double_t(c2[i]-c1[i])/Double_t(orbit2-orbit1);
	   cint[i]=ULong64_t(slope*(orbit-orbit1)) +c1[i] -base[i];
	}
	AliTriggerScalersESD* result = new AliTriggerScalersESD(classIndex,cint);
        return result;
      }
 }
 AliErrorClass(Form("Classindex %i not found.\n",classIndex));
 return 0;
}

//_____________________________________________________________________________
const AliTriggerScalersRecordESD* AliTriggerRunScalers::GetScalersDeltaForEvent(const AliTimeStamp* stamp) const
{
 // Find scalers for event for class in fScalersRecordESD
 // Assumes that fScalerRecord = fScalerRecordESD
 Int_t position = FindNearestScalersRecord(stamp);
 if ( position == -1 ) { 
  AliErrorClass("Event AliTimeStamp out of range!");
  return 0; 
 }
 // check also position=max
 if (fScalersRecordESD.GetEntriesFast()<2) return 0;
 
 AliTriggerScalersRecordESD* scalrec1 = (AliTriggerScalersRecordESD*)fScalersRecordESD.At(position);
 AliTriggerScalersRecordESD* scalrec2 = (AliTriggerScalersRecordESD*)fScalersRecordESD.At(position+1);
 TObjArray* scalers1 = (TObjArray*)scalrec1->GetTriggerScalers();
 TObjArray* scalers2 = (TObjArray*)scalrec2->GetTriggerScalers();

 if(scalers1->GetEntriesFast() != fnClasses){
  AliErrorClass("Internal error: #classes in RecordESD != fnClasses\n");
  return 0; 
 }
 AliTriggerScalersRecordESD *scalrec = new AliTriggerScalersRecordESD();
 AliTriggerScalersESD *s1,*s2;
 for ( Int_t ic=0; ic < (Int_t)fnClasses; ++ic ){

  s1 = (AliTriggerScalersESD*)scalers1->At(ic);
  s2 = (AliTriggerScalersESD*)scalers2->At(ic);

  ULong64_t c1[6],c2[6],cint[6];
  s1->GetAllScalers(c1);
  s2->GetAllScalers(c2);
  for(Int_t i=0;i<6;i++){
     cint[i]=c2[i]-c1[i];
  }
  AliTriggerScalersESD* result = new AliTriggerScalersESD(s1->GetClassIndex(),cint);
  scalrec->AddTriggerScalers(result);
 }
 UInt_t max = 16777216ul;
 UInt_t orbit, period;
 
 if (scalrec2->GetTimeStamp()->GetOrbit() > scalrec1->GetTimeStamp()->GetOrbit()) {
  orbit = scalrec2->GetTimeStamp()->GetOrbit() - scalrec1->GetTimeStamp()->GetOrbit();
  period = scalrec2->GetTimeStamp()->GetPeriod() - scalrec1->GetTimeStamp()->GetPeriod();
 }
 else {
  orbit = max - (scalrec1->GetTimeStamp()->GetOrbit() - scalrec2->GetTimeStamp()->GetOrbit());
  period = scalrec2->GetTimeStamp()->GetPeriod() - scalrec1->GetTimeStamp()->GetPeriod() - 1;
 }
 
 AliTimeStamp *timestamp = new AliTimeStamp(orbit, period, 0);
 scalrec->SetTimeStamp(timestamp);
 delete timestamp;
 return scalrec;
}
//_____________________________________________________________________________
const AliTriggerScalersRecordESD* AliTriggerRunScalers::GetScalersDeltaForRun() const
{
 // Find scalers for event for class in fScalersRecordESD
 // Assumes that fScalerRecord = fScalerRecordESD
 if (fScalersRecordESD.GetEntriesFast()<2) return 0;

 AliTriggerScalersRecordESD* scalrec1 = (AliTriggerScalersRecordESD*)fScalersRecordESD.At(0);
 AliTriggerScalersRecordESD* scalrec2 = (AliTriggerScalersRecordESD*)fScalersRecordESD.At(fScalersRecord.GetEntriesFast()-1);
 TObjArray* scalers1 = (TObjArray*)scalrec1->GetTriggerScalers();
 TObjArray* scalers2 = (TObjArray*)scalrec2->GetTriggerScalers();

 if(scalers1->GetEntriesFast() != fnClasses){
  AliErrorClass("Internal error: #classes in RecordESD != fnClasses\n");
  return 0; 
 }
 AliTriggerScalersESD *s1,*s2;
 AliTriggerScalersRecordESD *scalrec = new AliTriggerScalersRecordESD();
 for ( Int_t ic=0; ic < (Int_t)fnClasses; ++ic ){

  s1 = (AliTriggerScalersESD*)scalers1->At(ic);
  s2 = (AliTriggerScalersESD*)scalers2->At(ic);

  ULong64_t c1[6],c2[6],cint[6];
  s1->GetAllScalers(c1);
  s2->GetAllScalers(c2);
  for(Int_t i=0;i<6;i++){
     cint[i]=c2[i]-c1[i];
  }
  AliTriggerScalersESD* result = new AliTriggerScalersESD(s1->GetClassIndex(),cint);
  scalrec->AddTriggerScalers(result);
 }
 UInt_t max = 16777216ul;
 UInt_t orbit, period;
 
 if (scalrec2->GetTimeStamp()->GetOrbit() > scalrec1->GetTimeStamp()->GetOrbit()) {
  orbit = scalrec2->GetTimeStamp()->GetOrbit() - scalrec1->GetTimeStamp()->GetOrbit();
  period = scalrec2->GetTimeStamp()->GetPeriod() - scalrec1->GetTimeStamp()->GetPeriod();
 }
 else {
  orbit = max - (scalrec1->GetTimeStamp()->GetOrbit() - scalrec2->GetTimeStamp()->GetOrbit());
  period = scalrec2->GetTimeStamp()->GetPeriod() - scalrec1->GetTimeStamp()->GetPeriod() - 1;
 }
 
 AliTimeStamp *timestamp = new AliTimeStamp(orbit, period, 0);
 scalrec->SetTimeStamp(timestamp);
 delete timestamp;
 return scalrec;
}
//_____________________________________________________________________________
Bool_t AliTriggerRunScalers::CalculateMu(Double_t &mu, Double_t &errmu, ULong64_t countsB, ULong64_t countsAC, UShort_t nB, UShort_t nAC, UInt_t orbits, Bool_t bkgCorr, Double_t triggerEff, Double_t errorEff) 
{
 
   if (nB!=0 && orbits!=0)  {
      Double_t pB = (Double_t)countsB/((ULong64_t)nB*orbits); // probability for B trigger
      if (!bkgCorr || nAC==0 ) {
         mu = -log(1-pB)/triggerEff;
         errmu = TMath::Sqrt(pB/((1-pB)*nB*orbits) + mu*mu*errorEff*errorEff/(triggerEff*triggerEff)); //
         return kTRUE;
      }
      else 
      {
         Double_t pAC = (Double_t)countsAC/((ULong64_t)nAC*orbits); // probability for AC trigger (background)
         mu = ( log(1.-pAC) - log(1.-pB) )/triggerEff;
         // error
         errmu =  TMath::Sqrt(pB/((1.-pB)*nB*orbits) + pAC/((1.-pAC)*nAC*orbits) /*- 2*TMath::Sqrt(pB*pAC/(nB*nAC*(1.-pB)*(1.-pAC)))/orbits*/ + mu*mu*errorEff*errorEff/(triggerEff*triggerEff)); // assume no correlation between B and AC rates, hence no cov term in error
         return kTRUE;
      }
   }
   return kFALSE;
}
//_____________________________________________________________________________
Bool_t AliTriggerRunScalers::CalculateMu(Double_t &mu, Double_t &errmu, ULong64_t countsB, ULong64_t countsAC, ULong64_t beamB, UShort_t nB, UShort_t nAC, Bool_t bkgCorr, Double_t triggerEff, Double_t errorEff) 
{

   if (beamB!=0)  {
      Double_t pB = (Double_t)countsB/beamB; // probability for B trigger
      if (!bkgCorr || nAC==0 || nB==0) {
         mu = -log(1-pB)/triggerEff;
         errmu = TMath::Sqrt(pB/((1-pB)*beamB) + mu*mu*errorEff*errorEff/(triggerEff*triggerEff)); //
         return kTRUE;
      }
      else
      {
         Double_t pAC = (Double_t)countsAC/((ULong64_t)nAC*beamB/(Double_t)nB); // probability for AC trigger (background)
         mu = ( log(1-pAC) - log(1-pB) )/triggerEff;
         // error
         errmu =  TMath::Sqrt(pB/((1-pB)*beamB) + pAC/((1-pAC)*nAC*beamB/nB) + mu*mu*errorEff*errorEff/(triggerEff*triggerEff));
         return kTRUE;
      }
   }
   return kFALSE;
}
//_____________________________________________________________________________
ULong64_t AliTriggerRunScalers::GetDeltaScaler(const AliTriggerScalersRecordESD* scalRec1, const AliTriggerScalersRecordESD* scalRec2, Int_t classIndex, TString level)
{
  const AliTriggerScalersESD* scalers1 = scalRec1->GetTriggerScalersForClass(classIndex);
  const AliTriggerScalersESD* scalers2 = scalRec2->GetTriggerScalersForClass(classIndex);
  ULong64_t s1=0;
  ULong64_t s2=0;

  if (level == "l0b") {s1=scalers1->GetLOCB(); s2=scalers2->GetLOCB();}
  else if (level == "l0a") {s1=scalers1->GetLOCA(); s2=scalers2->GetLOCA();}
  else if (level == "l1b") {s1=scalers1->GetL1CB(); s2=scalers2->GetL1CB();}
  else if (level == "l1a") {s1=scalers1->GetL1CA(); s2=scalers2->GetL1CA();}
  else if (level == "l2b") {s1=scalers1->GetL2CB(); s2=scalers2->GetL2CB();}
  else if (level == "l2a") {s1=scalers1->GetL2CA(); s2=scalers2->GetL2CA();}
  else return 0;

  return s2-s1;
}
//_____________________________________________________________________________
Double_t AliTriggerRunScalers::GetDeltaTime(const AliTriggerScalersRecordESD* scalRec1, const AliTriggerScalersRecordESD* scalRec2)
{
  const AliTimeStamp* stamp1 = scalRec1->GetTimeStamp();
  const AliTimeStamp* stamp2 = scalRec2->GetTimeStamp();
  UInt_t orbit1 = stamp1->GetOrbit();
  UInt_t orbit2 = stamp2->GetOrbit();
  UInt_t period1 = stamp1->GetPeriod();
  UInt_t period2 = stamp2->GetPeriod();
  UInt_t max = 16777216;  // The period counter increases when 24 bit orbit counter overflow
  Double_t orbitSec = 89.1*1.e-6; // Length of 1 orbit in seconds
  return ((period2 - period1)*max + (orbit2-orbit1))*orbitSec;
}
//_____________________________________________________________________________
UInt_t AliTriggerRunScalers::GetDeltaOrbits(const AliTriggerScalersRecordESD* scalRec1, const AliTriggerScalersRecordESD* scalRec2)
{
  const AliTimeStamp* stamp1 = scalRec1->GetTimeStamp();
  const AliTimeStamp* stamp2 = scalRec2->GetTimeStamp();
  UInt_t orbit1 = stamp1->GetOrbit();
  UInt_t orbit2 = stamp2->GetOrbit();
  UInt_t period1 = stamp1->GetPeriod();
  UInt_t period2 = stamp2->GetPeriod();
  UInt_t max = 16777216;  // The period counter increases when 24 bit orbit counter overflow
  return (period2 - period1)*max + (orbit2-orbit1);
}
//_____________________________________________________________________________
Bool_t AliTriggerRunScalers::GetScalerRate(Double_t &rate, Double_t &error, const AliTriggerScalersRecordESD* scalRec1, const AliTriggerScalersRecordESD* scalRec2, Int_t classIndex, TString level)
{
  if (level != "l0b" && level != "l0a" && level != "l1b" && level != "l1a" && level != "l2b" && level != "l2a") return kFALSE;
 
  ULong64_t scaler = GetDeltaScaler(scalRec1, scalRec2, classIndex, level);
  Double_t time = GetDeltaTime(scalRec1, scalRec2 );
  if (time==0.) return kFALSE;
  rate = (Double_t)scaler/time;
  error = (Double_t)sqrt(scaler)/time;
  return kTRUE;
}
//_____________________________________________________________________________
Bool_t AliTriggerRunScalers::GetScalerRatePerBC(Double_t &rate, Double_t &error, const AliTriggerScalersRecordESD* scalRec1, const AliTriggerScalersRecordESD* scalRec2, AliTriggerConfiguration* cfg, Int_t classIndex, TString level)
{
  if (level != "l0b" && level != "l0a" && level != "l1b" && level != "l1a" && level != "l2b" && level != "l2a") return kFALSE;
  const AliTriggerClass* trgclass = cfg->GetTriggerClass(classIndex);
  AliTriggerBCMask* bcMask = new AliTriggerBCMask();
  if (trgclass) bcMask = (AliTriggerBCMask*)trgclass->GetBCMask();
  Int_t nBC=0;
  if (TString(bcMask->GetName()).CompareTo("NONE")!=0){
     nBC = (UShort_t)bcMask->GetNUnmaskedBCs();
  }
  if (nBC<1) return kFALSE;
  ULong64_t scaler = GetDeltaScaler(scalRec1, scalRec2, classIndex, level);
  Double_t time = GetDeltaTime(scalRec1, scalRec2 );
  if (time==0.) return kFALSE;
  rate = (Double_t)scaler/time/nBC;
  error = (Double_t)sqrt(scaler)/time/nBC;
  return kTRUE;
}
//_____________________________________________________________________________
Bool_t AliTriggerRunScalers::GetClassL2L0(Double_t &l2l0, Double_t &error, const AliTriggerScalersRecordESD* scalRec1, const AliTriggerScalersRecordESD* scalRec2, Int_t classIndex)
{
  ULong64_t l0 = GetDeltaScaler(scalRec1, scalRec2, classIndex, "l0b");
  ULong64_t l2 = GetDeltaScaler(scalRec1, scalRec2, classIndex, "l2a");
  if (l0!=0) {
    l2l0 = (Double_t)l2/l0;
    error = (Double_t)sqrt(l2-l2*l2/l0)/l0;
    return kTRUE;
  }
  return kFALSE;
}
//_____________________________________________________________________________
Bool_t AliTriggerRunScalers::GetMuFromClassScaler(Double_t &mu, Double_t &errmu, const char* className, const AliTriggerScalersRecordESD* scalRec1, const AliTriggerScalersRecordESD* scalRec2, AliTriggerConfiguration* cfg, Bool_t colBCsFromFillScheme, Bool_t bkgCorr, Double_t triggerEff, Double_t errorEff)
{
  // className = the first part of the class name. For example CINT1 for the class CINT1-B-NOPF-ALL
  // colBCsFromFillScheme=kTRUE - use filling scheme and orbit counter
  // colBCsFromFillScheme=kFALSE - use cbeamb scaler to get number of colliding BCs
  // One can also switch background correction ON or OFF with bkgCorr
  TObjArray cint1bNames;
  TObjArray cint1acNames;
  TObjArray cbeambNames;
  //
  cint1bNames.Add(new TNamed(Form("%s-ABCE",className),NULL));
  cint1bNames.Add(new TNamed(Form("%s-B",className),NULL));
  cint1acNames.Add(new TNamed(Form("%s-AC",className),NULL));
  cbeambNames.Add(new TNamed("CBEAMB",NULL));
  cbeambNames.Add(new TNamed("CTRUE",NULL));
  //
  Int_t cint1bIndex=-1, cint1acIndex=-1, cbeambIndex=-1;
  TString nameString;
  //
  for (Int_t i=0; i<cfg->GetClasses().GetEntriesFast(); i++ ) {
     nameString = cfg->GetClassNameFromIndex(i);
     for (Int_t j=0; j<cint1bNames.GetEntriesFast(); j++ ) {
       if (nameString.BeginsWith(cint1bNames.At(j)->GetName())) cint1bIndex = i;
     }
     for (Int_t j=0; j<cint1acNames.GetEntriesFast(); j++ ) {
       if (nameString.BeginsWith(cint1acNames.At(j)->GetName())) cint1acIndex = i;
     }
     for (Int_t j=0; j<cbeambNames.GetEntriesFast(); j++ ) {
       if (nameString.BeginsWith(cbeambNames.At(j)->GetName())) cbeambIndex = i;
     }
     nameString.Clear();
  }
  //
  ULong64_t cint1b = 0;
  UShort_t nB=0;
  if (cint1bIndex!=-1) {
     cint1b=GetDeltaScaler(scalRec1, scalRec2, cint1bIndex, "l0b");
     const AliTriggerClass* cint1bClass = cfg->GetTriggerClass(cint1bIndex);
     AliTriggerBCMask* cint1bBCMask = new AliTriggerBCMask();
     if (cint1bClass) cint1bBCMask = (AliTriggerBCMask*)cint1bClass->GetBCMask();
     if (TString(cint1bBCMask->GetName()).CompareTo("NONE")!=0){
        nB = (UShort_t)cint1bBCMask->GetNUnmaskedBCs();
     }
  } else return kFALSE;
  ULong64_t cint1ac = 0;
  UShort_t nAC=0;
  if (cint1acIndex!=-1) {
     cint1ac=GetDeltaScaler(scalRec1, scalRec2, cint1acIndex, "l0b");
     AliTriggerClass* cint1acClass = cfg->GetTriggerClass(cint1acIndex);
     AliTriggerBCMask* cint1acBCMask = new AliTriggerBCMask();
     if (cint1acClass) cint1acBCMask = (AliTriggerBCMask*)cint1acClass->GetBCMask();
     if (TString(cint1acBCMask->GetName()).CompareTo("NONE")!=0){
        nAC = (UShort_t)cint1acBCMask->GetNUnmaskedBCs();
     }
  }
  ULong64_t cbeamb = 0;
  if (cbeambIndex!=-1) cbeamb=GetDeltaScaler(scalRec1, scalRec2, cbeambIndex, "l0b");
  UInt_t orbits = GetDeltaOrbits(scalRec1, scalRec2);
  //
  if (bkgCorr && (nB==0 || nAC==0 )) return kFALSE;
  //
  if (colBCsFromFillScheme) {
     if (CalculateMu(mu, errmu, cint1b, cint1ac, nB, nAC, orbits, bkgCorr, triggerEff, errorEff)) return kTRUE;
  }
  else {
     if (cint1b!=0 && cbeamb==0) return kFALSE;
     if (CalculateMu(mu, errmu, cint1b, cint1ac, cbeamb, nB, nAC, bkgCorr, triggerEff, errorEff)) return kTRUE;
  }
  //
  return kFALSE;
}
//_____________________________________________________________________________
ULong64_t AliTriggerRunScalers::GetDeltaScalerForRun(Int_t classIndex, TString level)
{
  if (fScalersRecordESD.GetEntriesFast()==0) {
     if (CorrectScalersOverflow()==1) return 0;
  }
  if (level != "l0b" && level != "l0a" && level != "l1b" && level != "l1a" && level != "l2b" && level != "l2a") return 0;

  if (fScalersRecordESD.GetEntriesFast()>1) { 
     const AliTriggerScalersRecordESD* scalRec1 = (AliTriggerScalersRecordESD*)fScalersRecordESD.At(0);
     const AliTriggerScalersRecordESD* scalRec2 = (AliTriggerScalersRecordESD*)fScalersRecordESD.At(fScalersRecordESD.GetEntriesFast()-1);
     return GetDeltaScaler(scalRec1, scalRec2, classIndex, level);
  }
  return 0;
}
//_____________________________________________________________________________
Bool_t AliTriggerRunScalers::GetScalerRateForRun(Double_t &rate, Double_t &error, Int_t classIndex, TString level)
{
  if (level != "l0b" && level != "l0a" && level != "l1b" && level != "l1a" && level != "l2b" && level != "l2a") return 0;
  if (fScalersRecordESD.GetEntriesFast()==0) {
     if (CorrectScalersOverflow()==1) return kFALSE;
  }
  if (fScalersRecordESD.GetEntriesFast()>=4) {
     const AliTriggerScalersRecordESD* scalRec1 = (AliTriggerScalersRecordESD*)fScalersRecordESD.At(0);
     const AliTriggerScalersRecordESD* scalRec2 = (AliTriggerScalersRecordESD*)fScalersRecordESD.At(fScalersRecordESD.GetEntriesFast()-1);
     if (GetScalerRate(rate, error, scalRec1, scalRec2, classIndex, level)) return kTRUE;
  }
  return kFALSE;
}
//_____________________________________________________________________________
Bool_t AliTriggerRunScalers::GetClassL2L0ForRun(Double_t &l2l0, Double_t &error, Int_t classIndex)
{
  if (fScalersRecordESD.GetEntriesFast()==0) {
     if (CorrectScalersOverflow()==1) return kFALSE;
  }
  if (fScalersRecordESD.GetEntriesFast()>=4) {
     const AliTriggerScalersRecordESD* scalRec1 = (AliTriggerScalersRecordESD*)fScalersRecordESD.At(0);
     const AliTriggerScalersRecordESD* scalRec2 = (AliTriggerScalersRecordESD*)fScalersRecordESD.At(fScalersRecordESD.GetEntriesFast()-1);  // Skip the first and last scaler readings for timing reasons

     if ( GetClassL2L0(l2l0, error, scalRec1, scalRec2, classIndex)) return kTRUE;
  }
  return kFALSE;
}
//_____________________________________________________________________________
TGraphErrors* AliTriggerRunScalers::GetGraphScalerRate(const char* className, TString level, AliTriggerConfiguration* cfg)
{
  Int_t classIndex = cfg->GetClassIndexFromName(className);
  if (classIndex == -1) return 0;
  if (fScalersRecordESD.GetEntriesFast()==0) {
     if (CorrectScalersOverflow()==1) return 0;
  }
  if (level != "l0b" && level != "l0a" && level != "l1b" && level != "l1a" && level != "l2b" && level != "l2a") return 0;
  Int_t nent = fScalersRecordESD.GetEntriesFast();
  Double_t* time = new Double_t[nent];
  Double_t* etime = new Double_t[nent];
  Double_t* rate = new Double_t[nent];
  Double_t* erate = new Double_t[nent];
  for (Int_t i=0;i<nent-1;i++) {
     if (i>0) time[i] = time[i-1]+GetDeltaTime((AliTriggerScalersRecordESD*)fScalersRecordESD.At(i-1), (AliTriggerScalersRecordESD*)fScalersRecordESD.At(i))/2.+GetDeltaTime((AliTriggerScalersRecordESD*)fScalersRecordESD.At(i), (AliTriggerScalersRecordESD*)fScalersRecordESD.At(i+1))/2.;
     else time[0] = GetDeltaTime((AliTriggerScalersRecordESD*)fScalersRecordESD.At(0), (AliTriggerScalersRecordESD*)fScalersRecordESD.At(1))/2.;
     etime[i] = GetDeltaTime((AliTriggerScalersRecordESD*)fScalersRecordESD.At(i), (AliTriggerScalersRecordESD*)fScalersRecordESD.At(i+1))/2.;
     if (!GetScalerRate( rate[i], erate[i], (AliTriggerScalersRecordESD*)fScalersRecordESD.At(i), (AliTriggerScalersRecordESD*)fScalersRecordESD.At(i+1),  classIndex, level)) {rate[i]=-1; erate[i]=0.;}
  }
  TGraphErrors* graph = new TGraphErrors(nent-1, time, rate, etime, erate);  
  return graph;
}
//_____________________________________________________________________________
TGraphErrors* AliTriggerRunScalers::GetGraphScalerL2L0Ratio(const char* className, AliTriggerConfiguration* cfg)
{
  Int_t classIndex = cfg->GetClassIndexFromName(className);
  if (classIndex == -1) return 0;
  if (fScalersRecordESD.GetEntriesFast()==0) {
     if (CorrectScalersOverflow()==1) return 0;
  }
  Int_t nent = fScalersRecordESD.GetEntriesFast();
  Double_t* time = new Double_t[nent];
  Double_t* etime = new Double_t[nent];
  Double_t* ratio = new Double_t[nent];
  Double_t* eratio = new Double_t[nent];

  for (Int_t i=0;i<nent-1;i++) {
     if (i>0) time[i] = time[i-1]+GetDeltaTime((AliTriggerScalersRecordESD*)fScalersRecordESD.At(i-1), (AliTriggerScalersRecordESD*)fScalersRecordESD.At(i))/2.+GetDeltaTime((AliTriggerScalersRecordESD*)fScalersRecordESD.At(i), (AliTriggerScalersRecordESD*)fScalersRecordESD.At(i+1))/2.;
     else time[0] = GetDeltaTime((AliTriggerScalersRecordESD*)fScalersRecordESD.At(0), (AliTriggerScalersRecordESD*)fScalersRecordESD.At(1))/2.;
     etime[i] = GetDeltaTime((AliTriggerScalersRecordESD*)fScalersRecordESD.At(i), (AliTriggerScalersRecordESD*)fScalersRecordESD.At(i+1))/2.;
     if (!GetClassL2L0(ratio[i], eratio[i], (AliTriggerScalersRecordESD*)fScalersRecordESD.At(i), (AliTriggerScalersRecordESD*)fScalersRecordESD.At(i+1),  classIndex)) {ratio[i]=-1; eratio[i]=0.;}
  }
  TGraphErrors* graph = new TGraphErrors(nent-1, time, ratio, etime, eratio);  
  return graph;
}
//_____________________________________________________________________________
TGraphErrors* AliTriggerRunScalers::GetGraphMu(AliTriggerConfiguration* cfg, const char* className, Bool_t colBCsFromFillScheme, Bool_t bkgCorr, Double_t triggerEff, Double_t errorEff)
{
  if (fScalersRecordESD.GetEntriesFast()==0) {
     if (CorrectScalersOverflow()==1) return 0;
  }
  Int_t nent = fScalersRecordESD.GetEntriesFast();
  Double_t* time = new Double_t[nent];
  Double_t* etime = new Double_t[nent];
  Double_t* mu = new Double_t[nent];
  Double_t* emu = new Double_t[nent];
  
  for (Int_t i=0;i<nent-1;i++) {
     Double_t m=0.;
     Double_t err=0.;
     if (i!=0) time[i] = time[i-1]+GetDeltaTime((AliTriggerScalersRecordESD*)fScalersRecordESD.At(i-1), (AliTriggerScalersRecordESD*)fScalersRecordESD.At(i))/2.+GetDeltaTime((AliTriggerScalersRecordESD*)fScalersRecordESD.At(i), (AliTriggerScalersRecordESD*)fScalersRecordESD.At(i+1))/2.;
     else time[0] = GetDeltaTime((AliTriggerScalersRecordESD*)fScalersRecordESD.At(0), (AliTriggerScalersRecordESD*)fScalersRecordESD.At(1))/2.;
     etime[i] = GetDeltaTime((AliTriggerScalersRecordESD*)fScalersRecordESD.At(i), (AliTriggerScalersRecordESD*)fScalersRecordESD.At(i+1))/2.;
     if (GetMuFromClassScaler( m,err, className, (AliTriggerScalersRecordESD*)fScalersRecordESD.At(i), (AliTriggerScalersRecordESD*)fScalersRecordESD.At(i+1), cfg, colBCsFromFillScheme, bkgCorr, triggerEff, errorEff)!=kFALSE) {mu[i]=m; emu[i]=err;}
     else {mu[i]=-1.; emu[i]=0.;}
  }
  TGraphErrors* graph = new TGraphErrors(nent-1, time, mu, etime, emu);  
  return graph;
}
//_____________________________________________________________________________
void AliTriggerRunScalers::Print( const Option_t* ) const
{
   // Print
  cout << "Trigger Scalers Record per Run: " << endl;
  cout << "  File version :" <<  fVersion << endl;            
  cout << "  Run Number :" <<  fRunNumber << endl;          
  cout << "  Number of Classes :" <<  (Int_t)fnClasses << endl;          
  cout << "    Classes ID:";
  for( Int_t i=0; i<fnClasses; ++i ) 
    cout << "  " << (Int_t)fClassIndex[i];       
  cout << endl; 
    
  for( Int_t i=0; i<fScalersRecord.GetEntriesFast(); ++i ) 
     ((AliTriggerScalersRecord*)fScalersRecord.At(i))->Print();
}

 AliTriggerRunScalers.cxx:1
 AliTriggerRunScalers.cxx:2
 AliTriggerRunScalers.cxx:3
 AliTriggerRunScalers.cxx:4
 AliTriggerRunScalers.cxx:5
 AliTriggerRunScalers.cxx:6
 AliTriggerRunScalers.cxx:7
 AliTriggerRunScalers.cxx:8
 AliTriggerRunScalers.cxx:9
 AliTriggerRunScalers.cxx:10
 AliTriggerRunScalers.cxx:11
 AliTriggerRunScalers.cxx:12
 AliTriggerRunScalers.cxx:13
 AliTriggerRunScalers.cxx:14
 AliTriggerRunScalers.cxx:15
 AliTriggerRunScalers.cxx:16
 AliTriggerRunScalers.cxx:17
 AliTriggerRunScalers.cxx:18
 AliTriggerRunScalers.cxx:19
 AliTriggerRunScalers.cxx:20
 AliTriggerRunScalers.cxx:21
 AliTriggerRunScalers.cxx:22
 AliTriggerRunScalers.cxx:23
 AliTriggerRunScalers.cxx:24
 AliTriggerRunScalers.cxx:25
 AliTriggerRunScalers.cxx:26
 AliTriggerRunScalers.cxx:27
 AliTriggerRunScalers.cxx:28
 AliTriggerRunScalers.cxx:29
 AliTriggerRunScalers.cxx:30
 AliTriggerRunScalers.cxx:31
 AliTriggerRunScalers.cxx:32
 AliTriggerRunScalers.cxx:33
 AliTriggerRunScalers.cxx:34
 AliTriggerRunScalers.cxx:35
 AliTriggerRunScalers.cxx:36
 AliTriggerRunScalers.cxx:37
 AliTriggerRunScalers.cxx:38
 AliTriggerRunScalers.cxx:39
 AliTriggerRunScalers.cxx:40
 AliTriggerRunScalers.cxx:41
 AliTriggerRunScalers.cxx:42
 AliTriggerRunScalers.cxx:43
 AliTriggerRunScalers.cxx:44
 AliTriggerRunScalers.cxx:45
 AliTriggerRunScalers.cxx:46
 AliTriggerRunScalers.cxx:47
 AliTriggerRunScalers.cxx:48
 AliTriggerRunScalers.cxx:49
 AliTriggerRunScalers.cxx:50
 AliTriggerRunScalers.cxx:51
 AliTriggerRunScalers.cxx:52
 AliTriggerRunScalers.cxx:53
 AliTriggerRunScalers.cxx:54
 AliTriggerRunScalers.cxx:55
 AliTriggerRunScalers.cxx:56
 AliTriggerRunScalers.cxx:57
 AliTriggerRunScalers.cxx:58
 AliTriggerRunScalers.cxx:59
 AliTriggerRunScalers.cxx:60
 AliTriggerRunScalers.cxx:61
 AliTriggerRunScalers.cxx:62
 AliTriggerRunScalers.cxx:63
 AliTriggerRunScalers.cxx:64
 AliTriggerRunScalers.cxx:65
 AliTriggerRunScalers.cxx:66
 AliTriggerRunScalers.cxx:67
 AliTriggerRunScalers.cxx:68
 AliTriggerRunScalers.cxx:69
 AliTriggerRunScalers.cxx:70
 AliTriggerRunScalers.cxx:71
 AliTriggerRunScalers.cxx:72
 AliTriggerRunScalers.cxx:73
 AliTriggerRunScalers.cxx:74
 AliTriggerRunScalers.cxx:75
 AliTriggerRunScalers.cxx:76
 AliTriggerRunScalers.cxx:77
 AliTriggerRunScalers.cxx:78
 AliTriggerRunScalers.cxx:79
 AliTriggerRunScalers.cxx:80
 AliTriggerRunScalers.cxx:81
 AliTriggerRunScalers.cxx:82
 AliTriggerRunScalers.cxx:83
 AliTriggerRunScalers.cxx:84
 AliTriggerRunScalers.cxx:85
 AliTriggerRunScalers.cxx:86
 AliTriggerRunScalers.cxx:87
 AliTriggerRunScalers.cxx:88
 AliTriggerRunScalers.cxx:89
 AliTriggerRunScalers.cxx:90
 AliTriggerRunScalers.cxx:91
 AliTriggerRunScalers.cxx:92
 AliTriggerRunScalers.cxx:93
 AliTriggerRunScalers.cxx:94
 AliTriggerRunScalers.cxx:95
 AliTriggerRunScalers.cxx:96
 AliTriggerRunScalers.cxx:97
 AliTriggerRunScalers.cxx:98
 AliTriggerRunScalers.cxx:99
 AliTriggerRunScalers.cxx:100
 AliTriggerRunScalers.cxx:101
 AliTriggerRunScalers.cxx:102
 AliTriggerRunScalers.cxx:103
 AliTriggerRunScalers.cxx:104
 AliTriggerRunScalers.cxx:105
 AliTriggerRunScalers.cxx:106
 AliTriggerRunScalers.cxx:107
 AliTriggerRunScalers.cxx:108
 AliTriggerRunScalers.cxx:109
 AliTriggerRunScalers.cxx:110
 AliTriggerRunScalers.cxx:111
 AliTriggerRunScalers.cxx:112
 AliTriggerRunScalers.cxx:113
 AliTriggerRunScalers.cxx:114
 AliTriggerRunScalers.cxx:115
 AliTriggerRunScalers.cxx:116
 AliTriggerRunScalers.cxx:117
 AliTriggerRunScalers.cxx:118
 AliTriggerRunScalers.cxx:119
 AliTriggerRunScalers.cxx:120
 AliTriggerRunScalers.cxx:121
 AliTriggerRunScalers.cxx:122
 AliTriggerRunScalers.cxx:123
 AliTriggerRunScalers.cxx:124
 AliTriggerRunScalers.cxx:125
 AliTriggerRunScalers.cxx:126
 AliTriggerRunScalers.cxx:127
 AliTriggerRunScalers.cxx:128
 AliTriggerRunScalers.cxx:129
 AliTriggerRunScalers.cxx:130
 AliTriggerRunScalers.cxx:131
 AliTriggerRunScalers.cxx:132
 AliTriggerRunScalers.cxx:133
 AliTriggerRunScalers.cxx:134
 AliTriggerRunScalers.cxx:135
 AliTriggerRunScalers.cxx:136
 AliTriggerRunScalers.cxx:137
 AliTriggerRunScalers.cxx:138
 AliTriggerRunScalers.cxx:139
 AliTriggerRunScalers.cxx:140
 AliTriggerRunScalers.cxx:141
 AliTriggerRunScalers.cxx:142
 AliTriggerRunScalers.cxx:143
 AliTriggerRunScalers.cxx:144
 AliTriggerRunScalers.cxx:145
 AliTriggerRunScalers.cxx:146
 AliTriggerRunScalers.cxx:147
 AliTriggerRunScalers.cxx:148
 AliTriggerRunScalers.cxx:149
 AliTriggerRunScalers.cxx:150
 AliTriggerRunScalers.cxx:151
 AliTriggerRunScalers.cxx:152
 AliTriggerRunScalers.cxx:153
 AliTriggerRunScalers.cxx:154
 AliTriggerRunScalers.cxx:155
 AliTriggerRunScalers.cxx:156
 AliTriggerRunScalers.cxx:157
 AliTriggerRunScalers.cxx:158
 AliTriggerRunScalers.cxx:159
 AliTriggerRunScalers.cxx:160
 AliTriggerRunScalers.cxx:161
 AliTriggerRunScalers.cxx:162
 AliTriggerRunScalers.cxx:163
 AliTriggerRunScalers.cxx:164
 AliTriggerRunScalers.cxx:165
 AliTriggerRunScalers.cxx:166
 AliTriggerRunScalers.cxx:167
 AliTriggerRunScalers.cxx:168
 AliTriggerRunScalers.cxx:169
 AliTriggerRunScalers.cxx:170
 AliTriggerRunScalers.cxx:171
 AliTriggerRunScalers.cxx:172
 AliTriggerRunScalers.cxx:173
 AliTriggerRunScalers.cxx:174
 AliTriggerRunScalers.cxx:175
 AliTriggerRunScalers.cxx:176
 AliTriggerRunScalers.cxx:177
 AliTriggerRunScalers.cxx:178
 AliTriggerRunScalers.cxx:179
 AliTriggerRunScalers.cxx:180
 AliTriggerRunScalers.cxx:181
 AliTriggerRunScalers.cxx:182
 AliTriggerRunScalers.cxx:183
 AliTriggerRunScalers.cxx:184
 AliTriggerRunScalers.cxx:185
 AliTriggerRunScalers.cxx:186
 AliTriggerRunScalers.cxx:187
 AliTriggerRunScalers.cxx:188
 AliTriggerRunScalers.cxx:189
 AliTriggerRunScalers.cxx:190
 AliTriggerRunScalers.cxx:191
 AliTriggerRunScalers.cxx:192
 AliTriggerRunScalers.cxx:193
 AliTriggerRunScalers.cxx:194
 AliTriggerRunScalers.cxx:195
 AliTriggerRunScalers.cxx:196
 AliTriggerRunScalers.cxx:197
 AliTriggerRunScalers.cxx:198
 AliTriggerRunScalers.cxx:199
 AliTriggerRunScalers.cxx:200
 AliTriggerRunScalers.cxx:201
 AliTriggerRunScalers.cxx:202
 AliTriggerRunScalers.cxx:203
 AliTriggerRunScalers.cxx:204
 AliTriggerRunScalers.cxx:205
 AliTriggerRunScalers.cxx:206
 AliTriggerRunScalers.cxx:207
 AliTriggerRunScalers.cxx:208
 AliTriggerRunScalers.cxx:209
 AliTriggerRunScalers.cxx:210
 AliTriggerRunScalers.cxx:211
 AliTriggerRunScalers.cxx:212
 AliTriggerRunScalers.cxx:213
 AliTriggerRunScalers.cxx:214
 AliTriggerRunScalers.cxx:215
 AliTriggerRunScalers.cxx:216
 AliTriggerRunScalers.cxx:217
 AliTriggerRunScalers.cxx:218
 AliTriggerRunScalers.cxx:219
 AliTriggerRunScalers.cxx:220
 AliTriggerRunScalers.cxx:221
 AliTriggerRunScalers.cxx:222
 AliTriggerRunScalers.cxx:223
 AliTriggerRunScalers.cxx:224
 AliTriggerRunScalers.cxx:225
 AliTriggerRunScalers.cxx:226
 AliTriggerRunScalers.cxx:227
 AliTriggerRunScalers.cxx:228
 AliTriggerRunScalers.cxx:229
 AliTriggerRunScalers.cxx:230
 AliTriggerRunScalers.cxx:231
 AliTriggerRunScalers.cxx:232
 AliTriggerRunScalers.cxx:233
 AliTriggerRunScalers.cxx:234
 AliTriggerRunScalers.cxx:235
 AliTriggerRunScalers.cxx:236
 AliTriggerRunScalers.cxx:237
 AliTriggerRunScalers.cxx:238
 AliTriggerRunScalers.cxx:239
 AliTriggerRunScalers.cxx:240
 AliTriggerRunScalers.cxx:241
 AliTriggerRunScalers.cxx:242
 AliTriggerRunScalers.cxx:243
 AliTriggerRunScalers.cxx:244
 AliTriggerRunScalers.cxx:245
 AliTriggerRunScalers.cxx:246
 AliTriggerRunScalers.cxx:247
 AliTriggerRunScalers.cxx:248
 AliTriggerRunScalers.cxx:249
 AliTriggerRunScalers.cxx:250
 AliTriggerRunScalers.cxx:251
 AliTriggerRunScalers.cxx:252
 AliTriggerRunScalers.cxx:253
 AliTriggerRunScalers.cxx:254
 AliTriggerRunScalers.cxx:255
 AliTriggerRunScalers.cxx:256
 AliTriggerRunScalers.cxx:257
 AliTriggerRunScalers.cxx:258
 AliTriggerRunScalers.cxx:259
 AliTriggerRunScalers.cxx:260
 AliTriggerRunScalers.cxx:261
 AliTriggerRunScalers.cxx:262
 AliTriggerRunScalers.cxx:263
 AliTriggerRunScalers.cxx:264
 AliTriggerRunScalers.cxx:265
 AliTriggerRunScalers.cxx:266
 AliTriggerRunScalers.cxx:267
 AliTriggerRunScalers.cxx:268
 AliTriggerRunScalers.cxx:269
 AliTriggerRunScalers.cxx:270
 AliTriggerRunScalers.cxx:271
 AliTriggerRunScalers.cxx:272
 AliTriggerRunScalers.cxx:273
 AliTriggerRunScalers.cxx:274
 AliTriggerRunScalers.cxx:275
 AliTriggerRunScalers.cxx:276
 AliTriggerRunScalers.cxx:277
 AliTriggerRunScalers.cxx:278
 AliTriggerRunScalers.cxx:279
 AliTriggerRunScalers.cxx:280
 AliTriggerRunScalers.cxx:281
 AliTriggerRunScalers.cxx:282
 AliTriggerRunScalers.cxx:283
 AliTriggerRunScalers.cxx:284
 AliTriggerRunScalers.cxx:285
 AliTriggerRunScalers.cxx:286
 AliTriggerRunScalers.cxx:287
 AliTriggerRunScalers.cxx:288
 AliTriggerRunScalers.cxx:289
 AliTriggerRunScalers.cxx:290
 AliTriggerRunScalers.cxx:291
 AliTriggerRunScalers.cxx:292
 AliTriggerRunScalers.cxx:293
 AliTriggerRunScalers.cxx:294
 AliTriggerRunScalers.cxx:295
 AliTriggerRunScalers.cxx:296
 AliTriggerRunScalers.cxx:297
 AliTriggerRunScalers.cxx:298
 AliTriggerRunScalers.cxx:299
 AliTriggerRunScalers.cxx:300
 AliTriggerRunScalers.cxx:301
 AliTriggerRunScalers.cxx:302
 AliTriggerRunScalers.cxx:303
 AliTriggerRunScalers.cxx:304
 AliTriggerRunScalers.cxx:305
 AliTriggerRunScalers.cxx:306
 AliTriggerRunScalers.cxx:307
 AliTriggerRunScalers.cxx:308
 AliTriggerRunScalers.cxx:309
 AliTriggerRunScalers.cxx:310
 AliTriggerRunScalers.cxx:311
 AliTriggerRunScalers.cxx:312
 AliTriggerRunScalers.cxx:313
 AliTriggerRunScalers.cxx:314
 AliTriggerRunScalers.cxx:315
 AliTriggerRunScalers.cxx:316
 AliTriggerRunScalers.cxx:317
 AliTriggerRunScalers.cxx:318
 AliTriggerRunScalers.cxx:319
 AliTriggerRunScalers.cxx:320
 AliTriggerRunScalers.cxx:321
 AliTriggerRunScalers.cxx:322
 AliTriggerRunScalers.cxx:323
 AliTriggerRunScalers.cxx:324
 AliTriggerRunScalers.cxx:325
 AliTriggerRunScalers.cxx:326
 AliTriggerRunScalers.cxx:327
 AliTriggerRunScalers.cxx:328
 AliTriggerRunScalers.cxx:329
 AliTriggerRunScalers.cxx:330
 AliTriggerRunScalers.cxx:331
 AliTriggerRunScalers.cxx:332
 AliTriggerRunScalers.cxx:333
 AliTriggerRunScalers.cxx:334
 AliTriggerRunScalers.cxx:335
 AliTriggerRunScalers.cxx:336
 AliTriggerRunScalers.cxx:337
 AliTriggerRunScalers.cxx:338
 AliTriggerRunScalers.cxx:339
 AliTriggerRunScalers.cxx:340
 AliTriggerRunScalers.cxx:341
 AliTriggerRunScalers.cxx:342
 AliTriggerRunScalers.cxx:343
 AliTriggerRunScalers.cxx:344
 AliTriggerRunScalers.cxx:345
 AliTriggerRunScalers.cxx:346
 AliTriggerRunScalers.cxx:347
 AliTriggerRunScalers.cxx:348
 AliTriggerRunScalers.cxx:349
 AliTriggerRunScalers.cxx:350
 AliTriggerRunScalers.cxx:351
 AliTriggerRunScalers.cxx:352
 AliTriggerRunScalers.cxx:353
 AliTriggerRunScalers.cxx:354
 AliTriggerRunScalers.cxx:355
 AliTriggerRunScalers.cxx:356
 AliTriggerRunScalers.cxx:357
 AliTriggerRunScalers.cxx:358
 AliTriggerRunScalers.cxx:359
 AliTriggerRunScalers.cxx:360
 AliTriggerRunScalers.cxx:361
 AliTriggerRunScalers.cxx:362
 AliTriggerRunScalers.cxx:363
 AliTriggerRunScalers.cxx:364
 AliTriggerRunScalers.cxx:365
 AliTriggerRunScalers.cxx:366
 AliTriggerRunScalers.cxx:367
 AliTriggerRunScalers.cxx:368
 AliTriggerRunScalers.cxx:369
 AliTriggerRunScalers.cxx:370
 AliTriggerRunScalers.cxx:371
 AliTriggerRunScalers.cxx:372
 AliTriggerRunScalers.cxx:373
 AliTriggerRunScalers.cxx:374
 AliTriggerRunScalers.cxx:375
 AliTriggerRunScalers.cxx:376
 AliTriggerRunScalers.cxx:377
 AliTriggerRunScalers.cxx:378
 AliTriggerRunScalers.cxx:379
 AliTriggerRunScalers.cxx:380
 AliTriggerRunScalers.cxx:381
 AliTriggerRunScalers.cxx:382
 AliTriggerRunScalers.cxx:383
 AliTriggerRunScalers.cxx:384
 AliTriggerRunScalers.cxx:385
 AliTriggerRunScalers.cxx:386
 AliTriggerRunScalers.cxx:387
 AliTriggerRunScalers.cxx:388
 AliTriggerRunScalers.cxx:389
 AliTriggerRunScalers.cxx:390
 AliTriggerRunScalers.cxx:391
 AliTriggerRunScalers.cxx:392
 AliTriggerRunScalers.cxx:393
 AliTriggerRunScalers.cxx:394
 AliTriggerRunScalers.cxx:395
 AliTriggerRunScalers.cxx:396
 AliTriggerRunScalers.cxx:397
 AliTriggerRunScalers.cxx:398
 AliTriggerRunScalers.cxx:399
 AliTriggerRunScalers.cxx:400
 AliTriggerRunScalers.cxx:401
 AliTriggerRunScalers.cxx:402
 AliTriggerRunScalers.cxx:403
 AliTriggerRunScalers.cxx:404
 AliTriggerRunScalers.cxx:405
 AliTriggerRunScalers.cxx:406
 AliTriggerRunScalers.cxx:407
 AliTriggerRunScalers.cxx:408
 AliTriggerRunScalers.cxx:409
 AliTriggerRunScalers.cxx:410
 AliTriggerRunScalers.cxx:411
 AliTriggerRunScalers.cxx:412
 AliTriggerRunScalers.cxx:413
 AliTriggerRunScalers.cxx:414
 AliTriggerRunScalers.cxx:415
 AliTriggerRunScalers.cxx:416
 AliTriggerRunScalers.cxx:417
 AliTriggerRunScalers.cxx:418
 AliTriggerRunScalers.cxx:419
 AliTriggerRunScalers.cxx:420
 AliTriggerRunScalers.cxx:421
 AliTriggerRunScalers.cxx:422
 AliTriggerRunScalers.cxx:423
 AliTriggerRunScalers.cxx:424
 AliTriggerRunScalers.cxx:425
 AliTriggerRunScalers.cxx:426
 AliTriggerRunScalers.cxx:427
 AliTriggerRunScalers.cxx:428
 AliTriggerRunScalers.cxx:429
 AliTriggerRunScalers.cxx:430
 AliTriggerRunScalers.cxx:431
 AliTriggerRunScalers.cxx:432
 AliTriggerRunScalers.cxx:433
 AliTriggerRunScalers.cxx:434
 AliTriggerRunScalers.cxx:435
 AliTriggerRunScalers.cxx:436
 AliTriggerRunScalers.cxx:437
 AliTriggerRunScalers.cxx:438
 AliTriggerRunScalers.cxx:439
 AliTriggerRunScalers.cxx:440
 AliTriggerRunScalers.cxx:441
 AliTriggerRunScalers.cxx:442
 AliTriggerRunScalers.cxx:443
 AliTriggerRunScalers.cxx:444
 AliTriggerRunScalers.cxx:445
 AliTriggerRunScalers.cxx:446
 AliTriggerRunScalers.cxx:447
 AliTriggerRunScalers.cxx:448
 AliTriggerRunScalers.cxx:449
 AliTriggerRunScalers.cxx:450
 AliTriggerRunScalers.cxx:451
 AliTriggerRunScalers.cxx:452
 AliTriggerRunScalers.cxx:453
 AliTriggerRunScalers.cxx:454
 AliTriggerRunScalers.cxx:455
 AliTriggerRunScalers.cxx:456
 AliTriggerRunScalers.cxx:457
 AliTriggerRunScalers.cxx:458
 AliTriggerRunScalers.cxx:459
 AliTriggerRunScalers.cxx:460
 AliTriggerRunScalers.cxx:461
 AliTriggerRunScalers.cxx:462
 AliTriggerRunScalers.cxx:463
 AliTriggerRunScalers.cxx:464
 AliTriggerRunScalers.cxx:465
 AliTriggerRunScalers.cxx:466
 AliTriggerRunScalers.cxx:467
 AliTriggerRunScalers.cxx:468
 AliTriggerRunScalers.cxx:469
 AliTriggerRunScalers.cxx:470
 AliTriggerRunScalers.cxx:471
 AliTriggerRunScalers.cxx:472
 AliTriggerRunScalers.cxx:473
 AliTriggerRunScalers.cxx:474
 AliTriggerRunScalers.cxx:475
 AliTriggerRunScalers.cxx:476
 AliTriggerRunScalers.cxx:477
 AliTriggerRunScalers.cxx:478
 AliTriggerRunScalers.cxx:479
 AliTriggerRunScalers.cxx:480
 AliTriggerRunScalers.cxx:481
 AliTriggerRunScalers.cxx:482
 AliTriggerRunScalers.cxx:483
 AliTriggerRunScalers.cxx:484
 AliTriggerRunScalers.cxx:485
 AliTriggerRunScalers.cxx:486
 AliTriggerRunScalers.cxx:487
 AliTriggerRunScalers.cxx:488
 AliTriggerRunScalers.cxx:489
 AliTriggerRunScalers.cxx:490
 AliTriggerRunScalers.cxx:491
 AliTriggerRunScalers.cxx:492
 AliTriggerRunScalers.cxx:493
 AliTriggerRunScalers.cxx:494
 AliTriggerRunScalers.cxx:495
 AliTriggerRunScalers.cxx:496
 AliTriggerRunScalers.cxx:497
 AliTriggerRunScalers.cxx:498
 AliTriggerRunScalers.cxx:499
 AliTriggerRunScalers.cxx:500
 AliTriggerRunScalers.cxx:501
 AliTriggerRunScalers.cxx:502
 AliTriggerRunScalers.cxx:503
 AliTriggerRunScalers.cxx:504
 AliTriggerRunScalers.cxx:505
 AliTriggerRunScalers.cxx:506
 AliTriggerRunScalers.cxx:507
 AliTriggerRunScalers.cxx:508
 AliTriggerRunScalers.cxx:509
 AliTriggerRunScalers.cxx:510
 AliTriggerRunScalers.cxx:511
 AliTriggerRunScalers.cxx:512
 AliTriggerRunScalers.cxx:513
 AliTriggerRunScalers.cxx:514
 AliTriggerRunScalers.cxx:515
 AliTriggerRunScalers.cxx:516
 AliTriggerRunScalers.cxx:517
 AliTriggerRunScalers.cxx:518
 AliTriggerRunScalers.cxx:519
 AliTriggerRunScalers.cxx:520
 AliTriggerRunScalers.cxx:521
 AliTriggerRunScalers.cxx:522
 AliTriggerRunScalers.cxx:523
 AliTriggerRunScalers.cxx:524
 AliTriggerRunScalers.cxx:525
 AliTriggerRunScalers.cxx:526
 AliTriggerRunScalers.cxx:527
 AliTriggerRunScalers.cxx:528
 AliTriggerRunScalers.cxx:529
 AliTriggerRunScalers.cxx:530
 AliTriggerRunScalers.cxx:531
 AliTriggerRunScalers.cxx:532
 AliTriggerRunScalers.cxx:533
 AliTriggerRunScalers.cxx:534
 AliTriggerRunScalers.cxx:535
 AliTriggerRunScalers.cxx:536
 AliTriggerRunScalers.cxx:537
 AliTriggerRunScalers.cxx:538
 AliTriggerRunScalers.cxx:539
 AliTriggerRunScalers.cxx:540
 AliTriggerRunScalers.cxx:541
 AliTriggerRunScalers.cxx:542
 AliTriggerRunScalers.cxx:543
 AliTriggerRunScalers.cxx:544
 AliTriggerRunScalers.cxx:545
 AliTriggerRunScalers.cxx:546
 AliTriggerRunScalers.cxx:547
 AliTriggerRunScalers.cxx:548
 AliTriggerRunScalers.cxx:549
 AliTriggerRunScalers.cxx:550
 AliTriggerRunScalers.cxx:551
 AliTriggerRunScalers.cxx:552
 AliTriggerRunScalers.cxx:553
 AliTriggerRunScalers.cxx:554
 AliTriggerRunScalers.cxx:555
 AliTriggerRunScalers.cxx:556
 AliTriggerRunScalers.cxx:557
 AliTriggerRunScalers.cxx:558
 AliTriggerRunScalers.cxx:559
 AliTriggerRunScalers.cxx:560
 AliTriggerRunScalers.cxx:561
 AliTriggerRunScalers.cxx:562
 AliTriggerRunScalers.cxx:563
 AliTriggerRunScalers.cxx:564
 AliTriggerRunScalers.cxx:565
 AliTriggerRunScalers.cxx:566
 AliTriggerRunScalers.cxx:567
 AliTriggerRunScalers.cxx:568
 AliTriggerRunScalers.cxx:569
 AliTriggerRunScalers.cxx:570
 AliTriggerRunScalers.cxx:571
 AliTriggerRunScalers.cxx:572
 AliTriggerRunScalers.cxx:573
 AliTriggerRunScalers.cxx:574
 AliTriggerRunScalers.cxx:575
 AliTriggerRunScalers.cxx:576
 AliTriggerRunScalers.cxx:577
 AliTriggerRunScalers.cxx:578
 AliTriggerRunScalers.cxx:579
 AliTriggerRunScalers.cxx:580
 AliTriggerRunScalers.cxx:581
 AliTriggerRunScalers.cxx:582
 AliTriggerRunScalers.cxx:583
 AliTriggerRunScalers.cxx:584
 AliTriggerRunScalers.cxx:585
 AliTriggerRunScalers.cxx:586
 AliTriggerRunScalers.cxx:587
 AliTriggerRunScalers.cxx:588
 AliTriggerRunScalers.cxx:589
 AliTriggerRunScalers.cxx:590
 AliTriggerRunScalers.cxx:591
 AliTriggerRunScalers.cxx:592
 AliTriggerRunScalers.cxx:593
 AliTriggerRunScalers.cxx:594
 AliTriggerRunScalers.cxx:595
 AliTriggerRunScalers.cxx:596
 AliTriggerRunScalers.cxx:597
 AliTriggerRunScalers.cxx:598
 AliTriggerRunScalers.cxx:599
 AliTriggerRunScalers.cxx:600
 AliTriggerRunScalers.cxx:601
 AliTriggerRunScalers.cxx:602
 AliTriggerRunScalers.cxx:603
 AliTriggerRunScalers.cxx:604
 AliTriggerRunScalers.cxx:605
 AliTriggerRunScalers.cxx:606
 AliTriggerRunScalers.cxx:607
 AliTriggerRunScalers.cxx:608
 AliTriggerRunScalers.cxx:609
 AliTriggerRunScalers.cxx:610
 AliTriggerRunScalers.cxx:611
 AliTriggerRunScalers.cxx:612
 AliTriggerRunScalers.cxx:613
 AliTriggerRunScalers.cxx:614
 AliTriggerRunScalers.cxx:615
 AliTriggerRunScalers.cxx:616
 AliTriggerRunScalers.cxx:617
 AliTriggerRunScalers.cxx:618
 AliTriggerRunScalers.cxx:619
 AliTriggerRunScalers.cxx:620
 AliTriggerRunScalers.cxx:621
 AliTriggerRunScalers.cxx:622
 AliTriggerRunScalers.cxx:623
 AliTriggerRunScalers.cxx:624
 AliTriggerRunScalers.cxx:625
 AliTriggerRunScalers.cxx:626
 AliTriggerRunScalers.cxx:627
 AliTriggerRunScalers.cxx:628
 AliTriggerRunScalers.cxx:629
 AliTriggerRunScalers.cxx:630
 AliTriggerRunScalers.cxx:631
 AliTriggerRunScalers.cxx:632
 AliTriggerRunScalers.cxx:633
 AliTriggerRunScalers.cxx:634
 AliTriggerRunScalers.cxx:635
 AliTriggerRunScalers.cxx:636
 AliTriggerRunScalers.cxx:637
 AliTriggerRunScalers.cxx:638
 AliTriggerRunScalers.cxx:639
 AliTriggerRunScalers.cxx:640
 AliTriggerRunScalers.cxx:641
 AliTriggerRunScalers.cxx:642
 AliTriggerRunScalers.cxx:643
 AliTriggerRunScalers.cxx:644
 AliTriggerRunScalers.cxx:645
 AliTriggerRunScalers.cxx:646
 AliTriggerRunScalers.cxx:647
 AliTriggerRunScalers.cxx:648
 AliTriggerRunScalers.cxx:649
 AliTriggerRunScalers.cxx:650
 AliTriggerRunScalers.cxx:651
 AliTriggerRunScalers.cxx:652
 AliTriggerRunScalers.cxx:653
 AliTriggerRunScalers.cxx:654
 AliTriggerRunScalers.cxx:655
 AliTriggerRunScalers.cxx:656
 AliTriggerRunScalers.cxx:657
 AliTriggerRunScalers.cxx:658
 AliTriggerRunScalers.cxx:659
 AliTriggerRunScalers.cxx:660
 AliTriggerRunScalers.cxx:661
 AliTriggerRunScalers.cxx:662
 AliTriggerRunScalers.cxx:663
 AliTriggerRunScalers.cxx:664
 AliTriggerRunScalers.cxx:665
 AliTriggerRunScalers.cxx:666
 AliTriggerRunScalers.cxx:667
 AliTriggerRunScalers.cxx:668
 AliTriggerRunScalers.cxx:669
 AliTriggerRunScalers.cxx:670
 AliTriggerRunScalers.cxx:671
 AliTriggerRunScalers.cxx:672
 AliTriggerRunScalers.cxx:673
 AliTriggerRunScalers.cxx:674
 AliTriggerRunScalers.cxx:675
 AliTriggerRunScalers.cxx:676
 AliTriggerRunScalers.cxx:677
 AliTriggerRunScalers.cxx:678
 AliTriggerRunScalers.cxx:679
 AliTriggerRunScalers.cxx:680
 AliTriggerRunScalers.cxx:681
 AliTriggerRunScalers.cxx:682
 AliTriggerRunScalers.cxx:683
 AliTriggerRunScalers.cxx:684
 AliTriggerRunScalers.cxx:685
 AliTriggerRunScalers.cxx:686
 AliTriggerRunScalers.cxx:687
 AliTriggerRunScalers.cxx:688
 AliTriggerRunScalers.cxx:689
 AliTriggerRunScalers.cxx:690
 AliTriggerRunScalers.cxx:691
 AliTriggerRunScalers.cxx:692
 AliTriggerRunScalers.cxx:693
 AliTriggerRunScalers.cxx:694
 AliTriggerRunScalers.cxx:695
 AliTriggerRunScalers.cxx:696
 AliTriggerRunScalers.cxx:697
 AliTriggerRunScalers.cxx:698
 AliTriggerRunScalers.cxx:699
 AliTriggerRunScalers.cxx:700
 AliTriggerRunScalers.cxx:701
 AliTriggerRunScalers.cxx:702
 AliTriggerRunScalers.cxx:703
 AliTriggerRunScalers.cxx:704
 AliTriggerRunScalers.cxx:705
 AliTriggerRunScalers.cxx:706
 AliTriggerRunScalers.cxx:707
 AliTriggerRunScalers.cxx:708
 AliTriggerRunScalers.cxx:709
 AliTriggerRunScalers.cxx:710
 AliTriggerRunScalers.cxx:711
 AliTriggerRunScalers.cxx:712
 AliTriggerRunScalers.cxx:713
 AliTriggerRunScalers.cxx:714
 AliTriggerRunScalers.cxx:715
 AliTriggerRunScalers.cxx:716
 AliTriggerRunScalers.cxx:717
 AliTriggerRunScalers.cxx:718
 AliTriggerRunScalers.cxx:719
 AliTriggerRunScalers.cxx:720
 AliTriggerRunScalers.cxx:721
 AliTriggerRunScalers.cxx:722
 AliTriggerRunScalers.cxx:723
 AliTriggerRunScalers.cxx:724
 AliTriggerRunScalers.cxx:725
 AliTriggerRunScalers.cxx:726
 AliTriggerRunScalers.cxx:727
 AliTriggerRunScalers.cxx:728
 AliTriggerRunScalers.cxx:729
 AliTriggerRunScalers.cxx:730
 AliTriggerRunScalers.cxx:731
 AliTriggerRunScalers.cxx:732
 AliTriggerRunScalers.cxx:733
 AliTriggerRunScalers.cxx:734
 AliTriggerRunScalers.cxx:735
 AliTriggerRunScalers.cxx:736
 AliTriggerRunScalers.cxx:737
 AliTriggerRunScalers.cxx:738
 AliTriggerRunScalers.cxx:739
 AliTriggerRunScalers.cxx:740
 AliTriggerRunScalers.cxx:741
 AliTriggerRunScalers.cxx:742
 AliTriggerRunScalers.cxx:743
 AliTriggerRunScalers.cxx:744
 AliTriggerRunScalers.cxx:745
 AliTriggerRunScalers.cxx:746
 AliTriggerRunScalers.cxx:747
 AliTriggerRunScalers.cxx:748
 AliTriggerRunScalers.cxx:749
 AliTriggerRunScalers.cxx:750
 AliTriggerRunScalers.cxx:751
 AliTriggerRunScalers.cxx:752
 AliTriggerRunScalers.cxx:753
 AliTriggerRunScalers.cxx:754
 AliTriggerRunScalers.cxx:755
 AliTriggerRunScalers.cxx:756
 AliTriggerRunScalers.cxx:757
 AliTriggerRunScalers.cxx:758
 AliTriggerRunScalers.cxx:759
 AliTriggerRunScalers.cxx:760
 AliTriggerRunScalers.cxx:761
 AliTriggerRunScalers.cxx:762
 AliTriggerRunScalers.cxx:763
 AliTriggerRunScalers.cxx:764
 AliTriggerRunScalers.cxx:765
 AliTriggerRunScalers.cxx:766
 AliTriggerRunScalers.cxx:767
 AliTriggerRunScalers.cxx:768
 AliTriggerRunScalers.cxx:769
 AliTriggerRunScalers.cxx:770
 AliTriggerRunScalers.cxx:771
 AliTriggerRunScalers.cxx:772
 AliTriggerRunScalers.cxx:773
 AliTriggerRunScalers.cxx:774
 AliTriggerRunScalers.cxx:775
 AliTriggerRunScalers.cxx:776
 AliTriggerRunScalers.cxx:777
 AliTriggerRunScalers.cxx:778
 AliTriggerRunScalers.cxx:779
 AliTriggerRunScalers.cxx:780
 AliTriggerRunScalers.cxx:781
 AliTriggerRunScalers.cxx:782
 AliTriggerRunScalers.cxx:783
 AliTriggerRunScalers.cxx:784
 AliTriggerRunScalers.cxx:785
 AliTriggerRunScalers.cxx:786
 AliTriggerRunScalers.cxx:787
 AliTriggerRunScalers.cxx:788
 AliTriggerRunScalers.cxx:789
 AliTriggerRunScalers.cxx:790
 AliTriggerRunScalers.cxx:791
 AliTriggerRunScalers.cxx:792
 AliTriggerRunScalers.cxx:793
 AliTriggerRunScalers.cxx:794
 AliTriggerRunScalers.cxx:795
 AliTriggerRunScalers.cxx:796
 AliTriggerRunScalers.cxx:797
 AliTriggerRunScalers.cxx:798
 AliTriggerRunScalers.cxx:799
 AliTriggerRunScalers.cxx:800
 AliTriggerRunScalers.cxx:801
 AliTriggerRunScalers.cxx:802
 AliTriggerRunScalers.cxx:803
 AliTriggerRunScalers.cxx:804
 AliTriggerRunScalers.cxx:805
 AliTriggerRunScalers.cxx:806
 AliTriggerRunScalers.cxx:807
 AliTriggerRunScalers.cxx:808
 AliTriggerRunScalers.cxx:809
 AliTriggerRunScalers.cxx:810
 AliTriggerRunScalers.cxx:811
 AliTriggerRunScalers.cxx:812
 AliTriggerRunScalers.cxx:813
 AliTriggerRunScalers.cxx:814
 AliTriggerRunScalers.cxx:815
 AliTriggerRunScalers.cxx:816
 AliTriggerRunScalers.cxx:817
 AliTriggerRunScalers.cxx:818
 AliTriggerRunScalers.cxx:819
 AliTriggerRunScalers.cxx:820
 AliTriggerRunScalers.cxx:821
 AliTriggerRunScalers.cxx:822
 AliTriggerRunScalers.cxx:823
 AliTriggerRunScalers.cxx:824
 AliTriggerRunScalers.cxx:825
 AliTriggerRunScalers.cxx:826
 AliTriggerRunScalers.cxx:827
 AliTriggerRunScalers.cxx:828
 AliTriggerRunScalers.cxx:829
 AliTriggerRunScalers.cxx:830
 AliTriggerRunScalers.cxx:831
 AliTriggerRunScalers.cxx:832
 AliTriggerRunScalers.cxx:833
 AliTriggerRunScalers.cxx:834
 AliTriggerRunScalers.cxx:835
 AliTriggerRunScalers.cxx:836
 AliTriggerRunScalers.cxx:837
 AliTriggerRunScalers.cxx:838
 AliTriggerRunScalers.cxx:839
 AliTriggerRunScalers.cxx:840
 AliTriggerRunScalers.cxx:841
 AliTriggerRunScalers.cxx:842
 AliTriggerRunScalers.cxx:843
 AliTriggerRunScalers.cxx:844
 AliTriggerRunScalers.cxx:845
 AliTriggerRunScalers.cxx:846
 AliTriggerRunScalers.cxx:847
 AliTriggerRunScalers.cxx:848
 AliTriggerRunScalers.cxx:849
 AliTriggerRunScalers.cxx:850
 AliTriggerRunScalers.cxx:851
 AliTriggerRunScalers.cxx:852
 AliTriggerRunScalers.cxx:853
 AliTriggerRunScalers.cxx:854
 AliTriggerRunScalers.cxx:855
 AliTriggerRunScalers.cxx:856
 AliTriggerRunScalers.cxx:857
 AliTriggerRunScalers.cxx:858
 AliTriggerRunScalers.cxx:859
 AliTriggerRunScalers.cxx:860
 AliTriggerRunScalers.cxx:861
 AliTriggerRunScalers.cxx:862
 AliTriggerRunScalers.cxx:863
 AliTriggerRunScalers.cxx:864
 AliTriggerRunScalers.cxx:865
 AliTriggerRunScalers.cxx:866
 AliTriggerRunScalers.cxx:867
 AliTriggerRunScalers.cxx:868
 AliTriggerRunScalers.cxx:869
 AliTriggerRunScalers.cxx:870
 AliTriggerRunScalers.cxx:871
 AliTriggerRunScalers.cxx:872
 AliTriggerRunScalers.cxx:873
 AliTriggerRunScalers.cxx:874
 AliTriggerRunScalers.cxx:875
 AliTriggerRunScalers.cxx:876
 AliTriggerRunScalers.cxx:877
 AliTriggerRunScalers.cxx:878
 AliTriggerRunScalers.cxx:879
 AliTriggerRunScalers.cxx:880
 AliTriggerRunScalers.cxx:881
 AliTriggerRunScalers.cxx:882
 AliTriggerRunScalers.cxx:883
 AliTriggerRunScalers.cxx:884
 AliTriggerRunScalers.cxx:885
 AliTriggerRunScalers.cxx:886
 AliTriggerRunScalers.cxx:887
 AliTriggerRunScalers.cxx:888
 AliTriggerRunScalers.cxx:889
 AliTriggerRunScalers.cxx:890
 AliTriggerRunScalers.cxx:891
 AliTriggerRunScalers.cxx:892
 AliTriggerRunScalers.cxx:893
 AliTriggerRunScalers.cxx:894
 AliTriggerRunScalers.cxx:895
 AliTriggerRunScalers.cxx:896
 AliTriggerRunScalers.cxx:897
 AliTriggerRunScalers.cxx:898
 AliTriggerRunScalers.cxx:899
 AliTriggerRunScalers.cxx:900
 AliTriggerRunScalers.cxx:901
 AliTriggerRunScalers.cxx:902
 AliTriggerRunScalers.cxx:903
 AliTriggerRunScalers.cxx:904
 AliTriggerRunScalers.cxx:905
 AliTriggerRunScalers.cxx:906
 AliTriggerRunScalers.cxx:907
 AliTriggerRunScalers.cxx:908
 AliTriggerRunScalers.cxx:909
 AliTriggerRunScalers.cxx:910
 AliTriggerRunScalers.cxx:911
 AliTriggerRunScalers.cxx:912
 AliTriggerRunScalers.cxx:913
 AliTriggerRunScalers.cxx:914
 AliTriggerRunScalers.cxx:915
 AliTriggerRunScalers.cxx:916
 AliTriggerRunScalers.cxx:917
 AliTriggerRunScalers.cxx:918
 AliTriggerRunScalers.cxx:919
 AliTriggerRunScalers.cxx:920
 AliTriggerRunScalers.cxx:921
 AliTriggerRunScalers.cxx:922
 AliTriggerRunScalers.cxx:923
 AliTriggerRunScalers.cxx:924