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

/* $Id$ */

//-----------------------------------------------------------------------------
/// \class AliMUONVTrackReconstructor
/// Virtual MUON track reconstructor in ALICE (class renamed from AliMUONEventReconstructor)
///
/// This class contains as data a pointer to the array of reconstructed tracks
///
/// It contains as methods, among others:
/// * EventReconstruct to build the muon tracks
/// * EventReconstructTrigger to build the trigger tracks
/// * ValidateTracksWithTrigger to match tracker/trigger tracks
///
/// Several options and adjustable parameters are available for both KALMAN and ORIGINAL
/// tracking algorithms. They can be changed through the AliMUONRecoParam object
/// set in the reconstruction macro or read from the CDB
/// (see methods in AliMUONRecoParam.h file for details)
///
/// Main parameters and options are:
/// - *fgkSigmaToCutForTracking* : quality cut used to select new clusters to be
///   attached to the track candidate and to select good tracks.
/// - *fgkMakeTrackCandidatesFast* : if this flag is set to 'true', the track candidates
///   are made assuming linear propagation between stations 4 and 5.
/// - *fgkTrackAllTracks* : according to the value of this flag, in case that several
///   new clusters pass the quality cut, either we consider all the possibilities
///   (duplicating tracks) or we attach only the best cluster.
/// - *fgkRecoverTracks* : if this flag is set to 'true', we try to recover the tracks
///   lost during the tracking by removing the worst of the 2 clusters attached in the
///   previous station.
/// - *fgkImproveTracks* : if this flag is set to 'true', we try to improve the quality
///   of the tracks at the end of the tracking by removing clusters that do not pass
///   new quality cut (the track is removed is it does not contain enough cluster anymore).
/// - *fgkComplementTracks* : if this flag is set to 'true', we try to improve the quality
///   of the tracks at the end of the tracking by adding potentially missing clusters
///   (we may have 2 clusters in the same chamber because of the overlapping of detection
///   elements, which is not handle by the tracking algorithm).
/// - *fgkSigmaToCutForImprovement* : quality cut used when we try to improve the
///   quality of the tracks.
///
///  \author Philippe Pillot
//-----------------------------------------------------------------------------

#include "AliMUONVTrackReconstructor.h"

#include "AliMUONConstants.h"
#include "AliMUONObjectPair.h"
#include "AliMUONTriggerTrack.h"
#include "AliMUONTriggerCircuit.h"
#include "AliMUONLocalTrigger.h"
#include "AliMUONGlobalTrigger.h"
#include "AliMUONTrack.h"
#include "AliMUONTrackParam.h"
#include "AliMUONTrackExtrap.h"
#include "AliMUONTrackHitPattern.h"
#include "AliMUONVTrackStore.h"
#include "AliMUONVClusterStore.h"
#include "AliMUONVCluster.h"
#include "AliMUONVClusterServer.h"
#include "AliMUONVTriggerStore.h"
#include "AliMUONVTriggerTrackStore.h"
#include "AliMUONRecoParam.h"
#include "AliMUONGeometryTransformer.h"
#include "AliMUONVDigit.h"

#include "AliMpDEManager.h"
#include "AliMpArea.h"

#include "AliMpDDLStore.h"
#include "AliMpVSegmentation.h"
#include "AliMpSegmentation.h"
#include "AliMpPad.h"
#include "AliMpDetElement.h"
#include "AliMpCathodType.h"

#include "AliLog.h"
#include "AliCodeTimer.h"
#include "AliTracker.h"

#include <TClonesArray.h>
#include <TMath.h>
#include <TMatrixD.h>
#include <TVector2.h>

#include <Riostream.h>

using std::cout;
using std::endl;
/// \cond CLASSIMP
ClassImp(AliMUONVTrackReconstructor) // Class implementation in ROOT context
/// \endcond

  //__________________________________________________________________________
AliMUONVTrackReconstructor::AliMUONVTrackReconstructor(const AliMUONRecoParam* recoParam,
                                                       AliMUONVClusterServer* clusterServer,
						       const AliMUONGeometryTransformer* transformer)
: TObject(),
fRecTracksPtr(0x0),
fNRecTracks(0),
fClusterServer(clusterServer),
fkRecoParam(recoParam),
fkTransformer(transformer),
fMaxMCSAngle2(0x0)
{
  /// Constructor for class AliMUONVTrackReconstructor
  /// WARNING: if clusterServer=0x0, no clusterization will be possible at this level
  
  // Memory allocation for the TClonesArray of reconstructed tracks
  fRecTracksPtr = new TClonesArray("AliMUONTrack", 100);
  
  // set the magnetic field for track extrapolations
  AliMUONTrackExtrap::SetField();
  
  // set the maximum MCS angle in chamber from the minimum acceptable momentum
  AliMUONTrackParam param;
  Double_t inverseBendingP = (GetRecoParam()->GetMinBendingMomentum() > 0.) ? 1./GetRecoParam()->GetMinBendingMomentum() : 1.;
  param.SetInverseBendingMomentum(inverseBendingP);
  fMaxMCSAngle2 = new Double_t [AliMUONConstants::NTrackingCh()];
  for (Int_t iCh=0; iCh<AliMUONConstants::NTrackingCh(); iCh++)
    fMaxMCSAngle2[iCh] = AliMUONTrackExtrap::GetMCSAngle2(param, AliMUONConstants::ChamberThicknessInX0(iCh), 1.);
  
}

  //__________________________________________________________________________
AliMUONVTrackReconstructor::~AliMUONVTrackReconstructor()
{
  /// Destructor for class AliMUONVTrackReconstructor
  delete fRecTracksPtr;
  delete[] fMaxMCSAngle2;
}

  //__________________________________________________________________________
void AliMUONVTrackReconstructor::ResetTracks()
{
  /// To reset the TClonesArray of reconstructed tracks
  if (fRecTracksPtr) fRecTracksPtr->Clear("C");
  fNRecTracks = 0;
  return;
}

  //__________________________________________________________________________
void AliMUONVTrackReconstructor::EventReconstruct(AliMUONVClusterStore& clusterStore, AliMUONVTrackStore& trackStore)
{
  /// To reconstruct one event
  AliDebug(1,"");
  AliCodeTimerAuto("",0);
  
  // Reset array of tracks
  ResetTracks();
  
  // Look for candidates from clusters in stations(1..) 4 and 5 (abort in case of failure)
  if (!MakeTrackCandidates(clusterStore)) return;
  
  // Look for extra candidates from clusters in stations(1..) 4 and 5 (abort in case of failure)
  if (GetRecoParam()->MakeMoreTrackCandidates()) {
    if (!MakeMoreTrackCandidates(clusterStore)) return;
  }
  
  // Stop tracking if no candidate found
  if (fRecTracksPtr->GetEntriesFast() == 0) return;
  
  // Follow tracks in stations(1..) 3, 2 and 1 (abort in case of failure)
  if (!FollowTracks(clusterStore)) return;
  
  // Complement the reconstructed tracks
  if (GetRecoParam()->ComplementTracks()) {
    if (ComplementTracks(clusterStore)) RemoveIdenticalTracks();
  }
  
  // Improve the reconstructed tracks
  if (GetRecoParam()->ImproveTracks()) ImproveTracks();
  
  // Remove connected tracks
  RemoveConnectedTracks(3, 4, kFALSE);
  RemoveConnectedTracks(2, 2, kFALSE);
  if (GetRecoParam()->RemoveConnectedTracksInSt12()) RemoveConnectedTracks(0, 1, kFALSE);
  
  // Fill AliMUONTrack data members
  Finalize();
  if (!GetRecoParam()->RemoveConnectedTracksInSt12()) TagConnectedTracks(0, 1, kTRUE);
  
  // Make sure there is no bad track left
  RemoveBadTracks();
  
  // Refit the reconstructed tracks with a different resolution for mono-cathod clusters
  if (GetRecoParam()->DiscardMonoCathodClusters()) DiscardMonoCathodClusters();
  
  // Add tracks to MUON data container 
  for (Int_t i=0; i<fNRecTracks; ++i)
  {
    AliMUONTrack * track = (AliMUONTrack*) fRecTracksPtr->At(i);
    track->SetUniqueID(i+1);
    trackStore.Add(*track);
  }
  
}

