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

/*
  TO CHECK
  GetBz usage

 */

//-------------------------------------------------------------------------
//  Implementation of the ITS Upgrade tracker mother class.              
//-------------------------------------------------------------------------
#include <TTree.h>
#include <Riostream.h> 
#include <TMath.h>
#include <TFile.h>

#include "AliITSUTrackerGlo.h"
#include "AliESDEvent.h"
#include "AliESDtrack.h"
#include "AliITSURecoDet.h"
#include "AliITSURecoSens.h"
#include "AliITSUReconstructor.h"
#include "AliITSReconstructor.h"
#include "AliITSUSeed.h"
#include "AliITSUClusterPix.h"
#include "AliITSUGeomTGeo.h"
#include "AliCodeTimer.h"
#include "AliRefArray.h"
using namespace AliITSUAux;
using namespace TMath;

//----------------- tmp stuff -----------------

ClassImp(AliITSUTrackerGlo)

const Double_t AliITSUTrackerGlo::fgkToler =  1e-6;// tolerance for layer on-surface check


//_________________________________________________________________________
AliITSUTrackerGlo::AliITSUTrackerGlo(AliITSUReconstructor* rec)
:  fReconstructor(rec)
  ,fITS(0)
  ,fMatLUT(0)
  ,fCurrESDtrack(0)
  ,fCurrESDtrMClb(kDummyLabel)
  ,fCurrMass(kPionMass)
  ,fCountProlongationTrials(0)
  ,fCountITSin(0)
  ,fCountITSout(0)
  ,fCountITSrefit(0)
  ,fNTracksESD(0)
  ,fHypStore(100)
  ,fLayerMaxCandidates(1000)
  ,fLayerCandidates(0)
  ,fNBranchesAdded(0)
  ,fNCandidatesAdded(0)
  ,fCurrHyp(0)
  ,fWorkHyp(0)
  ,fSeedsPool("AliITSUSeed",0)
  ,fFreeSeedsID(0)
  ,fESDIndex(0)
  ,fNFreeSeeds(0)
  ,fLastSeedID(0)
  ,fNLrActive(0)
  ,fDefTrackConds(0)
  ,fCurrTrackCond(0)
  ,fCurrActLrID(-1)
  ,fCurrLayer(0)
  ,fTrackPhaseID(-1)
  ,fCurrPassID(-1)
  ,fUseMatLUT(kFALSE)
#ifdef  _ITSU_TUNING_MODE_
  ,fCHistoArrCorr(0)
  ,fCHistoArrFake(0)
#endif
{
  // Default constructor
  if (rec) Init(rec);
}

//_________________________________________________________________________
AliITSUTrackerGlo::~AliITSUTrackerGlo()
{
 // Default destructor
 //  
  delete[] fLayerCandidates;
  if (fWorkHyp) fWorkHyp->SetTPCSeed(0); // this hypothesis does not own the seed
  delete fWorkHyp;
  delete fMatLUT;
  //
#ifdef  _ITSU_TUNING_MODE_
  if (fCHistoArrCorr || fCHistoArrFake) {
    TFile* ctrOut = TFile::Open("itsuControlHistos.root","recreate");
    ctrOut->cd();
    AliInfo("Storing control histos");
    //    ctrOut->WriteObject(fCHistoArr,"controlH","kSingleKey");
    if (fCHistoArrCorr) {
      fCHistoArrCorr->Write();
      delete fCHistoArrCorr;
    }
    if (fCHistoArrFake) {
      fCHistoArrFake->Write();
      delete fCHistoArrFake;
    }
    ctrOut->Close();
    delete ctrOut;
    fCHistoArrCorr = 0;
    fCHistoArrFake = 0;    
  }
#endif 
  //
}

//_________________________________________________________________________
void AliITSUTrackerGlo::Init(AliITSUReconstructor* rec)
{
  // init with external reconstructor
  //
  fITS = rec->GetITSInterface();
  fNLrActive = fITS->GetNLayersActive();
  fWorkHyp = new AliITSUTrackHyp(fNLrActive);
  // create material lookup table
  const int kNTest = 1000;
  const double kStepsPerCM=5;
  fMatLUT  = new AliITSUMatLUT(fITS->GetRMin(),fITS->GetRMax(),Nint(kStepsPerCM*(fITS->GetRMax()-fITS->GetRMin())));
  double zmn = 1e6;
  for (int ilr=fITS->GetNLayers();ilr--;) {
    AliITSURecoLayer* lr = fITS->GetLayer(ilr);
    if (zmn>Abs(lr->GetZMin())) zmn = Abs(lr->GetZMin());
    if (zmn>Abs(lr->GetZMax())) zmn = Abs(lr->GetZMax());
  }
  fMatLUT->FillData(kNTest,-zmn,zmn);
  //
  if (fLayerMaxCandidates<1) fLayerMaxCandidates = 1000;
  fLayerCandidates = new AliITSUSeed*[fLayerMaxCandidates];
  fSeedsPool.ExpandCreateFast(1000); // RS TOCHECK
  fFreeSeedsID.Set(1000);
  fESDIndex.Set(1000);
  //
}

//_________________________________________________________________________
void AliITSUTrackerGlo::CreateDefaultTrackCond()
{
  // creates default tracking conditions to be used when recoparam does not provide them

  AliITSUTrackCond* cond = new AliITSUTrackCond();
  //
  cond->SetNLayers(fNLrActive);
  cond->AddNewCondition(fNLrActive);
  cond->AddGroupPattern( 0xffff, fNLrActive ); // require all layers hit
  cond->Init();
  //
  fDefTrackConds.AddLast(cond);
  //
  AliInfo("Created conditions: ");
  for (int i=0;i<fDefTrackConds.GetEntriesFast();i++) fDefTrackConds[i]->Print();
  //
}


//_________________________________________________________________________
Int_t AliITSUTrackerGlo::Clusters2Tracks(AliESDEvent *esdEv)
{
  //
  AliCodeTimerAuto("",0);
  SetTrackingPhase(kClus2Tracks);
  //
#ifdef  _ITSU_TUNING_MODE_
  if (!fCHistoArrCorr) fCHistoArrCorr = BookControlHistos("Corr");
  if (!fCHistoArrFake) fCHistoArrFake = BookControlHistos("Fake");
#endif
  static int evID = 0;
  static TArrayF esdTrPt(fESDIndex.GetSize()); 
  //
  TObjArray *trackConds = 0;
  //
  fCountProlongationTrials = 0;
  fCountITSin = 0;
  fCountITSout = 0;
  fCountITSrefit = 0;
  //
  ResetSeedsPool();
  fNTracksESD = esdEv->GetNumberOfTracks();
  AliInfo(Form("Will try to find prolongations for %d tracks",fNTracksESD));
  int nTrackCond = AliITSUReconstructor::GetRecoParam()->GetNTrackingConditions();
  fUseMatLUT = AliITSUReconstructor::GetRecoParam()->GetUseMatLUT(fTrackPhaseID);
  if (nTrackCond<1) {
    if (!fDefTrackConds.GetEntriesFast()) {
      AliInfo("No tracking conditions found in recoparams, creating default one requesting all layers hit");
      CreateDefaultTrackCond();
    }
    trackConds = &fDefTrackConds;
    nTrackCond = trackConds->GetEntriesFast();
  } 
  else trackConds = AliITSUReconstructor::GetRecoParam()->GetTrackingConditions();
  //
  static Bool_t first = kTRUE;
  if (first) {
    AliITSUReconstructor::GetRecoParam()->Print();
    first = kFALSE;
  }
  fHypStore.Delete();
  if (fHypStore.GetSize()<fNTracksESD) fHypStore.Expand(fNTracksESD+100);
  //
  fITS->SortClusters(AliITSUClusterPix::SortModeIdTrkYZ());
  fITS->ProcessClusters();
  //
#ifdef  _ITSU_TUNING_MODE_
  FlagSplitClusters(); // tmp RS
#endif
  //
  // the tracks will be reconstructed in decreasing pt order, sort them
  if (fESDIndex.GetSize()<fNTracksESD) {
    fESDIndex.Set(fNTracksESD+200);
    esdTrPt.Set(fNTracksESD+200);
  }
  int* esdTrackIndex = fESDIndex.GetArray();
  float* trPt = esdTrPt.GetArray();
  for (int itr=fNTracksESD;itr--;) {
    AliESDtrack* esdTr = esdEv->GetTrack(itr);
    esdTr->SetStatus(AliESDtrack::kITSupg);   // Flag all tracks as participating in the ITSup reco
    trPt[itr] = esdTr->Pt();
  }
  Sort(fNTracksESD,trPt,esdTrackIndex,kTRUE);    
  //
  for (int icnd=0;icnd<nTrackCond;icnd++) {
    fCurrPassID = icnd;
    fCurrTrackCond = (AliITSUTrackCond*)trackConds->UncheckedAt(icnd);
    if (!fCurrTrackCond->IsInitDone()) fCurrTrackCond->Init();
    // select ESD tracks to propagate
    for (int itr=0;itr<fNTracksESD;itr++) {
      int trID = esdTrackIndex[itr];
      fCurrESDtrack = esdEv->GetTrack(trID);
      fCurrESDtrMClb = fCurrESDtrack->GetLabel();
      //
      
      if (!NeedToProlong(fCurrESDtrack,trID)) continue;  // are we interested in this track?
      /*
	// if specific tracks should be checked, create a global array int wtc[] = {evITS*10000+trID, ...};
      Bool_t dbg = kFALSE;
      int nwtc = sizeof(wtc)/sizeof(int);

      for (int k=0;k<nwtc;k++) {
	if (wtc[k]==evID*10000+trID) {
	  dbg = kTRUE;
	  printf("At esdTr: %d MC: %d\n",wtc[k],fCurrESDtrMClb);
	  break;
	}
      }
      AliLog::SetClassDebugLevel("AliITSUTrackerGlo",dbg ? 10:0);
      */
#ifdef _ITSU_DEBUG_
      AliDebug(1,Form("Processing track %d(esd%d) | M=%.3f Pt=%.3f | MCLabel: %d",itr,trID,fCurrESDtrack->GetMassForTracking(kTRUE),fCurrESDtrack->Pt(),fCurrESDtrMClb));//RS
#endif
      FindTrack(fCurrESDtrack, trID);
    }   
    //
#ifdef _ITSU_DEBUG_
    if (AliDebugLevelClass()>+2) {
      AliInfo(Form("SeedsPool: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
      fHypStore.Print();
    }
#endif
    FinalizeHypotheses();
#ifdef  _ITSU_TUNING_MODE_
    CheckClusterUsage(); //!!RS
#endif
    //
    //AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0);     // in case wtc array was defined, uncomment this
  }
  //
  AliInfo(Form("%d ITSin for %d tried TPC seeds out of %d ESD tracks\n",fCountITSin,fCountProlongationTrials,fNTracksESD));
  evID++;
  return 0;
}

//_________________________________________________________________________
Int_t AliITSUTrackerGlo::PropagateBack(AliESDEvent *esdEv)
{
  //
  // Do outward fits in ITS
  AliCodeTimerAuto("",0);
  //
  SetTrackingPhase(kPropBack);
  fNTracksESD = esdEv->GetNumberOfTracks();
  fUseMatLUT = AliITSUReconstructor::GetRecoParam()->GetUseMatLUT(fTrackPhaseID);
#ifdef _ITSU_DEBUG_
  AliDebug(1,Form("Will propagate back %d tracks",fNTracksESD));
#endif
  //
  double bz0 = GetBz();
  Double_t xyzTrk[3],xyzVtx[3]={GetX(),GetY(),GetZ()};
  AliITSUTrackHyp dummyTr;
  const double kWatchStep=10.; // for larger steps watch arc vs segment difference
  Double_t times[AliPID::kSPECIESC];
  //
  for (int itr=0;itr<fNTracksESD;itr++) {
    fCurrESDtrack = esdEv->GetTrack(itr);
    fCurrESDtrMClb = fCurrESDtrack->GetLabel();
    // Start time integral and add distance from current position to vertex 
    if (fCurrESDtrack->IsOn(AliESDtrack::kITSout)) continue; // already done
    if (!fCurrESDtrack->IsOn(AliESDtrack::kTPCin)) continue; // skip ITS s.a.
    //
    fCurrESDtrack->GetXYZ(xyzTrk); 
    Double_t dst = 0.;     // set initial track lenght, tof
    {
      double dxs = xyzTrk[0] - xyzVtx[0];
      double dys = xyzTrk[1] - xyzVtx[1];
      //      double dzs = xyzTrk[2] - xyzVtx[2];
      // RS: for large segment steps d use approximation of cicrular arc s by
      // s = 2R*asin(d/2R) = d/p asin(p) \approx d/p (p + 1/6 p^3) = d (1+1/6 p^2)
      // where R is the track radius, p = d/2R = 2C*d (C is the abs curvature)
      // Hence s^2/d^2 = (1+1/6 p^2)^2
      dst = dxs*dxs + dys*dys;
      if (dst > kWatchStep*kWatchStep) { // correct circular part for arc/segment factor
	double crv = Abs(fCurrESDtrack->GetC(bz0));
	double fcarc = 1.+crv*crv*dst/6.;
	dst *= fcarc*fcarc;
      }
      // RS: temporary hack since we don't have SPD vertex:      
      //      dst += dzs*dzs;
      dst *= 1+fCurrESDtrack->GetTgl()*fCurrESDtrack->GetTgl();
      dst = Sqrt(dst); 
    }
    //    
    fCurrESDtrack->SetStatus(AliESDtrack::kTIME);
    //
    if (!fCurrESDtrack->IsOn(AliESDtrack::kITSin)) { // case of tracks w/o ITS prolongation: just set the time integration
      dummyTr.AliExternalTrackParam::operator=(*fCurrESDtrack);
      dummyTr.StartTimeIntegral();
      dummyTr.AddTimeStep(dst);
      dummyTr.GetIntegratedTimes(times,AliPID::kSPECIESC); 
      fCurrESDtrack->SetIntegratedTimes(times);
      fCurrESDtrack->SetIntegratedLength(dummyTr.GetIntegratedLength());
      continue;
    }
    //
    fCurrHyp = GetTrackHyp(itr);
    fCurrMass = fCurrHyp->GetMass();
    fCurrHyp->StartTimeIntegral();
    fCurrHyp->AddTimeStep(dst);
    fCurrHyp->ResetCovariance(10000);
    int nclFit = 0;
    double chi2 = RefitTrack(fCurrHyp,fITS->GetRMax(),nclFit);
    if (chi2>0) { // propagate to exit from the ITS/TPC screen
      int ndf = nclFit*2-5;
      if (ndf>0) chi2 /= ndf;
      fCurrHyp->SetChi2(chi2);
      UpdateESDTrack(fCurrHyp,AliESDtrack::kITSout);
    }
    else {
#ifdef _ITSU_DEBUG_
      AliDebug(2,Form("Refit Failed for track %d | ESDtrack#%d (MClb:%d)",itr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
#endif
      //fCurrHyp->AliExternalTrackParam::Print();
      //fCurrHyp->GetWinner()->Print();
    }
    //
  }
  //
  AliInfo(Form("%d ITSout in %d ITSin tracks for %d tried TPC seeds out of %d ESD tracks\n",
	       fCountITSout,fCountITSin,fCountProlongationTrials,fNTracksESD));
  //
  return 0;
}

//_________________________________________________________________________
Int_t AliITSUTrackerGlo::RefitInward(AliESDEvent *esdEv)
{
  //
  // refit the tracks inward, using current cov.matrix
  AliCodeTimerAuto("",0);
  //
  SetTrackingPhase(kRefitInw);
  fUseMatLUT = AliITSUReconstructor::GetRecoParam()->GetUseMatLUT(fTrackPhaseID);
  fNTracksESD = esdEv->GetNumberOfTracks();
  //  AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);

  AliDebug(1,Form("Will refit inward %d tracks",fNTracksESD));
  //  Bool_t uselogMS = AliExternalTrackParam::GetUseLogTermMS();
  //  AliExternalTrackParam::SetUseLogTermMS(kTRUE);
  //
  for (int itr=0;itr<fNTracksESD;itr++) {
    fCurrESDtrack = esdEv->GetTrack(itr);
    fCurrESDtrMClb = fCurrESDtrack->GetLabel();
    // Start time integral and add distance from current position to vertex 
    UInt_t trStat = fCurrESDtrack->GetStatus();
    if ( !(trStat & AliESDtrack::kITSout) ) continue;
    if (   trStat & AliESDtrack::kITSrefit ) continue; // already done
    if ( !(trStat & AliESDtrack::kTPCin)   ) continue; // skip ITS s.a.
    if (  (trStat & AliESDtrack::kTPCout) && !(trStat & AliESDtrack::kTPCrefit) ) continue;
    //
    fCurrHyp  = GetTrackHyp(itr);
    fCurrMass = fCurrHyp->GetMass();
    *fCurrHyp = *fCurrESDtrack;  // fetch current ESDtrack kinematics
    //
    int nclFit = 0;
    //    fCurrHyp->ResetCovariance(1000); // in case we want to do SA refit
    double chi2 = RefitTrack(fCurrHyp,fITS->GetRMin(),nclFit);
    if (chi2>0) { // propagate up to inside radius of the beam pipe      
      int ndf = nclFit*2-5;
      if (ndf>0) chi2 /= ndf;
      fCurrHyp->SetChi2(chi2);
      UpdateESDTrack(fCurrHyp,AliESDtrack::kITSrefit);
    }
    else {
#ifdef _ITSU_DEBUG_
      AliDebug(2,Form("Refit Failed for track %d |ESDtrack#%d (MClb:%d)",itr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
#endif
    }
  }    
  //  AliExternalTrackParam::SetUseLogTermMS(uselogMS);
  //
  AliInfo(Form("%d ITSrefit in %d ITSout in %d ITSin tracks for %d tried TPC seeds out of %d ESD tracks\n",
	       fCountITSrefit,fCountITSout,fCountITSin,fCountProlongationTrials,fNTracksESD));
  //
  //  AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0);
  return 0;
}

//_________________________________________________________________________
Int_t AliITSUTrackerGlo::LoadClusters(TTree * treeRP)
{
  // read from tree (if pointer provided) or directly from the ITS reco interface
  //
  return fReconstructor->LoadClusters(treeRP);
} 

//_________________________________________________________________________
void AliITSUTrackerGlo::UnloadClusters()
{
  //
  // Remove clusters from the memory 
  //
  /*
  AliITSURecoDet *det=fReconstructor->GetITSInterface();
  Int_t nlayers=det->GetNLayersActive();
  for (Int_t i=0; i<nlayers; i++) {
      TClonesArray *clusters=*(det->GetLayerActive(i)->GetClustersAddress());
      clusters->Delete();
  }
  */
} 

//_________________________________________________________________________
AliCluster * AliITSUTrackerGlo::GetCluster(Int_t /*index*/) const
{
  //
  // To be implemented 
  //
  
  Info("GetCluster","To be implemented");
  return 0x0;
} 

//_________________________________________________________________________
Bool_t AliITSUTrackerGlo::NeedToProlong(AliESDtrack* esdTr, Int_t esdID)
{
  // do we need to match this track to ITS?
  //
  static double bz = GetBz();
  //
  // is it already prolonged
  if (esdTr->IsOn(AliESDtrack::kITSin)) return kFALSE;
  AliITSUTrackHyp* hyp = GetTrackHyp(esdID); // in case there is unfinished hypothesis from prev.pass, clean it
  if (hyp) {
    CleanHypothesis(hyp);
    if (hyp->GetSkip()) return kFALSE; // need to skip this
    return kTRUE;
  }
  //
  if (!esdTr->IsOn(AliESDtrack::kTPCin) ||
      esdTr->IsOn(AliESDtrack::kTPCout) ||
      esdTr->IsOn(AliESDtrack::kITSin)  ||
      esdTr->GetKinkIndex(0)>0) return kFALSE;
  //
  if (esdTr->Pt()<AliITSUReconstructor::GetRecoParam()->GetMinPtForProlongation()) return kFALSE;
  //
  float dtz[2];
  esdTr->GetDZ(GetX(),GetY(),GetZ(),bz,dtz); 
  // if track is not V0 candidata but has large offset wrt IP, reject it. RS TOCHECK
  if ( !(esdTr->GetV0Index(0)>0 && dtz[0]>AliITSUReconstructor::GetRecoParam()->GetMaxDforV0dghtrForProlongation()) 
       && (Abs(dtz[0])>AliITSUReconstructor::GetRecoParam()->GetMaxDForProlongation() ||
	   Abs(dtz[1])>AliITSUReconstructor::GetRecoParam()->GetMaxDZForProlongation())) return kFALSE;
  //
  //  if (esdTr->Pt()<3) return kFALSE;//RS
  return kTRUE;
}

//_________________________________________________________________________
void AliITSUTrackerGlo::FindTrack(AliESDtrack* esdTr, Int_t esdID)
{
  // find prolongaion candidates finding for single seed
  //
  AliITSUSeed seedUC;  // copy of the seed from the upper layer
  AliITSUSeed seedT;   // transient seed between the seedUC and new prolongation hypothesis
  //
  //  AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
  AliITSUTrackHyp* hypTr = InitHypothesis(esdTr,esdID);  // initialize prolongations hypotheses tree
  if (!hypTr || hypTr->GetSkip()) return;
  //
  double bz = GetBz();
  fCurrHyp = fWorkHyp;
  fCurrHyp->InitFrom(hypTr);
  //
  AliITSURecoSens *hitSens[AliITSURecoLayer::kMaxSensMatching];
  //
  int ilaUp = fNLrActive;     // previous active layer really checked (some may be excluded!)
  //
  for (int ila=fNLrActive;ila--;ilaUp=fCurrActLrID) {
    if (fCurrTrackCond->IsLayerExcluded(ila)) continue;
    fCurrActLrID = ila;
    fCurrLayer = fITS->GetLayerActive(ila);
    Bool_t noClSharing = fCurrTrackCond->GetClSharing(ila)==0;
    //
    // for the outermost layer the seed is created from the ESD track
    int nSeedsUp = (ilaUp==fNLrActive) ? 1 : fCurrHyp->GetNSeeds(ilaUp);
    int maxNBranches   = fCurrTrackCond->GetMaxBranches(ila);
    int maxNCandidates = fCurrTrackCond->GetMaxCandidates(ila);
    //
    for (int isd=0;isd<nSeedsUp;isd++) {
      AliITSUSeed* seedU;
      if (ilaUp==fNLrActive) {
	seedU = 0;
	seedUC.InitFromSeed(fCurrHyp->GetTPCSeed()); // init with TPC seed from ref.R
      }
      else {
	seedU = fCurrHyp->GetSeed(ilaUp,isd);  // seed on prev.active layer to prolong	
	if (seedU->IsKilled()) continue;
	seedUC = *seedU;                       // its copy will be prolonged
	seedUC.SetParent(seedU);	
      }
      seedUC.ResetFMatrix();                    // reset the matrix for propagation to next layer
      // go till next active layer
#ifdef _ITSU_DEBUG_
      AliDebug(2,Form("working on Lr:%d Seed:%d of %d for esdID=%d (MClb:%d) | pT=%.3f",ila,isd,nSeedsUp,esdID,fCurrESDtrMClb,seedUC.Pt()));
#endif
      if (!TransportToLayer(&seedUC, fITS->GetLrIDActive(ilaUp), fITS->GetLrIDActive(ila)) ) { // external seed already prolonged
	//
#ifdef _ITSU_DEBUG_
	AliDebug(2,Form("Transport failed | esdID=%d (MClb:%d)",esdID,fCurrESDtrMClb));
#endif
	// Check if the seed satisfies to track definition
	if (NeedToKill(&seedUC,kTransportFailed) && seedU) KillSeed(seedU,kTRUE);
	continue; // RS TODO: decide what to do with tracks stopped on higher layers w/o killing
      }
      if (!GetRoadWidth(&seedUC, ila)) { // failed to find road width on the layer
	if (NeedToKill(&seedUC,kRWCheckFailed) && seedU) KillSeed(seedU,kTRUE);
	continue;
      }
      /*
      //RS toremove
      int mcquest = -1;
      if (!seedUC.ContainsFake() && AliDebugLevelClass()>2) {
	mcquest = fCurrESDtrMClb;
	seedUC.Print("etp");
	printf("FSParams: "); for (int ip=0;ip<kNTrImpData;ip++) printf("%+e ",fTrImpData[ip]); printf("\n");
      }
      //
      int nsens = fCurrLayer->FindSensors(&fTrImpData[kTrPhi0], hitSens, mcquest);  // find detectors which may be hit by the track
      */
      int nsens = fCurrLayer->FindSensors(&fTrImpData[kTrPhi0], hitSens);  // find detectors which may be hit by the track
#ifdef _ITSU_DEBUG_
      AliDebug(2,Form("Will check %d sensors on lr:%d | esdID=%d (MClb:%d)",nsens,ila,esdID,fCurrESDtrMClb));
#endif
      //
      seedUC.SetLr(ila);
      //
      for (int isn=nsens;isn--;) {
	AliITSURecoSens* sens = hitSens[isn];
	int ncl = sens->GetNClusters();
	if (!ncl) continue;
	seedT = seedUC;
	//
	// We need to propagate the seed to sensor on fCurrLayer staying in the frame of the sensor from prev.layer,
	// since the transport matrix should be defined in this frame.
	double xs; // X in the TF of current seed, corresponding to intersection with sensor plane
	if (!seedT.GetTrackingXAtXAlpha(sens->GetXTF(),sens->GetPhiTF(),bz, xs)) {
#ifdef _ITSU_DEBUG_
	  if (AliDebugLevelClass()>2) {
	    printf("Failed on GetTrackingXAtXAlpha: X=%.4f alp=%.3f\n",sens->GetXTF(),sens->GetPhiTF());
	    seedT.Print("etp");
	  }
#endif
	  continue;
	}
	if (xs<seedT.GetX()) {
	  if (!PropagateSeed(&seedT,xs,fCurrMass)) continue;
	}
	else { // some low precision tracks may hit the sensor plane outside of the layer radius
#ifdef _ITSU_DEBUG_
	  if (AliDebugLevelClass()>2) {
	    if (!seedT.ContainsFake()) {
	      printf("WRONG PROP on L%d, sens%d of %d: %.4f > %.4f\n",ila,isn,nsens,xs,seedT.GetX());
	      sens->Print();
	      seedT.Print("etp");	    
	    }
	  }
#endif
	  if (!seedT.PropagateParamOnlyTo(xs,bz)) continue;
	}
	//	if (!seedT.PropagateToX(xs,bz)) continue;
	//	if (!seedT.Rotate(sens->GetPhiTF())) continue;
	if (!seedT.RotateToAlpha(sens->GetPhiTF())) continue;
	//
	int clID0 = sens->GetFirstClusterId();
	for (int icl=ncl;icl--;) { // don't change the order, clusters are sorted
	  //AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
	  int clID = clID0 + icl;
	  if (noClSharing && fCurrLayer->GetCluster(clID)->IsClusterUsed()) continue;
	  int res = CheckCluster(&seedT,ila,clID0+icl);
	  //AliLog::SetClassDebugLevel("AliITSUTrackerGlo", 0);
	  //
	  if (res==kStopSearchOnSensor) break;     // stop looking on this sensor
	  if (res==kClusterNotMatching) continue;  // cluster does not match
	  // cluster is matching and it was added to the hypotheses tree
	}
      }
      // cluster search is done. Do we need to have a version of this seed skipping current layer
      if (!NeedToKill(&seedUC,kMissingCluster)) {
	AliITSUSeed* seedSkp = NewSeedFromPool(&seedUC);
	if (nsens) { // to do: make penalty to account for probability to miss the cluster for good reason
	  double penalty = -fCurrTrackCond->GetMissPenalty(ila);
	  seedSkp->SetChi2Cl(penalty);
	}
	if (seedSkp->GetChi2GloNrm()>fCurrTrackCond->GetMaxChi2GloNrm(ila)) {
	  MarkSeedFree(seedSkp);
	}
	else AddSeedBranch(seedSkp);
      }
      // transfer the new branches of the seed to the hypothesis container
      if (fNBranchesAdded) ValidateAllowedBranches(maxNBranches);
      //
    }
    if (fNCandidatesAdded) ValidateAllowedCandidates(ila,maxNCandidates);
    //    ((TObjArray*)fCurrHyp->GetLayerSeeds(ila))->Sort();
#ifdef _ITSU_DEBUG_
    if (AliDebugLevelClass()>2) { //RS
      printf(">>> All hypotheses on lr %d: \n",ila);
      for (int ih=0;ih<fCurrHyp->GetNSeeds(ila);ih++) {printf(" #%3d ",ih); fCurrHyp->GetSeed(ila,ih)->Print();}
    }
#endif
    //
    /*
    for (int ih=0;ih<fCurrHyp->GetNSeeds(ila);ih++) {
      if (ila!=0) continue;
      double vecL[5] = {0};
      double matL[15] = {0};
      AliITSUSeed* sp = fCurrHyp->GetSeed(ila,ih);
      while(sp->GetParent()) {
	sp->Smooth(vecL,matL);
	if (sp->GetLayerID()>=fNLrActive-1) break;
	sp = (AliITSUSeed*)sp->GetParent();
      }
      */
  }
  //
  SaveReducedHypothesesTree(hypTr);
#ifdef _ITSU_DEBUG_
  if (AliDebugLevelClass()>1) {
    printf("\nSaved hypotheses for esdTrack %d (MCLab:%d)\n",esdID,fCurrESDtrMClb);
    hypTr->Print("l");
  }
#endif
  fCurrHyp = 0;
  //  AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0);
  //
}

//_________________________________________________________________________
AliITSUTrackHyp* AliITSUTrackerGlo::InitHypothesis(AliESDtrack *esdTr, Int_t esdID)
{
  // init prolongaion candidates finding for single seed
  AliITSUTrackHyp* hyp = GetTrackHyp(esdID);
  if (hyp) return hyp;
  //
  fCountProlongationTrials++;
  //
  fCurrMass = esdTr->GetMassForTracking();
  //
  hyp = new AliITSUTrackHyp(fNLrActive);
  hyp->SetESDTrack(esdTr);
  hyp->SetUniqueID(esdID);
  hyp->SetMass(fCurrMass);
  hyp->SetTPCSeed( new AliExternalTrackParam(*esdTr) );
  SetTrackHyp(hyp,esdID);
  Bool_t lut = fUseMatLUT; // for propagation from TPC use TGeo material budget always
  fUseMatLUT = kFALSE;
  if (!TransportToLayer(hyp->GetTPCSeed(),fITS->GetNLayers(), fITS->GetLrIDActive(fNLrActive-1), fITS->GetRITSTPCRef())) hyp->SetSkip(); // propagate to outer R of ITS
  fUseMatLUT = lut;
  //
  return hyp;
  // TO DO
}

//_________________________________________________________________________
Bool_t AliITSUTrackerGlo::TransportToLayer(AliITSUSeed* seed, Int_t lFrom, Int_t lTo, Double_t rLim)
{
  // transport seed from layerFrom to the entrance (along the direction of the track) of layerTo or to rLim (if>0), wathever is closer
  //  
  //
  if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
  //
  int dir = lTo > lFrom ? 1:-1;
  AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
  Bool_t checkFirst = kTRUE;
  Bool_t limReached = kFALSE;
  while(lFrom!=lTo) {
    if (lrFr) {
      if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst))	return kFALSE; // go till the end of current layer
      checkFirst = kFALSE;
    }
    AliITSURecoLayer* lrTo =  fITS->GetLayer( (lFrom+=dir) );
    if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
    //
    // go the entrance of the layer, assuming no materials in between
    double xToGo = lrTo->GetR(-dir);
    if (rLim>0) {
      if (dir>0) {
	if (rLim<xToGo) {xToGo = rLim; limReached = kTRUE;}
      }
      else {
	if (rLim>xToGo) {xToGo = rLim; limReached = kTRUE;}
      }
    }
    //    double xts = xToGo;
    if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
      //printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
      //seed->Print("etp");
      //      
      return kFALSE;
    }
#ifdef _ITSU_DEBUG_
    AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
#endif
    if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
      //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
      //seed->Print("etp");
      return kFALSE;
    }
    lrFr = lrTo;
    if (limReached) break;
  }
  return kTRUE;
  //
}

//_________________________________________________________________________
Bool_t AliITSUTrackerGlo::TransportToLayer(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo, Double_t rLim)
{
  // transport track from layerFrom to the entrance of layerTo or to rLim (if>0), wathever is closer
  //  
  if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
  //
  int dir = lTo > lFrom ? 1:-1;
  AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
  Bool_t checkFirst = kTRUE;
  Bool_t limReached = kFALSE;
  while(lFrom!=lTo) {
    if (lrFr) {
      if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
      checkFirst = kFALSE;
    }
    AliITSURecoLayer* lrTo =  fITS->GetLayer( (lFrom+=dir) );
    if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
    //
    // go the entrance of the layer, assuming no materials in between
    double xToGo = lrTo->GetR(-dir);
    if (rLim>0) {
      if (dir>0) {
	if (rLim<xToGo) {xToGo = rLim; limReached = kTRUE;}
      }
      else {
	if (rLim>xToGo) {xToGo = rLim; limReached = kTRUE;}
      }
    }
    //    double xts = xToGo;
    if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
      //      printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
      //      seed->Print("etp");
      return kFALSE;
    }
#ifdef _ITSU_DEBUG_
    AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
#endif
    if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
      //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
      //seed->Print("etp");
      return kFALSE;
    }
    lrFr = lrTo;
    if (limReached) break;
  }
  return kTRUE;
  //
}

