]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - MUON/AliMUONTrack.cxx
Removing an AliInfo
[u/mrichter/AliRoot.git] / MUON / AliMUONTrack.cxx
index 1ac964eceef02e8e3e169dd674ac25187995f824..139e479afe2a86be4e0d70a8b4441b2b155a3e68 100644 (file)
 
 #include "AliMUONTrack.h"
 
-#include "AliMUONTrackParam.h" 
-#include "AliMUONVCluster.h" 
-#include "AliMUONObjectPair.h" 
+#include "AliMUONReconstructor.h"
+#include "AliMUONRecoParam.h"
+#include "AliMUONVCluster.h"
+#include "AliMUONVClusterStore.h"
+#include "AliMUONObjectPair.h"
 #include "AliMUONConstants.h"
-#include "AliMUONTrackExtrap.h" 
+#include "AliMUONTrackExtrap.h"
 
 #include "AliLog.h"
 
@@ -43,7 +45,7 @@ ClassImp(AliMUONTrack) // Class implementation in ROOT context
 //__________________________________________________________________________
 AliMUONTrack::AliMUONTrack()
   : TObject(),
-    fTrackParamAtCluster(new TClonesArray("AliMUONTrackParam",10)),
+    fTrackParamAtCluster(0x0),
     fFitWithVertex(kFALSE),
     fVertexErrXY2(),
     fFitWithMCS(kFALSE),
@@ -60,7 +62,6 @@ AliMUONTrack::AliMUONTrack()
     fLocalTrigger(0)
 {
   /// Default constructor
-  fTrackParamAtCluster->SetOwner(kTRUE);
   fVertexErrXY2[0] = 0.;
   fVertexErrXY2[1] = 0.;
 }