//__________________________________________________________________________
Bool_t AliMUONVTrackReconstructor::IsAcceptable(AliMUONTrackParam &trackParam)
{
  /// Return kTRUE if the track is within given limits on momentum/angle/origin
  
  const TMatrixD& kParamCov = trackParam.GetCovariances();
  Int_t chamber = trackParam.GetClusterPtr()->GetChamberId();
  Double_t z = trackParam.GetZ();
  Double_t sigmaCut = GetRecoParam()->GetSigmaCutForTracking();
  
  // MCS dipersion
  Double_t angleMCS2 = 0.;
  Double_t impactMCS2 = 0.;
  if (AliMUONTrackExtrap::IsFieldON() && chamber < 6) {
    
    // track momentum is known
    for (Int_t iCh=0; iCh<=chamber; iCh++) {
      Double_t localMCS2 = AliMUONTrackExtrap::GetMCSAngle2(trackParam, AliMUONConstants::ChamberThicknessInX0(iCh), 1.);
      angleMCS2 += localMCS2;
      impactMCS2 += AliMUONConstants::DefaultChamberZ(chamber) * AliMUONConstants::DefaultChamberZ(chamber) * localMCS2;
    }
    
  } else {
    
    // track momentum is unknown
    for (Int_t iCh=0; iCh<=chamber; iCh++) {
      angleMCS2 += fMaxMCSAngle2[iCh];
      impactMCS2 += AliMUONConstants::DefaultChamberZ(chamber) * AliMUONConstants::DefaultChamberZ(chamber) * fMaxMCSAngle2[iCh];
    }
    
  }
  
  // ------ track selection in non bending direction ------
  if (GetRecoParam()->SelectOnTrackSlope()) {
    
    // check if non bending slope is within tolerances
    Double_t nonBendingSlopeErr = TMath::Sqrt(kParamCov(1,1) + angleMCS2);
    if ((TMath::Abs(trackParam.GetNonBendingSlope()) - sigmaCut * nonBendingSlopeErr) > GetRecoParam()->GetMaxNonBendingSlope()) return kFALSE;
    
  } else {
    
    // or check if non bending impact parameter is within tolerances
    Double_t nonBendingImpactParam = TMath::Abs(trackParam.GetNonBendingCoor() - z * trackParam.GetNonBendingSlope());
    Double_t nonBendingImpactParamErr = TMath::Sqrt(kParamCov(0,0) + z * z * kParamCov(1,1) - 2. * z * kParamCov(0,1) + impactMCS2);
    if ((nonBendingImpactParam - sigmaCut * nonBendingImpactParamErr) > (3. * GetRecoParam()->GetNonBendingVertexDispersion())) return kFALSE;
    
  }
  
  // ------ track selection in bending direction ------
  if (AliMUONTrackExtrap::IsFieldON()) { // depending whether the field is ON or OFF
    
    // check if bending momentum is within tolerances
    Double_t bendingMomentum = TMath::Abs(1. / trackParam.GetInverseBendingMomentum());
    Double_t bendingMomentumErr = TMath::Sqrt(kParamCov(4,4)) * bendingMomentum * bendingMomentum;
    if (chamber < 6 && (bendingMomentum + sigmaCut * bendingMomentumErr) < GetRecoParam()->GetMinBendingMomentum()) return kFALSE;
    else if ((bendingMomentum + 3. * bendingMomentumErr) < GetRecoParam()->GetMinBendingMomentum()) return kFALSE;
    
  } else {
    
    if (GetRecoParam()->SelectOnTrackSlope()) {
      
      // check if bending slope is within tolerances
      Double_t bendingSlopeErr = TMath::Sqrt(kParamCov(3,3) + angleMCS2);
      if ((TMath::Abs(trackParam.GetBendingSlope()) - sigmaCut * bendingSlopeErr) > GetRecoParam()->GetMaxBendingSlope()) return kFALSE;
      
    } else {
      
      // or check if bending impact parameter is within tolerances
      Double_t bendingImpactParam = TMath::Abs(trackParam.GetBendingCoor() - z * trackParam.GetBendingSlope());
      Double_t bendingImpactParamErr = TMath::Sqrt(kParamCov(2,2) + z * z * kParamCov(3,3) - 2. * z * kParamCov(2,3) + impactMCS2);
      if ((bendingImpactParam - sigmaCut * bendingImpactParamErr) > (3. * GetRecoParam()->GetBendingVertexDispersion())) return kFALSE;
      
    }
    
  }
  
  return kTRUE;
  
}

//__________________________________________________________________________
TClonesArray* AliMUONVTrackReconstructor::MakeSegmentsBetweenChambers(const AliMUONVClusterStore& clusterStore, Int_t ch1, Int_t ch2)
{
  /// To make the list of segments from the list of clusters in the 2 given chambers.
  /// Return a TClonesArray of new segments (segments made in a previous call of this function are removed).
  AliDebug(1,Form("Enter MakeSegmentsBetweenChambers (1..) %d-%d", ch1+1, ch2+1));
  AliCodeTimerAuto("",0);
  
  AliMUONVCluster *cluster1, *cluster2;
  AliMUONObjectPair *segment;
  Double_t z1 = 0., z2 = 0., dZ = 0.;
  Double_t nonBendingSlope = 0., nonBendingSlopeErr = 0., nonBendingImpactParam = 0., nonBendingImpactParamErr = 0.;
  Double_t bendingSlope = 0., bendingSlopeErr = 0., bendingImpactParam = 0., bendingImpactParamErr = 0., bendingImpactParamErr2 = 0.;
  Double_t bendingMomentum = 0., bendingMomentumErr = 0.;
  Double_t bendingVertexDispersion2 = GetRecoParam()->GetBendingVertexDispersion() * GetRecoParam()->GetBendingVertexDispersion();
  Double_t angleMCS2 = 0.; // maximum angular dispersion**2 due to MCS in chamber
  Double_t impactMCS2 = 0.; // maximum impact parameter dispersion**2 due to MCS in chamber
  for (Int_t iCh=0; iCh<=ch1; iCh++) {
    angleMCS2 += fMaxMCSAngle2[iCh];
    impactMCS2 += AliMUONConstants::DefaultChamberZ(iCh) * AliMUONConstants::DefaultChamberZ(iCh) * fMaxMCSAngle2[iCh];
  }
  Double_t sigmaCut = GetRecoParam()->GetSigmaCutForTracking();
  
  // Create iterators to loop over clusters in both chambers
  TIter nextInCh1(clusterStore.CreateChamberIterator(ch1,ch1));
  TIter nextInCh2(clusterStore.CreateChamberIterator(ch2,ch2));
  
  // list of segments
  static TClonesArray *segments = new TClonesArray("AliMUONObjectPair", 100);
  segments->Clear("C");
  
  // Loop over clusters in the first chamber of the station
  while ( ( cluster1 = static_cast<AliMUONVCluster*>(nextInCh1()) ) ) {
    z1 = cluster1->GetZ();
    
    // reset cluster iterator of chamber 2
    nextInCh2.Reset();
    
    // Loop over clusters in the second chamber of the station
    while ( ( cluster2 = static_cast<AliMUONVCluster*>(nextInCh2()) ) ) {
      z2 = cluster2->GetZ();
      dZ = z1 - z2;
      
      // ------ track selection in non bending direction ------
      nonBendingSlope = (cluster1->GetX() - cluster2->GetX()) / dZ;
      if (GetRecoParam()->SelectOnTrackSlope()) {
	
	// check if non bending slope is within tolerances
	nonBendingSlopeErr = TMath::Sqrt((cluster1->GetErrX2() + cluster2->GetErrX2()) / dZ / dZ + angleMCS2);
	if ((TMath::Abs(nonBendingSlope) - sigmaCut * nonBendingSlopeErr) > GetRecoParam()->GetMaxNonBendingSlope()) continue;
	
      } else {
	
	// or check if non bending impact parameter is within tolerances
	nonBendingImpactParam = TMath::Abs(cluster1->GetX() - cluster1->GetZ() * nonBendingSlope);
	nonBendingImpactParamErr = TMath::Sqrt((z1 * z1 * cluster2->GetErrX2() + z2 * z2 * cluster1->GetErrX2()) / dZ / dZ + impactMCS2);
	if ((nonBendingImpactParam - sigmaCut * nonBendingImpactParamErr) > (3. * GetRecoParam()->GetNonBendingVertexDispersion())) continue;
	
      }
      
      // ------ track selection in bending direction ------
      bendingSlope = (cluster1->GetY() - cluster2->GetY()) / dZ;
      if (AliMUONTrackExtrap::IsFieldON()) { // depending whether the field is ON or OFF
	
	// check if bending momentum is within tolerances
	bendingImpactParam = cluster1->GetY() - cluster1->GetZ() * bendingSlope;
	bendingImpactParamErr2 = (z1 * z1 * cluster2->GetErrY2() + z2 * z2 * cluster1->GetErrY2()) / dZ / dZ + impactMCS2;
	bendingMomentum = TMath::Abs(AliMUONTrackExtrap::GetBendingMomentumFromImpactParam(bendingImpactParam));
	bendingMomentumErr = TMath::Sqrt((bendingVertexDispersion2 + bendingImpactParamErr2) /
					 bendingImpactParam / bendingImpactParam + 0.01) * bendingMomentum;
	if ((bendingMomentum + 3. * bendingMomentumErr) < GetRecoParam()->GetMinBendingMomentum()) continue;
	
      } else {
	
	if (GetRecoParam()->SelectOnTrackSlope()) {
	  
	  // check if bending slope is within tolerances
	  bendingSlopeErr = TMath::Sqrt((cluster1->GetErrY2() + cluster2->GetErrY2()) / dZ / dZ + angleMCS2);
	  if ((TMath::Abs(bendingSlope) - sigmaCut * bendingSlopeErr) > GetRecoParam()->GetMaxBendingSlope()) continue;
	  
	} else {
	  
	  // or check if bending impact parameter is within tolerances
	  bendingImpactParam = TMath::Abs(cluster1->GetY() - cluster1->GetZ() * bendingSlope);
	  bendingImpactParamErr = TMath::Sqrt((z1 * z1 * cluster2->GetErrY2() + z2 * z2 * cluster1->GetErrY2()) / dZ / dZ + impactMCS2);
	  if ((bendingImpactParam - sigmaCut * bendingImpactParamErr) > (3. * GetRecoParam()->GetBendingVertexDispersion())) continue;
	  
	}
	
      }
      
      // make new segment
      segment = new ((*segments)[segments->GetLast()+1]) AliMUONObjectPair(cluster1, cluster2, kFALSE, kFALSE);
      
      // Printout for debug
      if (AliLog::GetGlobalDebugLevel() > 1) {
	cout << "segmentIndex(0...): " << segments->GetLast() << endl;
	segment->Dump();
	cout << "Cluster in first chamber" << endl;
	cluster1->Print();
	cout << "Cluster in second chamber" << endl;
	cluster2->Print();
      }
      
    }
    
  }
  
  // Printout for debug
  AliDebug(1,Form("chambers%d-%d: NSegments =  %d ", ch1+1, ch2+1, segments->GetEntriesFast()));
  
  return segments;
}

  //__________________________________________________________________________