//_________________________________________________________________________
Bool_t AliITSUTrackerGlo::TransportToLayerX(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo, Double_t xStop)
{
  // transport track from layerFrom to the entrance of layerTo but do not pass control parameter X
  //  
  if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
  //
  int dir = lTo > lFrom ? 1:-1;
  AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
  Bool_t checkFirst = kTRUE;
  while(lFrom!=lTo) {
    if (lrFr) {
      if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
      checkFirst = kFALSE;
    }
    AliITSURecoLayer* lrTo =  fITS->GetLayer( (lFrom+=dir) );
    if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
    //
    // go the entrance of the layer, assuming no materials in between
    double xToGo = lrTo->GetR(-dir); // R of the entrance to layer
    //
    //    double xts = xToGo;
    if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
      //      printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
      //      seed->Print("etp");
      return kFALSE;
    }
    if ( (dir>0&&xToGo>xStop) || (dir<0&&xToGo<xStop) ) xToGo = xStop;
    //
#ifdef _ITSU_DEBUG_
    AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
#endif
    if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
      //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
      //seed->Print("etp");
      return kFALSE;
    }
    lrFr = lrTo;
  }
  return kTRUE;
  //
}

//_________________________________________________________________________
Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
{
  // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
  // If check is requested, do this only provided the track has not exited the layer already
  double xToGo = lr->GetR(dir);
  if (check) { // do we need to track till the surface of the current layer ?
    double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
#ifdef _ITSU_DEBUG_
    AliDebug(3,Form(" dir:%d Cur: %e Tgt: %e",dir,Sqrt(curR2),xToGo));
#endif
    if      (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // on the surface or outside of the layer
    else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // on the surface or outside of the layer
  }
  if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
  // go via layer to its boundary, applying material correction.
  if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
  //
  return kTRUE;
  //
}

//_________________________________________________________________________
Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
{
  // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
  // If check is requested, do this only provided the track has not exited the layer already
  double xToGo = lr->GetR(dir);
  if (check) { // do we need to track till the surface of the current layer ?
    double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
    if      (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // on the surface or outside of the layer
    else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // on the surface or outside of the layer
  }
#ifdef _ITSU_DEBUG_
  AliDebug(2,Form(" dir=%d : from R=%.4f -> R=%.4f",dir,Sqrt(seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY()), xToGo));
#endif
  if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
  // go via layer to its boundary, applying material correction.
  if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
  //
  return kTRUE;
  //
}


//_________________________________________________________________________
Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
{
  // go to the entrance of lr in direction dir, w/o applying material corrections.
  // If check is requested, do this only provided the track did not reach the layer already
  double xToGo = lr->GetR(-dir);
  if (check) { // do we need to track till the surface of the current layer ?
    double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
    if      (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // already passed
    else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // already passed
  }
  if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
  // go via layer to its boundary, applying material correction.
  if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
  return kTRUE;
  //
}

//_________________________________________________________________________
Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
{
  // go to the entrance of lr in direction dir, w/o applying material corrections.
  // If check is requested, do this only provided the track did not reach the layer already
  double xToGo = lr->GetR(-dir);
  if (check) { // do we need to track till the surface of the current layer ?
    double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
    if      (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // already passed
    else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // already passed
  }
  if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
  // go via layer to its boundary, applying material correction.
  if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
  return kTRUE;
  //
}

//_________________________________________________________________________
Bool_t AliITSUTrackerGlo::GetRoadWidth(AliITSUSeed* seed, int ilrA)
{
  // calculate road width in terms of phi and z for the track which MUST be on the external radius of the layer
  // as well as some aux info
  double bz = GetBz();
  AliITSURecoLayer* lrA = fITS->GetLayerActive(ilrA);
  seed->GetXYZ(&fTrImpData[kTrXIn]);    // lab position at the entrance from above
  static AliExternalTrackParam sc;   // seed copy for manipulations
  sc = *seed;
  //
  fTrImpData[kTrPhiIn] = ATan2(fTrImpData[kTrYIn],fTrImpData[kTrXIn]);
  if (!sc.Rotate(fTrImpData[kTrPhiIn])) return kFALSE; // go to the frame of the entry point into the layer
  BringTo02Pi(fTrImpData[kTrPhiIn]);
  double dr  = lrA->GetDR();                              // approximate X dist at the inner radius
  if (!sc.GetXYZAt(sc.GetX()-dr, bz, fTrImpData + kTrXOut)) {
    // special case: track does not reach inner radius, might be tangential
    double r = sc.GetD(0,0,bz);
    double x;
    if (!sc.GetXatLabR(r,x,bz,-1)) return kFALSE;
    dr = Abs(sc.GetX() - x);
    if (!sc.GetXYZAt(x, bz, fTrImpData + kTrXOut)) return kFALSE;
  }
  //
  fTrImpData[kTrPhiOut] = ATan2(fTrImpData[kTrYOut],fTrImpData[kTrXOut]);
  BringTo02Pi(fTrImpData[kTrPhiOut]);
  double sgy = sc.GetSigmaY2() + dr*dr*sc.GetSigmaSnp2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(ilrA);
  double sgz = sc.GetSigmaZ2() + dr*dr*sc.GetSigmaTgl2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(ilrA);
  sgy = Sqrt(sgy)*fCurrTrackCond->GetNSigmaRoadY(ilrA);
  sgz = Sqrt(sgz)*fCurrTrackCond->GetNSigmaRoadZ(ilrA);
  double phi0  = MeanPhiSmall(fTrImpData[kTrPhiOut],fTrImpData[kTrPhiIn]);
  double dphi0 = DeltaPhiSmall(fTrImpData[kTrPhiOut],fTrImpData[kTrPhiIn]);
  //
  fTrImpData[kTrPhi0] = phi0;
  fTrImpData[kTrZ0]   = 0.5*(fTrImpData[kTrZOut]+fTrImpData[kTrZIn]);
  dphi0 += sgy/lrA->GetR();
  fTrImpData[kTrDPhi] =  dphi0<PiOver2() ? dphi0 : PiOver2();
  fTrImpData[kTrDZ]   = 0.5*Abs(fTrImpData[kTrZOut]-fTrImpData[kTrZIn])   + sgz;
  //  
  return kTRUE;
}

//________________________________________
void AliITSUTrackerGlo::ResetSeedsPool()
{
  // mark all seeds in the pool as unused
  AliDebug(-1,Form("CurrentSize: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
  fNFreeSeeds = 0;
  fSeedsPool.Clear(); // seeds don't allocate memory
}


//________________________________________
void AliITSUTrackerGlo::MarkSeedFree(AliITSUSeed *sd) 
{
  // account that this seed is "deleted" 
  int id = sd->GetPoolID();
  if (id<0) {
    AliError(Form("Freeing of seed %p NOT from the pool is requested",sd)); 
    return;
  }
  //  AliInfo(Form("%d %p",id, seed));
  fSeedsPool.RemoveAt(id);
  if (fFreeSeedsID.GetSize()<=fNFreeSeeds) fFreeSeedsID.Set( 2*fNFreeSeeds + 100 );
  fFreeSeedsID.GetArray()[fNFreeSeeds++] = id;
}

//_________________________________________________________________________
Int_t AliITSUTrackerGlo::CheckCluster(AliITSUSeed* track, Int_t lr, Int_t clID) 
{
  // Check if the cluster (in tracking frame!) is matching to track. 
  // The track must be already propagated to sensor tracking frame.
  // Returns:  kStopSearchOnSensor if the search on given sensor should be stopped, 
  //           kClusterMatching    if the cluster is matching
  //           kClusterMatching    otherwise
  //
  const double kTolerX = 5e-4;
  //  AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
  AliCluster *cl = fCurrLayer->GetCluster(clID);
  //
  Bool_t goodCl = kFALSE;
  int currLabel = Abs(fCurrESDtrMClb);
  //
  if (cl->GetLabel(0)>=0) {
    for (int i=0;i<3;i++) if (cl->GetLabel(i)>=0 && cl->GetLabel(i)==currLabel) {goodCl = kTRUE; break;}
  }
  //
  if (TMath::Abs(cl->GetX())>kTolerX) { // if due to the misalingment X is large, propagate track only
    if (!track->PropagateParamOnlyTo(track->GetX()+cl->GetX(),GetBz())) {
      //
#ifdef _ITSU_DEBUG_
      if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
	AliDebug(2,Form("Lost good cl on L:%d failed propagation. |ESDtrack#%d (MClb:%d)",lr,fCurrESDtrack->GetID(),fCurrESDtrMClb)); 
	track->Print("etp");
	cl->Print();
      }
#endif
      return kStopSearchOnSensor; // propagation failed, seedT is intact
    }
  }
  double dy = cl->GetY()-track->GetY();
  double dz = cl->GetZ()-track->GetZ();
  //
  double dy2 = dy*dy;
  double tol2 = (track->GetSigmaY2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(lr))*
    fCurrTrackCond->GetNSigmaRoadY(lr)*fCurrTrackCond->GetNSigmaRoadY(lr); // RS TOOPTIMIZE
  if (dy2>tol2) {                          // the clusters are sorted in Z(col) then in Y(row). 
    //
#ifdef _ITSU_DEBUG_
    if (AliDebugLevelClass()>2 && goodCl &&  !track->ContainsFake()) {    
      AliDebug(2,Form("Lost good cl: dy2=%e > tol2=%e |ESDtrack#%d (MClb:%d)",dy2,tol2,fCurrESDtrack->GetID(),fCurrESDtrMClb)); 
      track->Print("etp");
      cl->Print();
      PrintSeedClusters(track);
    }
#endif
    // clusters are sorted in increasing Y and the loop goes from last (highers Y) to 1st.
    // if the track has too large y for this cluster (dy<0) it will be so for all other clusters of the sensor
    if (dy<0) return kStopSearchOnSensor;  // No chance that other cluster of this sensor will match (all Y's will be even larger)
    else      return kClusterNotMatching;   // Other clusters may match
  }
  double dz2 = dz*dz;
  tol2 = (track->GetSigmaZ2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(lr))*
    fCurrTrackCond->GetNSigmaRoadZ(lr)*fCurrTrackCond->GetNSigmaRoadZ(lr); // RS TOOPTIMIZE
  if (dz2>tol2) {
#ifdef _ITSU_DEBUG_
    if (AliDebugLevelClass()>2 && goodCl &&  !track->ContainsFake()) {
      AliDebug(2,Form("Lost good cl on L:%d : dz2=%e > tol2=%e |ESDtrack#%d (MClb:%d)",lr,dz2,tol2,fCurrESDtrack->GetID(),fCurrESDtrMClb)); 
      track->Print("etp");
      cl->Print();
      PrintSeedClusters(track);
    }
#endif
    return kClusterNotMatching; // Other clusters may match
  }
  //
  // check chi2
  Double_t p[2]={cl->GetY(), cl->GetZ()};
  Double_t cov[3]={cl->GetSigmaY2(), cl->GetSigmaYZ(), cl->GetSigmaZ2()};
  double chi2 = track->GetPredictedChi2(p,cov);
  //
  if (chi2>fCurrTrackCond->GetMaxTr2ClChi2(lr)) {
#ifdef _ITSU_DEBUG_
    if (AliDebugLevelClass()>2 && goodCl &&  !track->ContainsFake()) {
      AliDebug(2,Form("Lost good cl on L:%d : Chi2=%e > Chi2Max=%e |dy: %+.3e dz: %+.3e |ESDtrack#%d (MClb:%d)\n",
		      lr,chi2,fCurrTrackCond->GetMaxTr2ClChi2(lr),dy,dz,fCurrESDtrack->GetID(),fCurrESDtrMClb)); 
      track->Print("etp");
      cl->Print("");
      PrintSeedClusters(track);
    }
#endif
    return kClusterNotMatching;
  }
  //
  track = NewSeedFromPool(track);  // input track will be reused, use its clone for updates
  if (!track->Update()) {
#ifdef _ITSU_DEBUG_
    if (AliDebugLevelClass()>2 && goodCl &&  !track->ContainsFake()) {
      AliDebug(2,Form("Lost good cluster on L:%d : failed to update",lr));
      track->Print("etp");
      cl->Print("");
      PrintSeedClusters(track);
    }
#endif
    MarkSeedFree(track);
    return kClusterNotMatching;
  }
  track->SetChi2Cl(chi2);
  track->SetLrClusterID(lr,clID);
  //
  if (track->GetChi2GloNrm()>fCurrTrackCond->GetMaxChi2GloNrm(lr)) {
#ifdef _ITSU_DEBUG_
    if (AliDebugLevelClass()>2 && goodCl &&  !track->ContainsFake()) {
      AliDebug(2,Form("Lost good cl on L:%d : Chi2Glo=%e > Chi2Max=%e |dy: %+.3e dz: %+.3e |ESDtrack#%d (MClb:%d)\n",
		      lr,track->GetChi2GloNrm(),fCurrTrackCond->GetMaxChi2GloNrm(lr),dy,dz,fCurrESDtrack->GetID(),fCurrESDtrMClb)); 
      track->Print("etp");
      cl->Print("");
      PrintSeedClusters(track);
    }
#endif
    MarkSeedFree(track);
    return kClusterNotMatching;
  }
  //  cl->IncreaseClusterUsage(); // do this only for winners
  //
  track->SetFake(!goodCl);
  //
#ifdef _ITSU_DEBUG_
  AliDebug(2,Form("AddCl(%d) Cl%d lr:%d: dY:%+8.4f dZ:%+8.4f (MC: %5d %5d %5d) |Chi2=%f(%c)| ",
	       goodCl,clID,lr,dy,dz2,cl->GetLabel(0),cl->GetLabel(1),cl->GetLabel(2), chi2, track->IsFake() ? '-':'+'));
#endif
  //
  AddSeedBranch(track);
  //
  return kClusterMatching;
}

//_________________________________________________________________________
Bool_t AliITSUTrackerGlo::NeedToKill(AliITSUSeed *seed, Int_t flag)
{
  // check if the seed should not be discarded
  const UShort_t kMask = 0xffff;
  if (flag==kMissingCluster) {
    int lastChecked = seed->GetLayerID();
    UShort_t patt = seed->GetHitsPattern();
    if (lastChecked) patt |= (~(kMask<<lastChecked))&fCurrTrackCond->GetAllowedLayers(); // not all layers were checked, complete unchecked ones by potential hits
    Bool_t seedOK = fCurrTrackCond->CheckPattern(patt);
    return !seedOK;
  }
  return kTRUE;
}

//______________________________________________________________________________
Bool_t AliITSUTrackerGlo::PropagateSeed(AliITSUSeed *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr) 
{
  // propagate seed to given x applying material correction if requested
  const Double_t kEpsilon = 1e-5;
  Double_t xpos     = seed->GetX();
  Int_t dir         = (xpos<xToGo) ? 1:-1;
  Double_t xyz0[3],xyz1[3];
  //
  Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
  if (matCorr || updTime) seed->GetXYZ(xyz1);   //starting global position
  while ( (xToGo-xpos)*dir > kEpsilon){
    Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
    Double_t x    = xpos+step;
    Double_t bz=GetBz();   // getting the local Bz
    if (!seed->PropagateToX(x,bz)) return kFALSE;
    double ds = 0;
    if (matCorr || updTime) {
      xyz0[0]=xyz1[0]; // global pos at the beginning of step
      xyz0[1]=xyz1[1];
      xyz0[2]=xyz1[2];
      seed->GetXYZ(xyz1);    //  // global pos at the end of step
      if (matCorr) {
	Double_t xrho,xx0;
	ds = GetMaterialBudget(xyz0,xyz1,xx0,xrho);	
	if (dir>0) xrho = -xrho; // outward should be negative
	if (!seed->ApplyMaterialCorrection(xx0,xrho,mass,kFALSE)) {AliInfo("Fail2"); seed->Print("etp"); return kFALSE;}
      }
       else { // matCorr is not requested but time integral is
	double d0 = xyz1[0]-xyz0[0];
	double d1 = xyz1[1]-xyz0[1];
	double d2 = xyz1[2]-xyz0[2];	
	ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
      }     
    }
    if (updTime) seed->AddTimeStep(ds);
    xpos = seed->GetX();
  }
  return kTRUE;
}

//______________________________________________________________________________
Bool_t AliITSUTrackerGlo::PropagateSeed(AliExternalTrackParam *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr) 
{
  // propagate seed to given x applying material correction if requested
  const Double_t kEpsilon = 1e-5;
  Double_t xpos     = seed->GetX();
  Int_t dir         = (xpos<xToGo) ? 1:-1;
  Double_t xyz0[3],xyz1[3];
  //
#ifdef _ITSU_DEBUG_
  if (AliDebugLevelClass()>1) {
    AliDebug(2,Form("Before Propagate to X=%f with M=%.3f MaxStep=%.4f MatCorr=%d",xToGo,mass,maxStep,matCorr));
    seed->AliExternalTrackParam::Print();
  }
#endif
  Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
  if (matCorr || updTime) seed->GetXYZ(xyz1);   //starting global position
  while ( (xToGo-xpos)*dir > kEpsilon){
    Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
    Double_t x    = xpos+step;
    Double_t bz=GetBz();   // getting the local Bz
    if (!seed->PropagateTo(x,bz))  return kFALSE;
    double ds = 0;
    if (matCorr || updTime) {
      xyz0[0]=xyz1[0]; // global pos at the beginning of step
      xyz0[1]=xyz1[1];
      xyz0[2]=xyz1[2];
      seed->GetXYZ(xyz1);    //  // global pos at the end of step
      //
      if (matCorr) {
	Double_t xrho,xx0;
	ds = GetMaterialBudget(xyz0,xyz1,xx0,xrho);	
	if (dir>0) xrho = -xrho; // outward should be negative
	if (!seed->CorrectForMeanMaterial(xx0,xrho,mass)) return kFALSE;
      }
      else { // matCorr is not requested but time integral is
	double d0 = xyz1[0]-xyz0[0];
	double d1 = xyz1[1]-xyz0[1];
	double d2 = xyz1[2]-xyz0[2];	
	ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
      }
    }
    if (updTime) seed->AddTimeStep(ds);
    //
    xpos = seed->GetX();
  }
#ifdef _ITSU_DEBUG_
  if (AliDebugLevelClass()>1) {
    AliDebug(2,Form("After Propagate to X=%f",xToGo));
    seed->AliExternalTrackParam::Print();
  }
#endif
  return kTRUE;
}

//______________________________________________________________________________
Bool_t AliITSUTrackerGlo::FinalizeHypothesis(AliITSUTrackHyp* hyp)
{
  // finalize single TPC track hypothesis
  if (!hyp || hyp->GetESDTrack()->IsOn(AliESDtrack::kITSin)) return kFALSE;
  AliITSUSeed* winner = 0;
  fCurrHyp = hyp;
  fCurrMass = hyp->GetMass();
  if (!(winner=hyp->DefineWinner(fCurrTrackCond->GetActiveLrInner()))) return kFALSE;
  CookMCLabel(hyp);
  //
#ifdef  _ITSU_TUNING_MODE_  // fill tuning histos
  TObjArray* dest = hyp->GetLabel()>=0 ? fCHistoArrCorr : fCHistoArrFake;
  if (dest) {
    AliITSUSeed* sd = winner;
    double htrPt = hyp->Pt();
    do {
      int clID,lrID;
      if ( (clID=sd->GetLrCluster(lrID))<0 ) continue;
      ((TH2*)dest->At( GetHistoID(lrID,kHChi2Nrm     ,fCurrPassID) ))->Fill(htrPt,sd->GetChi2GloNrm());
      ((TH2*)dest->At( GetHistoID(lrID,kHBestInBranch,fCurrPassID) ))->Fill(htrPt,sd->GetOrdBranch());
      ((TH2*)dest->At( GetHistoID(lrID,kHBestInCand  ,fCurrPassID) ))->Fill(htrPt,sd->GetOrdCand());
      //
      if (dest==fCHistoArrFake && !sd->IsFake()) continue; // for the fake seeds fill only fake clusters part
      //
      ((TH2*)dest->At( GetHistoID(lrID,kHResY        ,fCurrPassID) ))->Fill(htrPt,sd->GetResidY());
      ((TH2*)dest->At( GetHistoID(lrID,kHResZ        ,fCurrPassID) ))->Fill(htrPt,sd->GetResidZ());
      ((TH2*)dest->At( GetHistoID(lrID,kHResYP       ,fCurrPassID) ))->Fill(htrPt,sd->GetPullY());
      ((TH2*)dest->At( GetHistoID(lrID,kHResZP       ,fCurrPassID) ))->Fill(htrPt,sd->GetPullZ());
      ((TH2*)dest->At( GetHistoID(lrID,kHChi2Cl      ,fCurrPassID) ))->Fill(htrPt,sd->GetChi2Cl());
      //
    } while((sd=(AliITSUSeed*)sd->GetParent()));
    //
    ((TH2*)dest->At( GetHistoID(-1,kHChiMatch,fCurrPassID) ))->Fill(htrPt,winner->GetChi2ITSTPC());
    ((TH2*)dest->At( GetHistoID(-1,kHChiITSSA,fCurrPassID) ))->Fill(htrPt,winner->GetChi2ITSSA());
  }
  //
#endif
  //
  CheckClusterSharingConflicts(hyp);
  return hyp->GetWinner();  // winner might change of disappear after resolving conflicts
  //
}

//______________________________________________________________________________
void AliITSUTrackerGlo::FinalizeHypotheses()
{
  // select winner for each hypothesis, remove cl. sharing conflicts
  AliCodeTimerAuto("",0); 
  //
  int* esdTrackIndex = fESDIndex.GetArray();
  for (int ih=0;ih<fNTracksESD;ih++) FinalizeHypothesis(GetTrackHyp(esdTrackIndex[ih])); // finlize in decreasing pt order
  //
  for (int ih=fNTracksESD;ih--;) UpdateESDTrack(GetTrackHyp(esdTrackIndex[ih]),AliESDtrack::kITSin); 
  //
}

//______________________________________________________________________________
void AliITSUTrackerGlo::CheckClusterSharingConflicts(AliITSUTrackHyp* hyp)
{
  // remove eventual cluster sharing conflicts
  AliITSUSeed* winner = 0;
  if (!(winner=hyp->GetWinner())) return;
  UShort_t idH = (UShort_t)hyp->GetUniqueID()+1;
  int lrID,clID;
  int inLr = fCurrTrackCond->GetActiveLrInner();
  AliITSUSeed* winSD=winner;
  do {
    if ( (clID=winSD->GetLrCluster(lrID))<0 ) continue;
    AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
    Int_t refID = cl->GetRecoInfo();                // was it already referred by some track (with refID-1 if >0)
    if (!refID) {cl->SetRecoInfo(idH); continue;}   // if not, refer from track to cluster
    if (refID==idH) continue;                          // ignore reference to itself
    // the cluster is already used by other track, need to resolve the conflict
    AliITSUTrackHyp* hypC = GetTrackHyp(refID-1);   // competitor
    //    printf("Ref to %d (%p) from %d Cl %d %d\n",refID-1,hypC, idH-1,clID,lrID);
    AliITSUSeed* winnerC = hypC->GetWinner();
    if (!winnerC) {
      printf("Missing winner, ERROR\n");
      continue;
    }
    //
#ifdef _ITSU_DEBUG_
    if (AliDebugLevelClass()>1) { //------------------------ DEBUG -------------------
      AliInfo(Form("Shared cl#%4d on Lr:%d: Hyp%3d/Q:%6.3f vs Hyp%3d/Q:%6.3f",
		   clID,lrID,idH-1,winner->GetQualityVar(),refID-1,winnerC->GetQualityVar()));
      // dump winner of hyp
      printf("#%4d Pt:%.3f Chi:%6.3f/%6.3f/%6.3f Ncl:%d MCits%+5d MCtpc:%+5d |",
	     idH-1,winner->Pt(),winner->GetChi2GloNrm(),winner->GetChi2ITSTPC(),winner->GetChi2ITSSA(),
	     winner->GetNLayersHit(),hyp->GetITSLabel(),hyp->GetESDTrack()->GetTPCLabel());
      int prevL=-1;
      AliITSUSeed* sd = winner;
      do {
	int lrs;
	int clIDt = sd->GetLrCluster(lrs);
	if (clIDt<0) continue;
	while( lrs>++prevL ) printf("%4s        ","----");
	printf("%4d (%5.1f)",clIDt,sd->GetChi2Cl());
      } while ((sd=(AliITSUSeed*)sd->GetParent()));
      printf("|\n");
      // dump winner of hypC
      printf("#%4d Pt:%.3f Chi:%6.3f/%6.3f/%6.3f Ncl:%d MCits%+5d MCtpc:%+5d |",
	     refID-1,winnerC->Pt(),winnerC->GetChi2GloNrm(),winnerC->GetChi2ITSTPC(),winnerC->GetChi2ITSSA(),
	     winnerC->GetNLayersHit(),hypC->GetITSLabel(),hypC->GetESDTrack()->GetTPCLabel());
      prevL=-1;
      sd = winnerC;
      do {
	int lrs;
	int clIDt = sd->GetLrCluster(lrs);
	if (clIDt<0) continue;
	while( lrs>++prevL ) printf("%4s        ","----");
	printf("%4d (%5.1f)",clIDt,sd->GetChi2Cl());
      } while ((sd=(AliITSUSeed*)sd->GetParent()));
      printf("|\n");
    }  //-------------------------------------------------- DEBUG -------------------
#endif
    //
    if (winner->GetQualityVar()<winnerC->GetQualityVar()) { // current tracks is better then competitor track
      FlagSeedClusters(winnerC,kFALSE,refID); // unflag cluster usage by loser
      cl->SetRecoInfo(idH);
      //
      //if ( hypC->GetLabel()>=0) AliInfo(Form("KILLING CORRECT: %d",refID-1));
      winnerC->Kill();
      if (hypC->DefineWinner(inLr)) {             // find new winner instead of suppressed candidate
	CookMCLabel(hypC);
	CheckClusterSharingConflicts(hypC);   // and check its sharing conflicts again
	if (winner->IsKilled()) break;        // the current winner might have been killed during check of new winner of hypC!
      }
    }    
    else { // competitor hypC is better than the hyp
      FlagSeedClusters(winner,kFALSE,idH); // unflag cluster usage by loser
      winner->Kill();
      //if ( hyp->GetLabel()>=0) AliInfo(Form("KILLING CORRECT: %d",idH-1));
      if (hyp->DefineWinner(inLr)) {
	CookMCLabel(hyp);
	CheckClusterSharingConflicts(hyp);
      }
      break;
    }
    //   
  } while ((winSD=(AliITSUSeed*)winSD->GetParent()));
  //
  return;
}

//______________________________________________________________________________
void AliITSUTrackerGlo::UpdateESDTrack(AliITSUTrackHyp* hyp,Int_t flag)
{
  // update ESD track with current best hypothesis
  if (!hyp) return;
  AliESDtrack* esdTr = hyp->GetESDTrack();
  if (!esdTr || esdTr->IsOn(flag)) return;
  AliITSUSeed* win = hyp->GetWinner();
  if (!win || win->IsKilled()) return;
  double chiSav;
  //
  switch (flag) {
  case AliESDtrack::kITSin: 
    esdTr->UpdateTrackParams(hyp,flag); // update kinematics
    fCountITSin++;
    // set fakes cluster info
    {
      UShort_t clfake = 0;
      Int_t    clSplit = 0;
      AliITSUSeed* sd = win;
      int ip=0;
      do {
	int lr, clID = sd->GetLrCluster(lr);
	if (sd->IsFake()) clfake |= 0x1<<lr;
	if (clID>=0) {
	  esdTr->SetITSModuleIndex(ip++, sd->GetLrClusterID());
	  AliITSUClusterPix *cl = (AliITSUClusterPix*)fITS->GetLayerActive(lr)->GetCluster(clID);
#ifdef  _ITSU_TUNING_MODE_
	  if (cl->IsSplit()) clSplit |= 0x1<<lr;
#endif
	}
      } while ((sd=(AliITSUSeed*)sd->GetParent()));
      //
      // RS: Temporary set special flag for tracks from the afterburner
      if (fCurrPassID>0) clfake |= 0x1<<7;
      //
      esdTr->SetITSSharedMap(clfake);
      esdTr->SetITSModuleIndex(10,clSplit);
    }
    // TEMPORARY: store iteration id
    esdTr->SetITSModuleIndex(11,fCurrPassID);
    break;
    //
  case AliESDtrack::kITSout: 
    // here the stored chi2 will correspond to backward ITS-SA fit
    esdTr->UpdateTrackParams(hyp,flag); // update kinematics
    fCountITSout++;
    // TODO: avoid setting friend
    break;
    //
  case AliESDtrack::kITSrefit: 

    // at the moment fill as a chi2 the TPC-ITS matching chi2
    chiSav = hyp->GetChi2();
    hyp->SetChi2(win->GetChi2ITSTPC());
    esdTr->UpdateTrackParams(hyp,flag); // update kinematics
    hyp->SetChi2(chiSav);
    fCountITSrefit++;
    // TODO: avoid setting cluster info
    break;
  default:
    AliFatal(Form("Unknown flag %d",flag));
  }
  //
  esdTr->SetITSLabel(hyp->GetITSLabel());
  // transfer chip indices
  // TODO
}

//______________________________________________________________________________
Double_t AliITSUTrackerGlo::RefitTrack(AliITSUTrackHyp* trc, Double_t rDest, Int_t &nclFit, Int_t stopCond)
{
  // refit track till radius rDest. The cluster,mass info is taken from fCurrHyp (and its winner seed)
  // if stopCond<0 : propagate till last cluster then stop
  // if stopCond==0: propagate till last cluster then try to go till limiting rDest, don't mind if fail
  // if stopCond>0 : rDest must be reached
  //
  double rCurr = Sqrt(trc->GetX()*trc->GetX() + trc->GetY()*trc->GetY());
  int dir,lrStart,lrStop;
  //
  dir = rCurr<rDest ? 1 : -1;
  lrStart = fITS->FindFirstLayerID(rCurr,dir);
  lrStop  = fITS->FindLastLayerID(rDest,dir); // lr id before which we have to stop
  //
#ifdef _ITSU_DEBUG_
  if (AliDebugLevelClass()>2) {
    printf("Refit %d: Lr:%d (%f) -> Lr:%d (%f)\n",dir,lrStart,rCurr, lrStop,rDest);
    printf("Before refit: "); trc->Print();
  }
#endif
  if (lrStop<0 || lrStart<0) AliFatal(Form("Failed to find start(%d) or last(%d) layers. Track from %.3f to %.3f",lrStart,lrStop,rCurr,rDest));
  //
  Int_t clInfo[2*AliITSUAux::kMaxLayers];
  Int_t nCl = fCurrHyp->FetchClusterInfo(clInfo);
  fCurrMass = fCurrHyp->GetMass();
  AliITSUTrackHyp tmpTr(*(AliKalmanTrack*)trc);
  double chi2 = 0;
  int iclLr[2],nclLr;
  nclFit = 0;
  //
  int lrStop1 = lrStop+dir;
  for (int ilr=lrStart;ilr!=lrStop1;ilr+=dir) {
    AliITSURecoLayer* lr = fITS->GetLayer(ilr);
    if ( dir*(rCurr-lr->GetR(dir))>0) continue; // this layer is already passed
    int ilrA2,ilrA = lr->GetActiveID();
    // passive layer or active w/o hits will be traversed on the way to next cluster
    if (!lr->IsActive() || clInfo[ilrA2=(ilrA<<1)]<0) continue; 
    //
    // select the order in which possible 2 clusters (in case of the overlap) will be traversed and fitted
    nclLr=0;
    if (dir>0) { // clusters are stored in increasing radius order
      iclLr[nclLr++]=clInfo[ilrA2++];
      if (clInfo[ilrA2]>=0) iclLr[nclLr++]=clInfo[ilrA2];
    }
    else {
      if ( clInfo[ilrA2+1]>=0 ) iclLr[nclLr++]=clInfo[ilrA2+1];
      iclLr[nclLr++]=clInfo[ilrA2];
    }
    //
    Bool_t transportedToLayer = kFALSE;
    for (int icl=0;icl<nclLr;icl++) {
      AliITSUClusterPix* clus =  (AliITSUClusterPix*)lr->GetCluster(iclLr[icl]);
      AliITSURecoSens* sens = lr->GetSensorFromID(clus->GetVolumeId());
      if (!tmpTr.Rotate(sens->GetPhiTF())) {
#ifdef _ITSU_DEBUG_
	AliESDtrack* trESD = fCurrHyp->GetESDTrack();
	AliDebug(2,Form("Failed on rotate to %f | ESDtrack#%d (MClb:%d)",sens->GetPhiTF(),trESD->GetID(),trESD->GetTPCLabel()));
#endif
	return -1;
      }
      //
      double xClus = sens->GetXTF()+clus->GetX();
      if (!transportedToLayer) {
	if (ilr!=lrStart && !TransportToLayerX(&tmpTr,lrStart,ilr,xClus)) {
#ifdef _ITSU_DEBUG_
	  AliESDtrack* trESD = fCurrHyp->GetESDTrack();
	  AliDebug(2,Form("Failed to transport %d -> %d | ESDtrack#%d (MClb:%d)\n",lrStart,ilr,trESD->GetID(),trESD->GetTPCLabel()));
#endif
	  return -1; // go to the entrance to the layer
	}
	lrStart = ilr;
	transportedToLayer = kTRUE;
      }
      //
#ifdef _ITSU_DEBUG_
      if (AliDebugLevelClass()>1) {
	AliDebug(2,Form("Propagate to cl:%d on lr %d Need to go %.4f -> %.4f",icl,ilrA,tmpTr.GetX(),xClus));
      }
#endif
      //
      if (!PropagateSeed(&tmpTr,xClus,fCurrMass)) {
#ifdef _ITSU_DEBUG_
	AliESDtrack* trESD = fCurrHyp->GetESDTrack();
	AliDebug(2,Form("Failed on propagate to %f (dir=%d) | ESDtrack#%d (MClb:%d)",xClus,dir,trESD->GetID(),trESD->GetTPCLabel()));
#endif
	return -1;
      }
      //
      Double_t p[2]={clus->GetY(), clus->GetZ()};
      Double_t cov[3]={clus->GetSigmaY2(), clus->GetSigmaYZ(), clus->GetSigmaZ2()};
      double chi2cl = tmpTr.GetPredictedChi2(p,cov);
      chi2 += chi2cl;
      //
#ifdef  _ITSU_TUNING_MODE_
      TObjArray* dest = trc->GetLabel()>=0 ? fCHistoArrCorr : fCHistoArrFake;
      if (dest && fTrackPhaseID>kClus2Tracks) {
	//
	double htrPt = tmpTr.Pt();
	double dy = p[0]-tmpTr.GetY();
	double dz = p[1]-tmpTr.GetZ();
	((TH2*)dest->At( GetHistoID(ilrA,kHResY,0,fTrackPhaseID) ))->Fill(htrPt,dy);
	((TH2*)dest->At( GetHistoID(ilrA,kHResZ,0,fTrackPhaseID) ))->Fill(htrPt,dz);
	double errY = tmpTr.GetSigmaY2();
	double errZ = tmpTr.GetSigmaZ2();
	if (errY>0) ((TH2*)dest->At( GetHistoID(ilrA,kHResYP,0,fTrackPhaseID) ))->Fill(htrPt,dy/Sqrt(errY));
	if (errZ>0) ((TH2*)dest->At( GetHistoID(ilrA,kHResZP,0,fTrackPhaseID) ))->Fill(htrPt,dz/Sqrt(errZ));
	((TH2*)dest->At( GetHistoID(ilrA,kHChi2Cl,0,fTrackPhaseID) ))->Fill(htrPt,chi2cl);
      }
#endif  
      //      
      if ( !tmpTr.Update(p,cov) ) {
#ifdef _ITSU_DEBUG_
	AliESDtrack* trESD = fCurrHyp->GetESDTrack();
	AliDebug(2,Form("Failed on Update | ESDtrack#%d (MClb:%d)",trESD->GetID(),trESD->GetTPCLabel()));
#endif
	return -1;
      }
#ifdef _ITSU_DEBUG_
      if (AliDebugLevelClass()>1) {
	printf("AfterRefit: "); tmpTr.Print();
      }
#endif
      if (++nclFit==nCl && stopCond<0) {
	*trc = (AliKalmanTrack&)tmpTr;
	return chi2; // it was requested to not propagate after last update
      }
    }
    //
  }
  // All clusters were succesfully fitted. Even if the track does not reach rDest, this is enough to validate it.
  // Still, try to go as close as possible to rDest.
  //
  if (lrStart!=lrStop) {
    if (!TransportToLayer(&tmpTr,lrStart,lrStop)) {
#ifdef _ITSU_DEBUG_
      AliDebug(1,Form("Failed on TransportToLayer %d->%d | ESDtrack#%d (MClb:%d), X=%f",lrStart,lrStop,fCurrESDtrack->GetID(),fCurrESDtrMClb,tmpTr.GetX()));
#endif
      return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
    }    
    if (!GoToExitFromLayer(&tmpTr,fITS->GetLayer(lrStop),dir)) {
#ifdef _ITSU_DEBUG_
      AliDebug(1,Form("Failed on GoToExitFromLayer %d | ESDtrack#%d (MClb:%d), X=%f",lrStop,fCurrESDtrack->GetID(),fCurrESDtrMClb,tmpTr.GetX()));
#endif
      return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
    }
  }
  // go to the destination radius. Note that here we don't select direction to avoid precision problems
  if (!tmpTr.GetXatLabR(rDest,rDest,GetBz(),0) || !PropagateSeed(&tmpTr,rDest,fCurrMass, 100, kFALSE)) {
    return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
  }
  *trc = (AliKalmanTrack&)tmpTr;
#ifdef _ITSU_DEBUG_
  if (AliDebugLevelClass()>2) {
    printf("After refit (now at lr %d): ",lrStop); trc->Print();
  }
#endif
  return chi2;
}

//__________________________________________________________________
void AliITSUTrackerGlo::CookMCLabel(AliITSUTrackHyp* hyp) 
{
  // build MC label
  //
  const int kMaxLbPerCl = 3;
  int lbID[kMaxLayers*6],lbStat[kMaxLayers*6];
  Int_t lr,nLab=0,nCl=0;
  AliITSUSeed *seed = hyp->GetWinner();
  while(seed) {
    int clID = seed->GetLrCluster(lr);
    if (clID>=0) {
      AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
      nCl++;
      for (int imc=0;imc<kMaxLbPerCl;imc++) { // labels within single cluster
	int trLb = cl->GetLabel(imc);
	if (trLb<0) break;
	// search this mc track in already accounted ones
	int iLab;
	for (iLab=0;iLab<nLab;iLab++) if (lbID[iLab]==trLb) break;
	if (iLab<nLab) lbStat[iLab]++;
	else {
	  lbID[nLab] = trLb;
	  lbStat[nLab++] = 1;
	}
      } // loop over given cluster's labels
    }
    seed = (AliITSUSeed*)seed->GetParent();
  } // loop over clusters
  // 
  AliESDtrack* esdTr = hyp->GetESDTrack();
  int tpcLab = esdTr ? Abs(esdTr->GetTPCLabel()) : -kDummyLabel;
  if (nCl && nLab) {
    int maxLab=0,nTPCok=0;
    for (int ilb=nLab;ilb--;) {
      int st = lbStat[ilb];
      if (lbStat[maxLab]<st) maxLab=ilb;
      if (tpcLab==lbID[ilb]) nTPCok += st;
    }
    hyp->SetFakeRatio(1.-float(nTPCok)/nCl);
    hyp->SetLabel( nTPCok==nCl ? tpcLab : -tpcLab);
    hyp->SetITSLabel( lbStat[maxLab]==nCl ? lbID[maxLab] : -lbID[maxLab]); // winner label
    return;
  }
  //
  hyp->SetFakeRatio(-1.);
  hyp->SetLabel( -tpcLab );
  hyp->SetITSLabel( kDummyLabel );
}

//__________________________________________________________________
Bool_t AliITSUTrackerGlo::AddSeedBranch(AliITSUSeed* seed)
{
  // add new prolongation branch for the seed to temporary array. The seeds are sorted according to Compare f-n
  if (fNCandidatesAdded+fNBranchesAdded>=fLayerMaxCandidates ) {
    AliInfo(Form("Number of candidates at layer %d for current seed reached %d, increasing buffer",fCurrActLrID,fLayerMaxCandidates));
    fLayerMaxCandidates = 2*(fLayerMaxCandidates+1);
    AliITSUSeed** tmpArr = fLayerCandidates;
    fLayerCandidates = new AliITSUSeed*[fLayerMaxCandidates];
    memcpy(fLayerCandidates,tmpArr,(fNCandidatesAdded+fNBranchesAdded)*sizeof(AliITSUSeed*));
    delete[] tmpArr; // delete only array, not objects
  }
  AliITSUSeed** branches = &fLayerCandidates[fNCandidatesAdded]; // note: fNCandidatesAdded is incremented after adding all branches of current seed
  int slot=fNBranchesAdded++;
  for (int slotF=slot;slotF--;) { // slotF is always slot-1
    AliITSUSeed* si = branches[slotF];
    if (si->Compare(seed)<0) break; // found the last seed with better quality
    // otherwise, shift the worse seed to the next slot
    branches[slot] = si;
    slot = slotF; // slot should be slotF+1
  }
  // if needed, move worse seeds
  branches[slot] = seed;
  return kTRUE;
  //
}

//__________________________________________________________________
void AliITSUTrackerGlo::ValidateAllowedBranches(Int_t acceptMax)
{
  // keep allowed number of branches for current seed and disable extras
  AliCodeTimerAuto("",0); 
  int nb = Min(fNBranchesAdded,acceptMax);
  //  if (nb<fNBranchesAdded) printf("ValidateAllowedBranches: %d of %d (%d) on lr %d\n",nb,fNBranchesAdded,fNCandidatesAdded,ilr);
  // disable unused branches
  AliITSUSeed** branches = &fLayerCandidates[fNCandidatesAdded];
  //
#ifdef  _ITSU_TUNING_MODE_
  for (int ib=0;ib<fNBranchesAdded;ib++) branches[ib]->SetOrdBranch(ib);
#endif
  //
  for (int ib=nb;ib<fNBranchesAdded;ib++) MarkSeedFree(branches[ib]);
  fNCandidatesAdded += nb; // update total candidates counter
  fNBranchesAdded = 0; // reset branches counter
  //
}

//__________________________________________________________________
void AliITSUTrackerGlo::ValidateAllowedCandidates(Int_t ilr, Int_t acceptMax)
{
  // transfer allowed number of branches to hypothesis container
  //
  AliCodeTimerAuto("",0); 
  // sort candidates in increasing order of chi2
  static int lastSize = 0;
  static int *index = 0;
  static float *chi2 = 0;
  if (fLayerMaxCandidates>lastSize) {
    lastSize = fLayerMaxCandidates;
    delete[] index;
    delete[] chi2;
    index = new int[lastSize];
    chi2  = new float[lastSize];
  }
  for (int i=fNCandidatesAdded;i--;) chi2[i] = fLayerCandidates[i]->GetChi2GloNrm();
  Sort(fNCandidatesAdded,chi2,index,kFALSE);
  //
  int nacc=0,nb=0;
  if (ilr>fCurrTrackCond->GetActiveLrInner()) { // just take 1st acceptMax candidates
    nb = Min(fNCandidatesAdded,acceptMax);
    for (int ib=0;ib<nb;ib++) AddProlongationHypothesis(fLayerCandidates[index[ib]],ilr);
  }
  else { // for innermost layer test ITS SA backward chi2 and matching to TPC
    AliITSUSeed* wn0 = fCurrHyp->GetWinner(); // in principle, should be NULL
    for (nb=0;nb<fNCandidatesAdded;nb++) {
      AliITSUSeed* sd = fLayerCandidates[index[nb]];    
      fCurrHyp->SetWinner(sd);
      if (!CheckBackwardMatching(sd)) MarkSeedFree(sd); // discard bad backward tracks
      else {
	AddProlongationHypothesis(sd,ilr);
	if (++nacc==acceptMax) {nb++; break;} // the rest will be discarded
      }
    }
    fCurrHyp->SetWinner(wn0); // restore original winner (NULL?)
  }
  // disable unused candidates
  for (int ib=nb;ib<fNCandidatesAdded;ib++) MarkSeedFree(fLayerCandidates[index[ib]]);    
  fNCandidatesAdded = 0; // reset candidates counter
  //
}

//__________________________________________________________________
Bool_t AliITSUTrackerGlo::CheckBackwardMatching(AliITSUSeed* seed)
{
  // check seed backward propagation chi2 and matching to TPC 
  double bz0 = GetBz();
  double rDest = fITS->GetRITSTPCRef(); // reference radius for comparison
  AliITSUTrackHyp trback;
  trback.AliExternalTrackParam::operator=(*seed);
  trback.ResetCovariance(10000);
  int nclFit = 0;
  double chi2sa = RefitTrack(&trback,rDest,nclFit,1);
  if (chi2sa<0) return kFALSE;
  int ndf = nclFit*2-5;
  if (ndf>0) chi2sa /= ndf;
  if (chi2sa>fCurrTrackCond->GetMaxITSSAChi2(nclFit)) return kFALSE;
  //
  // relate to TPC track at outer layer
  AliExternalTrackParam* tpcSeed = fCurrHyp->GetTPCSeed();
  if (!trback.Rotate(tpcSeed->GetAlpha()) || !trback.PropagateParamOnlyTo(tpcSeed->GetX(),bz0)) {
#ifdef _ITSU_DEBUG_
    if (AliDebugLevelClass()>+1 && !seed->ContainsFake()) {
      AliInfo(Form("Failed to propagate seed to TPC track @ X:%.3f Alpha:%.3f",tpcSeed->GetX(),tpcSeed->GetAlpha()));
      seed->Print("etp");
      trback.Print();
    }
#endif
    return kFALSE; 
  }
  double chi2Match = trback.GetPredictedChi2(tpcSeed)/5;
  if (chi2Match>fCurrTrackCond->GetMaxITSTPCMatchChi2()) return kFALSE;
  //
  seed->SetChi2ITSSA(chi2sa);
  seed->SetChi2ITSTPC(chi2Match);
  return kTRUE;
}

//__________________________________________________________________
void AliITSUTrackerGlo::SaveReducedHypothesesTree(AliITSUTrackHyp* dest)
{
  // remove those hypothesis seeds which dont lead to candidates at final layer
  //
  // 1st, flag the seeds to save
  int lr0 = fCurrTrackCond->GetActiveLrInner();
  for (int isd=0;isd<fCurrHyp->GetNSeeds(lr0);isd++) {
    AliITSUSeed* seed = fCurrHyp->RemoveSeed(lr0,isd);
    if (!seed) continue;
    seed->FlagTree(AliITSUSeed::kSave);
    dest->AddSeed(seed,lr0);
  }
  for (int ilr=0;ilr<fNLrActive;ilr++) {
    int nsd = fCurrHyp->GetNSeeds(ilr);
    for (int isd=0;isd<nsd;isd++) {
      AliITSUSeed* seed = fCurrHyp->RemoveSeed(ilr,isd);
      if (!seed) continue; // already discarded or saved
      if (seed->IsSaved()) dest->AddSeed(seed,ilr);
      else MarkSeedFree(seed);
    }
  }
  //
  //  AliInfo(Form("SeedsPool: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
}

//__________________________________________________________________
void AliITSUTrackerGlo::CleanHypothesis(AliITSUTrackHyp* hyp)
{
  // clean hypothesis
  hyp->SetWinner(0);
  for (int ilr=0;ilr<fNLrActive;ilr++) {
    int nsd = hyp->GetNSeeds(ilr);
    for (int isd=0;isd<nsd;isd++) {
      AliITSUSeed* seed = hyp->RemoveSeed(ilr,isd);
      if (!seed) continue; // already discarded or saved
      MarkSeedFree(seed);
    }
  }
}

//__________________________________________________________________
void AliITSUTrackerGlo::FlagSplitClusters()
{
  // set special bit on split clusters using MC info
  for (int ilr=fNLrActive;ilr--;) {
    int nsplit=0;
    AliITSURecoLayer* lr = fITS->GetLayerActive(ilr);
    for (int isn=lr->GetNSensors();isn--;) {
      AliITSURecoSens* sens = lr->GetSensor(isn);
      int nCl = sens->GetNClusters();
      if (!nCl) continue;
      int cl0 = sens->GetFirstClusterId();
      for (int icl=nCl;icl--;) {
	AliITSUClusterPix *cl = (AliITSUClusterPix*)lr->GetCluster(cl0+icl);
	for (int icl1=icl;icl1--;) {
	  AliITSUClusterPix *cl1 = (AliITSUClusterPix*)lr->GetCluster(cl0+icl1);
	  if (cl->HasCommonTrack(cl1)) {
	    if (!cl->IsSplit())  nsplit++;
	    if (!cl1->IsSplit()) nsplit++;
	    cl->SetSplit();
	    cl1->SetSplit();
	  }
	}
      }
    }
    AliInfo(Form("%4d out of %4d clusters are split",nsplit,lr->GetNClusters()));
  }
  //
}

//__________________________________________________________________
Bool_t AliITSUTrackerGlo::ContainsSplitCluster(const AliITSUSeed* seed, Int_t maxSize)
{
  // check if the seed contains split cluster with size < maxSize
  int lrID,clID;
  if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
    AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
    if (cl->IsSplit() && cl->GetNPix()<maxSize ) return kTRUE;
  }
  return seed->GetParent() ? ContainsSplitCluster((AliITSUSeed*)seed->GetParent(),maxSize) : kFALSE;
  //
}

//__________________________________________________________________
void AliITSUTrackerGlo::PrintSeedClusters(const AliITSUSeed* seed, Option_t* option)
{
  // print seeds clusters
  int lrID,clID;
  if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
    AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
    cl->Print(option);
  }
  if (seed->GetParent()) PrintSeedClusters((AliITSUSeed*)seed->GetParent(), option);
  //
}

//__________________________________________________________________
void AliITSUTrackerGlo::FlagSeedClusters(const AliITSUSeed* seed, Bool_t flg, UShort_t ref)
{
  // mark used clusters
  int lrID,clID;
  while (seed) {
    if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
      AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
      if (ref) { // do we need to set or delete cluster->track ref?
	if (flg) {
	  if (!cl->GetRecoInfo()) cl->SetRecoInfo(ref);  // set ref only if cluster already does not refer to other track, inc.counter
	}
	else {
	  if (cl->GetRecoInfo()==ref) cl->SetRecoInfo(0); // unset reference only if it refers to ref, decrease counter	    
	}
      }
    }
    seed = (AliITSUSeed*)seed->GetParent();
  }
  //
}

//__________________________________________________________________
void AliITSUTrackerGlo::CheckClusterUsage()
{
  // check cluster usage info
  printf("ClusterUsage at pass %d\n",fCurrPassID);
  for (int ilr=0;ilr<fNLrActive;ilr++) {
    AliITSURecoLayer* lr = fITS->GetLayerActive(ilr);
    int ncl = lr->GetNClusters();
    int nclUsed = 0;
    int nclUs0CntShr1 = 0;
    int nclUs1Cnt0 = 0;
    int nclUs1Shr1 = 0;
    int nclUseL = 0;
    for (int icl=0;icl<ncl;icl++) {
      AliITSUClusterPix *cl = (AliITSUClusterPix*)lr->GetCluster(icl);
      int usc = cl->GetClUsage();
      if (usc<1) {
	if (cl->IsClusterUsed() || cl->IsClusterShared()) {
	  nclUs0CntShr1++;
	  printf("Lr%d Cl%4d: UseCnt=0 but UseFlg=%d UseShr=%d\n",ilr,icl,cl->IsClusterUsed(),cl->IsClusterShared());
	}
	continue;
      }
      nclUsed++;
      if (!cl->IsClusterUsed()) {
	nclUs1Cnt0++;
	printf("Lr%d Cl%4d: UseCnt=%d but UseFlg=%d UseShr=%d\n",ilr,icl,usc, cl->IsClusterUsed(),cl->IsClusterShared());
      }
      if (cl->IsClusterShared()) {
	nclUs1Shr1++;
	printf("Lr%d Cl%4d: UseCnt=%d but UseFlg=%d UseShr=%d\n",ilr,icl,usc, cl->IsClusterUsed(),cl->IsClusterShared());
      }
      if (usc>1) {
	nclUseL++;
	printf("Lr%d Cl%4d: UseCnt=%d!! but UseFlg=%d UseShr=%d\n",ilr,icl,usc, cl->IsClusterUsed(),cl->IsClusterShared());
      }
    }
    printf("ClUsage Lr%d: Ncl=%5d Nusd=%5d (%5.2f) Use0CS1:%4d Use1C0:%4d Use1S1:%4d UseL:%d\n",
	   ilr,ncl,nclUsed,ncl? float(nclUsed)/ncl : 0, nclUs0CntShr1,nclUs1Cnt0,nclUs1Shr1,nclUseL);
  }
}


#ifdef  _ITSU_TUNING_MODE_
//__________________________________________________________________
TObjArray* AliITSUTrackerGlo::BookControlHistos(const char* pref)
{
  // book special control histos
  TString prefS = pref;
  prefS.ToLower();
  TObjArray* dest = new TObjArray;
  dest->SetOwner(kTRUE);
  //
  const int kNResDef=7;
  const double kResDef[kNResDef]={0.05,0.05,0.3, 0.05,1,0.5,1.5};
  const double ptMax=10;
  const double plMax=10;
  const double chiMax=100;
  const int nptbins=50;
  const int nresbins=400;
  const int nplbins=50;
  const int nchbins=200;
  const int maxBr  = 15;
  const int maxCand = 200;
  TString ttl;
  for (int phase=0;phase<kNTrackingPhases;phase++) {
    for (int pass=0;pass<AliITSUReconstructor::GetRecoParam()->GetNTrackingConditions();pass++) {
      //
      for (int ilr=0;ilr<fNLrActive;ilr++) {
	//
	// ----------------- These are histos to be filled in Cluster2Tracks of each pass. 
	// PropagateBack and RefitInward will be stored among the histos of 1st pass
	if (pass>0 && phase!=kClus2Tracks) continue;
	//
	double mxdf = ilr>=kNResDef ? kResDef[kNResDef-1] : kResDef[ilr];
	ttl = Form("Pass%d_S%d_residY%d_%s",pass,phase,ilr,pref);
	TH2F* hdy = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nresbins,-mxdf,mxdf);
	dest->AddAtAndExpand(hdy,GetHistoID(ilr,kHResY,pass,phase));
	hdy->SetDirectory(0);
	//
	ttl = Form("Pass%d_S%d_residYPull%d_%s",pass,phase,ilr,pref);	
	TH2F* hdyp = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nplbins,-plMax,plMax);
	dest->AddAtAndExpand(hdyp,GetHistoID(ilr,kHResYP,pass,phase));
	hdyp->SetDirectory(0);
	//
	ttl = Form("Pass%d_S%d_residZ%d_%s",pass,phase,ilr,pref);	
	TH2F* hdz = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nresbins,-mxdf,mxdf);
	dest->AddAtAndExpand(hdz,GetHistoID(ilr,kHResZ,pass,phase));
	hdz->SetDirectory(0);
	//
	ttl = Form("Pass%d_S%d_residZPull%d_%s",pass,phase,ilr,pref);		
	TH2F* hdzp = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nplbins,-plMax,plMax);
	hdzp->SetDirectory(0);
	dest->AddAtAndExpand(hdzp,GetHistoID(ilr,kHResZP,pass,phase));
	//
	ttl = Form("Pass%d_S%d_chi2Cl%d_%s",pass,phase,ilr,pref);		
	TH2F* hchi = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
	hchi->SetDirectory(0);
	dest->AddAtAndExpand(hchi,GetHistoID(ilr,kHChi2Cl,pass,phase));
	//
	// ------------------- These histos are filled for Clusters2Tracks only
	if (phase!=kClus2Tracks) continue;
	//
	ttl = Form("Pass%d_S%d_chi2Nrm%d_%s",pass,phase,ilr,pref);		
	TH2F* hchiN = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
	hchiN->SetDirectory(0);
	dest->AddAtAndExpand(hchiN,GetHistoID(ilr,kHChi2Nrm,pass,phase));
	//
	ttl = Form("Pass%d_S%d_bestInBranch%d_%s",pass,phase,ilr,pref);		
	TH2* hnbr = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, maxBr,-0.5,maxBr-0.5);
	hnbr->SetDirectory(0);
	dest->AddAtAndExpand(hnbr,GetHistoID(ilr,kHBestInBranch,pass,phase));
	//
	ttl = Form("Pass%d_S%d_bestInCands%d_%s",pass,phase,ilr,pref);		
	TH2* hncn = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, maxCand,-0.5,maxCand-0.5);
	hncn->SetDirectory(0);
	dest->AddAtAndExpand(hncn,GetHistoID(ilr,kHBestInCand,pass,phase));
	//
      } // loop over layers
      //
      // these are special histos, filled not per layer but in the end of track fit in Clusters2Tracks in EVERY pass
      //  
      if (phase==kClus2Tracks) {
	TH2* hchiMatch = 0; 
	ttl = Form("Pass%d_Chi2Match_%s",pass,pref);
	hchiMatch = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
	hchiMatch->SetDirectory(0);
	dest->AddAtAndExpand(hchiMatch,GetHistoID(-1,kHChiMatch,pass,phase));
	// 
	TH2* hchiSA = 0; 
	ttl = Form("Pass%d_Chi2ITSSA_%s",pass,pref);
	hchiSA = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins/2,0.,chiMax/2);
	hchiSA->SetDirectory(0);
	dest->AddAtAndExpand(hchiSA,GetHistoID(-1,kHChiITSSA,pass,phase));
	// 	
      }
    } // loop over tracking passes
  }// loop over tracking phases
  //
  return dest;
}