@@ -85,31 +86,13 @@ AliMUONTrack::AliMUONTrack(AliMUONObjectPair *segment)
     fLocalTrigger(0)
 {
   /// Constructor from two clusters
-  fTrackParamAtCluster->SetOwner(kTRUE);
   
   fVertexErrXY2[0] = 0.;
   fVertexErrXY2[1] = 0.;
   
   // Pointers to clusters from the segment
-  AliMUONVCluster* cluster1 = (AliMUONVCluster*) segment->First();
-  AliMUONVCluster* cluster2 = (AliMUONVCluster*) segment->Second();
-  
-  // check sorting in -Z (spectro z<0)
-  if (cluster1->GetZ() < cluster2->GetZ()) {
-    cluster1 = cluster2;
-    cluster2 = (AliMUONVCluster*) segment->First();
-  }
-  
-  // order the clusters into the track according to the station the segment belong to
-  //(the cluster first attached is the one from which we will start the tracking procedure)
-  AliMUONVCluster *firstCluster, *lastCluster;
-  if (cluster1->GetChamberId() == 8) {
-    firstCluster = cluster1;
-    lastCluster = cluster2;
-  } else {
-    firstCluster = cluster2;
-    lastCluster = cluster1;
-  }
+  AliMUONVCluster* firstCluster = (AliMUONVCluster*) segment->First();
+  AliMUONVCluster* lastCluster = (AliMUONVCluster*) segment->Second();
   
   // Compute track parameters
   Double_t z1 = firstCluster->GetZ();
@@ -127,7 +110,6 @@ AliMUONTrack::AliMUONTrack(AliMUONObjectPair *segment)
   Double_t bendingImpact = bendingCoor1 - z1 * bendingSlope;
   Double_t inverseBendingMomentum = 1. / AliMUONTrackExtrap::GetBendingMomentumFromImpactParam(bendingImpact);
   
-  
   // Set track parameters at first cluster
   AliMUONTrackParam trackParamAtFirstCluster;
   trackParamAtFirstCluster.SetZ(z1);
@@ -137,7 +119,6 @@ AliMUONTrack::AliMUONTrack(AliMUONObjectPair *segment)
   trackParamAtFirstCluster.SetBendingSlope(bendingSlope);
   trackParamAtFirstCluster.SetInverseBendingMomentum(inverseBendingMomentum);
   
-  
   // Set track parameters at last cluster
   AliMUONTrackParam trackParamAtLastCluster;
   trackParamAtLastCluster.SetZ(z2);
@@ -147,43 +128,40 @@ AliMUONTrack::AliMUONTrack(AliMUONObjectPair *segment)
   trackParamAtLastCluster.SetBendingSlope(bendingSlope);
   trackParamAtLastCluster.SetInverseBendingMomentum(inverseBendingMomentum);
   
-  
   // Compute and set track parameters covariances at first cluster
-  TMatrixD paramCov1(5,5);
-  paramCov1.Zero();
+  TMatrixD paramCov(5,5);
+  paramCov.Zero();
   // Non bending plane
-  paramCov1(0,0) = firstCluster->GetErrX2();
-  paramCov1(0,1) = firstCluster->GetErrX2() / dZ;
-  paramCov1(1,0) = paramCov1(0,1);
-  paramCov1(1,1) = ( firstCluster->GetErrX2() + lastCluster->GetErrX2() ) / dZ / dZ;
+  paramCov(0,0) = firstCluster->GetErrX2();
+  paramCov(0,1) = firstCluster->GetErrX2() / dZ;
+  paramCov(1,0) = paramCov(0,1);
+  paramCov(1,1) = ( firstCluster->GetErrX2() + lastCluster->GetErrX2() ) / dZ / dZ;
   // Bending plane
-  paramCov1(2,2) = firstCluster->GetErrY2();
-  paramCov1(2,3) = firstCluster->GetErrY2() / dZ;
-  paramCov1(3,2) = paramCov1(2,3);
-  paramCov1(3,3) = ( firstCluster->GetErrY2() + lastCluster->GetErrY2() ) / dZ / dZ;
-  // Inverse bending momentum (50% error)
-  paramCov1(4,4) = 0.5*inverseBendingMomentum * 0.5*inverseBendingMomentum;
-  // Set covariances
-  trackParamAtFirstCluster.SetCovariances(paramCov1);
-  
+  paramCov(2,2) = firstCluster->GetErrY2();
+  paramCov(2,3) = firstCluster->GetErrY2() / dZ;
+  paramCov(3,2) = paramCov(2,3);
+  paramCov(3,3) = ( firstCluster->GetErrY2() + lastCluster->GetErrY2() ) / dZ / dZ;
+  // Inverse bending momentum (vertex resolution + bending slope resolution + 10% error on dipole parameters+field)
+  paramCov(4,4) = ((AliMUONReconstructor::GetRecoParam()->GetBendingVertexDispersion() *
+                   AliMUONReconstructor::GetRecoParam()->GetBendingVertexDispersion() +
+                   (z1 * z1 * lastCluster->GetErrY2() + z2 * z2 * firstCluster->GetErrY2()) / dZ / dZ) /
+                  bendingImpact / bendingImpact + 0.1 * 0.1) * inverseBendingMomentum * inverseBendingMomentum;
+  paramCov(2,4) = - z2 * firstCluster->GetErrY2() * inverseBendingMomentum / bendingImpact / dZ;
+  paramCov(4,2) = paramCov(2,4);
+  paramCov(3,4) = - (z1 * lastCluster->GetErrY2() + z2 * firstCluster->GetErrY2()) * inverseBendingMomentum / bendingImpact / dZ / dZ;
+  paramCov(4,3) = paramCov(3,4);
   
-  // Compute and set track parameters covariances at last cluster (as if the first cluster did not exist)
-  TMatrixD paramCov2(5,5);
-  paramCov2.Zero();
-  // Non bending plane
-  paramCov2(0,0) = paramCov1(0,0);
-  paramCov2(1,1) = 100.*paramCov1(1,1);
-  // Bending plane
-  paramCov2(2,2) = paramCov1(2,2);
-  paramCov2(3,3) = 100.*paramCov1(3,3);
-  // Inverse bending momentum
-  paramCov2(4,4) = paramCov1(4,4);
   // Set covariances
-  trackParamAtLastCluster.SetCovariances(paramCov2);
+  trackParamAtFirstCluster.SetCovariances(paramCov);
   
-  // Flag clusters as being removable
-  trackParamAtFirstCluster.SetRemovable(kTRUE);
-  trackParamAtLastCluster.SetRemovable(kTRUE);
+  // Compute and set track parameters covariances at last cluster
+  paramCov(1,0) = - paramCov(1,0);
+  paramCov(0,1) = - paramCov(0,1);
+  paramCov(3,2) = - paramCov(3,2);
+  paramCov(2,3) = - paramCov(2,3);
+  paramCov(2,4) = z1 * lastCluster->GetErrY2() * inverseBendingMomentum / bendingImpact / dZ;
+  paramCov(4,2) = paramCov(2,4);
+  trackParamAtLastCluster.SetCovariances(paramCov);
   
   // Add track parameters at clusters
   AddTrackParamAtCluster(trackParamAtFirstCluster,*firstCluster);
@@ -191,10 +169,10 @@ AliMUONTrack::AliMUONTrack(AliMUONObjectPair *segment)
   
 }
 