void AliMUONVTrackReconstructor::RemoveUsedSegments(TClonesArray& segments)
{
  /// To remove pairs of clusters already attached to a track
  AliDebug(1,"Enter RemoveUsedSegments");
  Int_t nSegments = segments.GetEntriesFast();
  Int_t nTracks = fRecTracksPtr->GetEntriesFast();
  AliMUONObjectPair *segment;
  AliMUONTrack *track;
  AliMUONVCluster *cluster, *cluster1, *cluster2;
  Bool_t foundCluster1, foundCluster2, removeSegment;
  
  // Loop over segments
  for (Int_t iSegment=0; iSegment<nSegments; iSegment++) {
    segment = (AliMUONObjectPair*) segments.UncheckedAt(iSegment);
    
    cluster1 = (AliMUONVCluster*) segment->First();
    cluster2 = (AliMUONVCluster*) segment->Second();
    removeSegment = kFALSE;
    
    // Loop over tracks
    for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
      track = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iTrack);
      
      // skip empty slot
      if (!track) continue;
      
      foundCluster1 = kFALSE;
      foundCluster2 = kFALSE;
      
      // Loop over clusters
      Int_t nClusters = track->GetNClusters();
      for (Int_t iCluster = 0; iCluster < nClusters; iCluster++) {
        cluster = ((AliMUONTrackParam*) track->GetTrackParamAtCluster()->UncheckedAt(iCluster))->GetClusterPtr();
	
	// check if both clusters are in that track
	if (cluster == cluster1) foundCluster1 = kTRUE;
	else if (cluster == cluster2) foundCluster2 = kTRUE;
	
	if (foundCluster1 && foundCluster2) {
	  removeSegment = kTRUE;
	  break;
	}
	
      }
      
      if (removeSegment) break;
      
    }
    
    if (removeSegment) segments.RemoveAt(iSegment);
      
  }
  
  segments.Compress();
  
  // Printout for debug
  AliDebug(1,Form("NSegments =  %d ", segments.GetEntriesFast()));
}

  //__________________________________________________________________________
void AliMUONVTrackReconstructor::RemoveIdenticalTracks()
{
  /// To remove identical tracks:
  /// Tracks are considered identical if they have all their clusters in common.
  /// One keeps the track with the larger number of clusters if need be
  AliMUONTrack *track1, *track2;
  Int_t nTracks = fRecTracksPtr->GetEntriesFast();
  Int_t clustersInCommon, nClusters1, nClusters2;
  // Loop over first track of the pair
  for (Int_t iTrack1 = 0; iTrack1 < nTracks; iTrack1++) {
    track1 = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iTrack1);
    // skip empty slot
    if (!track1) continue;
    nClusters1 = track1->GetNClusters();
    // Loop over second track of the pair
    for (Int_t iTrack2 = iTrack1+1; iTrack2 < nTracks; iTrack2++) {
      track2 = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iTrack2);
      // skip empty slot
      if (!track2) continue;
      nClusters2 = track2->GetNClusters();
      // number of clusters in common between two tracks
      clustersInCommon = track1->ClustersInCommon(track2);
      // check for identical tracks
      if ((clustersInCommon == nClusters1) || (clustersInCommon == nClusters2)) {
        // decide which track to remove
        if (nClusters2 > nClusters1) {
	  // remove track1 and continue the first loop with the track next to track1
          fRecTracksPtr->RemoveAt(iTrack1);
	  fNRecTracks--;
	  break;
	} else {
	  // remove track2 and continue the second loop with the track next to track2
	  fRecTracksPtr->RemoveAt(iTrack2);
	  fNRecTracks--;
        }
      }
    } // track2
  } // track1
  fRecTracksPtr->Compress(); // this is essential to retrieve the TClonesArray afterwards
}

  //__________________________________________________________________________
void AliMUONVTrackReconstructor::RemoveDoubleTracks()
{
  /// To remove double tracks:
  /// Tracks are considered identical if more than half of the clusters of the track
  /// which has the smaller number of clusters are in common with the other track.
  /// Among two identical tracks, one keeps the track with the larger number of clusters
  /// or, if these numbers are equal, the track with the minimum chi2.
  AliMUONTrack *track1, *track2;
  Int_t nTracks = fRecTracksPtr->GetEntriesFast();
  Int_t clustersInCommon2, nClusters1, nClusters2;
  // Loop over first track of the pair
  for (Int_t iTrack1 = 0; iTrack1 < nTracks; iTrack1++) {
    track1 = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iTrack1);
    // skip empty slot
    if (!track1) continue;
    nClusters1 = track1->GetNClusters();
    // Loop over second track of the pair
    for (Int_t iTrack2 = iTrack1+1; iTrack2 < nTracks; iTrack2++) {
      track2 = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iTrack2);
      // skip empty slot
      if (!track2) continue;
      nClusters2 = track2->GetNClusters();
      // number of clusters in common between two tracks
      clustersInCommon2 = 2 * track1->ClustersInCommon(track2);
      // check for identical tracks
      if (clustersInCommon2 > nClusters1 || clustersInCommon2 > nClusters2) {
        // decide which track to remove
        if ((nClusters1 > nClusters2) || ((nClusters1 == nClusters2) && (track1->GetGlobalChi2() <= track2->GetGlobalChi2()))) {
	  // remove track2 and continue the second loop with the track next to track2
	  fRecTracksPtr->RemoveAt(iTrack2);
	  fNRecTracks--;
        } else {
	  // else remove track1 and continue the first loop with the track next to track1
          fRecTracksPtr->RemoveAt(iTrack1);
	  fNRecTracks--;
	  break;
        }
      }
    } // track2
  } // track1
  fRecTracksPtr->Compress(); // this is essential to retrieve the TClonesArray afterwards
}

  //__________________________________________________________________________
void AliMUONVTrackReconstructor::RemoveBadTracks()
{
  /// Remove tracks for which a problem occured somewhere during the tracking
  
  AliMUONTrack *track, *nextTrack;
  Bool_t trackRemoved = kFALSE;
  
  track = (AliMUONTrack*) fRecTracksPtr->First();
  while (track) {
    
    nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track);
    
    if (track->GetGlobalChi2() >= AliMUONTrack::MaxChi2()) {
      AliWarning("problem occured somewhere during the tracking --> discard track");
      fRecTracksPtr->Remove(track);
      fNRecTracks--;
      trackRemoved = kTRUE;
    }
    
    track = nextTrack;
    
  }
  
  // compress array of tracks if needed
  if (trackRemoved) fRecTracksPtr->Compress();
  
}

  //__________________________________________________________________________
