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$ */

//_________________________________________________________________________
//
// This is a TTask that made the calculation of the Time zero using TOF.
// Description: The algorithm used to calculate the time zero of
// interaction using TOF detector is the following.
// We select in the MonteCarlo some primary particles - or tracks in
// the following - that strike the TOF detector (the larger part are
// pions, kaons or protons).
// We choose a set of 10 selected tracks, for each track You have the
// length of the track when the TOF is reached (a standard TOF hit
// does not contain this additional information, this is the reason
// why we implemented a new time zero dedicated TOF hit class
// AliTOFhitT0; in order to store this type of hit You have to use the
// AliTOFv4T0 as TOF class in Your Config.C. In AliTOFv4T0 the
// StepManager was modified in order to fill the TOF hit branch with
// this type of hits; in fact the AliTOF::AddT0Hit is called rather
// that the usual AliTOF::AddHit), the momentum at generation (from
// TreeK) and the time of flight given by the TOF detector.
// (Observe that the ctor of the AliTOF class, when the AliTOFv4T0
// class is used, is called with the "tzero" option: it is in order
// create the fHits TClonesArray filled with AliTOFhitT0 objects,
// rather than with normal AliTOFhit)
// Then Momentum and time of flight for each track are smeared
// according to known experimental resolution (all sources of error
// have been token into account).
// Let consider now only one set of 10 tracks (the algorithm is the
// same for all sets).
// Assuming the (mass) hypothesis that each track can be AUT a pion,
// AUT a kaon, AUT a proton, we consider all the 3 at 10 possible
// cases.
// For each track in each (mass) configuration
// (a configuration can be
// e.g. pion/pion/kaon/proton/pion/proton/kaon/kaon/pion/pion)
// we calculate the time zero (we know in fact the velocity of the
// track after the assumption about its mass, the time of flight given
// by the TOF, and the corresponding path travelled till the TOF
// detector). Then for each mass configuration we have 10 time zero
// and we can calculate the ChiSquare for the current configuration
// using the weighted mean over all 10 time zero.
// We call the best assignment the mass configuration that gives the
// minimum value of the ChiSquare.
// We plot the weighted mean over all 10 time zero for the best
// assignment, the ChiSquare for the best assignment and the
// corresponding confidence level.
// The strong assumption is the MC selection of primary particles. It
// will be introduced in the future also some more realistic
// simulation about this point.
//
// Use case:
// root [0] AliTOFT0 * tzero = new AliTOFT0("galice.root")
// Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
// root [1] tzero->ExecuteTask()
// root [2] tzero->ExecuteTask("tim")
//             // available parameters:
//             tim - print benchmarking information
//             all - print usefull informations about the number of
//                   misidentified tracks and a comparison about the
//                   true configuration (known from MC) and the best
//                   assignment
//
//-- Author: F. Pierella
//
//_________________________________________________________________________

#include <TCanvas.h>
#include <TClonesArray.h>
#include <TFile.h>
//#include <TFolder.h>
#include <TFrame.h>
#include <TH1.h>
#include <TParticle.h>
#include <TBenchmark.h>
#include <TTask.h>
#include <TTree.h>
#include <TRandom.h>
#include <TROOT.h>

#include "AliLoader.h"
#include "AliMC.h"
#include "AliRun.h"

#include "AliTOFhitT0.h"
#include "AliTOFT0.h"
#include "AliTOF.h"

//extern TROOT *gROOT;
extern TRandom *gRandom;
extern TBenchmark *gBenchmark;

extern AliRun *gAlice;

ClassImp(AliTOFT0)

//____________________________________________________________________________ 
AliTOFT0::AliTOFT0():
  TTask("AliTOFT0",""),
  fNevents(0),
  fTimeResolution(0),
  fLowerMomBound(0),
  fUpperMomBound(0),
  fT0File(""),
  fHeadersFile("")
{
  // ctor
}
           
//____________________________________________________________________________ 
AliTOFT0::AliTOFT0(char* headerFile, Int_t nEvents):
  TTask("AliTOFT0",""), 
  fNevents(nEvents),
  fTimeResolution(1.2e-10),
  fLowerMomBound(1.5),
  fUpperMomBound(2.),
  fT0File(""),
  fHeadersFile(headerFile)
{
  //
  //
  //

  //fNevents=nEvents ; // Number of events for which calculate the T0, 
                     // default 0: it means all evens in current file
  //fLowerMomBound=1.5; // [GeV/c] default value
  //fUpperMomBound=2. ; // [GeV/c] default value
  //fTimeResolution   = 1.2e-10; // 120 ps by default	
  //fHeadersFile = headerFile ;

  TFile * file = (TFile*) gROOT->GetFile(fHeadersFile.Data() ) ;

  //File was not opened yet
  if(file == 0){
    if(fHeadersFile.Contains("rfio"))
      file =	TFile::Open(fHeadersFile,"update") ;
    else
      file = new TFile(fHeadersFile.Data(),"update") ;
    gAlice = (AliRun *) file->Get("gAlice") ;
  }

  // add Task to //root/Tasks folder
  TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
  roottasks->Add(this) ; 
}

//____________________________________________________________________________ 
AliTOFT0::AliTOFT0(const AliTOFT0 & tzero):
  TTask("AliTOFT0",""),
  fNevents(0),
  fTimeResolution(0),
  fLowerMomBound(0),
  fUpperMomBound(0),
  fT0File(""),
  fHeadersFile("")
{
  // copy ctr

( (AliTOFT0 &)tzero ).Copy(*this);
}