-  //__________________________________________________________________________
+//__________________________________________________________________________
 AliMUONTrack::AliMUONTrack(const AliMUONTrack& track)
   : TObject(track),
-    fTrackParamAtCluster(new TClonesArray("AliMUONTrackParam",10)),
+    fTrackParamAtCluster(0x0),
     fFitWithVertex(track.fFitWithVertex),
     fVertexErrXY2(),
     fFitWithMCS(track.fFitWithMCS),
@@ -213,10 +191,13 @@ AliMUONTrack::AliMUONTrack(const AliMUONTrack& track)
   ///copy constructor
   
   // necessary to make a copy of the objects and not only the pointers in TClonesArray.
-  AliMUONTrackParam *trackParamAtCluster = (AliMUONTrackParam*) track.fTrackParamAtCluster->First();
-  while (trackParamAtCluster) {
-    new ((*fTrackParamAtCluster)[GetNClusters()]) AliMUONTrackParam(*trackParamAtCluster);
-    trackParamAtCluster = (AliMUONTrackParam*) track.fTrackParamAtCluster->After(trackParamAtCluster);
+  if (track.fTrackParamAtCluster) {
+    fTrackParamAtCluster = new TClonesArray("AliMUONTrackParam",10);
+    AliMUONTrackParam *trackParamAtCluster = (AliMUONTrackParam*) track.fTrackParamAtCluster->First();
+    while (trackParamAtCluster) {
+      new ((*fTrackParamAtCluster)[GetNClusters()]) AliMUONTrackParam(*trackParamAtCluster);
+      trackParamAtCluster = (AliMUONTrackParam*) track.fTrackParamAtCluster->After(trackParamAtCluster);
+    }
   }
   
   // copy vertex resolution square used during the tracking procedure
@@ -242,13 +223,18 @@ AliMUONTrack & AliMUONTrack::operator=(const AliMUONTrack& track)
 
   // base class assignement
   TObject::operator=(track);
-
-  // necessary to make a copy of the objects and not only the pointers in TClonesArray.
-  fTrackParamAtCluster = new TClonesArray("AliMUONTrackParam",10);
-  AliMUONTrackParam *trackParamAtCluster = (AliMUONTrackParam*) track.fTrackParamAtCluster->First();
-  while (trackParamAtCluster) {
-    new ((*fTrackParamAtCluster)[GetNClusters()]) AliMUONTrackParam(*trackParamAtCluster);
-    trackParamAtCluster = (AliMUONTrackParam*) track.fTrackParamAtCluster->After(trackParamAtCluster);
+  
+  // clear memory
+  Clear();
+  
+  // necessary to make a copy of the objects and not only the pointers in TClonesArray
+  if (track.fTrackParamAtCluster) {
+    fTrackParamAtCluster = new TClonesArray("AliMUONTrackParam",10);
+    AliMUONTrackParam *trackParamAtCluster = (AliMUONTrackParam*) track.fTrackParamAtCluster->First();
+    while (trackParamAtCluster) {
+      new ((*fTrackParamAtCluster)[GetNClusters()]) AliMUONTrackParam(*trackParamAtCluster);
+      trackParamAtCluster = (AliMUONTrackParam*) track.fTrackParamAtCluster->After(trackParamAtCluster);
+    }
   }
   
   // copy cluster weights matrix if any
@@ -257,9 +243,6 @@ AliMUONTrack & AliMUONTrack::operator=(const AliMUONTrack& track)
       fClusterWeightsNonBending->ResizeTo(*(track.fClusterWeightsNonBending));
       *fClusterWeightsNonBending = *(track.fClusterWeightsNonBending);
     } else fClusterWeightsNonBending = new TMatrixD(*(track.fClusterWeightsNonBending));
-  } else if (fClusterWeightsNonBending) {
-    delete fClusterWeightsNonBending;
-    fClusterWeightsNonBending = 0x0;
   }
   
   // copy cluster weights matrix if any