void AliMUONVTrackReconstructor::RemoveConnectedTracks(Int_t stMin, Int_t stMax, Bool_t all)
{
  /// Find and remove tracks sharing 1 cluster or more in station(s) [stMin, stMax].
  /// For each couple of connected tracks, one removes the one with the smallest
  /// number of clusters or with the highest chi2 value in case of equality.
  /// If all=kTRUE: both tracks are removed.
  
  // tag the tracks to be removed
  TagConnectedTracks(stMin, stMax, all);
  
  // remove them
  Int_t nTracks = fRecTracksPtr->GetEntriesFast();
  for (Int_t i = 0; i < nTracks; i++) {
    if (((AliMUONTrack*) fRecTracksPtr->UncheckedAt(i))->IsConnected()) {
      fRecTracksPtr->RemoveAt(i);
      fNRecTracks--;
    }
  }
  
  // remove holes in the array if any
  fRecTracksPtr->Compress();
}

  //__________________________________________________________________________
void AliMUONVTrackReconstructor::TagConnectedTracks(Int_t stMin, Int_t stMax, Bool_t all)
{
  /// Find and tag tracks sharing 1 cluster or more in station(s) [stMin, stMax].
  /// For each couple of connected tracks, one tags the one with the smallest
  /// number of clusters or with the highest chi2 value in case of equality.
  /// If all=kTRUE: both tracks are tagged.
  
  AliMUONTrack *track1, *track2;
  Int_t nClusters1, nClusters2;
  Int_t nTracks = fRecTracksPtr->GetEntriesFast();
  
  // reset the tags
  for (Int_t i = 0; i < nTracks; i++) ((AliMUONTrack*) fRecTracksPtr->UncheckedAt(i))->Connected(kFALSE);
    
  // Loop over first track of the pair
  for (Int_t iTrack1 = 0; iTrack1 < nTracks; iTrack1++) {
    track1 = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iTrack1);
    
    // Loop over second track of the pair
    for (Int_t iTrack2 = iTrack1+1; iTrack2 < nTracks; iTrack2++) {
      track2 = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iTrack2);
      
      // check for connected tracks and tag them
      if (track1->ClustersInCommon(track2, stMin, stMax) > 0) {
        
	if (all) {
	  
	  // tag both tracks
	  track1->Connected();
	  track2->Connected();
	  
	} else {
	  
	  // tag only the worst track
	  nClusters1 = track1->GetNClusters();
	  nClusters2 = track2->GetNClusters();
	  if ((nClusters1 > nClusters2) || ((nClusters1 == nClusters2) && (track1->GetGlobalChi2() <= track2->GetGlobalChi2()))) {
	    track2->Connected();
	  } else {
	    track1->Connected();
	  }
	  
	}
	
      }
      
    }
    
  }
  
}

  //__________________________________________________________________________
void AliMUONVTrackReconstructor::AskForNewClustersInChamber(const AliMUONTrackParam &trackParam,
							    AliMUONVClusterStore& clusterStore, Int_t chamber)
{
  /// Ask the clustering to reconstruct new clusters around the track candidate position
  
  // check if the current chamber is useable
  if (!fClusterServer || !GetRecoParam()->UseChamber(chamber)) return;
  
  // maximum distance between the center of the chamber and a detection element
  // (accounting for the inclination of the chamber)
  static const Double_t kMaxDZ = 15.; // 15 cm
  
  // extrapolate track parameters to the chamber
  AliMUONTrackParam extrapTrackParam(trackParam);
  if (!AliMUONTrackExtrap::ExtrapToZCov(&extrapTrackParam, AliMUONConstants::DefaultChamberZ(chamber))) return;
  
  // build the searching area using the track and chamber resolutions and the maximum-distance-to-track value
  const TMatrixD& kParamCov = extrapTrackParam.GetCovariances();
  Double_t errX2 = kParamCov(0,0) + kMaxDZ * kMaxDZ * kParamCov(1,1) + 2. * kMaxDZ * TMath::Abs(kParamCov(0,1)) +
                   GetRecoParam()->GetDefaultNonBendingReso(chamber) * GetRecoParam()->GetDefaultNonBendingReso(chamber);
  Double_t errY2 = kParamCov(2,2) + kMaxDZ * kMaxDZ * kParamCov(3,3) + 2. * kMaxDZ * TMath::Abs(kParamCov(2,3)) +
		   GetRecoParam()->GetDefaultBendingReso(chamber) * GetRecoParam()->GetDefaultBendingReso(chamber);
  Double_t dX = TMath::Abs(trackParam.GetNonBendingSlope()) * kMaxDZ +
		GetRecoParam()->GetMaxNonBendingDistanceToTrack() +
		GetRecoParam()->GetSigmaCutForTracking() * TMath::Sqrt(2. * errX2);
  Double_t dY = TMath::Abs(trackParam.GetBendingSlope()) * kMaxDZ +
		GetRecoParam()->GetMaxBendingDistanceToTrack() +
		GetRecoParam()->GetSigmaCutForTracking() * TMath::Sqrt(2. * errY2);
  AliMpArea area(extrapTrackParam.GetNonBendingCoor(), 
                 extrapTrackParam.GetBendingCoor(),
                 dX, dY);
  
  // ask to cluterize in the given area of the given chamber
  fClusterServer->Clusterize(chamber, clusterStore, area, GetRecoParam());
  
}

  //__________________________________________________________________________
void AliMUONVTrackReconstructor::AskForNewClustersInStation(const AliMUONTrackParam &trackParam,
							    AliMUONVClusterStore& clusterStore, Int_t station)
{
  /// Ask the clustering to reconstruct new clusters around the track candidate position
  /// in the 2 chambers of the given station
  AskForNewClustersInChamber(trackParam, clusterStore, 2*station+1);
  AskForNewClustersInChamber(trackParam, clusterStore, 2*station);
}

  //__________________________________________________________________________
Double_t AliMUONVTrackReconstructor::TryOneCluster(const AliMUONTrackParam &trackParam, AliMUONVCluster* cluster,
						   AliMUONTrackParam &trackParamAtCluster, Bool_t updatePropagator)
{
/// Test the compatibility between the track and the cluster (using trackParam's covariance matrix):
/// return the corresponding Chi2
/// return trackParamAtCluster
  
  // extrapolate track parameters and covariances at the z position of the tested cluster
  // and set pointer to cluster into trackParamAtCluster
  trackParamAtCluster = trackParam;
  trackParamAtCluster.SetClusterPtr(cluster);
  if (!AliMUONTrackExtrap::ExtrapToZCov(&trackParamAtCluster, cluster->GetZ(), updatePropagator))
    return 2.*AliMUONTrack::MaxChi2();
  
  // Set differences between trackParam and cluster in the bending and non bending directions
  Double_t dX = cluster->GetX() - trackParamAtCluster.GetNonBendingCoor();
  Double_t dY = cluster->GetY() - trackParamAtCluster.GetBendingCoor();
  
  // Calculate errors and covariances
  const TMatrixD& kParamCov = trackParamAtCluster.GetCovariances();
  Double_t sigmaX2 = kParamCov(0,0) + cluster->GetErrX2();
  Double_t sigmaY2 = kParamCov(2,2) + cluster->GetErrY2();
  Double_t covXY   = kParamCov(0,2);
  Double_t det     = sigmaX2 * sigmaY2 - covXY * covXY;
  
  // Compute chi2
  if (det == 0.) return 2.*AliMUONTrack::MaxChi2();
  return (dX * dX * sigmaY2 + dY * dY * sigmaX2 - 2. * dX * dY * covXY) / det;
  
}

  //__________________________________________________________________________
Bool_t AliMUONVTrackReconstructor::TryOneClusterFast(const AliMUONTrackParam &trackParam, const AliMUONVCluster* cluster)
{
/// Test the compatibility between the track and the cluster
/// given the track and cluster resolutions + the maximum-distance-to-track value
/// and assuming linear propagation of the track:
/// return kTRUE if they are compatibles
  
  Double_t dZ = cluster->GetZ() - trackParam.GetZ();
  Double_t dX = cluster->GetX() - (trackParam.GetNonBendingCoor() + trackParam.GetNonBendingSlope() * dZ);
  Double_t dY = cluster->GetY() - (trackParam.GetBendingCoor() + trackParam.GetBendingSlope() * dZ);
  const TMatrixD& kParamCov = trackParam.GetCovariances();
  Double_t errX2 = kParamCov(0,0) + dZ * dZ * kParamCov(1,1) + 2. * dZ * kParamCov(0,1) + cluster->GetErrX2();
  Double_t errY2 = kParamCov(2,2) + dZ * dZ * kParamCov(3,3) + 2. * dZ * kParamCov(2,3) + cluster->GetErrY2();

  Double_t dXmax = GetRecoParam()->GetSigmaCutForTracking() * TMath::Sqrt(2. * errX2) +
                   GetRecoParam()->GetMaxNonBendingDistanceToTrack();
  Double_t dYmax = GetRecoParam()->GetSigmaCutForTracking() * TMath::Sqrt(2. * errY2) +
		   GetRecoParam()->GetMaxBendingDistanceToTrack();
  
  if (TMath::Abs(dX) > dXmax || TMath::Abs(dY) > dYmax) return kFALSE;
  
  return kTRUE;
  
}

  //__________________________________________________________________________