//____________________________________________________________________________ 
AliTOFT0& AliTOFT0::operator = (const AliTOFT0 & tzero)
{
  if (this==&tzero) return *this;

  fNevents=tzero.fNevents;
  fTimeResolution=tzero.fTimeResolution;
  fLowerMomBound=tzero.fLowerMomBound;
  fUpperMomBound=tzero.fUpperMomBound;
  fT0File=tzero.fT0File;
  fHeadersFile=tzero.fHeadersFile;
  return *this;

}

//____________________________________________________________________________ 
  AliTOFT0::~AliTOFT0()
{
  // dtor
}

//____________________________________________________________________________
void AliTOFT0::Exec(Option_t *option) 
{ 
  //
  // calculate T0 distribution for all events using chisquare 
  //
  Int_t ngood=0;
  Int_t nmisidentified=0;
  Int_t nmisidentified0=0;
  Int_t nmisidentified1=0;
  Int_t nmisidentified2=0;
  Int_t nmisidentified3=0;
  Int_t nmisidentified4=0;
  Int_t nmisidentified5=0;
  Int_t nmisidentified6=0;
  Int_t nmisidentified7=0;
  Int_t nmisidentified8=0;
  Int_t nmisidentified9=0;
  Int_t ipartold = -1;
  Int_t ipart;
  Int_t selected=0;
  Int_t istop=0;
  Float_t timeresolutioninns=fTimeResolution*(1.e+9); // convert in [ns]
  const Int_t kUPDATE = 5; // for visual option
  Int_t itimes=0;
  TCanvas* c1=0;
  TCanvas* c2=0;
  TCanvas* c3=0;

  if(strstr(option,"visual")){
    // Create a new canvas.
    //c1 = new TCanvas("c1","Dynamic Visual Filling of time zero histo",10,10,500,500);
    c1 = new TCanvas("c1","Dynamic Visual Filling of time zero histo",10,10,370,370);
    c1->SetFillColor(35);
    c1->GetFrame()->SetFillColor(21);
    c1->GetFrame()->SetBorderSize(6);
    c1->GetFrame()->SetBorderMode(-1);

    //c2 = new TCanvas("c2","Dynamic Visual Filling of chisquare histo",550,10,500,500);
    c2 = new TCanvas("c2","Dynamic Visual Filling of chisquare histo",380,10,370,370);
    c2->SetFillColor(35);
    c2->GetFrame()->SetFillColor(21);
    c2->GetFrame()->SetBorderSize(6);
    c2->GetFrame()->SetBorderMode(-1);

    //c3 = new TCanvas("c3","Dynamic Visual Filling of confidence level histo",280,550,500,500);
    c3 = new TCanvas("c3","Dynamic Visual Filling of confidence level histo",760,10,370,370);
    c3->SetFillColor(35);
    c3->GetFrame()->SetFillColor(21);
    c3->GetFrame()->SetBorderSize(6);
    c3->GetFrame()->SetBorderMode(-1);
  }

  if(strstr(option,"tim") || strstr(option,"all"))
    gBenchmark->Start("TOFT0");

  TH1F *htzerobest= new TH1F("htzerobest","T0 for best assignment",200,-1.,1.);
  TH1F* hchibest  = new TH1F("hchibest","ChiSquare Min Distribution",80,0.,40.);
  TH1F* hchibestconflevel  = new TH1F("hchibestconflevel","ChiSquare Min Confidence Level",10,0.,1.);

  // setting histo colors
  if(strstr(option,"visual")){
    htzerobest->SetFillColor(48);
    hchibest->SetFillColor(50);
    hchibestconflevel->SetFillColor(52);
  }

  Int_t   assparticle[10]={3,3,3,3,3,3,3,3,3,3};
  Int_t   truparticle[10]={3,3,3,3,3,3,3,3,3,3};
  Float_t t0best=999.;
  Float_t timeofflight[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
  Float_t momentum[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
  Float_t timezero[10];
  Float_t weightedtimezero[10];
  Float_t beta[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
  Float_t sqMomError[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
  Float_t sqTrackError[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
  Float_t massarray[3]={0.13957,0.493677,0.9382723};
  Float_t dummychisquare=0.;
  Float_t chisquare=999.;
  Float_t tracktoflen[10]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};

  AliTOF *detTOF = (AliTOF *) gAlice->GetDetector ("TOF");

  if (!detTOF) {
    Error("AliTOFT0","TOF not found");
    return;
  }

  if(strstr(option,"all")){
    AliInfo(Form("Selecting primary tracks with momentum between %f GeV/c and %f GeV/c",  fLowerMomBound, fUpperMomBound));
    AliInfo("Memorandum: 0 means PION | 1 means KAON | 2 means PROTON");
  }

  if (fNevents == 0) fNevents = (Int_t) AliRunLoader::Instance()->TreeE()->GetEntries();

  for (Int_t ievent = 0; ievent < fNevents; ievent++) {
    gAlice->GetEvent(ievent);
    TTree *hitTree = detTOF->GetLoader()->TreeH ();
    if (!hitTree)
      return;
    TParticle*    particle;
    AliTOFhitT0*  tofHit;
    TClonesArray* tofHits = detTOF->Hits();

    Int_t lasttrack=-1;
    Int_t nset=0;

    hitTree->SetBranchStatus("*",0); // switch off all branches
    hitTree->SetBranchStatus("TOF*",1); // switch on only TOF

    // Start loop on primary tracks in the hits containers

    Int_t ntracks = static_cast<Int_t>(hitTree->GetEntries());
    for (Int_t track = 0; track < ntracks; track++)
    {
      if(nset>=5) break; // check on the number of set analyzed
      
      gAlice->GetMCApp()->ResetHits();
      hitTree->GetEvent(track);

      //AliMC *mcApplication = (AliMC*)gAlice->GetMCApp();

      //particle = mcApplication->Particle(track);
      Int_t nhits = tofHits->GetEntriesFast();

      for (Int_t hit = 0; hit < nhits; hit++)
      {
	tofHit = (AliTOFhitT0 *) tofHits->UncheckedAt(hit);
	ipart    = tofHit->GetTrack();
	// check to discard the case when the same particle is selected more than one
	// time 

	if (ipart != ipartold){
	  
	  particle = (TParticle*)gAlice->GetMCApp()->Particle(ipart);
	  
	  Float_t idealtime=tofHit->GetTof();
	  //	   Float_t time=idealtime;
	  Float_t time   = gRandom->Gaus(idealtime, fTimeResolution);
	  Float_t toflen=tofHit->GetLen();
	  toflen=toflen/100.; // toflen given in m
	  Int_t pdg   = particle->GetPdgCode();
	  Int_t abspdg   =TMath::Abs(pdg);
	  Float_t idealmom  = particle->P();
	  Float_t momres=idealmom*0.025; // 2.5% res token into account for all momenta
	  Float_t mom =gRandom->Gaus(idealmom,momres);

	  Bool_t isgoodpart=(abspdg==211 || abspdg==2212 || abspdg==321);

	  time*=1.E+9; // tof given in nanoseconds	   
	  if (particle->GetFirstMother() < 0 && isgoodpart && mom<=fUpperMomBound && mom>=fLowerMomBound){
	    selected+=1;
	    istop=selected;
	    if(istop>10) break;
	    Int_t index=selected-1;
	    timeofflight[index]=time;
	    tracktoflen[index]=toflen;
	    momentum[index]=mom;
	    //AliInfo(Form(" %d  %d  %d ", timeofflight[index], tracktoflen[index], momentum[index]));
	    switch (abspdg) {
	    case 211:
	      truparticle[index]=0;
	      break ;
	    case 321:
	      truparticle[index]=1;
	      break ;
	    case 2212: 
	      truparticle[index]=2;
	      break ;
	    }
	    
	  }
	  ipartold = ipart;
	  
	  if(istop==10){ // start analysis on current set
	    nset+=1;
	    lasttrack=track;
	    istop=0;
	    selected=0;
	    //AliInfo("starting t0 calculation for current set");
	    for (Int_t i1=0; i1<3;i1++) {
	      beta[0]=momentum[0]/sqrt(massarray[i1]*massarray[i1]+momentum[0]*momentum[0]);
	      for (Int_t i2=0; i2<3;i2++) { 
		beta[1]=momentum[1]/sqrt(massarray[i2]*massarray[i2]+momentum[1]*momentum[1]);
		for (Int_t i3=0; i3<3;i3++) {
		  beta[2]=momentum[2]/sqrt(massarray[i3]*massarray[i3]+momentum[2]*momentum[2]);
		  for (Int_t i4=0; i4<3;i4++) {
		    beta[3]=momentum[3]/sqrt(massarray[i4]*massarray[i4]+momentum[3]*momentum[3]);
		    for (Int_t i5=0; i5<3;i5++) {
		      beta[4]=momentum[4]/sqrt(massarray[i5]*massarray[i5]+momentum[4]*momentum[4]);
		      for (Int_t i6=0; i6<3;i6++) {
			beta[5]=momentum[5]/sqrt(massarray[i6]*massarray[i6]+momentum[5]*momentum[5]);
			for (Int_t i7=0; i7<3;i7++) { 
			  beta[6]=momentum[6]/sqrt(massarray[i7]*massarray[i7]+momentum[6]*momentum[6]);
			  for (Int_t i8=0; i8<3;i8++) {
			    beta[7]=momentum[7]/sqrt(massarray[i8]*massarray[i8]+momentum[7]*momentum[7]);
			    for (Int_t i9=0; i9<3;i9++) {
			      beta[8]=momentum[8]/sqrt(massarray[i9]*massarray[i9]+momentum[8]*momentum[8]);
			      for (Int_t i10=0; i10<3;i10++) { 	
				beta[9]=momentum[9]/sqrt(massarray[i10]*massarray[i10]+momentum[9]*momentum[9]);
				
				Float_t meantzero=0.;
				Float_t sumAllweights=0.;
				for (Int_t itz=0; itz<10;itz++) {
				  sqMomError[itz]=((1.-beta[itz]*beta[itz])*0.025)*((1.-beta[itz]*beta[itz])*0.025)*(tracktoflen[itz]/(0.299792*beta[itz]))*(tracktoflen[itz]/(0.299792*beta[itz])); // this gives the square of the momentum error in nanoseconds
				  sqTrackError[itz]=(timeresolutioninns*timeresolutioninns+sqMomError[itz]); // total error for the current track
				  sumAllweights+=1./sqTrackError[itz];

				  timezero[itz]=(tracktoflen[itz]/(beta[itz]*0.299792))-timeofflight[itz];
				  weightedtimezero[itz]=((tracktoflen[itz]/(beta[itz]*0.299792))-timeofflight[itz])/sqTrackError[itz];// weighted time zero for current track
				  meantzero+=weightedtimezero[itz];
				} // end loop for (Int_t itz=0; itz<10;itz++)
				meantzero=meantzero/sumAllweights; // it is given in [ns]
				
				dummychisquare=0.;
				// calculate the chisquare for the current assignment
				for (Int_t icsq=0; icsq<10;icsq++) {
				  dummychisquare+=(timezero[icsq]-meantzero)*(timezero[icsq]-meantzero)/sqTrackError[icsq];
				} // end loop for (Int_t icsq=0; icsq<10;icsq++) 

				if(dummychisquare<=chisquare){
				  assparticle[0]=i1;
				  assparticle[1]=i2;
				  assparticle[2]=i3;
				  assparticle[3]=i4;
				  assparticle[4]=i5;
				  assparticle[5]=i6;
				  assparticle[6]=i7;
				  assparticle[7]=i8;
				  assparticle[8]=i9;
				  assparticle[9]=i10;
				  chisquare=dummychisquare;
				  t0best=meantzero;
				} // close if(dummychisquare<=chisquare)
				
			      } // end loop on i10
			    } // end loop on i9
			  } // end loop on i8
			} // end loop on i7
		      } // end loop on i6
		    } // end loop on i5
		  } // end loop on i4
		} // end loop on i3
	      } // end loop on i2
	    } // end loop on i1

	    if(truparticle[0]==assparticle[0] && truparticle[1]==assparticle[1] && truparticle[2]==assparticle[2]  && truparticle[3]==assparticle[3] && truparticle[4]==assparticle[4]&& truparticle[5]==assparticle[5] && truparticle[6]==assparticle[6] && truparticle[7]==assparticle[7]  && truparticle[8]==assparticle[8] && truparticle[9]==assparticle[9]) ngood+=1;
	    if(truparticle[0]!=assparticle[0]) nmisidentified0+=1;
	    if(truparticle[1]!=assparticle[1]) nmisidentified1+=1;
	    if(truparticle[2]!=assparticle[2]) nmisidentified2+=1;
	    if(truparticle[3]!=assparticle[3]) nmisidentified3+=1;
	    if(truparticle[4]!=assparticle[4]) nmisidentified4+=1;
	    if(truparticle[5]!=assparticle[5]) nmisidentified5+=1;
	    if(truparticle[6]!=assparticle[6]) nmisidentified6+=1;
	    if(truparticle[7]!=assparticle[7]) nmisidentified7+=1;
	    if(truparticle[8]!=assparticle[8]) nmisidentified8+=1;
	    if(truparticle[9]!=assparticle[9]) nmisidentified9+=1;
	    // filling histos
	    htzerobest->Fill(t0best);
	    hchibest->Fill(chisquare);
	    Double_t dblechisquare=(Double_t)chisquare;
	    Float_t confLevel=(Float_t)TMath::Prob(dblechisquare,9); // ndf 10-1=9
	    hchibestconflevel->Fill(confLevel);
	    itimes++;
	    if(strstr(option,"all")){
	      AliInfo(Form("True Assignment %d  %d  %d  %d  %d  %d  %d  %d  %d  %d", truparticle[0], truparticle[1], truparticle[2], truparticle[3], truparticle[4], truparticle[5], truparticle[6], truparticle[7], truparticle[8], truparticle[9]));
	      AliInfo(Form("Best Assignment %d  %d  %d  %d  %d  %d  %d  %d  %d  %d", assparticle[0], assparticle[1], assparticle[2], assparticle[3], assparticle[4], assparticle[5], assparticle[6], assparticle[7], assparticle[8], assparticle[9]));
	      AliInfo(Form("Minimum ChiSquare for current set   %f ", chisquare));
	      AliInfo(Form("Confidence Level (Minimum ChiSquare) %f", confLevel));
	    }
	    if (strstr(option,"visual") && itimes && (itimes%kUPDATE) == 0) {
	      if (itimes == kUPDATE){
		c1->cd();
		htzerobest->Draw();
		c2->cd();
		hchibest->Draw();
		c3->cd();
		hchibestconflevel->Draw();
	      }
	      c1->Modified();
	      c1->Update();
	      c2->Modified();
	      c2->Update();
	      c3->Modified();
	      c3->Update();
	      if (gSystem->ProcessEvents())
		break;
	    }
	    chisquare=999.;
	    t0best=999.;
	    
	  } // end for the current set. close if(istop==5)
	} // end condition on ipartold
      } // end loop on hits for the current track
      if(istop>=10) break;
    } // end loop on ntracks  
  } //event loop
  
  if(strstr(option,"all")){
    nmisidentified=(nmisidentified0+nmisidentified1+nmisidentified2+nmisidentified3+nmisidentified4+nmisidentified5+nmisidentified6+nmisidentified7+nmisidentified8+nmisidentified9);
    AliInfo(Form("total number of tracks token into account  %i", 10*5*fNevents));
    Float_t badPercentage=100.*(Float_t)nmisidentified/(10*5*fNevents);
    AliInfo(Form("total misidentified                       %i (%f %%) ", nmisidentified, badPercentage));
    AliInfo(Form("Total Number of set token into account     %i", 5*fNevents));
    Float_t goodSetPercentage=100.*(Float_t)ngood/(5*fNevents);
    AliInfo(Form("Number of set with no misidentified tracks %i (%f %%)", ngood, goodSetPercentage));
  }

  // free used memory for canvas
  delete c1; c1=0;
  delete c2; c2=0;
  delete c3; c3=0;

  // generating output filename only if not previously specified using SetTZeroFile

  /*const*/ Int_t kSize = 70+fHeadersFile.Length()+1;
  char outFileName[kSize];
  strncpy(outFileName,"ht010tr120ps",kSize); // global time resolution has to be converted from Int_t to char
                                      // in order to have in the output filename this parameter
  strncat(outFileName,fHeadersFile,kSize);

  if(fT0File.IsNull()) fT0File=outFileName;

  TFile* houtfile = new TFile(fT0File,"recreate");
  houtfile->cd();
  htzerobest->Write(0,TObject::kOverwrite);
  hchibest->Write(0,TObject::kOverwrite);
  hchibestconflevel->Write(0,TObject::kOverwrite);
  houtfile->Close();  
  
  
  if(strstr(option,"tim") || strstr(option,"all")){
    gBenchmark->Stop("TOFT0");
    AliInfo("AliTOFT0:");
    /*
    cout << "   took " << gBenchmark->GetCpuTime("TOFT0") << " seconds in order to calculate T0 " 
	 << gBenchmark->GetCpuTime("TOFT0")/fNevents << " seconds per event " << endl ;
    */
    gBenchmark->Print("TOFT0");
  }
}
 
//__________________________________________________________________
void AliTOFT0::SetTZeroFile(char * file )
{
  //
  // Set T0 file name
  //
  printf("Destination file : %s \n", file) ;
  fT0File=file;

}

//__________________________________________________________________
void AliTOFT0::Print(Option_t* /*option*/)const
{
  //
  // Print class content
  //
  printf("------------------- %s -------------\n", GetName()) ;
  if(!fT0File.IsNull())
    printf("  Writing T0 Distribution to file  %s \n",(char*) fT0File.Data());

}

//__________________________________________________________________
Bool_t AliTOFT0::operator==( AliTOFT0 const &tzero )const
{
  //
  // Equal operator
  // 

  if( (fTimeResolution==tzero.fTimeResolution)&&
      (fLowerMomBound==tzero.fLowerMomBound)&&
      (fUpperMomBound==tzero.fUpperMomBound))
    return kTRUE ;
  else
    return kFALSE ;
}
 AliTOFT0.cxx:1
 AliTOFT0.cxx:2
 AliTOFT0.cxx:3
 AliTOFT0.cxx:4
 AliTOFT0.cxx:5
 AliTOFT0.cxx:6
 AliTOFT0.cxx:7
 AliTOFT0.cxx:8
 AliTOFT0.cxx:9
 AliTOFT0.cxx:10
 AliTOFT0.cxx:11
 AliTOFT0.cxx:12
 AliTOFT0.cxx:13
 AliTOFT0.cxx:14
 AliTOFT0.cxx:15
 AliTOFT0.cxx:16
 AliTOFT0.cxx:17
 AliTOFT0.cxx:18
 AliTOFT0.cxx:19
 AliTOFT0.cxx:20
 AliTOFT0.cxx:21
 AliTOFT0.cxx:22
 AliTOFT0.cxx:23
 AliTOFT0.cxx:24
 AliTOFT0.cxx:25
 AliTOFT0.cxx:26
 AliTOFT0.cxx:27
 AliTOFT0.cxx:28
 AliTOFT0.cxx:29
 AliTOFT0.cxx:30
 AliTOFT0.cxx:31
 AliTOFT0.cxx:32
 AliTOFT0.cxx:33
 AliTOFT0.cxx:34
 AliTOFT0.cxx:35
 AliTOFT0.cxx:36
 AliTOFT0.cxx:37
 AliTOFT0.cxx:38
 AliTOFT0.cxx:39
 AliTOFT0.cxx:40
 AliTOFT0.cxx:41
 AliTOFT0.cxx:42
 AliTOFT0.cxx:43
 AliTOFT0.cxx:44
 AliTOFT0.cxx:45
 AliTOFT0.cxx:46
 AliTOFT0.cxx:47
 AliTOFT0.cxx:48
 AliTOFT0.cxx:49
 AliTOFT0.cxx:50
 AliTOFT0.cxx:51
 AliTOFT0.cxx:52
 AliTOFT0.cxx:53
 AliTOFT0.cxx:54
 AliTOFT0.cxx:55
 AliTOFT0.cxx:56
 AliTOFT0.cxx:57
 AliTOFT0.cxx:58
 AliTOFT0.cxx:59
 AliTOFT0.cxx:60
 AliTOFT0.cxx:61
 AliTOFT0.cxx:62
 AliTOFT0.cxx:63
 AliTOFT0.cxx:64
 AliTOFT0.cxx:65
 AliTOFT0.cxx:66
 AliTOFT0.cxx:67
 AliTOFT0.cxx:68
 AliTOFT0.cxx:69
 AliTOFT0.cxx:70
 AliTOFT0.cxx:71
 AliTOFT0.cxx:72
 AliTOFT0.cxx:73
 AliTOFT0.cxx:74
 AliTOFT0.cxx:75
 AliTOFT0.cxx:76
 AliTOFT0.cxx:77
 AliTOFT0.cxx:78
 AliTOFT0.cxx:79
 AliTOFT0.cxx:80
 AliTOFT0.cxx:81
 AliTOFT0.cxx:82
 AliTOFT0.cxx:83
 AliTOFT0.cxx:84
 AliTOFT0.cxx:85
 AliTOFT0.cxx:86
 AliTOFT0.cxx:87
 AliTOFT0.cxx:88
 AliTOFT0.cxx:89
 AliTOFT0.cxx:90
 AliTOFT0.cxx:91
 AliTOFT0.cxx:92
 AliTOFT0.cxx:93
 AliTOFT0.cxx:94
 AliTOFT0.cxx:95
 AliTOFT0.cxx:96
 AliTOFT0.cxx:97
 AliTOFT0.cxx:98
 AliTOFT0.cxx:99
 AliTOFT0.cxx:100
 AliTOFT0.cxx:101
 AliTOFT0.cxx:102
 AliTOFT0.cxx:103
 AliTOFT0.cxx:104
 AliTOFT0.cxx:105
 AliTOFT0.cxx:106
 AliTOFT0.cxx:107
 AliTOFT0.cxx:108
 AliTOFT0.cxx:109
 AliTOFT0.cxx:110
 AliTOFT0.cxx:111
 AliTOFT0.cxx:112
 AliTOFT0.cxx:113
 AliTOFT0.cxx:114
 AliTOFT0.cxx:115
 AliTOFT0.cxx:116
 AliTOFT0.cxx:117
 AliTOFT0.cxx:118
 AliTOFT0.cxx:119
 AliTOFT0.cxx:120
 AliTOFT0.cxx:121
 AliTOFT0.cxx:122
 AliTOFT0.cxx:123
 AliTOFT0.cxx:124
 AliTOFT0.cxx:125
 AliTOFT0.cxx:126
 AliTOFT0.cxx:127
 AliTOFT0.cxx:128
 AliTOFT0.cxx:129
 AliTOFT0.cxx:130
 AliTOFT0.cxx:131
 AliTOFT0.cxx:132
 AliTOFT0.cxx:133
 AliTOFT0.cxx:134
 AliTOFT0.cxx:135
 AliTOFT0.cxx:136
 AliTOFT0.cxx:137
 AliTOFT0.cxx:138
 AliTOFT0.cxx:139
 AliTOFT0.cxx:140
 AliTOFT0.cxx:141
 AliTOFT0.cxx:142
 AliTOFT0.cxx:143
 AliTOFT0.cxx:144
 AliTOFT0.cxx:145
 AliTOFT0.cxx:146
 AliTOFT0.cxx:147
 AliTOFT0.cxx:148
 AliTOFT0.cxx:149
 AliTOFT0.cxx:150
 AliTOFT0.cxx:151
 AliTOFT0.cxx:152
 AliTOFT0.cxx:153
 AliTOFT0.cxx:154
 AliTOFT0.cxx:155
 AliTOFT0.cxx:156
 AliTOFT0.cxx:157
 AliTOFT0.cxx:158
 AliTOFT0.cxx:159
 AliTOFT0.cxx:160
 AliTOFT0.cxx:161
 AliTOFT0.cxx:162
 AliTOFT0.cxx:163
 AliTOFT0.cxx:164
 AliTOFT0.cxx:165
 AliTOFT0.cxx:166
 AliTOFT0.cxx:167
 AliTOFT0.cxx:168
 AliTOFT0.cxx:169
 AliTOFT0.cxx:170
 AliTOFT0.cxx:171
 AliTOFT0.cxx:172
 AliTOFT0.cxx:173
 AliTOFT0.cxx:174
 AliTOFT0.cxx:175
 AliTOFT0.cxx:176
 AliTOFT0.cxx:177
 AliTOFT0.cxx:178
 AliTOFT0.cxx:179
 AliTOFT0.cxx:180
 AliTOFT0.cxx:181
 AliTOFT0.cxx:182
 AliTOFT0.cxx:183
 AliTOFT0.cxx:184
 AliTOFT0.cxx:185
 AliTOFT0.cxx:186
 AliTOFT0.cxx:187
 AliTOFT0.cxx:188
 AliTOFT0.cxx:189
 AliTOFT0.cxx:190
 AliTOFT0.cxx:191
 AliTOFT0.cxx:192
 AliTOFT0.cxx:193
 AliTOFT0.cxx:194
 AliTOFT0.cxx:195
 AliTOFT0.cxx:196
 AliTOFT0.cxx:197
 AliTOFT0.cxx:198
 AliTOFT0.cxx:199
 AliTOFT0.cxx:200
 AliTOFT0.cxx:201
 AliTOFT0.cxx:202
 AliTOFT0.cxx:203
 AliTOFT0.cxx:204
 AliTOFT0.cxx:205
 AliTOFT0.cxx:206
 AliTOFT0.cxx:207
 AliTOFT0.cxx:208
 AliTOFT0.cxx:209
 AliTOFT0.cxx:210
 AliTOFT0.cxx:211
 AliTOFT0.cxx:212
 AliTOFT0.cxx:213
 AliTOFT0.cxx:214
 AliTOFT0.cxx:215
 AliTOFT0.cxx:216
 AliTOFT0.cxx:217
 AliTOFT0.cxx:218
 AliTOFT0.cxx:219
 AliTOFT0.cxx:220
 AliTOFT0.cxx:221
 AliTOFT0.cxx:222
 AliTOFT0.cxx:223
 AliTOFT0.cxx:224
 AliTOFT0.cxx:225
 AliTOFT0.cxx:226
 AliTOFT0.cxx:227
 AliTOFT0.cxx:228
 AliTOFT0.cxx:229
 AliTOFT0.cxx:230
 AliTOFT0.cxx:231
 AliTOFT0.cxx:232
 AliTOFT0.cxx:233
 AliTOFT0.cxx:234
 AliTOFT0.cxx:235
 AliTOFT0.cxx:236
 AliTOFT0.cxx:237
 AliTOFT0.cxx:238
 AliTOFT0.cxx:239
 AliTOFT0.cxx:240
 AliTOFT0.cxx:241
 AliTOFT0.cxx:242
 AliTOFT0.cxx:243
 AliTOFT0.cxx:244
 AliTOFT0.cxx:245
 AliTOFT0.cxx:246
 AliTOFT0.cxx:247
 AliTOFT0.cxx:248
 AliTOFT0.cxx:249
 AliTOFT0.cxx:250
 AliTOFT0.cxx:251
 AliTOFT0.cxx:252
 AliTOFT0.cxx:253
 AliTOFT0.cxx:254
 AliTOFT0.cxx:255
 AliTOFT0.cxx:256
 AliTOFT0.cxx:257
 AliTOFT0.cxx:258
 AliTOFT0.cxx:259
 AliTOFT0.cxx:260
 AliTOFT0.cxx:261
 AliTOFT0.cxx:262
 AliTOFT0.cxx:263
 AliTOFT0.cxx:264
 AliTOFT0.cxx:265
 AliTOFT0.cxx:266
 AliTOFT0.cxx:267
 AliTOFT0.cxx:268
 AliTOFT0.cxx:269
 AliTOFT0.cxx:270
 AliTOFT0.cxx:271
 AliTOFT0.cxx:272
 AliTOFT0.cxx:273
 AliTOFT0.cxx:274
 AliTOFT0.cxx:275
 AliTOFT0.cxx:276
 AliTOFT0.cxx:277
 AliTOFT0.cxx:278
 AliTOFT0.cxx:279
 AliTOFT0.cxx:280
 AliTOFT0.cxx:281
 AliTOFT0.cxx:282
 AliTOFT0.cxx:283
 AliTOFT0.cxx:284
 AliTOFT0.cxx:285
 AliTOFT0.cxx:286
 AliTOFT0.cxx:287
 AliTOFT0.cxx:288
 AliTOFT0.cxx:289
 AliTOFT0.cxx:290
 AliTOFT0.cxx:291
 AliTOFT0.cxx:292
 AliTOFT0.cxx:293
 AliTOFT0.cxx:294
 AliTOFT0.cxx:295
 AliTOFT0.cxx:296
 AliTOFT0.cxx:297
 AliTOFT0.cxx:298
 AliTOFT0.cxx:299
 AliTOFT0.cxx:300
 AliTOFT0.cxx:301
 AliTOFT0.cxx:302
 AliTOFT0.cxx:303
 AliTOFT0.cxx:304
 AliTOFT0.cxx:305
 AliTOFT0.cxx:306
 AliTOFT0.cxx:307
 AliTOFT0.cxx:308
 AliTOFT0.cxx:309
 AliTOFT0.cxx:310
 AliTOFT0.cxx:311
 AliTOFT0.cxx:312
 AliTOFT0.cxx:313
 AliTOFT0.cxx:314
 AliTOFT0.cxx:315
 AliTOFT0.cxx:316
 AliTOFT0.cxx:317
 AliTOFT0.cxx:318
 AliTOFT0.cxx:319
 AliTOFT0.cxx:320
 AliTOFT0.cxx:321
 AliTOFT0.cxx:322
 AliTOFT0.cxx:323
 AliTOFT0.cxx:324
 AliTOFT0.cxx:325
 AliTOFT0.cxx:326
 AliTOFT0.cxx:327
 AliTOFT0.cxx:328
 AliTOFT0.cxx:329
 AliTOFT0.cxx:330
 AliTOFT0.cxx:331
 AliTOFT0.cxx:332
 AliTOFT0.cxx:333
 AliTOFT0.cxx:334
 AliTOFT0.cxx:335
 AliTOFT0.cxx:336
 AliTOFT0.cxx:337
 AliTOFT0.cxx:338
 AliTOFT0.cxx:339
 AliTOFT0.cxx:340
 AliTOFT0.cxx:341
 AliTOFT0.cxx:342
 AliTOFT0.cxx:343
 AliTOFT0.cxx:344
 AliTOFT0.cxx:345
 AliTOFT0.cxx:346
 AliTOFT0.cxx:347
 AliTOFT0.cxx:348
 AliTOFT0.cxx:349
 AliTOFT0.cxx:350
 AliTOFT0.cxx:351
 AliTOFT0.cxx:352
 AliTOFT0.cxx:353
 AliTOFT0.cxx:354
 AliTOFT0.cxx:355
 AliTOFT0.cxx:356
 AliTOFT0.cxx:357
 AliTOFT0.cxx:358
 AliTOFT0.cxx:359
 AliTOFT0.cxx:360
 AliTOFT0.cxx:361
 AliTOFT0.cxx:362
 AliTOFT0.cxx:363
 AliTOFT0.cxx:364
 AliTOFT0.cxx:365
 AliTOFT0.cxx:366
 AliTOFT0.cxx:367
 AliTOFT0.cxx:368
 AliTOFT0.cxx:369
 AliTOFT0.cxx:370
 AliTOFT0.cxx:371
 AliTOFT0.cxx:372
 AliTOFT0.cxx:373
 AliTOFT0.cxx:374
 AliTOFT0.cxx:375
 AliTOFT0.cxx:376
 AliTOFT0.cxx:377
 AliTOFT0.cxx:378
 AliTOFT0.cxx:379
 AliTOFT0.cxx:380
 AliTOFT0.cxx:381
 AliTOFT0.cxx:382
 AliTOFT0.cxx:383
 AliTOFT0.cxx:384
 AliTOFT0.cxx:385
 AliTOFT0.cxx:386
 AliTOFT0.cxx:387
 AliTOFT0.cxx:388
 AliTOFT0.cxx:389
 AliTOFT0.cxx:390
 AliTOFT0.cxx:391
 AliTOFT0.cxx:392
 AliTOFT0.cxx:393
 AliTOFT0.cxx:394
 AliTOFT0.cxx:395
 AliTOFT0.cxx:396
 AliTOFT0.cxx:397
 AliTOFT0.cxx:398
 AliTOFT0.cxx:399
 AliTOFT0.cxx:400
 AliTOFT0.cxx:401
 AliTOFT0.cxx:402
 AliTOFT0.cxx:403
 AliTOFT0.cxx:404
 AliTOFT0.cxx:405
 AliTOFT0.cxx:406
 AliTOFT0.cxx:407
 AliTOFT0.cxx:408
 AliTOFT0.cxx:409
 AliTOFT0.cxx:410
 AliTOFT0.cxx:411
 AliTOFT0.cxx:412
 AliTOFT0.cxx:413
 AliTOFT0.cxx:414
 AliTOFT0.cxx:415
 AliTOFT0.cxx:416
 AliTOFT0.cxx:417
 AliTOFT0.cxx:418
 AliTOFT0.cxx:419
 AliTOFT0.cxx:420
 AliTOFT0.cxx:421
 AliTOFT0.cxx:422
 AliTOFT0.cxx:423
 AliTOFT0.cxx:424
 AliTOFT0.cxx:425
 AliTOFT0.cxx:426
 AliTOFT0.cxx:427
 AliTOFT0.cxx:428
 AliTOFT0.cxx:429
 AliTOFT0.cxx:430
 AliTOFT0.cxx:431
 AliTOFT0.cxx:432
 AliTOFT0.cxx:433
 AliTOFT0.cxx:434
 AliTOFT0.cxx:435
 AliTOFT0.cxx:436
 AliTOFT0.cxx:437
 AliTOFT0.cxx:438
 AliTOFT0.cxx:439
 AliTOFT0.cxx:440
 AliTOFT0.cxx:441
 AliTOFT0.cxx:442
 AliTOFT0.cxx:443
 AliTOFT0.cxx:444
 AliTOFT0.cxx:445
 AliTOFT0.cxx:446
 AliTOFT0.cxx:447
 AliTOFT0.cxx:448
 AliTOFT0.cxx:449
 AliTOFT0.cxx:450
 AliTOFT0.cxx:451
 AliTOFT0.cxx:452
 AliTOFT0.cxx:453
 AliTOFT0.cxx:454
 AliTOFT0.cxx:455
 AliTOFT0.cxx:456
 AliTOFT0.cxx:457
 AliTOFT0.cxx:458
 AliTOFT0.cxx:459
 AliTOFT0.cxx:460
 AliTOFT0.cxx:461
 AliTOFT0.cxx:462
 AliTOFT0.cxx:463
 AliTOFT0.cxx:464
 AliTOFT0.cxx:465
 AliTOFT0.cxx:466
 AliTOFT0.cxx:467
 AliTOFT0.cxx:468
 AliTOFT0.cxx:469
 AliTOFT0.cxx:470
 AliTOFT0.cxx:471
 AliTOFT0.cxx:472
 AliTOFT0.cxx:473
 AliTOFT0.cxx:474
 AliTOFT0.cxx:475
 AliTOFT0.cxx:476
 AliTOFT0.cxx:477
 AliTOFT0.cxx:478
 AliTOFT0.cxx:479
 AliTOFT0.cxx:480
 AliTOFT0.cxx:481
 AliTOFT0.cxx:482
 AliTOFT0.cxx:483
 AliTOFT0.cxx:484
 AliTOFT0.cxx:485
 AliTOFT0.cxx:486
 AliTOFT0.cxx:487
 AliTOFT0.cxx:488
 AliTOFT0.cxx:489
 AliTOFT0.cxx:490
 AliTOFT0.cxx:491
 AliTOFT0.cxx:492
 AliTOFT0.cxx:493
 AliTOFT0.cxx:494
 AliTOFT0.cxx:495
 AliTOFT0.cxx:496
 AliTOFT0.cxx:497
 AliTOFT0.cxx:498
 AliTOFT0.cxx:499
 AliTOFT0.cxx:500
 AliTOFT0.cxx:501
 AliTOFT0.cxx:502
 AliTOFT0.cxx:503
 AliTOFT0.cxx:504
 AliTOFT0.cxx:505
 AliTOFT0.cxx:506
 AliTOFT0.cxx:507
 AliTOFT0.cxx:508
 AliTOFT0.cxx:509
 AliTOFT0.cxx:510
 AliTOFT0.cxx:511
 AliTOFT0.cxx:512
 AliTOFT0.cxx:513
 AliTOFT0.cxx:514
 AliTOFT0.cxx:515
 AliTOFT0.cxx:516
 AliTOFT0.cxx:517
 AliTOFT0.cxx:518
 AliTOFT0.cxx:519
 AliTOFT0.cxx:520
 AliTOFT0.cxx:521
 AliTOFT0.cxx:522
 AliTOFT0.cxx:523
 AliTOFT0.cxx:524
 AliTOFT0.cxx:525
 AliTOFT0.cxx:526
 AliTOFT0.cxx:527
 AliTOFT0.cxx:528
 AliTOFT0.cxx:529
 AliTOFT0.cxx:530
 AliTOFT0.cxx:531
 AliTOFT0.cxx:532
 AliTOFT0.cxx:533
 AliTOFT0.cxx:534
 AliTOFT0.cxx:535
 AliTOFT0.cxx:536
 AliTOFT0.cxx:537
 AliTOFT0.cxx:538
 AliTOFT0.cxx:539
 AliTOFT0.cxx:540
 AliTOFT0.cxx:541
 AliTOFT0.cxx:542
 AliTOFT0.cxx:543
 AliTOFT0.cxx:544
 AliTOFT0.cxx:545
 AliTOFT0.cxx:546
 AliTOFT0.cxx:547
 AliTOFT0.cxx:548
 AliTOFT0.cxx:549
 AliTOFT0.cxx:550
 AliTOFT0.cxx:551
 AliTOFT0.cxx:552
 AliTOFT0.cxx:553
 AliTOFT0.cxx:554
 AliTOFT0.cxx:555
 AliTOFT0.cxx:556
 AliTOFT0.cxx:557
 AliTOFT0.cxx:558
 AliTOFT0.cxx:559
 AliTOFT0.cxx:560
 AliTOFT0.cxx:561
 AliTOFT0.cxx:562
 AliTOFT0.cxx:563
 AliTOFT0.cxx:564
 AliTOFT0.cxx:565
 AliTOFT0.cxx:566
 AliTOFT0.cxx:567
 AliTOFT0.cxx:568
 AliTOFT0.cxx:569
 AliTOFT0.cxx:570
 AliTOFT0.cxx:571
 AliTOFT0.cxx:572
 AliTOFT0.cxx:573
 AliTOFT0.cxx:574
 AliTOFT0.cxx:575
 AliTOFT0.cxx:576