]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONVTrackReconstructor.cxx
A new method DrawPMDModule is added
[u/mrichter/AliRoot.git] / MUON / AliMUONVTrackReconstructor.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //-----------------------------------------------------------------------------
19 /// \class AliMUONVTrackReconstructor
20 /// Virtual MUON track reconstructor in ALICE (class renamed from AliMUONEventReconstructor)
21 ///
22 /// This class contains as data a pointer to the array of reconstructed tracks
23 ///
24 /// It contains as methods, among others:
25 /// * EventReconstruct to build the muon tracks
26 /// * EventReconstructTrigger to build the trigger tracks
27 /// * ValidateTracksWithTrigger to match tracker/trigger tracks
28 ///
29 /// Several options and adjustable parameters are available for both KALMAN and ORIGINAL
30 /// tracking algorithms. They can be changed by using:
31 /// AliMUONRecoParam *muonRecoParam = AliMUONRecoParam::GetLow(High)FluxParam();
32 /// muonRecoParam->Set...(); // see methods in AliMUONRecoParam.h for details
33 /// AliMUONReconstructor::SetRecoParam(muonRecoParam);
34 ///
35 /// Main parameters and options are:
36 /// - *fgkSigmaToCutForTracking* : quality cut used to select new clusters to be
37 ///   attached to the track candidate and to select good tracks.
38 /// - *fgkMakeTrackCandidatesFast* : if this flag is set to 'true', the track candidates
39 ///   are made assuming linear propagation between stations 4 and 5.
40 /// - *fgkTrackAllTracks* : according to the value of this flag, in case that several
41 ///   new clusters pass the quality cut, either we consider all the possibilities
42 ///   (duplicating tracks) or we attach only the best cluster.
43 /// - *fgkRecoverTracks* : if this flag is set to 'true', we try to recover the tracks
44 ///   lost during the tracking by removing the worst of the 2 clusters attached in the
45 ///   previous station.
46 /// - *fgkImproveTracks* : if this flag is set to 'true', we try to improve the quality
47 ///   of the tracks at the end of the tracking by removing clusters that do not pass
48 ///   new quality cut (the track is removed is it does not contain enough cluster anymore).
49 /// - *fgkComplementTracks* : if this flag is set to 'true', we try to improve the quality
50 ///   of the tracks at the end of the tracking by adding potentially missing clusters
51 ///   (we may have 2 clusters in the same chamber because of the overlapping of detection
52 ///   elements, which is not handle by the tracking algorithm).
53 /// - *fgkSigmaToCutForImprovement* : quality cut used when we try to improve the
54 ///   quality of the tracks.
55 ///
56 ///  \author Philippe Pillot
57 //-----------------------------------------------------------------------------
58
59 #include "AliMUONVTrackReconstructor.h"
60
61 #include "AliMUONConstants.h"
62 #include "AliMUONObjectPair.h"
63 #include "AliMUONTriggerTrack.h"
64 #include "AliMUONTriggerCircuit.h"
65 #include "AliMUONLocalTrigger.h"
66 #include "AliMUONGlobalTrigger.h"
67 #include "AliMUONTrack.h"
68 #include "AliMUONTrackParam.h"
69 #include "AliMUONTrackExtrap.h"
70 #include "AliMUONTrackHitPattern.h"
71 #include "AliMUONVTrackStore.h"
72 #include "AliMUONVClusterStore.h"
73 #include "AliMUONVCluster.h"
74 #include "AliMUONVClusterServer.h"
75 #include "AliMUONVTriggerStore.h"
76 #include "AliMUONVTriggerTrackStore.h"
77 #include "AliMpDEManager.h"
78 #include "AliMpArea.h"
79
80 #include "AliLog.h"
81 #include "AliCodeTimer.h"
82 #include "AliTracker.h"
83
84 #include <TClonesArray.h>
85 #include <TMath.h>
86 #include <TMatrixD.h>
87 #include <TVector2.h>
88
89 #include <Riostream.h>
90
91 /// \cond CLASSIMP
92 ClassImp(AliMUONVTrackReconstructor) // Class implementation in ROOT context
93 /// \endcond
94
95   //__________________________________________________________________________
96 AliMUONVTrackReconstructor::AliMUONVTrackReconstructor(AliMUONVClusterServer* clusterServer)
97   : TObject(),
98     fRecTracksPtr(0x0),
99     fNRecTracks(0),
100     fClusterServer(clusterServer)
101 {
102   /// Constructor for class AliMUONVTrackReconstructor
103   /// WARNING: if clusterServer=0x0, no clusterization will be possible at this level
104   
105   // Memory allocation for the TClonesArray of reconstructed tracks
106   fRecTracksPtr = new TClonesArray("AliMUONTrack", 100);
107   
108   // set the magnetic field for track extrapolations
109   const AliMagF* kField = AliTracker::GetFieldMap();
110   if (!kField) AliFatal("No field available");
111   AliMUONTrackExtrap::SetField(kField);
112 }
113
114   //__________________________________________________________________________
115 AliMUONVTrackReconstructor::~AliMUONVTrackReconstructor()
116 {
117   /// Destructor for class AliMUONVTrackReconstructor
118   delete fRecTracksPtr;
119 }
120
121   //__________________________________________________________________________
122 void AliMUONVTrackReconstructor::ResetTracks()
123 {
124   /// To reset the TClonesArray of reconstructed tracks
125   if (fRecTracksPtr) fRecTracksPtr->Clear("C");
126   fNRecTracks = 0;
127   return;
128 }
129
130   //__________________________________________________________________________
131 void AliMUONVTrackReconstructor::EventReconstruct(AliMUONVClusterStore& clusterStore, AliMUONVTrackStore& trackStore)
132 {
133   /// To reconstruct one event
134   AliDebug(1,"");
135   AliCodeTimerAuto("");
136   
137   // Reset array of tracks
138   ResetTracks();
139   
140   // Look for candidates from clusters in stations(1..) 4 and 5
141   MakeTrackCandidates(clusterStore);
142   
143   // Look for extra candidates from clusters in stations(1..) 4 and 5
144   if (AliMUONReconstructor::GetRecoParam()->MakeMoreTrackCandidates()) MakeMoreTrackCandidates(clusterStore);
145   
146   // Stop tracking if no candidate found
147   if (fRecTracksPtr->GetEntriesFast() == 0) return;
148   
149   // Follow tracks in stations(1..) 3, 2 and 1
150   FollowTracks(clusterStore);
151   
152   // Complement the reconstructed tracks
153   if (AliMUONReconstructor::GetRecoParam()->ComplementTracks()) ComplementTracks(clusterStore);
154   
155   // Improve the reconstructed tracks
156   if (AliMUONReconstructor::GetRecoParam()->ImproveTracks()) ImproveTracks();
157   
158   // Remove double tracks
159   RemoveDoubleTracks();
160   
161   // Fill AliMUONTrack data members
162   Finalize();
163   
164   // Add tracks to MUON data container 
165   for (Int_t i=0; i<fNRecTracks; ++i) 
166   {
167     AliMUONTrack * track = (AliMUONTrack*) fRecTracksPtr->At(i);
168     track->SetUniqueID(i+1);
169     trackStore.Add(*track);
170   }
171 }
172
173   //__________________________________________________________________________
174 TClonesArray* AliMUONVTrackReconstructor::MakeSegmentsBetweenChambers(const AliMUONVClusterStore& clusterStore, Int_t ch1, Int_t ch2)
175 {
176   /// To make the list of segments from the list of clusters in the 2 given chambers.
177   /// Return a new TClonesArray of segments.
178   /// It is the responsibility of the user to delete it afterward.
179   AliDebug(1,Form("Enter MakeSegmentsBetweenChambers (1..) %d-%d", ch1+1, ch2+1));
180   
181   AliMUONVCluster *cluster1, *cluster2;
182   AliMUONObjectPair *segment;
183   Double_t nonBendingSlope = 0, bendingSlope = 0, impactParam = 0., bendingMomentum = 0.; // to avoid compilation warning
184   
185   // Create iterators to loop over clusters in both chambers
186   TIter nextInCh1(clusterStore.CreateChamberIterator(ch1,ch1));
187   TIter nextInCh2(clusterStore.CreateChamberIterator(ch2,ch2));
188   
189   // list of segments
190   TClonesArray *segments = new TClonesArray("AliMUONObjectPair", 100);
191   
192   // Loop over clusters in the first chamber of the station
193   while ( ( cluster1 = static_cast<AliMUONVCluster*>(nextInCh1()) ) ) {
194     
195     // reset cluster iterator of chamber 2
196     nextInCh2.Reset();
197     
198     // Loop over clusters in the second chamber of the station
199     while ( ( cluster2 = static_cast<AliMUONVCluster*>(nextInCh2()) ) ) {
200       
201       // non bending slope
202       nonBendingSlope = (cluster1->GetX() - cluster2->GetX()) / (cluster1->GetZ() - cluster2->GetZ());
203       
204       // check if non bending slope is within tolerances
205       if (TMath::Abs(nonBendingSlope) > AliMUONReconstructor::GetRecoParam()->GetMaxNonBendingSlope()) continue;
206       
207       // bending slope
208       bendingSlope = (cluster1->GetY() - cluster2->GetY()) / (cluster1->GetZ() - cluster2->GetZ());
209       
210       // check the bending momentum of the bending slope depending if the field is ON or OFF
211       if (AliMUONTrackExtrap::IsFieldON()) {
212         
213         // impact parameter
214         impactParam = cluster1->GetY() - cluster1->GetZ() * bendingSlope;
215         
216         // absolute value of bending momentum
217         bendingMomentum = TMath::Abs(AliMUONTrackExtrap::GetBendingMomentumFromImpactParam(impactParam));
218         
219         // check if bending momentum is within tolerances
220         if (bendingMomentum < AliMUONReconstructor::GetRecoParam()->GetMinBendingMomentum() ||
221             bendingMomentum > AliMUONReconstructor::GetRecoParam()->GetMaxBendingMomentum()) continue;
222         
223       } else {
224         
225         // check if non bending slope is within tolerances
226         if (TMath::Abs(bendingSlope) > AliMUONReconstructor::GetRecoParam()->GetMaxBendingSlope()) continue;
227       
228       }
229       // make new segment
230       segment = new ((*segments)[segments->GetLast()+1]) AliMUONObjectPair(cluster1, cluster2, kFALSE, kFALSE);
231       
232       // Printout for debug
233       if (AliLog::GetGlobalDebugLevel() > 1) {
234         cout << "segmentIndex(0...): " << segments->GetLast() << endl;
235         segment->Dump();
236         cout << "Cluster in first chamber" << endl;
237         cluster1->Print();
238         cout << "Cluster in second chamber" << endl;
239         cluster2->Print();
240       }
241       
242     }
243     
244   }
245   
246   // Printout for debug
247   AliDebug(1,Form("chambers%d-%d: NSegments =  %d ", ch1+1, ch2+1, segments->GetEntriesFast()));
248   
249   return segments;
250 }
251
252   //__________________________________________________________________________
253 void AliMUONVTrackReconstructor::RemoveIdenticalTracks()
254 {
255   /// To remove identical tracks:
256   /// Tracks are considered identical if they have all their clusters in common.
257   /// One keeps the track with the larger number of clusters if need be
258   AliMUONTrack *track1, *track2, *trackToRemove;
259   Int_t clustersInCommon, nClusters1, nClusters2;
260   Bool_t removedTrack1;
261   // Loop over first track of the pair
262   track1 = (AliMUONTrack*) fRecTracksPtr->First();
263   while (track1) {
264     removedTrack1 = kFALSE;
265     nClusters1 = track1->GetNClusters();
266     // Loop over second track of the pair
267     track2 = (AliMUONTrack*) fRecTracksPtr->After(track1);
268     while (track2) {
269       nClusters2 = track2->GetNClusters();
270       // number of clusters in common between two tracks
271       clustersInCommon = track1->ClustersInCommon(track2);
272       // check for identical tracks
273       if ((clustersInCommon == nClusters1) || (clustersInCommon == nClusters2)) {
274         // decide which track to remove
275         if (nClusters2 > nClusters1) {
276           // remove track1 and continue the first loop with the track next to track1
277           trackToRemove = track1;
278           track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
279           fRecTracksPtr->Remove(trackToRemove);
280           fRecTracksPtr->Compress(); // this is essential to retrieve the TClonesArray afterwards
281           fNRecTracks--;
282           removedTrack1 = kTRUE;
283           break;
284         } else {
285           // remove track2 and continue the second loop with the track next to track2
286           trackToRemove = track2;
287           track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
288           fRecTracksPtr->Remove(trackToRemove);
289           fRecTracksPtr->Compress(); // this is essential to retrieve the TClonesArray afterwards
290           fNRecTracks--;
291         }
292       } else track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
293     } // track2
294     if (removedTrack1) continue;
295     track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
296   } // track1
297   return;
298 }
299
300   //__________________________________________________________________________
301 void AliMUONVTrackReconstructor::RemoveDoubleTracks()
302 {
303   /// To remove double tracks:
304   /// Tracks are considered identical if more than half of the clusters of the track
305   /// which has the smaller number of clusters are in common with the other track.
306   /// Among two identical tracks, one keeps the track with the larger number of clusters
307   /// or, if these numbers are equal, the track with the minimum chi2.
308   AliMUONTrack *track1, *track2, *trackToRemove;
309   Int_t clustersInCommon, nClusters1, nClusters2;
310   Bool_t removedTrack1;
311   // Loop over first track of the pair
312   track1 = (AliMUONTrack*) fRecTracksPtr->First();
313   while (track1) {
314     removedTrack1 = kFALSE;
315     nClusters1 = track1->GetNClusters();
316     // Loop over second track of the pair
317     track2 = (AliMUONTrack*) fRecTracksPtr->After(track1);
318     while (track2) {
319       nClusters2 = track2->GetNClusters();
320       // number of clusters in common between two tracks
321       clustersInCommon = track1->ClustersInCommon(track2);
322       // check for identical tracks
323       if (((nClusters1 < nClusters2) && (2 * clustersInCommon > nClusters1)) || (2 * clustersInCommon > nClusters2)) {
324         // decide which track to remove
325         if ((nClusters1 > nClusters2) || ((nClusters1 == nClusters2) && (track1->GetGlobalChi2() <= track2->GetGlobalChi2()))) {
326           // remove track2 and continue the second loop with the track next to track2
327           trackToRemove = track2;
328           track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
329           fRecTracksPtr->Remove(trackToRemove);
330           fRecTracksPtr->Compress(); // this is essential to retrieve the TClonesArray afterwards
331           fNRecTracks--;
332         } else {
333           // else remove track1 and continue the first loop with the track next to track1
334           trackToRemove = track1;
335           track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
336           fRecTracksPtr->Remove(trackToRemove);
337           fRecTracksPtr->Compress(); // this is essential to retrieve the TClonesArray afterwards
338           fNRecTracks--;
339           removedTrack1 = kTRUE;
340           break;
341         }
342       } else track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
343     } // track2
344     if (removedTrack1) continue;
345     track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
346   } // track1
347   return;
348 }
349
350   //__________________________________________________________________________
351 void AliMUONVTrackReconstructor::AskForNewClustersInChamber(const AliMUONTrackParam &trackParam,
352                                                             AliMUONVClusterStore& clusterStore, Int_t chamber)
353 {
354   /// Ask the clustering to reconstruct new clusters around the track candidate position
355   
356   // check if the current chamber is useable
357   if (!fClusterServer || !AliMUONReconstructor::GetRecoParam()->UseChamber(chamber)) return;
358   
359   // maximum distance between the center of the chamber and a detection element
360   // (accounting for the inclination of the chamber)
361   static const Double_t kMaxDZ = 15.; // 15 cm
362   
363   // extrapolate track parameters to the chamber
364   AliMUONTrackParam extrapTrackParam(trackParam);
365   AliMUONTrackExtrap::ExtrapToZCov(&extrapTrackParam, AliMUONConstants::DefaultChamberZ(chamber));
366   
367   // build the searching area using the track resolution and the maximum-distance-to-track value
368   const TMatrixD& kParamCov = extrapTrackParam.GetCovariances();
369   Double_t errX2 = kParamCov(0,0) + kMaxDZ * kMaxDZ * kParamCov(1,1) + 2. * kMaxDZ * TMath::Abs(kParamCov(0,1));
370   Double_t errY2 = kParamCov(2,2) + kMaxDZ * kMaxDZ * kParamCov(3,3) + 2. * kMaxDZ * TMath::Abs(kParamCov(2,3));
371   Double_t dX = TMath::Abs(trackParam.GetNonBendingSlope()) * kMaxDZ +
372                 AliMUONReconstructor::GetRecoParam()->GetMaxNonBendingDistanceToTrack() +
373                 AliMUONReconstructor::GetRecoParam()->GetSigmaCutForTracking() * TMath::Sqrt(errX2);
374   Double_t dY = TMath::Abs(trackParam.GetBendingSlope()) * kMaxDZ +
375                 AliMUONReconstructor::GetRecoParam()->GetMaxBendingDistanceToTrack() +
376                 AliMUONReconstructor::GetRecoParam()->GetSigmaCutForTracking() * TMath::Sqrt(errY2);
377   TVector2 dimensions(dX, dY);
378   TVector2 position(extrapTrackParam.GetNonBendingCoor(), extrapTrackParam.GetBendingCoor());
379   AliMpArea area(position, dimensions);
380   
381   // ask to cluterize in the given area of the given chamber
382   fClusterServer->Clusterize(chamber, clusterStore, area);
383   
384 }
385
386   //__________________________________________________________________________
387 void AliMUONVTrackReconstructor::AskForNewClustersInStation(const AliMUONTrackParam &trackParam,
388                                                             AliMUONVClusterStore& clusterStore, Int_t station)
389 {
390   /// Ask the clustering to reconstruct new clusters around the track candidate position
391   /// in the 2 chambers of the given station
392   AskForNewClustersInChamber(trackParam, clusterStore, 2*station+1);
393   AskForNewClustersInChamber(trackParam, clusterStore, 2*station);
394 }
395
396   //__________________________________________________________________________
397 Double_t AliMUONVTrackReconstructor::TryOneCluster(const AliMUONTrackParam &trackParam, AliMUONVCluster* cluster,
398                                                    AliMUONTrackParam &trackParamAtCluster, Bool_t updatePropagator)
399 {
400 /// Test the compatibility between the track and the cluster (using trackParam's covariance matrix):
401 /// return the corresponding Chi2
402 /// return trackParamAtCluster
403   
404   // extrapolate track parameters and covariances at the z position of the tested cluster
405   trackParamAtCluster = trackParam;
406   AliMUONTrackExtrap::ExtrapToZCov(&trackParamAtCluster, cluster->GetZ(), updatePropagator);
407   
408   // set pointer to cluster into trackParamAtCluster
409   trackParamAtCluster.SetClusterPtr(cluster);
410   
411   // Set differences between trackParam and cluster in the bending and non bending directions
412   Double_t dX = cluster->GetX() - trackParamAtCluster.GetNonBendingCoor();
413   Double_t dY = cluster->GetY() - trackParamAtCluster.GetBendingCoor();
414   
415   // Calculate errors and covariances
416   const TMatrixD& kParamCov = trackParamAtCluster.GetCovariances();
417   Double_t sigmaX2 = kParamCov(0,0) + cluster->GetErrX2();
418   Double_t sigmaY2 = kParamCov(2,2) + cluster->GetErrY2();
419   
420   // Compute chi2
421   return dX * dX / sigmaX2 + dY * dY / sigmaY2;
422   
423 }
424
425   //__________________________________________________________________________
426 Bool_t AliMUONVTrackReconstructor::TryOneClusterFast(const AliMUONTrackParam &trackParam, AliMUONVCluster* cluster)
427 {
428 /// Test the compatibility between the track and the cluster
429 /// given the track resolution + the maximum-distance-to-track value
430 /// and assuming linear propagation of the track:
431 /// return kTRUE if they are compatibles
432   
433   Double_t dZ = cluster->GetZ() - trackParam.GetZ();
434   Double_t dX = cluster->GetX() - (trackParam.GetNonBendingCoor() + trackParam.GetNonBendingSlope() * dZ);
435   Double_t dY = cluster->GetY() - (trackParam.GetBendingCoor() + trackParam.GetBendingSlope() * dZ);
436   const TMatrixD& kParamCov = trackParam.GetCovariances();
437   Double_t errX2 = kParamCov(0,0) + dZ * dZ * kParamCov(1,1) + 2. * dZ * kParamCov(0,1);
438   Double_t errY2 = kParamCov(2,2) + dZ * dZ * kParamCov(3,3) + 2. * dZ * kParamCov(2,3);
439
440   Double_t dXmax = AliMUONReconstructor::GetRecoParam()->GetSigmaCutForTracking() * TMath::Sqrt(errX2) +
441                    AliMUONReconstructor::GetRecoParam()->GetMaxNonBendingDistanceToTrack();
442   Double_t dYmax = AliMUONReconstructor::GetRecoParam()->GetSigmaCutForTracking() * TMath::Sqrt(errY2) +
443                    AliMUONReconstructor::GetRecoParam()->GetMaxBendingDistanceToTrack();
444   
445   if (TMath::Abs(dX) > dXmax || TMath::Abs(dY) > dYmax) return kFALSE;
446   
447   return kTRUE;
448   
449 }
450
451   //__________________________________________________________________________
452 Double_t AliMUONVTrackReconstructor::TryTwoClustersFast(const AliMUONTrackParam &trackParamAtCluster1, AliMUONVCluster* cluster2,
453                                                         AliMUONTrackParam &trackParamAtCluster2)
454 {
455 /// Test the compatibility between the track and the 2 clusters together (using trackParam's covariance matrix)
456 /// assuming linear propagation between the two clusters:
457 /// return the corresponding Chi2 accounting for covariances between the 2 clusters
458 /// return trackParamAtCluster2
459   
460   // extrapolate linearly track parameters and covariances at the z position of the second cluster
461   trackParamAtCluster2 = trackParamAtCluster1;
462   AliMUONTrackExtrap::LinearExtrapToZ(&trackParamAtCluster2, cluster2->GetZ());
463   
464   // set pointer to cluster2 into trackParamAtCluster2
465   trackParamAtCluster2.SetClusterPtr(cluster2);
466   
467   // Set differences between track and clusters in the bending and non bending directions
468   AliMUONVCluster* cluster1 = trackParamAtCluster1.GetClusterPtr();
469   Double_t dX1 = cluster1->GetX() - trackParamAtCluster1.GetNonBendingCoor();
470   Double_t dX2 = cluster2->GetX() - trackParamAtCluster2.GetNonBendingCoor();
471   Double_t dY1 = cluster1->GetY() - trackParamAtCluster1.GetBendingCoor();
472   Double_t dY2 = cluster2->GetY() - trackParamAtCluster2.GetBendingCoor();
473   
474   // Calculate errors and covariances
475   const TMatrixD& kParamCov1 = trackParamAtCluster1.GetCovariances();
476   const TMatrixD& kParamCov2 = trackParamAtCluster2.GetCovariances();
477   Double_t dZ = trackParamAtCluster2.GetZ() - trackParamAtCluster1.GetZ();
478   Double_t sigma2X1 = kParamCov1(0,0) + cluster1->GetErrX2();
479   Double_t sigma2X2 = kParamCov2(0,0) + cluster2->GetErrX2();
480   Double_t covX1X2  = kParamCov1(0,0) + dZ * kParamCov1(0,1);
481   Double_t sigma2Y1 = kParamCov1(2,2) + cluster1->GetErrY2();
482   Double_t sigma2Y2 = kParamCov2(2,2) + cluster2->GetErrY2();
483   Double_t covY1Y2  = kParamCov1(2,2) + dZ * kParamCov1(2,3);
484   
485   // Compute chi2
486   Double_t detX = sigma2X1 * sigma2X2 - covX1X2 * covX1X2;
487   Double_t detY = sigma2Y1 * sigma2Y2 - covY1Y2 * covY1Y2;
488   if (detX == 0. || detY == 0.) return 1.e10;
489   return   (dX1 * dX1 * sigma2X2 + dX2 * dX2 * sigma2X1 - 2. * dX1 * dX2 * covX1X2) / detX
490          + (dY1 * dY1 * sigma2Y2 + dY2 * dY2 * sigma2Y1 - 2. * dY1 * dY2 * covY1Y2) / detY;
491   
492 }
493
494   //__________________________________________________________________________
495 Bool_t AliMUONVTrackReconstructor::FollowLinearTrackInChamber(AliMUONTrack &trackCandidate, const AliMUONVClusterStore& clusterStore,
496                                                               Int_t nextChamber)
497 {
498   /// Follow trackCandidate in chamber(0..) nextChamber assuming linear propagation, and search for compatible cluster(s)
499   /// Keep all possibilities or only the best one(s) according to the flag fgkTrackAllTracks:
500   /// kTRUE:  duplicate "trackCandidate" if there are several possibilities and add the new tracks at the end of
501   ///         fRecTracksPtr to avoid conficts with other track candidates at this current stage of the tracking procedure.
502   ///         Remove the obsolete "trackCandidate" at the end.
503   /// kFALSE: add only the best cluster(s) to the "trackCandidate". Try to add a couple of clusters in priority.
504   /// return kTRUE if new cluster(s) have been found (otherwise return kFALSE)
505   AliDebug(1,Form("Enter FollowLinearTrackInChamber(1..) %d", nextChamber+1));
506   
507   Double_t chi2WithOneCluster = 1.e10;
508   Double_t maxChi2WithOneCluster = 2. * AliMUONReconstructor::GetRecoParam()->GetSigmaCutForTracking() *
509                                         AliMUONReconstructor::GetRecoParam()->GetSigmaCutForTracking(); // 2 because 2 quantities in chi2
510   Double_t bestChi2WithOneCluster = maxChi2WithOneCluster;
511   Bool_t foundOneCluster = kFALSE;
512   AliMUONTrack *newTrack = 0x0;
513   AliMUONVCluster *cluster;
514   AliMUONTrackParam trackParam;
515   AliMUONTrackParam extrapTrackParamAtCluster;
516   AliMUONTrackParam bestTrackParamAtCluster;
517   
518   // Get track parameters according to the propagation direction
519   if (nextChamber > 7) trackParam = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->Last();
520   else trackParam = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->First();
521   
522   // Printout for debuging
523   if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 2) || (AliLog::GetGlobalDebugLevel() >= 2)) {
524     cout<<endl<<"Track parameters and covariances at first cluster:"<<endl;
525     trackParam.GetParameters().Print();
526     trackParam.GetCovariances().Print();
527   }
528   
529   // Add MCS effect
530   AliMUONTrackExtrap::AddMCSEffect(&trackParam,AliMUONConstants::ChamberThicknessInX0(),1.);
531   
532   // Printout for debuging
533   if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
534     cout << "FollowLinearTrackInChamber: look for cluster in chamber(1..): " << nextChamber+1 << endl;
535   }
536   
537   // Create iterators to loop over clusters in chamber
538   TIter next(clusterStore.CreateChamberIterator(nextChamber,nextChamber));
539   
540   // look for candidates in chamber
541   while ( ( cluster = static_cast<AliMUONVCluster*>(next()) ) ) {
542     
543     // try to add the current cluster fast
544     if (!TryOneClusterFast(trackParam, cluster)) continue;
545     
546     // try to add the current cluster accuratly
547     extrapTrackParamAtCluster = trackParam;
548     AliMUONTrackExtrap::LinearExtrapToZ(&extrapTrackParamAtCluster, cluster->GetZ());
549     chi2WithOneCluster = TryOneCluster(extrapTrackParamAtCluster, cluster, extrapTrackParamAtCluster);
550     
551     // if good chi2 then consider to add cluster
552     if (chi2WithOneCluster < maxChi2WithOneCluster) {
553       foundOneCluster = kTRUE;
554       
555       // Printout for debuging
556       if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
557         cout << "FollowLinearTrackInChamber: found one cluster in chamber(1..): " << nextChamber+1
558         << " (Chi2 = " << chi2WithOneCluster << ")" << endl;
559         cluster->Print();
560       }
561       
562       if (AliMUONReconstructor::GetRecoParam()->TrackAllTracks()) {
563         // copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new cluster
564         newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
565         if (AliMUONReconstructor::GetRecoParam()->RequestStation(nextChamber/2))
566           extrapTrackParamAtCluster.SetRemovable(kFALSE);
567         else extrapTrackParamAtCluster.SetRemovable(kTRUE);
568         newTrack->AddTrackParamAtCluster(extrapTrackParamAtCluster,*cluster);
569         fNRecTracks++;
570         
571         // Printout for debuging
572         if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
573           cout << "FollowLinearTrackInChamber: added one cluster in chamber(1..): " << nextChamber+1 << endl;
574           if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
575         }
576         
577       } else if (chi2WithOneCluster < bestChi2WithOneCluster) {
578         // keep track of the best cluster
579         bestChi2WithOneCluster = chi2WithOneCluster;
580         bestTrackParamAtCluster = extrapTrackParamAtCluster;
581       }
582       
583     }
584     
585   }
586   
587   // fill out the best track if required else clean up the fRecTracksPtr array
588   if (!AliMUONReconstructor::GetRecoParam()->TrackAllTracks()) {
589     if (foundOneCluster) {
590       if (AliMUONReconstructor::GetRecoParam()->RequestStation(nextChamber/2))
591         bestTrackParamAtCluster.SetRemovable(kFALSE);
592       else bestTrackParamAtCluster.SetRemovable(kTRUE);
593       trackCandidate.AddTrackParamAtCluster(bestTrackParamAtCluster,*(bestTrackParamAtCluster.GetClusterPtr()));
594       
595       // Printout for debuging
596       if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
597         cout << "FollowLinearTrackInChamber: added the best cluster in chamber(1..): " << bestTrackParamAtCluster.GetClusterPtr()->GetChamberId()+1 << endl;
598         if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
599       }
600       
601     } else return kFALSE;
602     
603   } else if (foundOneCluster) {
604     
605     // remove obsolete track
606     fRecTracksPtr->Remove(&trackCandidate);
607     fNRecTracks--;
608     
609   } else return kFALSE;
610   
611   return kTRUE;
612   
613 }
614
615 //__________________________________________________________________________
616 Bool_t AliMUONVTrackReconstructor::FollowLinearTrackInStation(AliMUONTrack &trackCandidate, const AliMUONVClusterStore& clusterStore,
617                                                               Int_t nextStation)
618 {
619   /// Follow trackCandidate in station(0..) nextStation assuming linear propagation, and search for compatible cluster(s)
620   /// Keep all possibilities or only the best one(s) according to the flag fgkTrackAllTracks:
621   /// kTRUE:  duplicate "trackCandidate" if there are several possibilities and add the new tracks at the end of
622   ///         fRecTracksPtr to avoid conficts with other track candidates at this current stage of the tracking procedure.
623   ///         Remove the obsolete "trackCandidate" at the end.
624   /// kFALSE: add only the best cluster(s) to the "trackCandidate". Try to add a couple of clusters in priority.
625   /// return kTRUE if new cluster(s) have been found (otherwise return kFALSE)
626   AliDebug(1,Form("Enter FollowLinearTrackInStation(1..) %d", nextStation+1));
627   
628   // Order the chamber according to the propagation direction (tracking starts with chamber 2):
629   // - nextStation == station(1...) 5 => forward propagation
630   // - nextStation < station(1...) 5 => backward propagation
631   Int_t ch1, ch2;
632   if (nextStation==4) {
633     ch1 = 2*nextStation+1;
634     ch2 = 2*nextStation;
635   } else {
636     ch1 = 2*nextStation;
637     ch2 = 2*nextStation+1;
638   }
639   
640   Double_t chi2WithOneCluster = 1.e10;
641   Double_t chi2WithTwoClusters = 1.e10;
642   Double_t maxChi2WithOneCluster = 2. * AliMUONReconstructor::GetRecoParam()->GetSigmaCutForTracking() *
643                                         AliMUONReconstructor::GetRecoParam()->GetSigmaCutForTracking(); // 2 because 2 quantities in chi2
644   Double_t maxChi2WithTwoClusters = 4. * AliMUONReconstructor::GetRecoParam()->GetSigmaCutForTracking() *
645                                          AliMUONReconstructor::GetRecoParam()->GetSigmaCutForTracking(); // 4 because 4 quantities in chi2
646   Double_t bestChi2WithOneCluster = maxChi2WithOneCluster;
647   Double_t bestChi2WithTwoClusters = maxChi2WithTwoClusters;
648   Bool_t foundOneCluster = kFALSE;
649   Bool_t foundTwoClusters = kFALSE;
650   AliMUONTrack *newTrack = 0x0;
651   AliMUONVCluster *clusterCh1, *clusterCh2;
652   AliMUONTrackParam trackParam;
653   AliMUONTrackParam extrapTrackParamAtCluster1;
654   AliMUONTrackParam extrapTrackParamAtCluster2;
655   AliMUONTrackParam bestTrackParamAtCluster1;
656   AliMUONTrackParam bestTrackParamAtCluster2;
657   
658   Int_t nClusters = clusterStore.GetSize();
659   Bool_t *clusterCh1Used = new Bool_t[nClusters];
660   for (Int_t i = 0; i < nClusters; i++) clusterCh1Used[i] = kFALSE;
661   Int_t iCluster1;
662   
663   // Get track parameters according to the propagation direction
664   if (nextStation==4) trackParam = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->Last();
665   else trackParam = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->First();
666   
667   // Printout for debuging
668   if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 2) || (AliLog::GetGlobalDebugLevel() >= 2)) {
669     cout<<endl<<"Track parameters and covariances at first cluster:"<<endl;
670     trackParam.GetParameters().Print();
671     trackParam.GetCovariances().Print();
672   }
673   
674   // Add MCS effect
675   AliMUONTrackExtrap::AddMCSEffect(&trackParam,AliMUONConstants::ChamberThicknessInX0(),1.);
676   
677   // Printout for debuging
678   if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
679     cout << "FollowLinearTrackInStation: look for clusters in chamber(1..): " << ch2+1 << endl;
680   }
681   
682   // Create iterators to loop over clusters in both chambers
683   TIter nextInCh1(clusterStore.CreateChamberIterator(ch1,ch1));
684   TIter nextInCh2(clusterStore.CreateChamberIterator(ch2,ch2));
685   
686   // look for candidates in chamber 2
687   while ( ( clusterCh2 = static_cast<AliMUONVCluster*>(nextInCh2()) ) ) {
688     
689     // try to add the current cluster fast
690     if (!TryOneClusterFast(trackParam, clusterCh2)) continue;
691     
692     // try to add the current cluster accuratly
693     extrapTrackParamAtCluster2 = trackParam;
694     AliMUONTrackExtrap::LinearExtrapToZ(&extrapTrackParamAtCluster2, clusterCh2->GetZ());
695     chi2WithOneCluster = TryOneCluster(extrapTrackParamAtCluster2, clusterCh2, extrapTrackParamAtCluster2);
696     
697     // if good chi2 then try to attach a cluster in the other chamber too
698     if (chi2WithOneCluster < maxChi2WithOneCluster) {
699       Bool_t foundSecondCluster = kFALSE;
700       
701       // Printout for debuging
702       if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
703         cout << "FollowLinearTrackInStation: found one cluster in chamber(1..): " << ch2+1
704              << " (Chi2 = " << chi2WithOneCluster << ")" << endl;
705         clusterCh2->Print();
706         cout << "                      look for second clusters in chamber(1..): " << ch1+1 << " ..." << endl;
707       }
708       
709       // add MCS effect
710       AliMUONTrackExtrap::AddMCSEffect(&extrapTrackParamAtCluster2,AliMUONConstants::ChamberThicknessInX0(),1.);
711       
712       // reset cluster iterator of chamber 1
713       nextInCh1.Reset();
714       iCluster1 = -1;
715       
716       // look for second candidates in chamber 1
717       while ( ( clusterCh1 = static_cast<AliMUONVCluster*>(nextInCh1()) ) ) {
718         iCluster1++;
719         
720         // try to add the current cluster fast
721         if (!TryOneClusterFast(extrapTrackParamAtCluster2, clusterCh1)) continue;
722         
723         // try to add the current cluster in addition to the one found in the previous chamber
724         chi2WithTwoClusters = TryTwoClustersFast(extrapTrackParamAtCluster2, clusterCh1, extrapTrackParamAtCluster1);
725         
726         // if good chi2 then consider to add the 2 clusters to the "trackCandidate"
727         if (chi2WithTwoClusters < maxChi2WithTwoClusters) {
728           foundSecondCluster = kTRUE;
729           foundTwoClusters = kTRUE;
730           
731           // Printout for debuging
732           if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
733             cout << "FollowLinearTrackInStation: found one cluster in chamber(1..): " << ch1+1
734                  << " (Chi2 = " << chi2WithTwoClusters << ")" << endl;
735             clusterCh1->Print();
736           }
737           
738           if (AliMUONReconstructor::GetRecoParam()->TrackAllTracks()) {
739             // copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new clusters
740             newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
741             extrapTrackParamAtCluster1.SetRemovable(kTRUE);
742             newTrack->AddTrackParamAtCluster(extrapTrackParamAtCluster1,*clusterCh1);
743             extrapTrackParamAtCluster2.SetRemovable(kTRUE);
744             newTrack->AddTrackParamAtCluster(extrapTrackParamAtCluster2,*clusterCh2);
745             fNRecTracks++;
746             
747             // Tag clusterCh1 as used
748             clusterCh1Used[iCluster1] = kTRUE;
749             
750             // Printout for debuging
751             if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
752               cout << "FollowLinearTrackInStation: added two clusters in station(1..): " << nextStation+1 << endl;
753               if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
754             }
755             
756           } else if (chi2WithTwoClusters < bestChi2WithTwoClusters) {
757             // keep track of the best couple of clusters
758             bestChi2WithTwoClusters = chi2WithTwoClusters;
759             bestTrackParamAtCluster1 = extrapTrackParamAtCluster1;
760             bestTrackParamAtCluster2 = extrapTrackParamAtCluster2;
761           }
762           
763         }
764         
765       }
766       
767       // if no cluster found in chamber1 then consider to add clusterCh2 only
768       if (!foundSecondCluster) {
769         foundOneCluster = kTRUE;
770         
771         if (AliMUONReconstructor::GetRecoParam()->TrackAllTracks()) {
772           // copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new cluster
773           newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
774           if (AliMUONReconstructor::GetRecoParam()->RequestStation(nextStation))
775             extrapTrackParamAtCluster2.SetRemovable(kFALSE);
776           else extrapTrackParamAtCluster2.SetRemovable(kTRUE);
777           newTrack->AddTrackParamAtCluster(extrapTrackParamAtCluster2,*clusterCh2);
778           fNRecTracks++;
779           
780           // Printout for debuging
781           if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
782             cout << "FollowLinearTrackInStation: added one cluster in chamber(1..): " << ch2+1 << endl;
783             if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
784           }
785           
786         } else if (!foundTwoClusters && chi2WithOneCluster < bestChi2WithOneCluster) {
787           // keep track of the best cluster except if a couple of clusters has already been found
788           bestChi2WithOneCluster = chi2WithOneCluster;
789           bestTrackParamAtCluster1 = extrapTrackParamAtCluster2;
790         }
791         
792       }
793       
794     }
795     
796   }
797   
798   // look for candidates in chamber 1 not already attached to a track
799   // if we want to keep all possible tracks or if no good couple of clusters has been found
800   if (AliMUONReconstructor::GetRecoParam()->TrackAllTracks() || !foundTwoClusters) {
801     
802     // Printout for debuging
803     if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
804       cout << "FollowLinearTrackInStation: look for single cluster in chamber(1..): " << ch1+1 << endl;
805     }
806     
807     //Extrapolate trackCandidate to chamber "ch2"
808     AliMUONTrackExtrap::LinearExtrapToZ(&trackParam, AliMUONConstants::DefaultChamberZ(ch2));
809     
810     // add MCS effect for next step
811     AliMUONTrackExtrap::AddMCSEffect(&trackParam,AliMUONConstants::ChamberThicknessInX0(),1.);
812     
813     // reset cluster iterator of chamber 1
814     nextInCh1.Reset();
815     iCluster1 = -1;
816     
817     // look for second candidates in chamber 1
818     while ( ( clusterCh1 = static_cast<AliMUONVCluster*>(nextInCh1()) ) ) {
819       iCluster1++;
820       
821       if (clusterCh1Used[iCluster1]) continue; // Skip clusters already used
822       
823       // try to add the current cluster fast
824       if (!TryOneClusterFast(trackParam, clusterCh1)) continue;
825       
826       // try to add the current cluster accuratly
827       extrapTrackParamAtCluster1 = trackParam;
828       AliMUONTrackExtrap::LinearExtrapToZ(&extrapTrackParamAtCluster1, clusterCh1->GetZ());
829       chi2WithOneCluster = TryOneCluster(extrapTrackParamAtCluster1, clusterCh1, extrapTrackParamAtCluster1);
830       
831       // if good chi2 then consider to add clusterCh1
832       // We do not try to attach a cluster in the other chamber too since it has already been done above
833       if (chi2WithOneCluster < maxChi2WithOneCluster) {
834         foundOneCluster = kTRUE;
835         
836         // Printout for debuging
837         if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
838           cout << "FollowLinearTrackInStation: found one cluster in chamber(1..): " << ch1+1
839                << " (Chi2 = " << chi2WithOneCluster << ")" << endl;
840           clusterCh1->Print();
841         }
842         
843         if (AliMUONReconstructor::GetRecoParam()->TrackAllTracks()) {
844           // copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new cluster
845           newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
846           if (AliMUONReconstructor::GetRecoParam()->RequestStation(nextStation))
847             extrapTrackParamAtCluster1.SetRemovable(kFALSE);
848           else extrapTrackParamAtCluster1.SetRemovable(kTRUE);
849           newTrack->AddTrackParamAtCluster(extrapTrackParamAtCluster1,*clusterCh1);
850           fNRecTracks++;
851           
852           // Printout for debuging
853           if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
854             cout << "FollowLinearTrackInStation: added one cluster in chamber(1..): " << ch1+1 << endl;
855             if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
856           }
857           
858         } else if (chi2WithOneCluster < bestChi2WithOneCluster) {
859           // keep track of the best cluster except if a couple of clusters has already been found
860           bestChi2WithOneCluster = chi2WithOneCluster;
861           bestTrackParamAtCluster1 = extrapTrackParamAtCluster1;
862         }
863         
864       }
865       
866     }
867     
868   }
869   
870   // fill out the best track if required else clean up the fRecTracksPtr array
871   if (!AliMUONReconstructor::GetRecoParam()->TrackAllTracks()) {
872     if (foundTwoClusters) {
873       bestTrackParamAtCluster1.SetRemovable(kTRUE);
874       trackCandidate.AddTrackParamAtCluster(bestTrackParamAtCluster1,*(bestTrackParamAtCluster1.GetClusterPtr()));
875       bestTrackParamAtCluster2.SetRemovable(kTRUE);
876       trackCandidate.AddTrackParamAtCluster(bestTrackParamAtCluster2,*(bestTrackParamAtCluster2.GetClusterPtr()));
877       
878       // Printout for debuging
879       if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
880         cout << "FollowLinearTrackInStation: added the two best clusters in station(1..): " << nextStation+1 << endl;
881         if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
882       }
883       
884     } else if (foundOneCluster) {
885       if (AliMUONReconstructor::GetRecoParam()->RequestStation(nextStation))
886         bestTrackParamAtCluster1.SetRemovable(kFALSE);
887       else bestTrackParamAtCluster1.SetRemovable(kTRUE);
888       trackCandidate.AddTrackParamAtCluster(bestTrackParamAtCluster1,*(bestTrackParamAtCluster1.GetClusterPtr()));
889       
890       // Printout for debuging
891       if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
892         cout << "FollowLinearTrackInStation: added the best cluster in chamber(1..): " << bestTrackParamAtCluster1.GetClusterPtr()->GetChamberId()+1 << endl;
893         if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
894       }
895       
896     } else {
897       delete [] clusterCh1Used;
898       return kFALSE;
899     }
900     
901   } else if (foundOneCluster || foundTwoClusters) {
902     
903     // remove obsolete track
904     fRecTracksPtr->Remove(&trackCandidate);
905     fNRecTracks--;
906     
907   } else {
908     delete [] clusterCh1Used;  
909     return kFALSE;
910   }
911   
912   delete [] clusterCh1Used;
913   return kTRUE;
914   
915 }
916
917 //__________________________________________________________________________
918 void AliMUONVTrackReconstructor::ImproveTracks()
919 {
920   /// Improve tracks by removing clusters with local chi2 highter than the defined cut
921   /// Recompute track parameters and covariances at the remaining clusters
922   AliDebug(1,"Enter ImproveTracks");
923   
924   AliMUONTrack *track, *nextTrack;
925   
926   track = (AliMUONTrack*) fRecTracksPtr->First();
927   while (track) {
928     
929     // prepare next track in case the actual track is suppressed
930     nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track);
931     
932     ImproveTrack(*track);
933     
934     // remove track if improvement failed
935     if (!track->IsImproved()) {
936       fRecTracksPtr->Remove(track);
937       fNRecTracks--;
938     }
939     
940     track = nextTrack;
941   }
942   
943   // compress the array in case of some tracks have been removed
944   fRecTracksPtr->Compress();
945   
946 }
947
948 //__________________________________________________________________________
949 void AliMUONVTrackReconstructor::Finalize()
950 {
951   /// Recompute track parameters and covariances at each attached cluster from those at the first one
952   
953   AliMUONTrack *track;
954   
955   track = (AliMUONTrack*) fRecTracksPtr->First();
956   while (track) {
957     
958     FinalizeTrack(*track);
959     
960     track = (AliMUONTrack*) fRecTracksPtr->After(track);
961     
962   }
963   
964 }
965
966 //__________________________________________________________________________
967 void AliMUONVTrackReconstructor::ValidateTracksWithTrigger(AliMUONVTrackStore& trackStore,
968                                                            const AliMUONVTriggerTrackStore& triggerTrackStore,
969                                                            const AliMUONVTriggerStore& triggerStore,
970                                                            const AliMUONTrackHitPattern& trackHitPattern)
971 {
972   /// Try to match track from tracking system with trigger track
973   AliCodeTimerAuto("");
974
975   trackHitPattern.ExecuteValidation(trackStore, triggerTrackStore, triggerStore);
976 }
977
978   //__________________________________________________________________________
979 void AliMUONVTrackReconstructor::EventReconstructTrigger(const AliMUONTriggerCircuit& circuit,
980                                                          const AliMUONVTriggerStore& triggerStore,
981                                                          AliMUONVTriggerTrackStore& triggerTrackStore)
982 {
983   /// To make the trigger tracks from Local Trigger
984   AliDebug(1, "");
985   AliCodeTimerAuto("");
986   
987   AliMUONGlobalTrigger* globalTrigger = triggerStore.Global();
988   
989   UChar_t gloTrigPat = 0;
990
991   if (globalTrigger)
992   {
993     gloTrigPat = globalTrigger->GetGlobalResponse();
994   }
995   
996   TIter next(triggerStore.CreateIterator());
997   AliMUONLocalTrigger* locTrg(0x0);
998
999   Float_t z11 = AliMUONConstants::DefaultChamberZ(10);
1000   Float_t z21 = AliMUONConstants::DefaultChamberZ(12);
1001       
1002   AliMUONTriggerTrack triggerTrack;
1003   
1004   while ( ( locTrg = static_cast<AliMUONLocalTrigger*>(next()) ) )
1005   {
1006     Bool_t xTrig=locTrg->IsTrigX();
1007     Bool_t yTrig=locTrg->IsTrigY();
1008     
1009     Int_t localBoardId = locTrg->LoCircuit();
1010     
1011     if (xTrig && yTrig) 
1012     { // make Trigger Track if trigger in X and Y
1013       
1014       Float_t y11 = circuit.GetY11Pos(localBoardId, locTrg->LoStripX()); 
1015       // need first to convert deviation to [0-30] 
1016       // (see AliMUONLocalTriggerBoard::LocalTrigger)
1017       Int_t deviation = locTrg->GetDeviation(); 
1018       Int_t stripX21 = locTrg->LoStripX()+deviation+1;
1019       Float_t y21 = circuit.GetY21Pos(localBoardId, stripX21);       
1020       Float_t x11 = circuit.GetX11Pos(localBoardId, locTrg->LoStripY());
1021       
1022       AliDebug(1, Form(" MakeTriggerTrack %d %d %d %d %f %f %f \n",locTrg->LoCircuit(),
1023                        locTrg->LoStripX(),locTrg->LoStripX()+locTrg->LoDev()+1,locTrg->LoStripY(),y11, y21, x11));
1024       
1025       Float_t thetax = TMath::ATan2( x11 , z11 );
1026       Float_t thetay = TMath::ATan2( (y21-y11) , (z21-z11) );
1027       
1028       triggerTrack.SetX11(x11);
1029       triggerTrack.SetY11(y11);
1030       triggerTrack.SetThetax(thetax);
1031       triggerTrack.SetThetay(thetay);
1032       triggerTrack.SetGTPattern(gloTrigPat);
1033       triggerTrack.SetLoTrgNum(localBoardId);
1034       
1035       triggerTrackStore.Add(triggerTrack);
1036     } // board is fired 
1037   } // end of loop on Local Trigger
1038 }
1039