Double_t AliMUONVTrackReconstructor::TryTwoClustersFast(const AliMUONTrackParam &trackParamAtCluster1, AliMUONVCluster* cluster2,
						        AliMUONTrackParam &trackParamAtCluster2)
{
/// Test the compatibility between the track and the 2 clusters together (using trackParam's covariance matrix)
/// assuming linear propagation between the two clusters:
/// return the corresponding Chi2 accounting for covariances between the 2 clusters
/// return trackParamAtCluster2
  
  // extrapolate linearly track parameters and covariances at the z position of the second cluster
  trackParamAtCluster2 = trackParamAtCluster1;
  AliMUONTrackExtrap::LinearExtrapToZCov(&trackParamAtCluster2, cluster2->GetZ());
  
  // set pointer to cluster2 into trackParamAtCluster2
  trackParamAtCluster2.SetClusterPtr(cluster2);
  
  // Set differences between track and clusters in the bending and non bending directions
  AliMUONVCluster* cluster1 = trackParamAtCluster1.GetClusterPtr();
  Double_t dX1 = cluster1->GetX() - trackParamAtCluster1.GetNonBendingCoor();
  Double_t dX2 = cluster2->GetX() - trackParamAtCluster2.GetNonBendingCoor();
  Double_t dY1 = cluster1->GetY() - trackParamAtCluster1.GetBendingCoor();
  Double_t dY2 = cluster2->GetY() - trackParamAtCluster2.GetBendingCoor();
  
  // Calculate errors and covariances
  const TMatrixD& kParamCov1 = trackParamAtCluster1.GetCovariances();
  const TMatrixD& kParamCov2 = trackParamAtCluster2.GetCovariances();
  Double_t dZ = trackParamAtCluster2.GetZ() - trackParamAtCluster1.GetZ();
  Double_t sigma2X1 = kParamCov1(0,0) + cluster1->GetErrX2();
  Double_t sigma2X2 = kParamCov2(0,0) + cluster2->GetErrX2();
  Double_t covX1X2  = kParamCov1(0,0) + dZ * kParamCov1(0,1);
  Double_t sigma2Y1 = kParamCov1(2,2) + cluster1->GetErrY2();
  Double_t sigma2Y2 = kParamCov2(2,2) + cluster2->GetErrY2();
  Double_t covY1Y2  = kParamCov1(2,2) + dZ * kParamCov1(2,3);
  
  // Compute chi2
  Double_t detX = sigma2X1 * sigma2X2 - covX1X2 * covX1X2;
  Double_t detY = sigma2Y1 * sigma2Y2 - covY1Y2 * covY1Y2;
  if (detX == 0. || detY == 0.) return 2.*AliMUONTrack::MaxChi2();
  return   (dX1 * dX1 * sigma2X2 + dX2 * dX2 * sigma2X1 - 2. * dX1 * dX2 * covX1X2) / detX
  	 + (dY1 * dY1 * sigma2Y2 + dY2 * dY2 * sigma2Y1 - 2. * dY1 * dY2 * covY1Y2) / detY;
  
}

  //__________________________________________________________________________
Bool_t AliMUONVTrackReconstructor::FollowLinearTrackInChamber(AliMUONTrack &trackCandidate, const AliMUONVClusterStore& clusterStore,
							      Int_t nextChamber)
{
  /// Follow trackCandidate in chamber(0..) nextChamber assuming linear propagation, and search for compatible cluster(s)
  /// Keep all possibilities or only the best one(s) according to the flag fgkTrackAllTracks:
  /// kTRUE:  duplicate "trackCandidate" if there are several possibilities and add the new tracks at the end of
  ///         fRecTracksPtr to avoid conficts with other track candidates at this current stage of the tracking procedure.
  ///         Remove the obsolete "trackCandidate" at the end.
  /// kFALSE: add only the best cluster(s) to the "trackCandidate". Try to add a couple of clusters in priority.
  /// return kTRUE if new cluster(s) have been found (otherwise return kFALSE)
  AliDebug(1,Form("Enter FollowLinearTrackInChamber(1..) %d", nextChamber+1));
  
  Double_t chi2WithOneCluster = AliMUONTrack::MaxChi2();
  Double_t maxChi2WithOneCluster = 2. * GetRecoParam()->GetSigmaCutForTracking() *
					GetRecoParam()->GetSigmaCutForTracking(); // 2 because 2 quantities in chi2
  Double_t bestChi2WithOneCluster = maxChi2WithOneCluster;
  Bool_t foundOneCluster = kFALSE;
  AliMUONTrack *newTrack = 0x0;
  AliMUONVCluster *cluster;
  AliMUONTrackParam trackParam;
  AliMUONTrackParam extrapTrackParamAtCluster;
  AliMUONTrackParam bestTrackParamAtCluster;
  
  // Get track parameters according to the propagation direction
  if (nextChamber > 7) trackParam = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->Last();
  else trackParam = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->First();
  
  // Printout for debuging
  if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 2) || (AliLog::GetGlobalDebugLevel() >= 2)) {
    cout<<endl<<"Track parameters and covariances at first cluster:"<<endl;
    trackParam.GetParameters().Print();
    trackParam.GetCovariances().Print();
  }
  
  // Add MCS effect
  AliMUONTrackExtrap::AddMCSEffect(&trackParam,AliMUONConstants::ChamberThicknessInX0(trackParam.GetClusterPtr()->GetChamberId()),-1.);
  
  // Printout for debuging
  if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
    cout << "FollowLinearTrackInChamber: look for cluster in chamber(1..): " << nextChamber+1 << endl;
  }
  
  // Create iterators to loop over clusters in chamber
  TIter next(clusterStore.CreateChamberIterator(nextChamber,nextChamber));
  
  // look for candidates in chamber
  while ( ( cluster = static_cast<AliMUONVCluster*>(next()) ) ) {
    
    // try to add the current cluster fast
    if (!TryOneClusterFast(trackParam, cluster)) continue;
    
    // try to add the current cluster accuratly
    extrapTrackParamAtCluster = trackParam;
    AliMUONTrackExtrap::LinearExtrapToZCov(&extrapTrackParamAtCluster, cluster->GetZ());
    chi2WithOneCluster = TryOneCluster(extrapTrackParamAtCluster, cluster, extrapTrackParamAtCluster);
    
    // if good chi2 then consider to add cluster
    if (chi2WithOneCluster < maxChi2WithOneCluster) {
      foundOneCluster = kTRUE;
      
      // Printout for debuging
      if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
	cout << "FollowLinearTrackInChamber: found one cluster in chamber(1..): " << nextChamber+1
	<< " (Chi2 = " << chi2WithOneCluster << ")" << endl;
	cluster->Print();
      }
      
      if (GetRecoParam()->TrackAllTracks()) {
	// copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new cluster
	newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
	if (GetRecoParam()->RequestStation(nextChamber/2))
	  extrapTrackParamAtCluster.SetRemovable(kFALSE);
	else extrapTrackParamAtCluster.SetRemovable(kTRUE);
	newTrack->AddTrackParamAtCluster(extrapTrackParamAtCluster,*cluster);
	newTrack->SetGlobalChi2(trackCandidate.GetGlobalChi2()+chi2WithOneCluster);
	fNRecTracks++;
	
	// Printout for debuging
	if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
	  cout << "FollowLinearTrackInChamber: added one cluster in chamber(1..): " << nextChamber+1 << endl;
	  if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
	}
	
      } else if (chi2WithOneCluster < bestChi2WithOneCluster) {
	// keep track of the best cluster
	bestChi2WithOneCluster = chi2WithOneCluster;
	bestTrackParamAtCluster = extrapTrackParamAtCluster;
      }
      
    }
    
  }
  
  // fill out the best track if required else clean up the fRecTracksPtr array
  if (!GetRecoParam()->TrackAllTracks()) {
    if (foundOneCluster) {
      if (GetRecoParam()->RequestStation(nextChamber/2))
	bestTrackParamAtCluster.SetRemovable(kFALSE);
      else bestTrackParamAtCluster.SetRemovable(kTRUE);
      trackCandidate.AddTrackParamAtCluster(bestTrackParamAtCluster,*(bestTrackParamAtCluster.GetClusterPtr()));
      trackCandidate.SetGlobalChi2(trackCandidate.GetGlobalChi2()+bestChi2WithOneCluster);
      
      // Printout for debuging
      if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
        cout << "FollowLinearTrackInChamber: added the best cluster in chamber(1..): " << bestTrackParamAtCluster.GetClusterPtr()->GetChamberId()+1 << endl;
        if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
      }
      
    } else return kFALSE;
    
  } else if (foundOneCluster) {
    
    // remove obsolete track
    fRecTracksPtr->Remove(&trackCandidate);
    fNRecTracks--;
    
  } else return kFALSE;
  
  return kTRUE;
  
}