//__________________________________________________________________
Int_t AliITSUTrackerGlo::GetHistoID(Int_t lr, Int_t hid, Int_t pass, Int_t phase)
{
  // get id for the requested histo
  if (lr<0) lr=-1;
  lr++;
  return pass*kHistosPass + phase*kHistosPhase + lr*kMaxHID + hid;
  //
}

#endif
 AliITSUTrackerGlo.cxx:1
 AliITSUTrackerGlo.cxx:2
 AliITSUTrackerGlo.cxx:3
 AliITSUTrackerGlo.cxx:4
 AliITSUTrackerGlo.cxx:5
 AliITSUTrackerGlo.cxx:6
 AliITSUTrackerGlo.cxx:7
 AliITSUTrackerGlo.cxx:8
 AliITSUTrackerGlo.cxx:9
 AliITSUTrackerGlo.cxx:10
 AliITSUTrackerGlo.cxx:11
 AliITSUTrackerGlo.cxx:12
 AliITSUTrackerGlo.cxx:13
 AliITSUTrackerGlo.cxx:14
 AliITSUTrackerGlo.cxx:15
 AliITSUTrackerGlo.cxx:16
 AliITSUTrackerGlo.cxx:17
 AliITSUTrackerGlo.cxx:18
 AliITSUTrackerGlo.cxx:19
 AliITSUTrackerGlo.cxx:20
 AliITSUTrackerGlo.cxx:21
 AliITSUTrackerGlo.cxx:22
 AliITSUTrackerGlo.cxx:23
 AliITSUTrackerGlo.cxx:24
 AliITSUTrackerGlo.cxx:25
 AliITSUTrackerGlo.cxx:26
 AliITSUTrackerGlo.cxx:27
 AliITSUTrackerGlo.cxx:28
 AliITSUTrackerGlo.cxx:29
 AliITSUTrackerGlo.cxx:30
 AliITSUTrackerGlo.cxx:31
 AliITSUTrackerGlo.cxx:32
 AliITSUTrackerGlo.cxx:33
 AliITSUTrackerGlo.cxx:34
 AliITSUTrackerGlo.cxx:35
 AliITSUTrackerGlo.cxx:36
 AliITSUTrackerGlo.cxx:37
 AliITSUTrackerGlo.cxx:38
 AliITSUTrackerGlo.cxx:39
 AliITSUTrackerGlo.cxx:40
 AliITSUTrackerGlo.cxx:41
 AliITSUTrackerGlo.cxx:42
 AliITSUTrackerGlo.cxx:43
 AliITSUTrackerGlo.cxx:44
 AliITSUTrackerGlo.cxx:45
 AliITSUTrackerGlo.cxx:46
 AliITSUTrackerGlo.cxx:47
 AliITSUTrackerGlo.cxx:48
 AliITSUTrackerGlo.cxx:49
 AliITSUTrackerGlo.cxx:50
 AliITSUTrackerGlo.cxx:51
 AliITSUTrackerGlo.cxx:52
 AliITSUTrackerGlo.cxx:53
 AliITSUTrackerGlo.cxx:54
 AliITSUTrackerGlo.cxx:55
 AliITSUTrackerGlo.cxx:56
 AliITSUTrackerGlo.cxx:57
 AliITSUTrackerGlo.cxx:58
 AliITSUTrackerGlo.cxx:59
 AliITSUTrackerGlo.cxx:60
 AliITSUTrackerGlo.cxx:61
 AliITSUTrackerGlo.cxx:62
 AliITSUTrackerGlo.cxx:63
 AliITSUTrackerGlo.cxx:64
 AliITSUTrackerGlo.cxx:65
 AliITSUTrackerGlo.cxx:66
 AliITSUTrackerGlo.cxx:67
 AliITSUTrackerGlo.cxx:68
 AliITSUTrackerGlo.cxx:69
 AliITSUTrackerGlo.cxx:70
 AliITSUTrackerGlo.cxx:71
 AliITSUTrackerGlo.cxx:72
 AliITSUTrackerGlo.cxx:73
 AliITSUTrackerGlo.cxx:74
 AliITSUTrackerGlo.cxx:75
 AliITSUTrackerGlo.cxx:76
 AliITSUTrackerGlo.cxx:77
 AliITSUTrackerGlo.cxx:78
 AliITSUTrackerGlo.cxx:79
 AliITSUTrackerGlo.cxx:80
 AliITSUTrackerGlo.cxx:81
 AliITSUTrackerGlo.cxx:82
 AliITSUTrackerGlo.cxx:83
 AliITSUTrackerGlo.cxx:84
 AliITSUTrackerGlo.cxx:85
 AliITSUTrackerGlo.cxx:86
 AliITSUTrackerGlo.cxx:87
 AliITSUTrackerGlo.cxx:88
 AliITSUTrackerGlo.cxx:89
 AliITSUTrackerGlo.cxx:90
 AliITSUTrackerGlo.cxx:91
 AliITSUTrackerGlo.cxx:92
 AliITSUTrackerGlo.cxx:93
 AliITSUTrackerGlo.cxx:94
 AliITSUTrackerGlo.cxx:95
 AliITSUTrackerGlo.cxx:96
 AliITSUTrackerGlo.cxx:97
 AliITSUTrackerGlo.cxx:98
 AliITSUTrackerGlo.cxx:99
 AliITSUTrackerGlo.cxx:100
 AliITSUTrackerGlo.cxx:101
 AliITSUTrackerGlo.cxx:102
 AliITSUTrackerGlo.cxx:103
 AliITSUTrackerGlo.cxx:104
 AliITSUTrackerGlo.cxx:105
 AliITSUTrackerGlo.cxx:106
 AliITSUTrackerGlo.cxx:107
 AliITSUTrackerGlo.cxx:108
 AliITSUTrackerGlo.cxx:109
 AliITSUTrackerGlo.cxx:110
 AliITSUTrackerGlo.cxx:111
 AliITSUTrackerGlo.cxx:112
 AliITSUTrackerGlo.cxx:113
 AliITSUTrackerGlo.cxx:114
 AliITSUTrackerGlo.cxx:115
 AliITSUTrackerGlo.cxx:116
 AliITSUTrackerGlo.cxx:117
 AliITSUTrackerGlo.cxx:118
 AliITSUTrackerGlo.cxx:119
 AliITSUTrackerGlo.cxx:120
 AliITSUTrackerGlo.cxx:121
 AliITSUTrackerGlo.cxx:122
 AliITSUTrackerGlo.cxx:123
 AliITSUTrackerGlo.cxx:124
 AliITSUTrackerGlo.cxx:125
 AliITSUTrackerGlo.cxx:126
 AliITSUTrackerGlo.cxx:127
 AliITSUTrackerGlo.cxx:128
 AliITSUTrackerGlo.cxx:129
 AliITSUTrackerGlo.cxx:130
 AliITSUTrackerGlo.cxx:131
 AliITSUTrackerGlo.cxx:132
 AliITSUTrackerGlo.cxx:133
 AliITSUTrackerGlo.cxx:134
 AliITSUTrackerGlo.cxx:135
 AliITSUTrackerGlo.cxx:136
 AliITSUTrackerGlo.cxx:137
 AliITSUTrackerGlo.cxx:138
 AliITSUTrackerGlo.cxx:139
 AliITSUTrackerGlo.cxx:140
 AliITSUTrackerGlo.cxx:141
 AliITSUTrackerGlo.cxx:142
 AliITSUTrackerGlo.cxx:143
 AliITSUTrackerGlo.cxx:144
 AliITSUTrackerGlo.cxx:145
 AliITSUTrackerGlo.cxx:146
 AliITSUTrackerGlo.cxx:147
 AliITSUTrackerGlo.cxx:148
 AliITSUTrackerGlo.cxx:149
 AliITSUTrackerGlo.cxx:150
 AliITSUTrackerGlo.cxx:151
 AliITSUTrackerGlo.cxx:152
 AliITSUTrackerGlo.cxx:153
 AliITSUTrackerGlo.cxx:154
 AliITSUTrackerGlo.cxx:155
 AliITSUTrackerGlo.cxx:156
 AliITSUTrackerGlo.cxx:157
 AliITSUTrackerGlo.cxx:158
 AliITSUTrackerGlo.cxx:159
 AliITSUTrackerGlo.cxx:160
 AliITSUTrackerGlo.cxx:161
 AliITSUTrackerGlo.cxx:162
 AliITSUTrackerGlo.cxx:163
 AliITSUTrackerGlo.cxx:164
 AliITSUTrackerGlo.cxx:165
 AliITSUTrackerGlo.cxx:166
 AliITSUTrackerGlo.cxx:167
 AliITSUTrackerGlo.cxx:168
 AliITSUTrackerGlo.cxx:169
 AliITSUTrackerGlo.cxx:170
 AliITSUTrackerGlo.cxx:171
 AliITSUTrackerGlo.cxx:172
 AliITSUTrackerGlo.cxx:173
 AliITSUTrackerGlo.cxx:174
 AliITSUTrackerGlo.cxx:175
 AliITSUTrackerGlo.cxx:176
 AliITSUTrackerGlo.cxx:177
 AliITSUTrackerGlo.cxx:178
 AliITSUTrackerGlo.cxx:179
 AliITSUTrackerGlo.cxx:180
 AliITSUTrackerGlo.cxx:181
 AliITSUTrackerGlo.cxx:182
 AliITSUTrackerGlo.cxx:183
 AliITSUTrackerGlo.cxx:184
 AliITSUTrackerGlo.cxx:185
 AliITSUTrackerGlo.cxx:186
 AliITSUTrackerGlo.cxx:187
 AliITSUTrackerGlo.cxx:188
 AliITSUTrackerGlo.cxx:189
 AliITSUTrackerGlo.cxx:190
 AliITSUTrackerGlo.cxx:191
 AliITSUTrackerGlo.cxx:192
 AliITSUTrackerGlo.cxx:193
 AliITSUTrackerGlo.cxx:194
 AliITSUTrackerGlo.cxx:195
 AliITSUTrackerGlo.cxx:196
 AliITSUTrackerGlo.cxx:197
 AliITSUTrackerGlo.cxx:198
 AliITSUTrackerGlo.cxx:199
 AliITSUTrackerGlo.cxx:200
 AliITSUTrackerGlo.cxx:201
 AliITSUTrackerGlo.cxx:202
 AliITSUTrackerGlo.cxx:203
 AliITSUTrackerGlo.cxx:204
 AliITSUTrackerGlo.cxx:205
 AliITSUTrackerGlo.cxx:206
 AliITSUTrackerGlo.cxx:207
 AliITSUTrackerGlo.cxx:208
 AliITSUTrackerGlo.cxx:209
 AliITSUTrackerGlo.cxx:210
 AliITSUTrackerGlo.cxx:211
 AliITSUTrackerGlo.cxx:212
 AliITSUTrackerGlo.cxx:213
 AliITSUTrackerGlo.cxx:214
 AliITSUTrackerGlo.cxx:215
 AliITSUTrackerGlo.cxx:216
 AliITSUTrackerGlo.cxx:217
 AliITSUTrackerGlo.cxx:218
 AliITSUTrackerGlo.cxx:219
 AliITSUTrackerGlo.cxx:220
 AliITSUTrackerGlo.cxx:221
 AliITSUTrackerGlo.cxx:222
 AliITSUTrackerGlo.cxx:223
 AliITSUTrackerGlo.cxx:224
 AliITSUTrackerGlo.cxx:225
 AliITSUTrackerGlo.cxx:226
 AliITSUTrackerGlo.cxx:227
 AliITSUTrackerGlo.cxx:228
 AliITSUTrackerGlo.cxx:229
 AliITSUTrackerGlo.cxx:230
 AliITSUTrackerGlo.cxx:231
 AliITSUTrackerGlo.cxx:232
 AliITSUTrackerGlo.cxx:233
 AliITSUTrackerGlo.cxx:234
 AliITSUTrackerGlo.cxx:235
 AliITSUTrackerGlo.cxx:236
 AliITSUTrackerGlo.cxx:237
 AliITSUTrackerGlo.cxx:238
 AliITSUTrackerGlo.cxx:239
 AliITSUTrackerGlo.cxx:240
 AliITSUTrackerGlo.cxx:241
 AliITSUTrackerGlo.cxx:242
 AliITSUTrackerGlo.cxx:243
 AliITSUTrackerGlo.cxx:244
 AliITSUTrackerGlo.cxx:245
 AliITSUTrackerGlo.cxx:246
 AliITSUTrackerGlo.cxx:247
 AliITSUTrackerGlo.cxx:248
 AliITSUTrackerGlo.cxx:249
 AliITSUTrackerGlo.cxx:250
 AliITSUTrackerGlo.cxx:251
 AliITSUTrackerGlo.cxx:252
 AliITSUTrackerGlo.cxx:253
 AliITSUTrackerGlo.cxx:254
 AliITSUTrackerGlo.cxx:255
 AliITSUTrackerGlo.cxx:256
 AliITSUTrackerGlo.cxx:257
 AliITSUTrackerGlo.cxx:258
 AliITSUTrackerGlo.cxx:259
 AliITSUTrackerGlo.cxx:260
 AliITSUTrackerGlo.cxx:261
 AliITSUTrackerGlo.cxx:262
 AliITSUTrackerGlo.cxx:263
 AliITSUTrackerGlo.cxx:264
 AliITSUTrackerGlo.cxx:265
 AliITSUTrackerGlo.cxx:266
 AliITSUTrackerGlo.cxx:267
 AliITSUTrackerGlo.cxx:268
 AliITSUTrackerGlo.cxx:269
 AliITSUTrackerGlo.cxx:270
 AliITSUTrackerGlo.cxx:271
 AliITSUTrackerGlo.cxx:272
 AliITSUTrackerGlo.cxx:273
 AliITSUTrackerGlo.cxx:274
 AliITSUTrackerGlo.cxx:275
 AliITSUTrackerGlo.cxx:276
 AliITSUTrackerGlo.cxx:277
 AliITSUTrackerGlo.cxx:278
 AliITSUTrackerGlo.cxx:279
 AliITSUTrackerGlo.cxx:280
 AliITSUTrackerGlo.cxx:281
 AliITSUTrackerGlo.cxx:282
 AliITSUTrackerGlo.cxx:283
 AliITSUTrackerGlo.cxx:284
 AliITSUTrackerGlo.cxx:285
 AliITSUTrackerGlo.cxx:286
 AliITSUTrackerGlo.cxx:287
 AliITSUTrackerGlo.cxx:288
 AliITSUTrackerGlo.cxx:289
 AliITSUTrackerGlo.cxx:290
 AliITSUTrackerGlo.cxx:291
 AliITSUTrackerGlo.cxx:292
 AliITSUTrackerGlo.cxx:293
 AliITSUTrackerGlo.cxx:294
 AliITSUTrackerGlo.cxx:295
 AliITSUTrackerGlo.cxx:296
 AliITSUTrackerGlo.cxx:297
 AliITSUTrackerGlo.cxx:298
 AliITSUTrackerGlo.cxx:299
 AliITSUTrackerGlo.cxx:300
 AliITSUTrackerGlo.cxx:301
 AliITSUTrackerGlo.cxx:302
 AliITSUTrackerGlo.cxx:303
 AliITSUTrackerGlo.cxx:304
 AliITSUTrackerGlo.cxx:305
 AliITSUTrackerGlo.cxx:306
 AliITSUTrackerGlo.cxx:307
 AliITSUTrackerGlo.cxx:308
 AliITSUTrackerGlo.cxx:309
 AliITSUTrackerGlo.cxx:310
 AliITSUTrackerGlo.cxx:311
 AliITSUTrackerGlo.cxx:312
 AliITSUTrackerGlo.cxx:313
 AliITSUTrackerGlo.cxx:314
 AliITSUTrackerGlo.cxx:315
 AliITSUTrackerGlo.cxx:316
 AliITSUTrackerGlo.cxx:317
 AliITSUTrackerGlo.cxx:318
 AliITSUTrackerGlo.cxx:319
 AliITSUTrackerGlo.cxx:320
 AliITSUTrackerGlo.cxx:321
 AliITSUTrackerGlo.cxx:322
 AliITSUTrackerGlo.cxx:323
 AliITSUTrackerGlo.cxx:324
 AliITSUTrackerGlo.cxx:325
 AliITSUTrackerGlo.cxx:326
 AliITSUTrackerGlo.cxx:327
 AliITSUTrackerGlo.cxx:328
 AliITSUTrackerGlo.cxx:329
 AliITSUTrackerGlo.cxx:330
 AliITSUTrackerGlo.cxx:331
 AliITSUTrackerGlo.cxx:332
 AliITSUTrackerGlo.cxx:333
 AliITSUTrackerGlo.cxx:334
 AliITSUTrackerGlo.cxx:335
 AliITSUTrackerGlo.cxx:336
 AliITSUTrackerGlo.cxx:337
 AliITSUTrackerGlo.cxx:338
 AliITSUTrackerGlo.cxx:339
 AliITSUTrackerGlo.cxx:340
 AliITSUTrackerGlo.cxx:341
 AliITSUTrackerGlo.cxx:342
 AliITSUTrackerGlo.cxx:343
 AliITSUTrackerGlo.cxx:344
 AliITSUTrackerGlo.cxx:345
 AliITSUTrackerGlo.cxx:346
 AliITSUTrackerGlo.cxx:347
 AliITSUTrackerGlo.cxx:348
 AliITSUTrackerGlo.cxx:349
 AliITSUTrackerGlo.cxx:350
 AliITSUTrackerGlo.cxx:351
 AliITSUTrackerGlo.cxx:352
 AliITSUTrackerGlo.cxx:353
 AliITSUTrackerGlo.cxx:354
 AliITSUTrackerGlo.cxx:355
 AliITSUTrackerGlo.cxx:356
 AliITSUTrackerGlo.cxx:357
 AliITSUTrackerGlo.cxx:358
 AliITSUTrackerGlo.cxx:359
 AliITSUTrackerGlo.cxx:360
 AliITSUTrackerGlo.cxx:361
 AliITSUTrackerGlo.cxx:362
 AliITSUTrackerGlo.cxx:363
 AliITSUTrackerGlo.cxx:364
 AliITSUTrackerGlo.cxx:365
 AliITSUTrackerGlo.cxx:366
 AliITSUTrackerGlo.cxx:367
 AliITSUTrackerGlo.cxx:368
 AliITSUTrackerGlo.cxx:369
 AliITSUTrackerGlo.cxx:370
 AliITSUTrackerGlo.cxx:371
 AliITSUTrackerGlo.cxx:372
 AliITSUTrackerGlo.cxx:373
 AliITSUTrackerGlo.cxx:374
 AliITSUTrackerGlo.cxx:375
 AliITSUTrackerGlo.cxx:376
 AliITSUTrackerGlo.cxx:377
 AliITSUTrackerGlo.cxx:378
 AliITSUTrackerGlo.cxx:379
 AliITSUTrackerGlo.cxx:380
 AliITSUTrackerGlo.cxx:381
 AliITSUTrackerGlo.cxx:382
 AliITSUTrackerGlo.cxx:383
 AliITSUTrackerGlo.cxx:384
 AliITSUTrackerGlo.cxx:385
 AliITSUTrackerGlo.cxx:386
 AliITSUTrackerGlo.cxx:387
 AliITSUTrackerGlo.cxx:388
 AliITSUTrackerGlo.cxx:389
 AliITSUTrackerGlo.cxx:390
 AliITSUTrackerGlo.cxx:391
 AliITSUTrackerGlo.cxx:392
 AliITSUTrackerGlo.cxx:393
 AliITSUTrackerGlo.cxx:394
 AliITSUTrackerGlo.cxx:395
 AliITSUTrackerGlo.cxx:396
 AliITSUTrackerGlo.cxx:397
 AliITSUTrackerGlo.cxx:398
 AliITSUTrackerGlo.cxx:399
 AliITSUTrackerGlo.cxx:400
 AliITSUTrackerGlo.cxx:401
 AliITSUTrackerGlo.cxx:402
 AliITSUTrackerGlo.cxx:403
 AliITSUTrackerGlo.cxx:404
 AliITSUTrackerGlo.cxx:405
 AliITSUTrackerGlo.cxx:406
 AliITSUTrackerGlo.cxx:407
 AliITSUTrackerGlo.cxx:408
 AliITSUTrackerGlo.cxx:409
 AliITSUTrackerGlo.cxx:410
 AliITSUTrackerGlo.cxx:411
 AliITSUTrackerGlo.cxx:412
 AliITSUTrackerGlo.cxx:413
 AliITSUTrackerGlo.cxx:414
 AliITSUTrackerGlo.cxx:415
 AliITSUTrackerGlo.cxx:416
 AliITSUTrackerGlo.cxx:417
 AliITSUTrackerGlo.cxx:418
 AliITSUTrackerGlo.cxx:419
 AliITSUTrackerGlo.cxx:420
 AliITSUTrackerGlo.cxx:421
 AliITSUTrackerGlo.cxx:422
 AliITSUTrackerGlo.cxx:423
 AliITSUTrackerGlo.cxx:424
 AliITSUTrackerGlo.cxx:425
 AliITSUTrackerGlo.cxx:426
 AliITSUTrackerGlo.cxx:427
 AliITSUTrackerGlo.cxx:428
 AliITSUTrackerGlo.cxx:429
 AliITSUTrackerGlo.cxx:430
 AliITSUTrackerGlo.cxx:431
 AliITSUTrackerGlo.cxx:432
 AliITSUTrackerGlo.cxx:433
 AliITSUTrackerGlo.cxx:434
 AliITSUTrackerGlo.cxx:435
 AliITSUTrackerGlo.cxx:436
 AliITSUTrackerGlo.cxx:437
 AliITSUTrackerGlo.cxx:438
 AliITSUTrackerGlo.cxx:439
 AliITSUTrackerGlo.cxx:440
 AliITSUTrackerGlo.cxx:441
 AliITSUTrackerGlo.cxx:442
 AliITSUTrackerGlo.cxx:443
 AliITSUTrackerGlo.cxx:444
 AliITSUTrackerGlo.cxx:445
 AliITSUTrackerGlo.cxx:446
 AliITSUTrackerGlo.cxx:447
 AliITSUTrackerGlo.cxx:448
 AliITSUTrackerGlo.cxx:449
 AliITSUTrackerGlo.cxx:450
 AliITSUTrackerGlo.cxx:451
 AliITSUTrackerGlo.cxx:452
 AliITSUTrackerGlo.cxx:453
 AliITSUTrackerGlo.cxx:454
 AliITSUTrackerGlo.cxx:455
 AliITSUTrackerGlo.cxx:456
 AliITSUTrackerGlo.cxx:457
 AliITSUTrackerGlo.cxx:458
 AliITSUTrackerGlo.cxx:459
 AliITSUTrackerGlo.cxx:460
 AliITSUTrackerGlo.cxx:461
 AliITSUTrackerGlo.cxx:462
 AliITSUTrackerGlo.cxx:463
 AliITSUTrackerGlo.cxx:464
 AliITSUTrackerGlo.cxx:465
 AliITSUTrackerGlo.cxx:466
 AliITSUTrackerGlo.cxx:467
 AliITSUTrackerGlo.cxx:468
 AliITSUTrackerGlo.cxx:469
 AliITSUTrackerGlo.cxx:470
 AliITSUTrackerGlo.cxx:471
 AliITSUTrackerGlo.cxx:472
 AliITSUTrackerGlo.cxx:473
 AliITSUTrackerGlo.cxx:474
 AliITSUTrackerGlo.cxx:475
 AliITSUTrackerGlo.cxx:476
 AliITSUTrackerGlo.cxx:477
 AliITSUTrackerGlo.cxx:478
 AliITSUTrackerGlo.cxx:479
 AliITSUTrackerGlo.cxx:480
 AliITSUTrackerGlo.cxx:481
 AliITSUTrackerGlo.cxx:482
 AliITSUTrackerGlo.cxx:483
 AliITSUTrackerGlo.cxx:484
 AliITSUTrackerGlo.cxx:485
 AliITSUTrackerGlo.cxx:486
 AliITSUTrackerGlo.cxx:487
 AliITSUTrackerGlo.cxx:488
 AliITSUTrackerGlo.cxx:489
 AliITSUTrackerGlo.cxx:490
 AliITSUTrackerGlo.cxx:491
 AliITSUTrackerGlo.cxx:492
 AliITSUTrackerGlo.cxx:493
 AliITSUTrackerGlo.cxx:494
 AliITSUTrackerGlo.cxx:495
 AliITSUTrackerGlo.cxx:496
 AliITSUTrackerGlo.cxx:497
 AliITSUTrackerGlo.cxx:498
 AliITSUTrackerGlo.cxx:499
 AliITSUTrackerGlo.cxx:500
 AliITSUTrackerGlo.cxx:501
 AliITSUTrackerGlo.cxx:502
 AliITSUTrackerGlo.cxx:503
 AliITSUTrackerGlo.cxx:504
 AliITSUTrackerGlo.cxx:505
 AliITSUTrackerGlo.cxx:506
 AliITSUTrackerGlo.cxx:507
 AliITSUTrackerGlo.cxx:508
 AliITSUTrackerGlo.cxx:509
 AliITSUTrackerGlo.cxx:510
 AliITSUTrackerGlo.cxx:511
 AliITSUTrackerGlo.cxx:512
 AliITSUTrackerGlo.cxx:513
 AliITSUTrackerGlo.cxx:514
 AliITSUTrackerGlo.cxx:515
 AliITSUTrackerGlo.cxx:516
 AliITSUTrackerGlo.cxx:517
 AliITSUTrackerGlo.cxx:518
 AliITSUTrackerGlo.cxx:519
 AliITSUTrackerGlo.cxx:520
 AliITSUTrackerGlo.cxx:521
 AliITSUTrackerGlo.cxx:522
 AliITSUTrackerGlo.cxx:523
 AliITSUTrackerGlo.cxx:524
 AliITSUTrackerGlo.cxx:525
 AliITSUTrackerGlo.cxx:526
 AliITSUTrackerGlo.cxx:527
 AliITSUTrackerGlo.cxx:528
 AliITSUTrackerGlo.cxx:529
 AliITSUTrackerGlo.cxx:530
 AliITSUTrackerGlo.cxx:531
 AliITSUTrackerGlo.cxx:532
 AliITSUTrackerGlo.cxx:533
 AliITSUTrackerGlo.cxx:534
 AliITSUTrackerGlo.cxx:535
 AliITSUTrackerGlo.cxx:536
 AliITSUTrackerGlo.cxx:537
 AliITSUTrackerGlo.cxx:538
 AliITSUTrackerGlo.cxx:539
 AliITSUTrackerGlo.cxx:540
 AliITSUTrackerGlo.cxx:541
 AliITSUTrackerGlo.cxx:542
 AliITSUTrackerGlo.cxx:543
 AliITSUTrackerGlo.cxx:544
 AliITSUTrackerGlo.cxx:545
 AliITSUTrackerGlo.cxx:546
 AliITSUTrackerGlo.cxx:547
 AliITSUTrackerGlo.cxx:548
 AliITSUTrackerGlo.cxx:549
 AliITSUTrackerGlo.cxx:550
 AliITSUTrackerGlo.cxx:551
 AliITSUTrackerGlo.cxx:552
 AliITSUTrackerGlo.cxx:553
 AliITSUTrackerGlo.cxx:554
 AliITSUTrackerGlo.cxx:555
 AliITSUTrackerGlo.cxx:556
 AliITSUTrackerGlo.cxx:557
 AliITSUTrackerGlo.cxx:558
 AliITSUTrackerGlo.cxx:559
 AliITSUTrackerGlo.cxx:560
 AliITSUTrackerGlo.cxx:561
 AliITSUTrackerGlo.cxx:562
 AliITSUTrackerGlo.cxx:563
 AliITSUTrackerGlo.cxx:564
 AliITSUTrackerGlo.cxx:565
 AliITSUTrackerGlo.cxx:566
 AliITSUTrackerGlo.cxx:567
 AliITSUTrackerGlo.cxx:568
 AliITSUTrackerGlo.cxx:569
 AliITSUTrackerGlo.cxx:570
 AliITSUTrackerGlo.cxx:571
 AliITSUTrackerGlo.cxx:572
 AliITSUTrackerGlo.cxx:573
 AliITSUTrackerGlo.cxx:574
 AliITSUTrackerGlo.cxx:575
 AliITSUTrackerGlo.cxx:576
 AliITSUTrackerGlo.cxx:577
 AliITSUTrackerGlo.cxx:578
 AliITSUTrackerGlo.cxx:579
 AliITSUTrackerGlo.cxx:580
 AliITSUTrackerGlo.cxx:581
 AliITSUTrackerGlo.cxx:582
 AliITSUTrackerGlo.cxx:583
 AliITSUTrackerGlo.cxx:584
 AliITSUTrackerGlo.cxx:585
 AliITSUTrackerGlo.cxx:586
 AliITSUTrackerGlo.cxx:587
 AliITSUTrackerGlo.cxx:588
 AliITSUTrackerGlo.cxx:589
 AliITSUTrackerGlo.cxx:590
 AliITSUTrackerGlo.cxx:591
 AliITSUTrackerGlo.cxx:592
 AliITSUTrackerGlo.cxx:593
 AliITSUTrackerGlo.cxx:594
 AliITSUTrackerGlo.cxx:595
 AliITSUTrackerGlo.cxx:596
 AliITSUTrackerGlo.cxx:597
 AliITSUTrackerGlo.cxx:598
 AliITSUTrackerGlo.cxx:599
 AliITSUTrackerGlo.cxx:600
 AliITSUTrackerGlo.cxx:601
 AliITSUTrackerGlo.cxx:602
 AliITSUTrackerGlo.cxx:603
 AliITSUTrackerGlo.cxx:604
 AliITSUTrackerGlo.cxx:605
 AliITSUTrackerGlo.cxx:606
 AliITSUTrackerGlo.cxx:607
 AliITSUTrackerGlo.cxx:608
 AliITSUTrackerGlo.cxx:609
 AliITSUTrackerGlo.cxx:610
 AliITSUTrackerGlo.cxx:611
 AliITSUTrackerGlo.cxx:612
 AliITSUTrackerGlo.cxx:613
 AliITSUTrackerGlo.cxx:614
 AliITSUTrackerGlo.cxx:615
 AliITSUTrackerGlo.cxx:616
 AliITSUTrackerGlo.cxx:617
 AliITSUTrackerGlo.cxx:618
 AliITSUTrackerGlo.cxx:619
 AliITSUTrackerGlo.cxx:620
 AliITSUTrackerGlo.cxx:621
 AliITSUTrackerGlo.cxx:622
 AliITSUTrackerGlo.cxx:623
 AliITSUTrackerGlo.cxx:624
 AliITSUTrackerGlo.cxx:625
 AliITSUTrackerGlo.cxx:626
 AliITSUTrackerGlo.cxx:627
 AliITSUTrackerGlo.cxx:628
 AliITSUTrackerGlo.cxx:629
 AliITSUTrackerGlo.cxx:630
 AliITSUTrackerGlo.cxx:631
 AliITSUTrackerGlo.cxx:632
 AliITSUTrackerGlo.cxx:633
 AliITSUTrackerGlo.cxx:634
 AliITSUTrackerGlo.cxx:635
 AliITSUTrackerGlo.cxx:636
 AliITSUTrackerGlo.cxx:637
 AliITSUTrackerGlo.cxx:638
 AliITSUTrackerGlo.cxx:639
 AliITSUTrackerGlo.cxx:640
 AliITSUTrackerGlo.cxx:641
 AliITSUTrackerGlo.cxx:642
 AliITSUTrackerGlo.cxx:643
 AliITSUTrackerGlo.cxx:644
 AliITSUTrackerGlo.cxx:645
 AliITSUTrackerGlo.cxx:646
 AliITSUTrackerGlo.cxx:647
 AliITSUTrackerGlo.cxx:648
 AliITSUTrackerGlo.cxx:649
 AliITSUTrackerGlo.cxx:650
 AliITSUTrackerGlo.cxx:651
 AliITSUTrackerGlo.cxx:652
 AliITSUTrackerGlo.cxx:653
 AliITSUTrackerGlo.cxx:654
 AliITSUTrackerGlo.cxx:655
 AliITSUTrackerGlo.cxx:656
 AliITSUTrackerGlo.cxx:657
 AliITSUTrackerGlo.cxx:658
 AliITSUTrackerGlo.cxx:659
 AliITSUTrackerGlo.cxx:660
 AliITSUTrackerGlo.cxx:661
 AliITSUTrackerGlo.cxx:662
 AliITSUTrackerGlo.cxx:663
 AliITSUTrackerGlo.cxx:664
 AliITSUTrackerGlo.cxx:665
 AliITSUTrackerGlo.cxx:666
 AliITSUTrackerGlo.cxx:667
 AliITSUTrackerGlo.cxx:668
 AliITSUTrackerGlo.cxx:669
 AliITSUTrackerGlo.cxx:670
 AliITSUTrackerGlo.cxx:671
 AliITSUTrackerGlo.cxx:672
 AliITSUTrackerGlo.cxx:673
 AliITSUTrackerGlo.cxx:674
 AliITSUTrackerGlo.cxx:675
 AliITSUTrackerGlo.cxx:676
 AliITSUTrackerGlo.cxx:677
 AliITSUTrackerGlo.cxx:678
 AliITSUTrackerGlo.cxx:679
 AliITSUTrackerGlo.cxx:680
 AliITSUTrackerGlo.cxx:681
 AliITSUTrackerGlo.cxx:682
 AliITSUTrackerGlo.cxx:683
 AliITSUTrackerGlo.cxx:684
 AliITSUTrackerGlo.cxx:685
 AliITSUTrackerGlo.cxx:686
 AliITSUTrackerGlo.cxx:687
 AliITSUTrackerGlo.cxx:688
 AliITSUTrackerGlo.cxx:689
 AliITSUTrackerGlo.cxx:690
 AliITSUTrackerGlo.cxx:691
 AliITSUTrackerGlo.cxx:692
 AliITSUTrackerGlo.cxx:693
 AliITSUTrackerGlo.cxx:694
 AliITSUTrackerGlo.cxx:695
 AliITSUTrackerGlo.cxx:696
 AliITSUTrackerGlo.cxx:697
 AliITSUTrackerGlo.cxx:698
 AliITSUTrackerGlo.cxx:699
 AliITSUTrackerGlo.cxx:700
 AliITSUTrackerGlo.cxx:701
 AliITSUTrackerGlo.cxx:702
 AliITSUTrackerGlo.cxx:703
 AliITSUTrackerGlo.cxx:704
 AliITSUTrackerGlo.cxx:705
 AliITSUTrackerGlo.cxx:706
 AliITSUTrackerGlo.cxx:707
 AliITSUTrackerGlo.cxx:708
 AliITSUTrackerGlo.cxx:709
 AliITSUTrackerGlo.cxx:710
 AliITSUTrackerGlo.cxx:711
 AliITSUTrackerGlo.cxx:712
 AliITSUTrackerGlo.cxx:713
 AliITSUTrackerGlo.cxx:714
 AliITSUTrackerGlo.cxx:715
 AliITSUTrackerGlo.cxx:716
 AliITSUTrackerGlo.cxx:717
 AliITSUTrackerGlo.cxx:718
 AliITSUTrackerGlo.cxx:719
 AliITSUTrackerGlo.cxx:720
 AliITSUTrackerGlo.cxx:721
 AliITSUTrackerGlo.cxx:722
 AliITSUTrackerGlo.cxx:723
 AliITSUTrackerGlo.cxx:724
 AliITSUTrackerGlo.cxx:725
 AliITSUTrackerGlo.cxx:726
 AliITSUTrackerGlo.cxx:727
 AliITSUTrackerGlo.cxx:728
 AliITSUTrackerGlo.cxx:729
 AliITSUTrackerGlo.cxx:730
 AliITSUTrackerGlo.cxx:731
 AliITSUTrackerGlo.cxx:732
 AliITSUTrackerGlo.cxx:733
 AliITSUTrackerGlo.cxx:734
 AliITSUTrackerGlo.cxx:735
 AliITSUTrackerGlo.cxx:736
 AliITSUTrackerGlo.cxx:737
 AliITSUTrackerGlo.cxx:738
 AliITSUTrackerGlo.cxx:739
 AliITSUTrackerGlo.cxx:740
 AliITSUTrackerGlo.cxx:741
 AliITSUTrackerGlo.cxx:742
 AliITSUTrackerGlo.cxx:743
 AliITSUTrackerGlo.cxx:744
 AliITSUTrackerGlo.cxx:745
 AliITSUTrackerGlo.cxx:746
 AliITSUTrackerGlo.cxx:747
 AliITSUTrackerGlo.cxx:748
 AliITSUTrackerGlo.cxx:749
 AliITSUTrackerGlo.cxx:750
 AliITSUTrackerGlo.cxx:751
 AliITSUTrackerGlo.cxx:752
 AliITSUTrackerGlo.cxx:753
 AliITSUTrackerGlo.cxx:754
 AliITSUTrackerGlo.cxx:755
 AliITSUTrackerGlo.cxx:756
 AliITSUTrackerGlo.cxx:757
 AliITSUTrackerGlo.cxx:758
 AliITSUTrackerGlo.cxx:759
 AliITSUTrackerGlo.cxx:760
 AliITSUTrackerGlo.cxx:761
 AliITSUTrackerGlo.cxx:762
 AliITSUTrackerGlo.cxx:763
 AliITSUTrackerGlo.cxx:764
 AliITSUTrackerGlo.cxx:765
 AliITSUTrackerGlo.cxx:766
 AliITSUTrackerGlo.cxx:767
 AliITSUTrackerGlo.cxx:768
 AliITSUTrackerGlo.cxx:769
 AliITSUTrackerGlo.cxx:770
 AliITSUTrackerGlo.cxx:771
 AliITSUTrackerGlo.cxx:772
 AliITSUTrackerGlo.cxx:773
 AliITSUTrackerGlo.cxx:774
 AliITSUTrackerGlo.cxx:775
 AliITSUTrackerGlo.cxx:776
 AliITSUTrackerGlo.cxx:777
 AliITSUTrackerGlo.cxx:778
 AliITSUTrackerGlo.cxx:779
 AliITSUTrackerGlo.cxx:780
 AliITSUTrackerGlo.cxx:781
 AliITSUTrackerGlo.cxx:782
 AliITSUTrackerGlo.cxx:783
 AliITSUTrackerGlo.cxx:784
 AliITSUTrackerGlo.cxx:785
 AliITSUTrackerGlo.cxx:786
 AliITSUTrackerGlo.cxx:787
 AliITSUTrackerGlo.cxx:788
 AliITSUTrackerGlo.cxx:789
 AliITSUTrackerGlo.cxx:790
 AliITSUTrackerGlo.cxx:791
 AliITSUTrackerGlo.cxx:792
 AliITSUTrackerGlo.cxx:793
 AliITSUTrackerGlo.cxx:794
 AliITSUTrackerGlo.cxx:795
 AliITSUTrackerGlo.cxx:796
 AliITSUTrackerGlo.cxx:797
 AliITSUTrackerGlo.cxx:798
 AliITSUTrackerGlo.cxx:799
 AliITSUTrackerGlo.cxx:800
 AliITSUTrackerGlo.cxx:801
 AliITSUTrackerGlo.cxx:802
 AliITSUTrackerGlo.cxx:803
 AliITSUTrackerGlo.cxx:804
 AliITSUTrackerGlo.cxx:805
 AliITSUTrackerGlo.cxx:806
 AliITSUTrackerGlo.cxx:807
 AliITSUTrackerGlo.cxx:808
 AliITSUTrackerGlo.cxx:809
 AliITSUTrackerGlo.cxx:810
 AliITSUTrackerGlo.cxx:811
 AliITSUTrackerGlo.cxx:812
 AliITSUTrackerGlo.cxx:813
 AliITSUTrackerGlo.cxx:814
 AliITSUTrackerGlo.cxx:815
 AliITSUTrackerGlo.cxx:816
 AliITSUTrackerGlo.cxx:817
 AliITSUTrackerGlo.cxx:818
 AliITSUTrackerGlo.cxx:819
 AliITSUTrackerGlo.cxx:820
 AliITSUTrackerGlo.cxx:821
 AliITSUTrackerGlo.cxx:822
 AliITSUTrackerGlo.cxx:823
 AliITSUTrackerGlo.cxx:824
 AliITSUTrackerGlo.cxx:825
 AliITSUTrackerGlo.cxx:826
 AliITSUTrackerGlo.cxx:827
 AliITSUTrackerGlo.cxx:828
 AliITSUTrackerGlo.cxx:829
 AliITSUTrackerGlo.cxx:830
 AliITSUTrackerGlo.cxx:831
 AliITSUTrackerGlo.cxx:832
 AliITSUTrackerGlo.cxx:833
 AliITSUTrackerGlo.cxx:834
 AliITSUTrackerGlo.cxx:835
 AliITSUTrackerGlo.cxx:836
 AliITSUTrackerGlo.cxx:837
 AliITSUTrackerGlo.cxx:838
 AliITSUTrackerGlo.cxx:839
 AliITSUTrackerGlo.cxx:840
 AliITSUTrackerGlo.cxx:841
 AliITSUTrackerGlo.cxx:842
 AliITSUTrackerGlo.cxx:843
 AliITSUTrackerGlo.cxx:844
 AliITSUTrackerGlo.cxx:845
 AliITSUTrackerGlo.cxx:846
 AliITSUTrackerGlo.cxx:847
 AliITSUTrackerGlo.cxx:848
 AliITSUTrackerGlo.cxx:849
 AliITSUTrackerGlo.cxx:850
 AliITSUTrackerGlo.cxx:851
 AliITSUTrackerGlo.cxx:852
 AliITSUTrackerGlo.cxx:853
 AliITSUTrackerGlo.cxx:854
 AliITSUTrackerGlo.cxx:855
 AliITSUTrackerGlo.cxx:856
 AliITSUTrackerGlo.cxx:857
 AliITSUTrackerGlo.cxx:858
 AliITSUTrackerGlo.cxx:859
 AliITSUTrackerGlo.cxx:860
 AliITSUTrackerGlo.cxx:861
 AliITSUTrackerGlo.cxx:862
 AliITSUTrackerGlo.cxx:863
 AliITSUTrackerGlo.cxx:864
 AliITSUTrackerGlo.cxx:865
 AliITSUTrackerGlo.cxx:866
 AliITSUTrackerGlo.cxx:867
 AliITSUTrackerGlo.cxx:868
 AliITSUTrackerGlo.cxx:869
 AliITSUTrackerGlo.cxx:870
 AliITSUTrackerGlo.cxx:871
 AliITSUTrackerGlo.cxx:872
 AliITSUTrackerGlo.cxx:873
 AliITSUTrackerGlo.cxx:874
 AliITSUTrackerGlo.cxx:875
 AliITSUTrackerGlo.cxx:876
 AliITSUTrackerGlo.cxx:877
 AliITSUTrackerGlo.cxx:878
 AliITSUTrackerGlo.cxx:879
 AliITSUTrackerGlo.cxx:880
 AliITSUTrackerGlo.cxx:881
 AliITSUTrackerGlo.cxx:882
 AliITSUTrackerGlo.cxx:883
 AliITSUTrackerGlo.cxx:884
 AliITSUTrackerGlo.cxx:885
 AliITSUTrackerGlo.cxx:886
 AliITSUTrackerGlo.cxx:887
 AliITSUTrackerGlo.cxx:888
 AliITSUTrackerGlo.cxx:889
 AliITSUTrackerGlo.cxx:890
 AliITSUTrackerGlo.cxx:891
 AliITSUTrackerGlo.cxx:892
 AliITSUTrackerGlo.cxx:893
 AliITSUTrackerGlo.cxx:894
 AliITSUTrackerGlo.cxx:895
 AliITSUTrackerGlo.cxx:896
 AliITSUTrackerGlo.cxx:897
 AliITSUTrackerGlo.cxx:898
 AliITSUTrackerGlo.cxx:899
 AliITSUTrackerGlo.cxx:900
 AliITSUTrackerGlo.cxx:901
 AliITSUTrackerGlo.cxx:902
 AliITSUTrackerGlo.cxx:903
 AliITSUTrackerGlo.cxx:904
 AliITSUTrackerGlo.cxx:905
 AliITSUTrackerGlo.cxx:906
 AliITSUTrackerGlo.cxx:907
 AliITSUTrackerGlo.cxx:908
 AliITSUTrackerGlo.cxx:909
 AliITSUTrackerGlo.cxx:910
 AliITSUTrackerGlo.cxx:911
 AliITSUTrackerGlo.cxx:912
 AliITSUTrackerGlo.cxx:913
 AliITSUTrackerGlo.cxx:914
 AliITSUTrackerGlo.cxx:915
 AliITSUTrackerGlo.cxx:916
 AliITSUTrackerGlo.cxx:917
 AliITSUTrackerGlo.cxx:918
 AliITSUTrackerGlo.cxx:919
 AliITSUTrackerGlo.cxx:920
 AliITSUTrackerGlo.cxx:921
 AliITSUTrackerGlo.cxx:922
 AliITSUTrackerGlo.cxx:923
 AliITSUTrackerGlo.cxx:924
 AliITSUTrackerGlo.cxx:925
 AliITSUTrackerGlo.cxx:926
 AliITSUTrackerGlo.cxx:927
 AliITSUTrackerGlo.cxx:928
 AliITSUTrackerGlo.cxx:929
 AliITSUTrackerGlo.cxx:930
 AliITSUTrackerGlo.cxx:931
 AliITSUTrackerGlo.cxx:932
 AliITSUTrackerGlo.cxx:933
 AliITSUTrackerGlo.cxx:934
 AliITSUTrackerGlo.cxx:935
 AliITSUTrackerGlo.cxx:936
 AliITSUTrackerGlo.cxx:937
 AliITSUTrackerGlo.cxx:938
 AliITSUTrackerGlo.cxx:939
 AliITSUTrackerGlo.cxx:940
 AliITSUTrackerGlo.cxx:941
 AliITSUTrackerGlo.cxx:942
 AliITSUTrackerGlo.cxx:943
 AliITSUTrackerGlo.cxx:944
 AliITSUTrackerGlo.cxx:945
 AliITSUTrackerGlo.cxx:946
 AliITSUTrackerGlo.cxx:947
 AliITSUTrackerGlo.cxx:948
 AliITSUTrackerGlo.cxx:949
 AliITSUTrackerGlo.cxx:950
 AliITSUTrackerGlo.cxx:951
 AliITSUTrackerGlo.cxx:952
 AliITSUTrackerGlo.cxx:953
 AliITSUTrackerGlo.cxx:954
 AliITSUTrackerGlo.cxx:955
 AliITSUTrackerGlo.cxx:956
 AliITSUTrackerGlo.cxx:957
 AliITSUTrackerGlo.cxx:958
 AliITSUTrackerGlo.cxx:959
 AliITSUTrackerGlo.cxx:960
 AliITSUTrackerGlo.cxx:961
 AliITSUTrackerGlo.cxx:962
 AliITSUTrackerGlo.cxx:963
 AliITSUTrackerGlo.cxx:964
 AliITSUTrackerGlo.cxx:965
 AliITSUTrackerGlo.cxx:966
 AliITSUTrackerGlo.cxx:967
 AliITSUTrackerGlo.cxx:968
 AliITSUTrackerGlo.cxx:969
 AliITSUTrackerGlo.cxx:970
 AliITSUTrackerGlo.cxx:971
 AliITSUTrackerGlo.cxx:972
 AliITSUTrackerGlo.cxx:973
 AliITSUTrackerGlo.cxx:974
 AliITSUTrackerGlo.cxx:975
 AliITSUTrackerGlo.cxx:976
 AliITSUTrackerGlo.cxx:977
 AliITSUTrackerGlo.cxx:978
 AliITSUTrackerGlo.cxx:979
 AliITSUTrackerGlo.cxx:980
 AliITSUTrackerGlo.cxx:981
 AliITSUTrackerGlo.cxx:982
 AliITSUTrackerGlo.cxx:983
 AliITSUTrackerGlo.cxx:984
 AliITSUTrackerGlo.cxx:985
 AliITSUTrackerGlo.cxx:986
 AliITSUTrackerGlo.cxx:987
 AliITSUTrackerGlo.cxx:988
 AliITSUTrackerGlo.cxx:989
 AliITSUTrackerGlo.cxx:990
 AliITSUTrackerGlo.cxx:991
 AliITSUTrackerGlo.cxx:992
 AliITSUTrackerGlo.cxx:993
 AliITSUTrackerGlo.cxx:994
 AliITSUTrackerGlo.cxx:995
 AliITSUTrackerGlo.cxx:996
 AliITSUTrackerGlo.cxx:997
 AliITSUTrackerGlo.cxx:998
 AliITSUTrackerGlo.cxx:999
 AliITSUTrackerGlo.cxx:1000
 AliITSUTrackerGlo.cxx:1001
 AliITSUTrackerGlo.cxx:1002
 AliITSUTrackerGlo.cxx:1003
 AliITSUTrackerGlo.cxx:1004
 AliITSUTrackerGlo.cxx:1005
 AliITSUTrackerGlo.cxx:1006
 AliITSUTrackerGlo.cxx:1007
 AliITSUTrackerGlo.cxx:1008
 AliITSUTrackerGlo.cxx:1009
 AliITSUTrackerGlo.cxx:1010
 AliITSUTrackerGlo.cxx:1011
 AliITSUTrackerGlo.cxx:1012
 AliITSUTrackerGlo.cxx:1013
 AliITSUTrackerGlo.cxx:1014
 AliITSUTrackerGlo.cxx:1015
 AliITSUTrackerGlo.cxx:1016
 AliITSUTrackerGlo.cxx:1017
 AliITSUTrackerGlo.cxx:1018
 AliITSUTrackerGlo.cxx:1019
 AliITSUTrackerGlo.cxx:1020
 AliITSUTrackerGlo.cxx:1021
 AliITSUTrackerGlo.cxx:1022
 AliITSUTrackerGlo.cxx:1023
 AliITSUTrackerGlo.cxx:1024
 AliITSUTrackerGlo.cxx:1025
 AliITSUTrackerGlo.cxx:1026
 AliITSUTrackerGlo.cxx:1027
 AliITSUTrackerGlo.cxx:1028
 AliITSUTrackerGlo.cxx:1029
 AliITSUTrackerGlo.cxx:1030
 AliITSUTrackerGlo.cxx:1031
 AliITSUTrackerGlo.cxx:1032
 AliITSUTrackerGlo.cxx:1033
 AliITSUTrackerGlo.cxx:1034
 AliITSUTrackerGlo.cxx:1035
 AliITSUTrackerGlo.cxx:1036
 AliITSUTrackerGlo.cxx:1037
 AliITSUTrackerGlo.cxx:1038
 AliITSUTrackerGlo.cxx:1039
 AliITSUTrackerGlo.cxx:1040
 AliITSUTrackerGlo.cxx:1041
 AliITSUTrackerGlo.cxx:1042
 AliITSUTrackerGlo.cxx:1043
 AliITSUTrackerGlo.cxx:1044
 AliITSUTrackerGlo.cxx:1045
 AliITSUTrackerGlo.cxx:1046
 AliITSUTrackerGlo.cxx:1047
 AliITSUTrackerGlo.cxx:1048
 AliITSUTrackerGlo.cxx:1049
 AliITSUTrackerGlo.cxx:1050
 AliITSUTrackerGlo.cxx:1051
 AliITSUTrackerGlo.cxx:1052
 AliITSUTrackerGlo.cxx:1053
 AliITSUTrackerGlo.cxx:1054
 AliITSUTrackerGlo.cxx:1055
 AliITSUTrackerGlo.cxx:1056
 AliITSUTrackerGlo.cxx:1057
 AliITSUTrackerGlo.cxx:1058
 AliITSUTrackerGlo.cxx:1059
 AliITSUTrackerGlo.cxx:1060
 AliITSUTrackerGlo.cxx:1061
 AliITSUTrackerGlo.cxx:1062
 AliITSUTrackerGlo.cxx:1063
 AliITSUTrackerGlo.cxx:1064
 AliITSUTrackerGlo.cxx:1065
 AliITSUTrackerGlo.cxx:1066
 AliITSUTrackerGlo.cxx:1067
 AliITSUTrackerGlo.cxx:1068
 AliITSUTrackerGlo.cxx:1069
 AliITSUTrackerGlo.cxx:1070
 AliITSUTrackerGlo.cxx:1071
 AliITSUTrackerGlo.cxx:1072
 AliITSUTrackerGlo.cxx:1073
 AliITSUTrackerGlo.cxx:1074
 AliITSUTrackerGlo.cxx:1075
 AliITSUTrackerGlo.cxx:1076
 AliITSUTrackerGlo.cxx:1077
 AliITSUTrackerGlo.cxx:1078
 AliITSUTrackerGlo.cxx:1079
 AliITSUTrackerGlo.cxx:1080
 AliITSUTrackerGlo.cxx:1081
 AliITSUTrackerGlo.cxx:1082
 AliITSUTrackerGlo.cxx:1083
 AliITSUTrackerGlo.cxx:1084
 AliITSUTrackerGlo.cxx:1085
 AliITSUTrackerGlo.cxx:1086
 AliITSUTrackerGlo.cxx:1087
 AliITSUTrackerGlo.cxx:1088
 AliITSUTrackerGlo.cxx:1089
 AliITSUTrackerGlo.cxx:1090
 AliITSUTrackerGlo.cxx:1091
 AliITSUTrackerGlo.cxx:1092
 AliITSUTrackerGlo.cxx:1093
 AliITSUTrackerGlo.cxx:1094
 AliITSUTrackerGlo.cxx:1095
 AliITSUTrackerGlo.cxx:1096
 AliITSUTrackerGlo.cxx:1097
 AliITSUTrackerGlo.cxx:1098
 AliITSUTrackerGlo.cxx:1099
 AliITSUTrackerGlo.cxx:1100
 AliITSUTrackerGlo.cxx:1101
 AliITSUTrackerGlo.cxx:1102
 AliITSUTrackerGlo.cxx:1103
 AliITSUTrackerGlo.cxx:1104
 AliITSUTrackerGlo.cxx:1105
 AliITSUTrackerGlo.cxx:1106
 AliITSUTrackerGlo.cxx:1107
 AliITSUTrackerGlo.cxx:1108
 AliITSUTrackerGlo.cxx:1109
 AliITSUTrackerGlo.cxx:1110
 AliITSUTrackerGlo.cxx:1111
 AliITSUTrackerGlo.cxx:1112
 AliITSUTrackerGlo.cxx:1113
 AliITSUTrackerGlo.cxx:1114
 AliITSUTrackerGlo.cxx:1115
 AliITSUTrackerGlo.cxx:1116
 AliITSUTrackerGlo.cxx:1117
 AliITSUTrackerGlo.cxx:1118
 AliITSUTrackerGlo.cxx:1119
 AliITSUTrackerGlo.cxx:1120
 AliITSUTrackerGlo.cxx:1121
 AliITSUTrackerGlo.cxx:1122
 AliITSUTrackerGlo.cxx:1123
 AliITSUTrackerGlo.cxx:1124
 AliITSUTrackerGlo.cxx:1125
 AliITSUTrackerGlo.cxx:1126
 AliITSUTrackerGlo.cxx:1127
 AliITSUTrackerGlo.cxx:1128
 AliITSUTrackerGlo.cxx:1129
 AliITSUTrackerGlo.cxx:1130
 AliITSUTrackerGlo.cxx:1131
 AliITSUTrackerGlo.cxx:1132
 AliITSUTrackerGlo.cxx:1133
 AliITSUTrackerGlo.cxx:1134
 AliITSUTrackerGlo.cxx:1135
 AliITSUTrackerGlo.cxx:1136
 AliITSUTrackerGlo.cxx:1137
 AliITSUTrackerGlo.cxx:1138
 AliITSUTrackerGlo.cxx:1139
 AliITSUTrackerGlo.cxx:1140
 AliITSUTrackerGlo.cxx:1141
 AliITSUTrackerGlo.cxx:1142
 AliITSUTrackerGlo.cxx:1143
 AliITSUTrackerGlo.cxx:1144
 AliITSUTrackerGlo.cxx:1145
 AliITSUTrackerGlo.cxx:1146
 AliITSUTrackerGlo.cxx:1147
 AliITSUTrackerGlo.cxx:1148
 AliITSUTrackerGlo.cxx:1149
 AliITSUTrackerGlo.cxx:1150
 AliITSUTrackerGlo.cxx:1151
 AliITSUTrackerGlo.cxx:1152
 AliITSUTrackerGlo.cxx:1153
 AliITSUTrackerGlo.cxx:1154
 AliITSUTrackerGlo.cxx:1155
 AliITSUTrackerGlo.cxx:1156
 AliITSUTrackerGlo.cxx:1157
 AliITSUTrackerGlo.cxx:1158
 AliITSUTrackerGlo.cxx:1159
 AliITSUTrackerGlo.cxx:1160
 AliITSUTrackerGlo.cxx:1161
 AliITSUTrackerGlo.cxx:1162
 AliITSUTrackerGlo.cxx:1163
 AliITSUTrackerGlo.cxx:1164
 AliITSUTrackerGlo.cxx:1165
 AliITSUTrackerGlo.cxx:1166
 AliITSUTrackerGlo.cxx:1167
 AliITSUTrackerGlo.cxx:1168
 AliITSUTrackerGlo.cxx:1169
 AliITSUTrackerGlo.cxx:1170
 AliITSUTrackerGlo.cxx:1171
 AliITSUTrackerGlo.cxx:1172
 AliITSUTrackerGlo.cxx:1173
 AliITSUTrackerGlo.cxx:1174
 AliITSUTrackerGlo.cxx:1175
 AliITSUTrackerGlo.cxx:1176
 AliITSUTrackerGlo.cxx:1177
 AliITSUTrackerGlo.cxx:1178
 AliITSUTrackerGlo.cxx:1179
 AliITSUTrackerGlo.cxx:1180
 AliITSUTrackerGlo.cxx:1181
 AliITSUTrackerGlo.cxx:1182
 AliITSUTrackerGlo.cxx:1183
 AliITSUTrackerGlo.cxx:1184
 AliITSUTrackerGlo.cxx:1185
 AliITSUTrackerGlo.cxx:1186
 AliITSUTrackerGlo.cxx:1187
 AliITSUTrackerGlo.cxx:1188
 AliITSUTrackerGlo.cxx:1189
 AliITSUTrackerGlo.cxx:1190
 AliITSUTrackerGlo.cxx:1191
 AliITSUTrackerGlo.cxx:1192
 AliITSUTrackerGlo.cxx:1193
 AliITSUTrackerGlo.cxx:1194
 AliITSUTrackerGlo.cxx:1195
 AliITSUTrackerGlo.cxx:1196
 AliITSUTrackerGlo.cxx:1197
 AliITSUTrackerGlo.cxx:1198
 AliITSUTrackerGlo.cxx:1199
 AliITSUTrackerGlo.cxx:1200
 AliITSUTrackerGlo.cxx:1201
 AliITSUTrackerGlo.cxx:1202
 AliITSUTrackerGlo.cxx:1203
 AliITSUTrackerGlo.cxx:1204
 AliITSUTrackerGlo.cxx:1205
 AliITSUTrackerGlo.cxx:1206
 AliITSUTrackerGlo.cxx:1207
 AliITSUTrackerGlo.cxx:1208
 AliITSUTrackerGlo.cxx:1209
 AliITSUTrackerGlo.cxx:1210
 AliITSUTrackerGlo.cxx:1211
 AliITSUTrackerGlo.cxx:1212
 AliITSUTrackerGlo.cxx:1213
 AliITSUTrackerGlo.cxx:1214
 AliITSUTrackerGlo.cxx:1215
 AliITSUTrackerGlo.cxx:1216
 AliITSUTrackerGlo.cxx:1217
 AliITSUTrackerGlo.cxx:1218
 AliITSUTrackerGlo.cxx:1219
 AliITSUTrackerGlo.cxx:1220
 AliITSUTrackerGlo.cxx:1221
 AliITSUTrackerGlo.cxx:1222
 AliITSUTrackerGlo.cxx:1223
 AliITSUTrackerGlo.cxx:1224
 AliITSUTrackerGlo.cxx:1225
 AliITSUTrackerGlo.cxx:1226
 AliITSUTrackerGlo.cxx:1227
 AliITSUTrackerGlo.cxx:1228
 AliITSUTrackerGlo.cxx:1229
 AliITSUTrackerGlo.cxx:1230
 AliITSUTrackerGlo.cxx:1231
 AliITSUTrackerGlo.cxx:1232
 AliITSUTrackerGlo.cxx:1233
 AliITSUTrackerGlo.cxx:1234
 AliITSUTrackerGlo.cxx:1235
 AliITSUTrackerGlo.cxx:1236
 AliITSUTrackerGlo.cxx:1237
 AliITSUTrackerGlo.cxx:1238
 AliITSUTrackerGlo.cxx:1239
 AliITSUTrackerGlo.cxx:1240
 AliITSUTrackerGlo.cxx:1241
 AliITSUTrackerGlo.cxx:1242
 AliITSUTrackerGlo.cxx:1243
 AliITSUTrackerGlo.cxx:1244
 AliITSUTrackerGlo.cxx:1245
 AliITSUTrackerGlo.cxx:1246
 AliITSUTrackerGlo.cxx:1247
 AliITSUTrackerGlo.cxx:1248
 AliITSUTrackerGlo.cxx:1249
 AliITSUTrackerGlo.cxx:1250
 AliITSUTrackerGlo.cxx:1251
 AliITSUTrackerGlo.cxx:1252
 AliITSUTrackerGlo.cxx:1253
 AliITSUTrackerGlo.cxx:1254
 AliITSUTrackerGlo.cxx:1255
 AliITSUTrackerGlo.cxx:1256
 AliITSUTrackerGlo.cxx:1257
 AliITSUTrackerGlo.cxx:1258
 AliITSUTrackerGlo.cxx:1259
 AliITSUTrackerGlo.cxx:1260
 AliITSUTrackerGlo.cxx:1261
 AliITSUTrackerGlo.cxx:1262
 AliITSUTrackerGlo.cxx:1263
 AliITSUTrackerGlo.cxx:1264
 AliITSUTrackerGlo.cxx:1265
 AliITSUTrackerGlo.cxx:1266
 AliITSUTrackerGlo.cxx:1267
 AliITSUTrackerGlo.cxx:1268
 AliITSUTrackerGlo.cxx:1269
 AliITSUTrackerGlo.cxx:1270
 AliITSUTrackerGlo.cxx:1271
 AliITSUTrackerGlo.cxx:1272
 AliITSUTrackerGlo.cxx:1273
 AliITSUTrackerGlo.cxx:1274
 AliITSUTrackerGlo.cxx:1275
 AliITSUTrackerGlo.cxx:1276
 AliITSUTrackerGlo.cxx:1277
 AliITSUTrackerGlo.cxx:1278
 AliITSUTrackerGlo.cxx:1279
 AliITSUTrackerGlo.cxx:1280
 AliITSUTrackerGlo.cxx:1281
 AliITSUTrackerGlo.cxx:1282
 AliITSUTrackerGlo.cxx:1283
 AliITSUTrackerGlo.cxx:1284
 AliITSUTrackerGlo.cxx:1285
 AliITSUTrackerGlo.cxx:1286
 AliITSUTrackerGlo.cxx:1287
 AliITSUTrackerGlo.cxx:1288
 AliITSUTrackerGlo.cxx:1289
 AliITSUTrackerGlo.cxx:1290
 AliITSUTrackerGlo.cxx:1291
 AliITSUTrackerGlo.cxx:1292
 AliITSUTrackerGlo.cxx:1293
 AliITSUTrackerGlo.cxx:1294
 AliITSUTrackerGlo.cxx:1295
 AliITSUTrackerGlo.cxx:1296
 AliITSUTrackerGlo.cxx:1297
 AliITSUTrackerGlo.cxx:1298
 AliITSUTrackerGlo.cxx:1299
 AliITSUTrackerGlo.cxx:1300
 AliITSUTrackerGlo.cxx:1301
 AliITSUTrackerGlo.cxx:1302
 AliITSUTrackerGlo.cxx:1303
 AliITSUTrackerGlo.cxx:1304
 AliITSUTrackerGlo.cxx:1305
 AliITSUTrackerGlo.cxx:1306
 AliITSUTrackerGlo.cxx:1307
 AliITSUTrackerGlo.cxx:1308
 AliITSUTrackerGlo.cxx:1309
 AliITSUTrackerGlo.cxx:1310
 AliITSUTrackerGlo.cxx:1311
 AliITSUTrackerGlo.cxx:1312
 AliITSUTrackerGlo.cxx:1313
 AliITSUTrackerGlo.cxx:1314
 AliITSUTrackerGlo.cxx:1315
 AliITSUTrackerGlo.cxx:1316
 AliITSUTrackerGlo.cxx:1317
 AliITSUTrackerGlo.cxx:1318
 AliITSUTrackerGlo.cxx:1319
 AliITSUTrackerGlo.cxx:1320
 AliITSUTrackerGlo.cxx:1321
 AliITSUTrackerGlo.cxx:1322
 AliITSUTrackerGlo.cxx:1323
 AliITSUTrackerGlo.cxx:1324
 AliITSUTrackerGlo.cxx:1325
 AliITSUTrackerGlo.cxx:1326
 AliITSUTrackerGlo.cxx:1327
 AliITSUTrackerGlo.cxx:1328
 AliITSUTrackerGlo.cxx:1329
 AliITSUTrackerGlo.cxx:1330
 AliITSUTrackerGlo.cxx:1331
 AliITSUTrackerGlo.cxx:1332
 AliITSUTrackerGlo.cxx:1333
 AliITSUTrackerGlo.cxx:1334
 AliITSUTrackerGlo.cxx:1335
 AliITSUTrackerGlo.cxx:1336
 AliITSUTrackerGlo.cxx:1337
 AliITSUTrackerGlo.cxx:1338
 AliITSUTrackerGlo.cxx:1339
 AliITSUTrackerGlo.cxx:1340
 AliITSUTrackerGlo.cxx:1341
 AliITSUTrackerGlo.cxx:1342
 AliITSUTrackerGlo.cxx:1343
 AliITSUTrackerGlo.cxx:1344
 AliITSUTrackerGlo.cxx:1345
 AliITSUTrackerGlo.cxx:1346
 AliITSUTrackerGlo.cxx:1347
 AliITSUTrackerGlo.cxx:1348
 AliITSUTrackerGlo.cxx:1349
 AliITSUTrackerGlo.cxx:1350
 AliITSUTrackerGlo.cxx:1351
 AliITSUTrackerGlo.cxx:1352
 AliITSUTrackerGlo.cxx:1353
 AliITSUTrackerGlo.cxx:1354
 AliITSUTrackerGlo.cxx:1355
 AliITSUTrackerGlo.cxx:1356
 AliITSUTrackerGlo.cxx:1357
 AliITSUTrackerGlo.cxx:1358
 AliITSUTrackerGlo.cxx:1359
 AliITSUTrackerGlo.cxx:1360
 AliITSUTrackerGlo.cxx:1361
 AliITSUTrackerGlo.cxx:1362
 AliITSUTrackerGlo.cxx:1363
 AliITSUTrackerGlo.cxx:1364
 AliITSUTrackerGlo.cxx:1365
 AliITSUTrackerGlo.cxx:1366
 AliITSUTrackerGlo.cxx:1367
 AliITSUTrackerGlo.cxx:1368
 AliITSUTrackerGlo.cxx:1369
 AliITSUTrackerGlo.cxx:1370
 AliITSUTrackerGlo.cxx:1371
 AliITSUTrackerGlo.cxx:1372
 AliITSUTrackerGlo.cxx:1373
 AliITSUTrackerGlo.cxx:1374
 AliITSUTrackerGlo.cxx:1375
 AliITSUTrackerGlo.cxx:1376
 AliITSUTrackerGlo.cxx:1377
 AliITSUTrackerGlo.cxx:1378
 AliITSUTrackerGlo.cxx:1379
 AliITSUTrackerGlo.cxx:1380
 AliITSUTrackerGlo.cxx:1381
 AliITSUTrackerGlo.cxx:1382
 AliITSUTrackerGlo.cxx:1383
 AliITSUTrackerGlo.cxx:1384
 AliITSUTrackerGlo.cxx:1385
 AliITSUTrackerGlo.cxx:1386
 AliITSUTrackerGlo.cxx:1387
 AliITSUTrackerGlo.cxx:1388
 AliITSUTrackerGlo.cxx:1389
 AliITSUTrackerGlo.cxx:1390
 AliITSUTrackerGlo.cxx:1391
 AliITSUTrackerGlo.cxx:1392
 AliITSUTrackerGlo.cxx:1393
 AliITSUTrackerGlo.cxx:1394
 AliITSUTrackerGlo.cxx:1395
 AliITSUTrackerGlo.cxx:1396
 AliITSUTrackerGlo.cxx:1397
 AliITSUTrackerGlo.cxx:1398
 AliITSUTrackerGlo.cxx:1399
 AliITSUTrackerGlo.cxx:1400
 AliITSUTrackerGlo.cxx:1401
 AliITSUTrackerGlo.cxx:1402
 AliITSUTrackerGlo.cxx:1403
 AliITSUTrackerGlo.cxx:1404
 AliITSUTrackerGlo.cxx:1405
 AliITSUTrackerGlo.cxx:1406
 AliITSUTrackerGlo.cxx:1407
 AliITSUTrackerGlo.cxx:1408
 AliITSUTrackerGlo.cxx:1409
 AliITSUTrackerGlo.cxx:1410
 AliITSUTrackerGlo.cxx:1411
 AliITSUTrackerGlo.cxx:1412
 AliITSUTrackerGlo.cxx:1413
 AliITSUTrackerGlo.cxx:1414
 AliITSUTrackerGlo.cxx:1415
 AliITSUTrackerGlo.cxx:1416
 AliITSUTrackerGlo.cxx:1417
 AliITSUTrackerGlo.cxx:1418
 AliITSUTrackerGlo.cxx:1419
 AliITSUTrackerGlo.cxx:1420
 AliITSUTrackerGlo.cxx:1421
 AliITSUTrackerGlo.cxx:1422
 AliITSUTrackerGlo.cxx:1423
 AliITSUTrackerGlo.cxx:1424
 AliITSUTrackerGlo.cxx:1425
 AliITSUTrackerGlo.cxx:1426
 AliITSUTrackerGlo.cxx:1427
 AliITSUTrackerGlo.cxx:1428
 AliITSUTrackerGlo.cxx:1429
 AliITSUTrackerGlo.cxx:1430
 AliITSUTrackerGlo.cxx:1431
 AliITSUTrackerGlo.cxx:1432
 AliITSUTrackerGlo.cxx:1433
 AliITSUTrackerGlo.cxx:1434
 AliITSUTrackerGlo.cxx:1435
 AliITSUTrackerGlo.cxx:1436
 AliITSUTrackerGlo.cxx:1437
 AliITSUTrackerGlo.cxx:1438
 AliITSUTrackerGlo.cxx:1439
 AliITSUTrackerGlo.cxx:1440
 AliITSUTrackerGlo.cxx:1441
 AliITSUTrackerGlo.cxx:1442
 AliITSUTrackerGlo.cxx:1443
 AliITSUTrackerGlo.cxx:1444
 AliITSUTrackerGlo.cxx:1445
 AliITSUTrackerGlo.cxx:1446
 AliITSUTrackerGlo.cxx:1447
 AliITSUTrackerGlo.cxx:1448
 AliITSUTrackerGlo.cxx:1449
 AliITSUTrackerGlo.cxx:1450
 AliITSUTrackerGlo.cxx:1451
 AliITSUTrackerGlo.cxx:1452
 AliITSUTrackerGlo.cxx:1453
 AliITSUTrackerGlo.cxx:1454
 AliITSUTrackerGlo.cxx:1455
 AliITSUTrackerGlo.cxx:1456
 AliITSUTrackerGlo.cxx:1457
 AliITSUTrackerGlo.cxx:1458
 AliITSUTrackerGlo.cxx:1459
 AliITSUTrackerGlo.cxx:1460
 AliITSUTrackerGlo.cxx:1461
 AliITSUTrackerGlo.cxx:1462
 AliITSUTrackerGlo.cxx:1463
 AliITSUTrackerGlo.cxx:1464
 AliITSUTrackerGlo.cxx:1465
 AliITSUTrackerGlo.cxx:1466
 AliITSUTrackerGlo.cxx:1467
 AliITSUTrackerGlo.cxx:1468
 AliITSUTrackerGlo.cxx:1469
 AliITSUTrackerGlo.cxx:1470
 AliITSUTrackerGlo.cxx:1471
 AliITSUTrackerGlo.cxx:1472
 AliITSUTrackerGlo.cxx:1473
 AliITSUTrackerGlo.cxx:1474
 AliITSUTrackerGlo.cxx:1475
 AliITSUTrackerGlo.cxx:1476
 AliITSUTrackerGlo.cxx:1477
 AliITSUTrackerGlo.cxx:1478
 AliITSUTrackerGlo.cxx:1479
 AliITSUTrackerGlo.cxx:1480
 AliITSUTrackerGlo.cxx:1481
 AliITSUTrackerGlo.cxx:1482
 AliITSUTrackerGlo.cxx:1483
 AliITSUTrackerGlo.cxx:1484
 AliITSUTrackerGlo.cxx:1485
 AliITSUTrackerGlo.cxx:1486
 AliITSUTrackerGlo.cxx:1487
 AliITSUTrackerGlo.cxx:1488
 AliITSUTrackerGlo.cxx:1489
 AliITSUTrackerGlo.cxx:1490
 AliITSUTrackerGlo.cxx:1491
 AliITSUTrackerGlo.cxx:1492
 AliITSUTrackerGlo.cxx:1493
 AliITSUTrackerGlo.cxx:1494
 AliITSUTrackerGlo.cxx:1495
 AliITSUTrackerGlo.cxx:1496
 AliITSUTrackerGlo.cxx:1497
 AliITSUTrackerGlo.cxx:1498
 AliITSUTrackerGlo.cxx:1499
 AliITSUTrackerGlo.cxx:1500
 AliITSUTrackerGlo.cxx:1501
 AliITSUTrackerGlo.cxx:1502
 AliITSUTrackerGlo.cxx:1503
 AliITSUTrackerGlo.cxx:1504
 AliITSUTrackerGlo.cxx:1505
 AliITSUTrackerGlo.cxx:1506
 AliITSUTrackerGlo.cxx:1507
 AliITSUTrackerGlo.cxx:1508
 AliITSUTrackerGlo.cxx:1509
 AliITSUTrackerGlo.cxx:1510
 AliITSUTrackerGlo.cxx:1511
 AliITSUTrackerGlo.cxx:1512
 AliITSUTrackerGlo.cxx:1513
 AliITSUTrackerGlo.cxx:1514
 AliITSUTrackerGlo.cxx:1515
 AliITSUTrackerGlo.cxx:1516
 AliITSUTrackerGlo.cxx:1517
 AliITSUTrackerGlo.cxx:1518
 AliITSUTrackerGlo.cxx:1519
 AliITSUTrackerGlo.cxx:1520
 AliITSUTrackerGlo.cxx:1521
 AliITSUTrackerGlo.cxx:1522
 AliITSUTrackerGlo.cxx:1523
 AliITSUTrackerGlo.cxx:1524
 AliITSUTrackerGlo.cxx:1525
 AliITSUTrackerGlo.cxx:1526
 AliITSUTrackerGlo.cxx:1527
 AliITSUTrackerGlo.cxx:1528
 AliITSUTrackerGlo.cxx:1529
 AliITSUTrackerGlo.cxx:1530
 AliITSUTrackerGlo.cxx:1531
 AliITSUTrackerGlo.cxx:1532
 AliITSUTrackerGlo.cxx:1533
 AliITSUTrackerGlo.cxx:1534
 AliITSUTrackerGlo.cxx:1535
 AliITSUTrackerGlo.cxx:1536
 AliITSUTrackerGlo.cxx:1537
 AliITSUTrackerGlo.cxx:1538
 AliITSUTrackerGlo.cxx:1539
 AliITSUTrackerGlo.cxx:1540
 AliITSUTrackerGlo.cxx:1541
 AliITSUTrackerGlo.cxx:1542
 AliITSUTrackerGlo.cxx:1543
 AliITSUTrackerGlo.cxx:1544
 AliITSUTrackerGlo.cxx:1545
 AliITSUTrackerGlo.cxx:1546
 AliITSUTrackerGlo.cxx:1547
 AliITSUTrackerGlo.cxx:1548
 AliITSUTrackerGlo.cxx:1549
 AliITSUTrackerGlo.cxx:1550
 AliITSUTrackerGlo.cxx:1551
 AliITSUTrackerGlo.cxx:1552
 AliITSUTrackerGlo.cxx:1553
 AliITSUTrackerGlo.cxx:1554
 AliITSUTrackerGlo.cxx:1555
 AliITSUTrackerGlo.cxx:1556
 AliITSUTrackerGlo.cxx:1557
 AliITSUTrackerGlo.cxx:1558
 AliITSUTrackerGlo.cxx:1559
 AliITSUTrackerGlo.cxx:1560
 AliITSUTrackerGlo.cxx:1561
 AliITSUTrackerGlo.cxx:1562
 AliITSUTrackerGlo.cxx:1563
 AliITSUTrackerGlo.cxx:1564
 AliITSUTrackerGlo.cxx:1565
 AliITSUTrackerGlo.cxx:1566
 AliITSUTrackerGlo.cxx:1567
 AliITSUTrackerGlo.cxx:1568
 AliITSUTrackerGlo.cxx:1569
 AliITSUTrackerGlo.cxx:1570
 AliITSUTrackerGlo.cxx:1571
 AliITSUTrackerGlo.cxx:1572
 AliITSUTrackerGlo.cxx:1573
 AliITSUTrackerGlo.cxx:1574
 AliITSUTrackerGlo.cxx:1575
 AliITSUTrackerGlo.cxx:1576
 AliITSUTrackerGlo.cxx:1577
 AliITSUTrackerGlo.cxx:1578
 AliITSUTrackerGlo.cxx:1579
 AliITSUTrackerGlo.cxx:1580
 AliITSUTrackerGlo.cxx:1581
 AliITSUTrackerGlo.cxx:1582
 AliITSUTrackerGlo.cxx:1583
 AliITSUTrackerGlo.cxx:1584
 AliITSUTrackerGlo.cxx:1585
 AliITSUTrackerGlo.cxx:1586
 AliITSUTrackerGlo.cxx:1587
 AliITSUTrackerGlo.cxx:1588
 AliITSUTrackerGlo.cxx:1589
 AliITSUTrackerGlo.cxx:1590
 AliITSUTrackerGlo.cxx:1591
 AliITSUTrackerGlo.cxx:1592
 AliITSUTrackerGlo.cxx:1593
 AliITSUTrackerGlo.cxx:1594
 AliITSUTrackerGlo.cxx:1595
 AliITSUTrackerGlo.cxx:1596
 AliITSUTrackerGlo.cxx:1597
 AliITSUTrackerGlo.cxx:1598
 AliITSUTrackerGlo.cxx:1599
 AliITSUTrackerGlo.cxx:1600
 AliITSUTrackerGlo.cxx:1601
 AliITSUTrackerGlo.cxx:1602
 AliITSUTrackerGlo.cxx:1603
 AliITSUTrackerGlo.cxx:1604
 AliITSUTrackerGlo.cxx:1605
 AliITSUTrackerGlo.cxx:1606
 AliITSUTrackerGlo.cxx:1607
 AliITSUTrackerGlo.cxx:1608
 AliITSUTrackerGlo.cxx:1609
 AliITSUTrackerGlo.cxx:1610
 AliITSUTrackerGlo.cxx:1611
 AliITSUTrackerGlo.cxx:1612
 AliITSUTrackerGlo.cxx:1613
 AliITSUTrackerGlo.cxx:1614
 AliITSUTrackerGlo.cxx:1615
 AliITSUTrackerGlo.cxx:1616
 AliITSUTrackerGlo.cxx:1617
 AliITSUTrackerGlo.cxx:1618
 AliITSUTrackerGlo.cxx:1619
 AliITSUTrackerGlo.cxx:1620
 AliITSUTrackerGlo.cxx:1621
 AliITSUTrackerGlo.cxx:1622
 AliITSUTrackerGlo.cxx:1623
 AliITSUTrackerGlo.cxx:1624
 AliITSUTrackerGlo.cxx:1625
 AliITSUTrackerGlo.cxx:1626
 AliITSUTrackerGlo.cxx:1627
 AliITSUTrackerGlo.cxx:1628
 AliITSUTrackerGlo.cxx:1629
 AliITSUTrackerGlo.cxx:1630
 AliITSUTrackerGlo.cxx:1631
 AliITSUTrackerGlo.cxx:1632
 AliITSUTrackerGlo.cxx:1633
 AliITSUTrackerGlo.cxx:1634
 AliITSUTrackerGlo.cxx:1635
 AliITSUTrackerGlo.cxx:1636
 AliITSUTrackerGlo.cxx:1637
 AliITSUTrackerGlo.cxx:1638
 AliITSUTrackerGlo.cxx:1639
 AliITSUTrackerGlo.cxx:1640
 AliITSUTrackerGlo.cxx:1641
 AliITSUTrackerGlo.cxx:1642
 AliITSUTrackerGlo.cxx:1643
 AliITSUTrackerGlo.cxx:1644
 AliITSUTrackerGlo.cxx:1645
 AliITSUTrackerGlo.cxx:1646
 AliITSUTrackerGlo.cxx:1647
 AliITSUTrackerGlo.cxx:1648
 AliITSUTrackerGlo.cxx:1649
 AliITSUTrackerGlo.cxx:1650
 AliITSUTrackerGlo.cxx:1651
 AliITSUTrackerGlo.cxx:1652
 AliITSUTrackerGlo.cxx:1653
 AliITSUTrackerGlo.cxx:1654
 AliITSUTrackerGlo.cxx:1655
 AliITSUTrackerGlo.cxx:1656
 AliITSUTrackerGlo.cxx:1657
 AliITSUTrackerGlo.cxx:1658
 AliITSUTrackerGlo.cxx:1659
 AliITSUTrackerGlo.cxx:1660
 AliITSUTrackerGlo.cxx:1661
 AliITSUTrackerGlo.cxx:1662
 AliITSUTrackerGlo.cxx:1663
 AliITSUTrackerGlo.cxx:1664
 AliITSUTrackerGlo.cxx:1665
 AliITSUTrackerGlo.cxx:1666
 AliITSUTrackerGlo.cxx:1667
 AliITSUTrackerGlo.cxx:1668
 AliITSUTrackerGlo.cxx:1669
 AliITSUTrackerGlo.cxx:1670
 AliITSUTrackerGlo.cxx:1671
 AliITSUTrackerGlo.cxx:1672
 AliITSUTrackerGlo.cxx:1673
 AliITSUTrackerGlo.cxx:1674
 AliITSUTrackerGlo.cxx:1675
 AliITSUTrackerGlo.cxx:1676
 AliITSUTrackerGlo.cxx:1677
 AliITSUTrackerGlo.cxx:1678
 AliITSUTrackerGlo.cxx:1679
 AliITSUTrackerGlo.cxx:1680
 AliITSUTrackerGlo.cxx:1681
 AliITSUTrackerGlo.cxx:1682
 AliITSUTrackerGlo.cxx:1683
 AliITSUTrackerGlo.cxx:1684
 AliITSUTrackerGlo.cxx:1685
 AliITSUTrackerGlo.cxx:1686
 AliITSUTrackerGlo.cxx:1687
 AliITSUTrackerGlo.cxx:1688
 AliITSUTrackerGlo.cxx:1689
 AliITSUTrackerGlo.cxx:1690
 AliITSUTrackerGlo.cxx:1691
 AliITSUTrackerGlo.cxx:1692
 AliITSUTrackerGlo.cxx:1693
 AliITSUTrackerGlo.cxx:1694
 AliITSUTrackerGlo.cxx:1695
 AliITSUTrackerGlo.cxx:1696
 AliITSUTrackerGlo.cxx:1697
 AliITSUTrackerGlo.cxx:1698
 AliITSUTrackerGlo.cxx:1699
 AliITSUTrackerGlo.cxx:1700
 AliITSUTrackerGlo.cxx:1701
 AliITSUTrackerGlo.cxx:1702
 AliITSUTrackerGlo.cxx:1703
 AliITSUTrackerGlo.cxx:1704
 AliITSUTrackerGlo.cxx:1705
 AliITSUTrackerGlo.cxx:1706
 AliITSUTrackerGlo.cxx:1707
 AliITSUTrackerGlo.cxx:1708
 AliITSUTrackerGlo.cxx:1709
 AliITSUTrackerGlo.cxx:1710
 AliITSUTrackerGlo.cxx:1711
 AliITSUTrackerGlo.cxx:1712
 AliITSUTrackerGlo.cxx:1713
 AliITSUTrackerGlo.cxx:1714
 AliITSUTrackerGlo.cxx:1715
 AliITSUTrackerGlo.cxx:1716
 AliITSUTrackerGlo.cxx:1717
 AliITSUTrackerGlo.cxx:1718
 AliITSUTrackerGlo.cxx:1719
 AliITSUTrackerGlo.cxx:1720
 AliITSUTrackerGlo.cxx:1721
 AliITSUTrackerGlo.cxx:1722
 AliITSUTrackerGlo.cxx:1723
 AliITSUTrackerGlo.cxx:1724
 AliITSUTrackerGlo.cxx:1725
 AliITSUTrackerGlo.cxx:1726
 AliITSUTrackerGlo.cxx:1727
 AliITSUTrackerGlo.cxx:1728
 AliITSUTrackerGlo.cxx:1729
 AliITSUTrackerGlo.cxx:1730
 AliITSUTrackerGlo.cxx:1731
 AliITSUTrackerGlo.cxx:1732
 AliITSUTrackerGlo.cxx:1733
 AliITSUTrackerGlo.cxx:1734
 AliITSUTrackerGlo.cxx:1735
 AliITSUTrackerGlo.cxx:1736
 AliITSUTrackerGlo.cxx:1737
 AliITSUTrackerGlo.cxx:1738
 AliITSUTrackerGlo.cxx:1739
 AliITSUTrackerGlo.cxx:1740
 AliITSUTrackerGlo.cxx:1741
 AliITSUTrackerGlo.cxx:1742
 AliITSUTrackerGlo.cxx:1743
 AliITSUTrackerGlo.cxx:1744
 AliITSUTrackerGlo.cxx:1745
 AliITSUTrackerGlo.cxx:1746
 AliITSUTrackerGlo.cxx:1747
 AliITSUTrackerGlo.cxx:1748
 AliITSUTrackerGlo.cxx:1749
 AliITSUTrackerGlo.cxx:1750
 AliITSUTrackerGlo.cxx:1751
 AliITSUTrackerGlo.cxx:1752
 AliITSUTrackerGlo.cxx:1753
 AliITSUTrackerGlo.cxx:1754
 AliITSUTrackerGlo.cxx:1755
 AliITSUTrackerGlo.cxx:1756
 AliITSUTrackerGlo.cxx:1757
 AliITSUTrackerGlo.cxx:1758
 AliITSUTrackerGlo.cxx:1759
 AliITSUTrackerGlo.cxx:1760
 AliITSUTrackerGlo.cxx:1761
 AliITSUTrackerGlo.cxx:1762
 AliITSUTrackerGlo.cxx:1763
 AliITSUTrackerGlo.cxx:1764
 AliITSUTrackerGlo.cxx:1765
 AliITSUTrackerGlo.cxx:1766
 AliITSUTrackerGlo.cxx:1767
 AliITSUTrackerGlo.cxx:1768
 AliITSUTrackerGlo.cxx:1769
 AliITSUTrackerGlo.cxx:1770
 AliITSUTrackerGlo.cxx:1771
 AliITSUTrackerGlo.cxx:1772
 AliITSUTrackerGlo.cxx:1773
 AliITSUTrackerGlo.cxx:1774
 AliITSUTrackerGlo.cxx:1775
 AliITSUTrackerGlo.cxx:1776
 AliITSUTrackerGlo.cxx:1777
 AliITSUTrackerGlo.cxx:1778
 AliITSUTrackerGlo.cxx:1779
 AliITSUTrackerGlo.cxx:1780
 AliITSUTrackerGlo.cxx:1781
 AliITSUTrackerGlo.cxx:1782
 AliITSUTrackerGlo.cxx:1783
 AliITSUTrackerGlo.cxx:1784
 AliITSUTrackerGlo.cxx:1785
 AliITSUTrackerGlo.cxx:1786
 AliITSUTrackerGlo.cxx:1787
 AliITSUTrackerGlo.cxx:1788
 AliITSUTrackerGlo.cxx:1789
 AliITSUTrackerGlo.cxx:1790
 AliITSUTrackerGlo.cxx:1791
 AliITSUTrackerGlo.cxx:1792
 AliITSUTrackerGlo.cxx:1793
 AliITSUTrackerGlo.cxx:1794
 AliITSUTrackerGlo.cxx:1795
 AliITSUTrackerGlo.cxx:1796
 AliITSUTrackerGlo.cxx:1797
 AliITSUTrackerGlo.cxx:1798
 AliITSUTrackerGlo.cxx:1799
 AliITSUTrackerGlo.cxx:1800
 AliITSUTrackerGlo.cxx:1801
 AliITSUTrackerGlo.cxx:1802
 AliITSUTrackerGlo.cxx:1803
 AliITSUTrackerGlo.cxx:1804
 AliITSUTrackerGlo.cxx:1805
 AliITSUTrackerGlo.cxx:1806
 AliITSUTrackerGlo.cxx:1807
 AliITSUTrackerGlo.cxx:1808
 AliITSUTrackerGlo.cxx:1809
 AliITSUTrackerGlo.cxx:1810
 AliITSUTrackerGlo.cxx:1811
 AliITSUTrackerGlo.cxx:1812
 AliITSUTrackerGlo.cxx:1813
 AliITSUTrackerGlo.cxx:1814
 AliITSUTrackerGlo.cxx:1815
 AliITSUTrackerGlo.cxx:1816
 AliITSUTrackerGlo.cxx:1817
 AliITSUTrackerGlo.cxx:1818
 AliITSUTrackerGlo.cxx:1819
 AliITSUTrackerGlo.cxx:1820
 AliITSUTrackerGlo.cxx:1821
 AliITSUTrackerGlo.cxx:1822
 AliITSUTrackerGlo.cxx:1823
 AliITSUTrackerGlo.cxx:1824
 AliITSUTrackerGlo.cxx:1825
 AliITSUTrackerGlo.cxx:1826
 AliITSUTrackerGlo.cxx:1827
 AliITSUTrackerGlo.cxx:1828
 AliITSUTrackerGlo.cxx:1829
 AliITSUTrackerGlo.cxx:1830
 AliITSUTrackerGlo.cxx:1831
 AliITSUTrackerGlo.cxx:1832
 AliITSUTrackerGlo.cxx:1833
 AliITSUTrackerGlo.cxx:1834
 AliITSUTrackerGlo.cxx:1835
 AliITSUTrackerGlo.cxx:1836
 AliITSUTrackerGlo.cxx:1837
 AliITSUTrackerGlo.cxx:1838
 AliITSUTrackerGlo.cxx:1839
 AliITSUTrackerGlo.cxx:1840
 AliITSUTrackerGlo.cxx:1841
 AliITSUTrackerGlo.cxx:1842
 AliITSUTrackerGlo.cxx:1843
 AliITSUTrackerGlo.cxx:1844
 AliITSUTrackerGlo.cxx:1845
 AliITSUTrackerGlo.cxx:1846
 AliITSUTrackerGlo.cxx:1847
 AliITSUTrackerGlo.cxx:1848
 AliITSUTrackerGlo.cxx:1849
 AliITSUTrackerGlo.cxx:1850
 AliITSUTrackerGlo.cxx:1851
 AliITSUTrackerGlo.cxx:1852
 AliITSUTrackerGlo.cxx:1853
 AliITSUTrackerGlo.cxx:1854
 AliITSUTrackerGlo.cxx:1855
 AliITSUTrackerGlo.cxx:1856
 AliITSUTrackerGlo.cxx:1857
 AliITSUTrackerGlo.cxx:1858
 AliITSUTrackerGlo.cxx:1859
 AliITSUTrackerGlo.cxx:1860
 AliITSUTrackerGlo.cxx:1861
 AliITSUTrackerGlo.cxx:1862
 AliITSUTrackerGlo.cxx:1863
 AliITSUTrackerGlo.cxx:1864
 AliITSUTrackerGlo.cxx:1865
 AliITSUTrackerGlo.cxx:1866
 AliITSUTrackerGlo.cxx:1867
 AliITSUTrackerGlo.cxx:1868
 AliITSUTrackerGlo.cxx:1869
 AliITSUTrackerGlo.cxx:1870
 AliITSUTrackerGlo.cxx:1871
 AliITSUTrackerGlo.cxx:1872
 AliITSUTrackerGlo.cxx:1873
 AliITSUTrackerGlo.cxx:1874
 AliITSUTrackerGlo.cxx:1875
 AliITSUTrackerGlo.cxx:1876
 AliITSUTrackerGlo.cxx:1877
 AliITSUTrackerGlo.cxx:1878
 AliITSUTrackerGlo.cxx:1879
 AliITSUTrackerGlo.cxx:1880
 AliITSUTrackerGlo.cxx:1881
 AliITSUTrackerGlo.cxx:1882
 AliITSUTrackerGlo.cxx:1883
 AliITSUTrackerGlo.cxx:1884
 AliITSUTrackerGlo.cxx:1885
 AliITSUTrackerGlo.cxx:1886
 AliITSUTrackerGlo.cxx:1887
 AliITSUTrackerGlo.cxx:1888
 AliITSUTrackerGlo.cxx:1889
 AliITSUTrackerGlo.cxx:1890
 AliITSUTrackerGlo.cxx:1891
 AliITSUTrackerGlo.cxx:1892
 AliITSUTrackerGlo.cxx:1893
 AliITSUTrackerGlo.cxx:1894
 AliITSUTrackerGlo.cxx:1895
 AliITSUTrackerGlo.cxx:1896
 AliITSUTrackerGlo.cxx:1897
 AliITSUTrackerGlo.cxx:1898
 AliITSUTrackerGlo.cxx:1899
 AliITSUTrackerGlo.cxx:1900
 AliITSUTrackerGlo.cxx:1901
 AliITSUTrackerGlo.cxx:1902
 AliITSUTrackerGlo.cxx:1903
 AliITSUTrackerGlo.cxx:1904
 AliITSUTrackerGlo.cxx:1905
 AliITSUTrackerGlo.cxx:1906
 AliITSUTrackerGlo.cxx:1907
 AliITSUTrackerGlo.cxx:1908
 AliITSUTrackerGlo.cxx:1909
 AliITSUTrackerGlo.cxx:1910
 AliITSUTrackerGlo.cxx:1911
 AliITSUTrackerGlo.cxx:1912
 AliITSUTrackerGlo.cxx:1913
 AliITSUTrackerGlo.cxx:1914
 AliITSUTrackerGlo.cxx:1915
 AliITSUTrackerGlo.cxx:1916
 AliITSUTrackerGlo.cxx:1917
 AliITSUTrackerGlo.cxx:1918
 AliITSUTrackerGlo.cxx:1919
 AliITSUTrackerGlo.cxx:1920
 AliITSUTrackerGlo.cxx:1921
 AliITSUTrackerGlo.cxx:1922
 AliITSUTrackerGlo.cxx:1923
 AliITSUTrackerGlo.cxx:1924
 AliITSUTrackerGlo.cxx:1925
 AliITSUTrackerGlo.cxx:1926
 AliITSUTrackerGlo.cxx:1927
 AliITSUTrackerGlo.cxx:1928
 AliITSUTrackerGlo.cxx:1929
 AliITSUTrackerGlo.cxx:1930
 AliITSUTrackerGlo.cxx:1931
 AliITSUTrackerGlo.cxx:1932
 AliITSUTrackerGlo.cxx:1933
 AliITSUTrackerGlo.cxx:1934
 AliITSUTrackerGlo.cxx:1935
 AliITSUTrackerGlo.cxx:1936
 AliITSUTrackerGlo.cxx:1937
 AliITSUTrackerGlo.cxx:1938
 AliITSUTrackerGlo.cxx:1939
 AliITSUTrackerGlo.cxx:1940
 AliITSUTrackerGlo.cxx:1941
 AliITSUTrackerGlo.cxx:1942
 AliITSUTrackerGlo.cxx:1943
 AliITSUTrackerGlo.cxx:1944
 AliITSUTrackerGlo.cxx:1945
 AliITSUTrackerGlo.cxx:1946
 AliITSUTrackerGlo.cxx:1947
 AliITSUTrackerGlo.cxx:1948
 AliITSUTrackerGlo.cxx:1949
 AliITSUTrackerGlo.cxx:1950
 AliITSUTrackerGlo.cxx:1951
 AliITSUTrackerGlo.cxx:1952
 AliITSUTrackerGlo.cxx:1953
 AliITSUTrackerGlo.cxx:1954
 AliITSUTrackerGlo.cxx:1955
 AliITSUTrackerGlo.cxx:1956
 AliITSUTrackerGlo.cxx:1957
 AliITSUTrackerGlo.cxx:1958
 AliITSUTrackerGlo.cxx:1959
 AliITSUTrackerGlo.cxx:1960
 AliITSUTrackerGlo.cxx:1961
 AliITSUTrackerGlo.cxx:1962
 AliITSUTrackerGlo.cxx:1963
 AliITSUTrackerGlo.cxx:1964
 AliITSUTrackerGlo.cxx:1965
 AliITSUTrackerGlo.cxx:1966
 AliITSUTrackerGlo.cxx:1967
 AliITSUTrackerGlo.cxx:1968
 AliITSUTrackerGlo.cxx:1969
 AliITSUTrackerGlo.cxx:1970
 AliITSUTrackerGlo.cxx:1971
 AliITSUTrackerGlo.cxx:1972
 AliITSUTrackerGlo.cxx:1973
 AliITSUTrackerGlo.cxx:1974
 AliITSUTrackerGlo.cxx:1975
 AliITSUTrackerGlo.cxx:1976
 AliITSUTrackerGlo.cxx:1977
 AliITSUTrackerGlo.cxx:1978
 AliITSUTrackerGlo.cxx:1979
 AliITSUTrackerGlo.cxx:1980
 AliITSUTrackerGlo.cxx:1981
 AliITSUTrackerGlo.cxx:1982
 AliITSUTrackerGlo.cxx:1983
 AliITSUTrackerGlo.cxx:1984
 AliITSUTrackerGlo.cxx:1985
 AliITSUTrackerGlo.cxx:1986
 AliITSUTrackerGlo.cxx:1987
 AliITSUTrackerGlo.cxx:1988
 AliITSUTrackerGlo.cxx:1989
 AliITSUTrackerGlo.cxx:1990
 AliITSUTrackerGlo.cxx:1991
 AliITSUTrackerGlo.cxx:1992
 AliITSUTrackerGlo.cxx:1993
 AliITSUTrackerGlo.cxx:1994
 AliITSUTrackerGlo.cxx:1995
 AliITSUTrackerGlo.cxx:1996
 AliITSUTrackerGlo.cxx:1997
 AliITSUTrackerGlo.cxx:1998
 AliITSUTrackerGlo.cxx:1999
 AliITSUTrackerGlo.cxx:2000
 AliITSUTrackerGlo.cxx:2001
 AliITSUTrackerGlo.cxx:2002
 AliITSUTrackerGlo.cxx:2003
 AliITSUTrackerGlo.cxx:2004
 AliITSUTrackerGlo.cxx:2005
 AliITSUTrackerGlo.cxx:2006
 AliITSUTrackerGlo.cxx:2007
 AliITSUTrackerGlo.cxx:2008
 AliITSUTrackerGlo.cxx:2009
 AliITSUTrackerGlo.cxx:2010
 AliITSUTrackerGlo.cxx:2011
 AliITSUTrackerGlo.cxx:2012
 AliITSUTrackerGlo.cxx:2013
 AliITSUTrackerGlo.cxx:2014
 AliITSUTrackerGlo.cxx:2015
 AliITSUTrackerGlo.cxx:2016
 AliITSUTrackerGlo.cxx:2017
 AliITSUTrackerGlo.cxx:2018
 AliITSUTrackerGlo.cxx:2019
 AliITSUTrackerGlo.cxx:2020
 AliITSUTrackerGlo.cxx:2021
 AliITSUTrackerGlo.cxx:2022
 AliITSUTrackerGlo.cxx:2023
 AliITSUTrackerGlo.cxx:2024
 AliITSUTrackerGlo.cxx:2025
 AliITSUTrackerGlo.cxx:2026
 AliITSUTrackerGlo.cxx:2027
 AliITSUTrackerGlo.cxx:2028
 AliITSUTrackerGlo.cxx:2029
 AliITSUTrackerGlo.cxx:2030
 AliITSUTrackerGlo.cxx:2031
 AliITSUTrackerGlo.cxx:2032
 AliITSUTrackerGlo.cxx:2033
 AliITSUTrackerGlo.cxx:2034
 AliITSUTrackerGlo.cxx:2035
 AliITSUTrackerGlo.cxx:2036
 AliITSUTrackerGlo.cxx:2037
 AliITSUTrackerGlo.cxx:2038
 AliITSUTrackerGlo.cxx:2039
 AliITSUTrackerGlo.cxx:2040
 AliITSUTrackerGlo.cxx:2041
 AliITSUTrackerGlo.cxx:2042
 AliITSUTrackerGlo.cxx:2043
 AliITSUTrackerGlo.cxx:2044
 AliITSUTrackerGlo.cxx:2045
 AliITSUTrackerGlo.cxx:2046
 AliITSUTrackerGlo.cxx:2047
 AliITSUTrackerGlo.cxx:2048
 AliITSUTrackerGlo.cxx:2049
 AliITSUTrackerGlo.cxx:2050
 AliITSUTrackerGlo.cxx:2051
 AliITSUTrackerGlo.cxx:2052
 AliITSUTrackerGlo.cxx:2053
 AliITSUTrackerGlo.cxx:2054
 AliITSUTrackerGlo.cxx:2055
 AliITSUTrackerGlo.cxx:2056
 AliITSUTrackerGlo.cxx:2057
 AliITSUTrackerGlo.cxx:2058
 AliITSUTrackerGlo.cxx:2059
 AliITSUTrackerGlo.cxx:2060
 AliITSUTrackerGlo.cxx:2061
 AliITSUTrackerGlo.cxx:2062
 AliITSUTrackerGlo.cxx:2063
 AliITSUTrackerGlo.cxx:2064
 AliITSUTrackerGlo.cxx:2065
 AliITSUTrackerGlo.cxx:2066