@@ -268,18 +251,12 @@ AliMUONTrack & AliMUONTrack::operator=(const AliMUONTrack& track)
       fClusterWeightsBending->ResizeTo(*(track.fClusterWeightsBending));
       *fClusterWeightsBending = *(track.fClusterWeightsBending);
     } else fClusterWeightsBending = new TMatrixD(*(track.fClusterWeightsBending));
-  } else if (fClusterWeightsBending) {
-    delete fClusterWeightsBending;
-    fClusterWeightsBending = 0x0;
   }
   
   // copy track parameters at vertex if any
   if (track.fTrackParamAtVertex) {
     if (fTrackParamAtVertex) *fTrackParamAtVertex = *(track.fTrackParamAtVertex);
     else fTrackParamAtVertex = new AliMUONTrackParam(*(track.fTrackParamAtVertex));
-  } else if (fTrackParamAtVertex) {
-    delete fTrackParamAtVertex;
-    fTrackParamAtVertex = 0x0;
   }
   
   fFitWithVertex      =  track.fFitWithVertex;
@@ -312,12 +289,47 @@ AliMUONTrack::~AliMUONTrack()
 void AliMUONTrack::Clear(Option_t* opt)
 {
   /// Clear arrays
-  fTrackParamAtCluster->Clear(opt);
+  if (opt && opt[0] == 'C' && fTrackParamAtCluster) fTrackParamAtCluster->Clear("C");
+  else {
+    delete fTrackParamAtCluster;
+    fTrackParamAtCluster = 0x0;
+  }
+  delete fClusterWeightsNonBending; fClusterWeightsNonBending = 0x0;
+  delete fClusterWeightsBending; fClusterWeightsBending = 0x0;
+  delete fTrackParamAtVertex; fTrackParamAtVertex = 0x0;
+}
+
+  //__________________________________________________________________________
+void AliMUONTrack::Reset()
+{
+  /// Reset to default values
+  SetUniqueID(0);
+  fFitWithVertex = kFALSE;
+  fVertexErrXY2[0] = 0.;
+  fVertexErrXY2[1] = 0.;
+  fFitWithMCS = kFALSE;
+  fGlobalChi2 = -1.;
+  fImproved = kFALSE;
+  fMatchTrigger = -1;
+  floTrgNum = -1;
+  fChi2MatchTrigger = 0.;
+  fTrackID = 0;
+  fHitsPatternInTrigCh = 0;
+  fLocalTrigger = 0;
+  delete fTrackParamAtCluster; fTrackParamAtCluster = 0x0;
   delete fClusterWeightsNonBending; fClusterWeightsNonBending = 0x0;
   delete fClusterWeightsBending; fClusterWeightsBending = 0x0;
   delete fTrackParamAtVertex; fTrackParamAtVertex = 0x0;
 }
 