//__________________________________________________________________________
Bool_t AliMUONVTrackReconstructor::FollowLinearTrackInStation(AliMUONTrack &trackCandidate, const AliMUONVClusterStore& clusterStore,
							      Int_t nextStation)
{
  /// Follow trackCandidate in station(0..) nextStation assuming linear propagation, and search for compatible cluster(s)
  /// Keep all possibilities or only the best one(s) according to the flag fgkTrackAllTracks:
  /// kTRUE:  duplicate "trackCandidate" if there are several possibilities and add the new tracks at the end of
  ///         fRecTracksPtr to avoid conficts with other track candidates at this current stage of the tracking procedure.
  ///         Remove the obsolete "trackCandidate" at the end.
  /// kFALSE: add only the best cluster(s) to the "trackCandidate". Try to add a couple of clusters in priority.
  /// return kTRUE if new cluster(s) have been found (otherwise return kFALSE)
  AliDebug(1,Form("Enter FollowLinearTrackInStation(1..) %d", nextStation+1));
  
  // Order the chamber according to the propagation direction (tracking starts with chamber 2):
  // - nextStation == station(1...) 5 => forward propagation
  // - nextStation < station(1...) 5 => backward propagation
  Int_t ch1, ch2;
  if (nextStation==4) {
    ch1 = 2*nextStation+1;
    ch2 = 2*nextStation;
  } else {
    ch1 = 2*nextStation;
    ch2 = 2*nextStation+1;
  }
  
  Double_t chi2WithOneCluster = AliMUONTrack::MaxChi2();
  Double_t chi2WithTwoClusters = AliMUONTrack::MaxChi2();
  Double_t maxChi2WithOneCluster = 2. * GetRecoParam()->GetSigmaCutForTracking() *
					GetRecoParam()->GetSigmaCutForTracking(); // 2 because 2 quantities in chi2
  Double_t maxChi2WithTwoClusters = 4. * GetRecoParam()->GetSigmaCutForTracking() *
					 GetRecoParam()->GetSigmaCutForTracking(); // 4 because 4 quantities in chi2
  Double_t bestChi2WithOneCluster = maxChi2WithOneCluster;
  Double_t bestChi2WithTwoClusters = maxChi2WithTwoClusters;
  Bool_t foundOneCluster = kFALSE;
  Bool_t foundTwoClusters = kFALSE;
  AliMUONTrack *newTrack = 0x0;
  AliMUONVCluster *clusterCh1, *clusterCh2;
  AliMUONTrackParam trackParam;
  AliMUONTrackParam extrapTrackParamAtCluster1;
  AliMUONTrackParam extrapTrackParamAtCluster2;
  AliMUONTrackParam bestTrackParamAtCluster1;
  AliMUONTrackParam bestTrackParamAtCluster2;
  
  Int_t nClusters = clusterStore.GetSize();
  Bool_t *clusterCh1Used = new Bool_t[nClusters];
  for (Int_t i = 0; i < nClusters; i++) clusterCh1Used[i] = kFALSE;
  Int_t iCluster1;
  
  // Get track parameters according to the propagation direction
  if (nextStation==4) trackParam = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->Last();
  else trackParam = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->First();
  
  // Printout for debuging
  if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 2) || (AliLog::GetGlobalDebugLevel() >= 2)) {
    cout<<endl<<"Track parameters and covariances at first cluster:"<<endl;
    trackParam.GetParameters().Print();
    trackParam.GetCovariances().Print();
  }
  
  // Add MCS effect
  AliMUONTrackExtrap::AddMCSEffect(&trackParam,AliMUONConstants::ChamberThicknessInX0(trackParam.GetClusterPtr()->GetChamberId()),-1.);
  
  // Printout for debuging
  if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
    cout << "FollowLinearTrackInStation: look for clusters in chamber(1..): " << ch2+1 << endl;
  }
  
  // Create iterators to loop over clusters in both chambers
  TIter nextInCh1(clusterStore.CreateChamberIterator(ch1,ch1));
  TIter nextInCh2(clusterStore.CreateChamberIterator(ch2,ch2));
  
  // look for candidates in chamber 2
  while ( ( clusterCh2 = static_cast<AliMUONVCluster*>(nextInCh2()) ) ) {
    
    // try to add the current cluster fast
    if (!TryOneClusterFast(trackParam, clusterCh2)) continue;
    
    // try to add the current cluster accuratly
    extrapTrackParamAtCluster2 = trackParam;
    AliMUONTrackExtrap::LinearExtrapToZCov(&extrapTrackParamAtCluster2, clusterCh2->GetZ());
    chi2WithOneCluster = TryOneCluster(extrapTrackParamAtCluster2, clusterCh2, extrapTrackParamAtCluster2);
    
    // if good chi2 then try to attach a cluster in the other chamber too
    if (chi2WithOneCluster < maxChi2WithOneCluster) {
      Bool_t foundSecondCluster = kFALSE;
      
      // Printout for debuging
      if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
        cout << "FollowLinearTrackInStation: found one cluster in chamber(1..): " << ch2+1
	     << " (Chi2 = " << chi2WithOneCluster << ")" << endl;
	clusterCh2->Print();
        cout << "                      look for second clusters in chamber(1..): " << ch1+1 << " ..." << endl;
      }
      
      // add MCS effect
      AliMUONTrackExtrap::AddMCSEffect(&extrapTrackParamAtCluster2,AliMUONConstants::ChamberThicknessInX0(ch2),-1.);
      
      // reset cluster iterator of chamber 1
      nextInCh1.Reset();
      iCluster1 = -1;
      
      // look for second candidates in chamber 1
      while ( ( clusterCh1 = static_cast<AliMUONVCluster*>(nextInCh1()) ) ) {
        iCluster1++;
	
        // try to add the current cluster fast
        if (!TryOneClusterFast(extrapTrackParamAtCluster2, clusterCh1)) continue;
	
	// try to add the current cluster in addition to the one found in the previous chamber
	chi2WithTwoClusters = TryTwoClustersFast(extrapTrackParamAtCluster2, clusterCh1, extrapTrackParamAtCluster1);
        
	// if good chi2 then consider to add the 2 clusters to the "trackCandidate"
	if (chi2WithTwoClusters < maxChi2WithTwoClusters) {
	  foundSecondCluster = kTRUE;
	  foundTwoClusters = kTRUE;
          
	  // Printout for debuging
	  if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
	    cout << "FollowLinearTrackInStation: found one cluster in chamber(1..): " << ch1+1
		 << " (Chi2 = " << chi2WithTwoClusters << ")" << endl;
	    clusterCh1->Print();
	  }
	  
	  if (GetRecoParam()->TrackAllTracks()) {
	    // copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new clusters
            newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
	    extrapTrackParamAtCluster1.SetRemovable(kTRUE);
	    newTrack->AddTrackParamAtCluster(extrapTrackParamAtCluster1,*clusterCh1);
	    extrapTrackParamAtCluster2.SetRemovable(kTRUE);
	    newTrack->AddTrackParamAtCluster(extrapTrackParamAtCluster2,*clusterCh2);
	    newTrack->SetGlobalChi2(newTrack->GetGlobalChi2()+chi2WithTwoClusters);
	    fNRecTracks++;
	    
	    // Tag clusterCh1 as used
	    clusterCh1Used[iCluster1] = kTRUE;
	    
	    // Printout for debuging
	    if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
	      cout << "FollowLinearTrackInStation: added two clusters in station(1..): " << nextStation+1 << endl;
	      if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
	    }
	    
          } else if (chi2WithTwoClusters < bestChi2WithTwoClusters) {
	    // keep track of the best couple of clusters
	    bestChi2WithTwoClusters = chi2WithTwoClusters;
	    bestTrackParamAtCluster1 = extrapTrackParamAtCluster1;
	    bestTrackParamAtCluster2 = extrapTrackParamAtCluster2;
          }
	  
	}
	
      }
      
      // if no cluster found in chamber1 then consider to add clusterCh2 only
      if (!foundSecondCluster) {
	foundOneCluster = kTRUE;
        
	if (GetRecoParam()->TrackAllTracks()) {
	  // copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new cluster
          newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
	  if (GetRecoParam()->RequestStation(nextStation))
	    extrapTrackParamAtCluster2.SetRemovable(kFALSE);
	  else extrapTrackParamAtCluster2.SetRemovable(kTRUE);
	  newTrack->AddTrackParamAtCluster(extrapTrackParamAtCluster2,*clusterCh2);
	  newTrack->SetGlobalChi2(newTrack->GetGlobalChi2()+chi2WithOneCluster);
	  fNRecTracks++;
	  
	  // Printout for debuging
	  if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
	    cout << "FollowLinearTrackInStation: added one cluster in chamber(1..): " << ch2+1 << endl;
	    if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
	  }
	  
	} else if (!foundTwoClusters && chi2WithOneCluster < bestChi2WithOneCluster) {
	  // keep track of the best cluster except if a couple of clusters has already been found
	  bestChi2WithOneCluster = chi2WithOneCluster;
	  bestTrackParamAtCluster1 = extrapTrackParamAtCluster2;
        }
	
      }
      
    }
    
  }
  
  // look for candidates in chamber 1 not already attached to a track
  // if we want to keep all possible tracks or if no good couple of clusters has been found
  if (GetRecoParam()->TrackAllTracks() || !foundTwoClusters) {
    
    // Printout for debuging
    if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
      cout << "FollowLinearTrackInStation: look for single cluster in chamber(1..): " << ch1+1 << endl;
    }
    
    //Extrapolate trackCandidate to chamber "ch2"
    AliMUONTrackExtrap::LinearExtrapToZCov(&trackParam, AliMUONConstants::DefaultChamberZ(ch2));
    
    // add MCS effect for next step
    AliMUONTrackExtrap::AddMCSEffect(&trackParam,AliMUONConstants::ChamberThicknessInX0(ch2),-1.);
    
    // reset cluster iterator of chamber 1
    nextInCh1.Reset();
    iCluster1 = -1;
    
    // look for second candidates in chamber 1
    while ( ( clusterCh1 = static_cast<AliMUONVCluster*>(nextInCh1()) ) ) {
      iCluster1++;
      
      if (clusterCh1Used[iCluster1]) continue; // Skip clusters already used
      
      // try to add the current cluster fast
      if (!TryOneClusterFast(trackParam, clusterCh1)) continue;
      
      // try to add the current cluster accuratly
      extrapTrackParamAtCluster1 = trackParam;
      AliMUONTrackExtrap::LinearExtrapToZCov(&extrapTrackParamAtCluster1, clusterCh1->GetZ());
      chi2WithOneCluster = TryOneCluster(extrapTrackParamAtCluster1, clusterCh1, extrapTrackParamAtCluster1);
      
      // if good chi2 then consider to add clusterCh1
      // We do not try to attach a cluster in the other chamber too since it has already been done above
      if (chi2WithOneCluster < maxChi2WithOneCluster) {
	foundOneCluster = kTRUE;
  	
	// Printout for debuging
  	if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
  	  cout << "FollowLinearTrackInStation: found one cluster in chamber(1..): " << ch1+1
  	       << " (Chi2 = " << chi2WithOneCluster << ")" << endl;
	  clusterCh1->Print();
  	}
	
	if (GetRecoParam()->TrackAllTracks()) {
	  // copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new cluster
  	  newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
	  if (GetRecoParam()->RequestStation(nextStation))
	    extrapTrackParamAtCluster1.SetRemovable(kFALSE);
	  else extrapTrackParamAtCluster1.SetRemovable(kTRUE);
	  newTrack->AddTrackParamAtCluster(extrapTrackParamAtCluster1,*clusterCh1);
	  newTrack->SetGlobalChi2(newTrack->GetGlobalChi2()+chi2WithOneCluster);
	  fNRecTracks++;
  	  
	  // Printout for debuging
  	  if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
  	    cout << "FollowLinearTrackInStation: added one cluster in chamber(1..): " << ch1+1 << endl;
  	    if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
  	  }
	  
	} else if (chi2WithOneCluster < bestChi2WithOneCluster) {
	  // keep track of the best cluster except if a couple of clusters has already been found
	  bestChi2WithOneCluster = chi2WithOneCluster;
	  bestTrackParamAtCluster1 = extrapTrackParamAtCluster1;
  	}
	
      }
      
    }
    
  }
  
  // fill out the best track if required else clean up the fRecTracksPtr array
  if (!GetRecoParam()->TrackAllTracks()) {
    if (foundTwoClusters) {
      bestTrackParamAtCluster1.SetRemovable(kTRUE);
      trackCandidate.AddTrackParamAtCluster(bestTrackParamAtCluster1,*(bestTrackParamAtCluster1.GetClusterPtr()));
      bestTrackParamAtCluster2.SetRemovable(kTRUE);
      trackCandidate.AddTrackParamAtCluster(bestTrackParamAtCluster2,*(bestTrackParamAtCluster2.GetClusterPtr()));
      trackCandidate.SetGlobalChi2(trackCandidate.GetGlobalChi2()+bestChi2WithTwoClusters);
      
      // Printout for debuging
      if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
        cout << "FollowLinearTrackInStation: added the two best clusters in station(1..): " << nextStation+1 << endl;
        if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
      }
      
    } else if (foundOneCluster) {
      if (GetRecoParam()->RequestStation(nextStation))
	bestTrackParamAtCluster1.SetRemovable(kFALSE);
      else bestTrackParamAtCluster1.SetRemovable(kTRUE);
      trackCandidate.AddTrackParamAtCluster(bestTrackParamAtCluster1,*(bestTrackParamAtCluster1.GetClusterPtr()));
      trackCandidate.SetGlobalChi2(trackCandidate.GetGlobalChi2()+bestChi2WithOneCluster);
      
      // Printout for debuging
      if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
        cout << "FollowLinearTrackInStation: added the best cluster in chamber(1..): " << bestTrackParamAtCluster1.GetClusterPtr()->GetChamberId()+1 << endl;
        if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
      }
      
    } else {
      delete [] clusterCh1Used;
      return kFALSE;
    }
    
  } else if (foundOneCluster || foundTwoClusters) {
    
    // remove obsolete track
    fRecTracksPtr->Remove(&trackCandidate);
    fNRecTracks--;
    
  } else {
    delete [] clusterCh1Used;  
    return kFALSE;
  }
  
  delete [] clusterCh1Used;
  return kTRUE;
  
}

//__________________________________________________________________________
void AliMUONVTrackReconstructor::ImproveTracks()
{
  /// Improve tracks by removing clusters with local chi2 highter than the defined cut
  /// Recompute track parameters and covariances at the remaining clusters
  AliDebug(1,"Enter ImproveTracks");
  
  AliMUONTrack *track, *nextTrack;
  
  track = (AliMUONTrack*) fRecTracksPtr->First();
  while (track) {
    
    // prepare next track in case the actual track is suppressed
    nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track);
    
    ImproveTrack(*track);
    
    // remove track if improvement failed
    if (!track->IsImproved()) {
      fRecTracksPtr->Remove(track);
      fNRecTracks--;
    }
    
    track = nextTrack;
  }
  
  // compress the array in case of some tracks have been removed
  fRecTracksPtr->Compress();
  
}

//__________________________________________________________________________
void AliMUONVTrackReconstructor::Finalize()
{
  /// Recompute track parameters and covariances at each attached cluster
  /// Set the label pointing to the corresponding MC track
  /// Remove the track if finalization failed
  
  AliMUONTrack *track, *nextTrack;
  Bool_t trackRemoved = kFALSE;
  
  track = (AliMUONTrack*) fRecTracksPtr->First();
  while (track) {
    
    nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track);
    
    if (FinalizeTrack(*track)) track->FindMCLabel();
    else {
      fRecTracksPtr->Remove(track);
      fNRecTracks--;
      trackRemoved = kTRUE;
    }
    
    track = nextTrack;
    
  }
  
  // compress array of tracks if needed
  if (trackRemoved) fRecTracksPtr->Compress();
  
}

//__________________________________________________________________________
void AliMUONVTrackReconstructor::DiscardMonoCathodClusters()
{
  /// Assign a different resolution to the mono-cathod clusters
  /// in the direction of the missing plane and refit the track
  /// Remove the track in case of failure
  
  if (!fkTransformer) AliFatal("missing geometry transformer");
  
  AliMUONTrack *track, *nextTrack;
  Bool_t trackRemoved = kFALSE;
  
  track = (AliMUONTrack*) fRecTracksPtr->First();
  while (track) {
    
    nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track);
    
    ChangeMonoCathodClusterRes(*track);
    
    if (!RefitTrack(*track) || (GetRecoParam()->ImproveTracks() && !track->IsImproved())) {
      fRecTracksPtr->Remove(track);
      fNRecTracks--;
      trackRemoved = kTRUE;
    }
    
    track = nextTrack;
    
  }
  
  // compress array of tracks if needed
  if (trackRemoved) fRecTracksPtr->Compress();
  
}