+  //__________________________________________________________________________
+TClonesArray* AliMUONTrack::GetTrackParamAtCluster() const
+{
+  /// return array of track parameters at cluster (create it if needed)
+  if (!fTrackParamAtCluster) fTrackParamAtCluster = new TClonesArray("AliMUONTrackParam",10);
+  return fTrackParamAtCluster;
+}
+
   //__________________________________________________________________________
 void AliMUONTrack::AddTrackParamAtCluster(const AliMUONTrackParam &trackParam, AliMUONVCluster &cluster, Bool_t copy)
 {
@@ -339,6 +351,7 @@ void AliMUONTrack::AddTrackParamAtCluster(const AliMUONTrackParam &trackParam, A
   }
   
   // add parameters to the array of track parameters
+  if (!fTrackParamAtCluster) fTrackParamAtCluster = new TClonesArray("AliMUONTrackParam",10);
   AliMUONTrackParam* trackParamAtCluster = new ((*fTrackParamAtCluster)[GetNClusters()]) AliMUONTrackParam(trackParam);
   
   // link parameters with the associated cluster or its copy
@@ -346,13 +359,16 @@ void AliMUONTrack::AddTrackParamAtCluster(const AliMUONTrackParam &trackParam, A
     AliMUONVCluster *clusterCopy = static_cast<AliMUONVCluster*>(cluster.Clone());
     trackParamAtCluster->SetClusterPtr(clusterCopy, kTRUE);
   } else trackParamAtCluster->SetClusterPtr(&cluster);
+  
+  // sort the array of track parameters
+  fTrackParamAtCluster->Sort();
 }
 
   //__________________________________________________________________________
 void AliMUONTrack::RemoveTrackParamAtCluster(AliMUONTrackParam *trackParam)
 {
   /// Remove trackParam from the array of TrackParamAtCluster
-  if (!fTrackParamAtCluster->Remove(trackParam)) {
+  if (!fTrackParamAtCluster || !fTrackParamAtCluster->Remove(trackParam)) {
     AliWarning("object to remove does not exist in array fTrackParamAtCluster");
     return;
   }
@@ -432,6 +448,78 @@ void AliMUONTrack::UpdateCovTrackParamAtCluster()
     trackParamAtCluster = (AliMUONTrackParam*) (fTrackParamAtCluster->After(trackParamAtCluster));
   }
   
+}
+
+  //__________________________________________________________________________
+Bool_t AliMUONTrack::IsValid()
+{
+  /// check the validity of the current track (at least one cluster per requested station)
+  
+  Int_t nClusters = GetNClusters();
+  AliMUONTrackParam *trackParam;
+  Int_t currentStation = 0, expectedStation = 0;
+  
+  for (Int_t i = 0; i < nClusters; i++) {
+    trackParam = (AliMUONTrackParam*) fTrackParamAtCluster->UncheckedAt(i);
+    
+    // skip unrequested stations
+    while (expectedStation < AliMUONConstants::NTrackingSt() &&
+          !AliMUONReconstructor::GetRecoParam()->RequestStation(expectedStation)) expectedStation++;
+    
+    currentStation = trackParam->GetClusterPtr()->GetChamberId()/2;
+    
+    // missing station
+    if (currentStation > expectedStation) return kFALSE;
+    
+    // found station --> look for next one
+    if (currentStation == expectedStation) expectedStation++;
+    
+  }
+  
+  return expectedStation == AliMUONConstants::NTrackingSt();
+  
+}
+
+  //__________________________________________________________________________
+void AliMUONTrack::TagRemovableClusters() {
+  /// Identify clusters that can be removed from the track,
+  /// with the only requirement to have at least 1 cluster per requested station
+  
+  Int_t nClusters = GetNClusters();
+  AliMUONTrackParam *trackParam, *nextTrackParam;
+  Int_t currentCh, nextCh;
+  
+  // reset flags to kFALSE for all clusters in required station
+  for (Int_t i = 0; i < nClusters; i++) {
+    trackParam = (AliMUONTrackParam*) fTrackParamAtCluster->UncheckedAt(i);
+    if (AliMUONReconstructor::GetRecoParam()->RequestStation(trackParam->GetClusterPtr()->GetChamberId()/2))
+      trackParam->SetRemovable(kFALSE);
+    else trackParam->SetRemovable(kTRUE);
+  }
+  
+  // loop over track parameters
+  for (Int_t i = 0; i < nClusters; i++) {
+    trackParam = (AliMUONTrackParam*) fTrackParamAtCluster->UncheckedAt(i);
+    
+    currentCh = trackParam->GetClusterPtr()->GetChamberId();
+    
+    // loop over next track parameters
+    for (Int_t j = i+1; j < nClusters; j++) {
+      nextTrackParam = (AliMUONTrackParam*) fTrackParamAtCluster->UncheckedAt(j);
+      
+      nextCh = nextTrackParam->GetClusterPtr()->GetChamberId();
+      
+      // check if the 2 clusters are on the same station
+      if (nextCh/2 != currentCh/2) break;
+      
+      // set clusters in the same station as being removable
+      trackParam->SetRemovable(kTRUE);
+      nextTrackParam->SetRemovable(kTRUE);
+      
+    }
+    
+  }
+    
 }
 
   //__________________________________________________________________________