//__________________________________________________________________________
void AliMUONVTrackReconstructor::ChangeMonoCathodClusterRes(AliMUONTrack &track)
{
  /// Assign a different resolution to the mono-cathod clusters
  /// in the direction of the missing plane and refit the track
  
  // Loop over clusters
  AliMUONVCluster *cluster;
  Int_t nClusters = track.GetNClusters();
  for (Int_t iCluster = 0; iCluster < nClusters; iCluster++) {
    cluster = ((AliMUONTrackParam*) track.GetTrackParamAtCluster()->UncheckedAt(iCluster))->GetClusterPtr();
    
    // do it only for stations 3, 4 & 5
    if (cluster->GetChamberId() < 4) continue;
    
    // get the cathod corresponding to the bending/non-bending plane
    Int_t deId = cluster->GetDetElemId();
    AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(deId, kFALSE);
    if (!de) continue;
    AliMp::CathodType cath1 = de->GetCathodType(AliMp::kBendingPlane); 
    AliMp::CathodType cath2 = de->GetCathodType(AliMp::kNonBendingPlane); 
    
    // get the corresponding segmentation
    const AliMpVSegmentation* seg1 = AliMpSegmentation::Instance()->GetMpSegmentation(deId, cath1);
    const AliMpVSegmentation* seg2 = AliMpSegmentation::Instance()->GetMpSegmentation(deId, cath2);
    if (!seg1 || !seg2) continue;
    
    // get local coordinate of the cluster
    Double_t lX,lY,lZ;
    Double_t gX = cluster->GetX();
    Double_t gY = cluster->GetY();
    Double_t gZ = cluster->GetZ();
    fkTransformer->Global2Local(deId,gX,gY,gZ,lX,lY,lZ);
    
    // find pads below the cluster
    AliMpPad pad1 = seg1->PadByPosition(lX, lY, kFALSE);
    AliMpPad pad2 = seg2->PadByPosition(lX, lY, kFALSE);
    
    // build their ID if pads are valid
    UInt_t padId1 = (pad1.IsValid()) ? AliMUONVDigit::BuildUniqueID(deId, pad1.GetManuId(), pad1.GetManuChannel(), cath1) : 0;
    UInt_t padId2 = (pad2.IsValid()) ? AliMUONVDigit::BuildUniqueID(deId, pad2.GetManuId(), pad2.GetManuChannel(), cath2) : 0;
    
    // check if the cluster contains these pads 
    Bool_t hasNonBending = kFALSE;
    Bool_t hasBending = kFALSE;
    for (Int_t iDigit = 0; iDigit < cluster->GetNDigits(); iDigit++) {
      
      UInt_t digitId = cluster->GetDigitId(iDigit);
      
      if (digitId == padId1) {
	
	hasBending = kTRUE;
	if (hasNonBending) break;
	
      } else if (digitId == padId2) {
	
	hasNonBending = kTRUE;
	if (hasBending) break;
	
      }
      
    }
    
    // modify the cluster resolution if needed
    if (!hasNonBending) cluster->SetErrXY(GetRecoParam()->GetMonoCathodClNonBendingRes(), cluster->GetErrY());
    if (!hasBending) cluster->SetErrXY(cluster->GetErrX(), GetRecoParam()->GetMonoCathodClBendingRes());
    
  }
  
}

//__________________________________________________________________________
void AliMUONVTrackReconstructor::ValidateTracksWithTrigger(AliMUONVTrackStore& trackStore,
                                                           const AliMUONVTriggerTrackStore& triggerTrackStore,
                                                           const AliMUONVTriggerStore& triggerStore,
                                                           const AliMUONTrackHitPattern& trackHitPattern)
{
  /// Try to match track from tracking system with trigger track
  AliCodeTimerAuto("",0);

  trackHitPattern.ExecuteValidation(trackStore, triggerTrackStore, triggerStore);
}


//__________________________________________________________________________
void AliMUONVTrackReconstructor::EventReconstructTrigger(const AliMUONTriggerCircuit& circuit,
                                                         const AliMUONVTriggerStore& triggerStore,
                                                         AliMUONVTriggerTrackStore& triggerTrackStore)
{
  /// Fill trigger track store from local trigger
  AliDebug(1, "");
  AliCodeTimerAuto("",0);

  AliMUONGlobalTrigger* globalTrigger = triggerStore.Global();
  
  UChar_t gloTrigPat = 0;

  if (globalTrigger)
  {
    gloTrigPat = globalTrigger->GetGlobalResponse();
  }
  
  AliMUONTriggerTrack triggerTrack;
  
  TIter next(triggerStore.CreateIterator());
  AliMUONLocalTrigger* locTrg(0x0);
  
  while ( ( locTrg = static_cast<AliMUONLocalTrigger*>(next()) ) )
  {
    if ( locTrg->IsTrigX() && locTrg->IsTrigY() ) 
    { // make Trigger Track if trigger in X and Y
      
      if (TriggerToTrack(circuit, *locTrg, triggerTrack, gloTrigPat))
	triggerTrackStore.Add(triggerTrack);

    } // board is fired 
  } // end of loop on Local Trigger
}

//__________________________________________________________________________
Bool_t AliMUONVTrackReconstructor::TriggerToTrack(const AliMUONTriggerCircuit& circuit,
                                                const AliMUONLocalTrigger& locTrg,
                                                AliMUONTriggerTrack& triggerTrack,
                                                UChar_t globalTriggerPattern)
{
  /// To make the trigger tracks from Local Trigger
  const Double_t kTrigNonBendReso = AliMUONConstants::TriggerNonBendingReso();
  const Double_t kTrigBendReso = AliMUONConstants::TriggerBendingReso();
  const Double_t kSqrt12 = TMath::Sqrt(12.);
  
  TMatrixD trigCov(3,3);

  Int_t localBoardId = locTrg.LoCircuit();
      
  Float_t y11 = circuit.GetY11Pos(localBoardId, locTrg.LoStripX()); 
  Float_t z11 = circuit.GetZ11Pos(localBoardId, locTrg.LoStripX());
  // need first to convert deviation to [0-30] 
  // (see AliMUONLocalTriggerBoard::LocalTrigger)
  Int_t deviation = locTrg.GetDeviation(); 
  Int_t stripX21 = locTrg.LoStripX()+deviation+1;
  Float_t y21 = circuit.GetY21Pos(localBoardId, stripX21);       
  Float_t z21 = circuit.GetZ21Pos(localBoardId, stripX21);
  Float_t x11 = circuit.GetX11Pos(localBoardId, locTrg.LoStripY());
      
  AliDebug(1, Form(" MakeTriggerTrack %3d %2d %2d %2d (%f %f %f) (%f %f)\n",locTrg.LoCircuit(),
                   locTrg.LoStripX(),locTrg.LoStripX()+deviation+1,locTrg.LoStripY(),x11, y11, z11, y21, z21));
      
  if (TMath::Abs(z11) < 0.00001) return kFALSE;

  Double_t deltaZ = z11 - z21;
      
  Float_t slopeX = x11/z11;
  Float_t slopeY = (y11-y21) / deltaZ;
      
  Float_t sigmaX = circuit.GetX11Width(localBoardId, locTrg.LoStripY()) / kSqrt12;
  Float_t sigmaY = circuit.GetY11Width(localBoardId, locTrg.LoStripX()) / kSqrt12;
  Float_t sigmaY21 = circuit.GetY21Width(localBoardId, locTrg.LoStripX()) / kSqrt12;
      
  trigCov.Zero();
  trigCov(0,0) = kTrigNonBendReso * kTrigNonBendReso + sigmaX * sigmaX;
  trigCov(1,1) = kTrigBendReso * kTrigBendReso + sigmaY * sigmaY;
  trigCov(2,2) = 
    (2. * kTrigBendReso * kTrigBendReso + sigmaY * sigmaY + sigmaY21 * sigmaY21 ) / deltaZ / deltaZ;
    trigCov(1,2) = trigCov(2,1) = trigCov(1,1) / deltaZ;
      
  triggerTrack.SetX11(x11);
  triggerTrack.SetY11(y11);
  triggerTrack.SetZ11(z11);
  triggerTrack.SetZ21(z21);
  triggerTrack.SetSlopeX(slopeX);
  triggerTrack.SetSlopeY(slopeY);
  triggerTrack.SetGTPattern(globalTriggerPattern);
  triggerTrack.SetLoTrgNum(localBoardId);
  triggerTrack.SetCovariances(trigCov);
  triggerTrack.SetUniqueID(locTrg.GetUniqueID());

  return kTRUE;

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