@@ -442,6 +530,12 @@ Bool_t AliMUONTrack::ComputeLocalChi2(Bool_t accountForMCS)
   /// - Also recompute the weight matrices of the attached clusters if accountForMCS=kTRUE
   /// - Assume that track parameters at each cluster are corrects
   /// - Return kFALSE if computation failed
+  AliDebug(1,"Enter ComputeLocalChi2");
+  
+  if (!fTrackParamAtCluster) {
+    AliWarning("no cluster attached to this track");
+    return kFALSE;
+  }
   
   if (accountForMCS) { // Compute local chi2 taking into account multiple scattering effects
       
@@ -558,6 +652,12 @@ Double_t AliMUONTrack::ComputeGlobalChi2(Bool_t accountForMCS)
   /// - Assume that track parameters at each cluster are corrects
   /// - Assume the cluster weights matrices are corrects
   /// - Return negative value if chi2 computation failed
+  AliDebug(1,"Enter ComputeGlobalChi2");
+  
+  if (!fTrackParamAtCluster) {
+    AliWarning("no cluster attached to this track");
+    return 1.e10;
+  }
   
   Double_t chi2 = 0.;
   
@@ -622,6 +722,12 @@ Bool_t AliMUONTrack::ComputeClusterWeights(TMatrixD* mcsCovariances)
   /// - Use the provided MCS covariance matrix if any (otherwise build it temporarily)
   /// - Assume that track parameters at each cluster are corrects
   /// - Return kFALSE if computation failed
+  AliDebug(1,"Enter ComputeClusterWeights1");
+  
+  if (!fTrackParamAtCluster) {
+    AliWarning("no cluster attached to this track");
+    return kFALSE;
+  }
   
   // Alocate memory
   Int_t nClusters = GetNClusters();
@@ -644,6 +750,7 @@ Bool_t AliMUONTrack::ComputeClusterWeights(TMatrixD& clusterWeightsNB, TMatrixD&
   /// accounting for multiple scattering correlations and cluster resolution
   /// - Use the provided MCS covariance matrix if any (otherwise build it temporarily)
   /// - Return kFALSE if computation failed
+  AliDebug(1,"Enter ComputeClusterWeights2");
   
   // Check MCS covariance matrix and recompute it if need
   Int_t nClusters = GetNClusters();
@@ -733,6 +840,7 @@ void AliMUONTrack::ComputeMCSCovariances(TMatrixD& mcsCovariances) const
 {
   /// Compute the multiple scattering covariance matrix
   /// (assume that track parameters at each cluster are corrects)
+  AliDebug(1,"Enter ComputeMCSCovariances");
   
   // Reset the size of the covariance matrix if needed
   Int_t nClusters = GetNClusters();
@@ -821,6 +929,7 @@ Int_t AliMUONTrack::ClustersInCommon(AliMUONTrack* track) const
 {
   /// Returns the number of clusters in common between the current track ("this")
   /// and the track pointed to by "track".
+  if (!fTrackParamAtCluster || !this->fTrackParamAtCluster) return 0;
   Int_t clustersInCommon = 0;
   AliMUONTrackParam *trackParamAtCluster1, *trackParamAtCluster2;
   // Loop over clusters of first track
@@ -848,22 +957,24 @@ Double_t AliMUONTrack::GetNormalizedChi2() const
   
   Double_t numberOfDegFree = (2. * GetNClusters() - 5.);
   if (numberOfDegFree > 0.) return fGlobalChi2 / numberOfDegFree;
-  else return 1.e10;
+  else return fGlobalChi2; // system is under-constraint
 }
 
   //__________________________________________________________________________
-Bool_t* AliMUONTrack::CompatibleTrack(AliMUONTrack *track, Double_t sigma2Cut) const
+Bool_t* AliMUONTrack::CompatibleTrack(AliMUONTrack *track, Double_t sigmaCut) const
 {
   /// for each chamber: return kTRUE (kFALSE) if clusters are compatible (not compatible)
   AliMUONTrackParam *trackParamAtCluster1, *trackParamAtCluster2;
   AliMUONVCluster *cluster1, *cluster2;
   Double_t chi2, dX, dY, dZ;
-  Double_t chi2Max = sigma2Cut * sigma2Cut;
+  Double_t chi2Max = sigmaCut * sigmaCut;
   Double_t dZMax = 1.; // 1 cm
   
   Bool_t *compatibleCluster = new Bool_t[AliMUONConstants::NTrackingCh()]; 
   for ( Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) compatibleCluster[ch] = kFALSE;
 
+  if (!fTrackParamAtCluster || !this->fTrackParamAtCluster) return compatibleCluster;
+  
   // Loop over clusters of first track
   trackParamAtCluster1 = (AliMUONTrackParam*) this->fTrackParamAtCluster->First();
   while (trackParamAtCluster1) {
@@ -903,14 +1014,6 @@ Bool_t* AliMUONTrack::CompatibleTrack(AliMUONTrack *track, Double_t sigma2Cut) c
   return compatibleCluster;
 }
 
-//__________________________________________________________________________
-AliMUONTrackParam* AliMUONTrack::GetTrackParamAtVertex()
-{
-  /// return reference to track parameters at vertex (create it before if needed)
-  if (!fTrackParamAtVertex) fTrackParamAtVertex = new AliMUONTrackParam();
-  return fTrackParamAtVertex;
-}
-
 //__________________________________________________________________________
 void AliMUONTrack::SetTrackParamAtVertex(const AliMUONTrackParam* trackParam)
 {
@@ -952,7 +1055,7 @@ void AliMUONTrack::Print(Option_t*) const
       ", LoTrgNum=" << setw(3) << GetLoTrgNum()  << 
     ", Chi2-tracking-trigger=" << setw(8) << setprecision(5) <<  GetChi2MatchTrigger();
   cout << Form(" HitTriggerPattern %x",fHitsPatternInTrigCh) << endl;
-  fTrackParamAtCluster->First()->Print("FULL");
+  if (fTrackParamAtCluster) fTrackParamAtCluster->First()->Print("FULL");
 }
 
 //__________________________________________________________________________