]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONTrackReconstructorK.cxx
Writing back ITS alignement data (Ruben).
[u/mrichter/AliRoot.git] / MUON / AliMUONTrackReconstructorK.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 AliMUONTrackReconstructorK
20 ///
21 /// MUON track reconstructor using the kalman method
22 ///
23 /// This class contains as data:
24 /// - the parameters for the track reconstruction
25 ///
26 /// It contains as methods, among others:
27 /// - MakeTracks to build the tracks
28 ///
29 //-----------------------------------------------------------------------------
30
31 #include "AliMUONTrackReconstructorK.h"
32
33 #include "AliMUONConstants.h"
34 #include "AliMUONVCluster.h"
35 #include "AliMUONVClusterServer.h"
36 #include "AliMUONVClusterStore.h"
37 #include "AliMUONTrack.h"
38 #include "AliMUONTrackParam.h"
39 #include "AliMUONTrackExtrap.h"
40 #include "AliMUONRecoParam.h"
41
42 #include "AliMpArea.h"
43
44 #include "AliLog.h"
45
46 #include <Riostream.h>
47 #include <TMath.h>
48 #include <TMatrixD.h>
49 #include <TClonesArray.h>
50
51 using std::endl;
52 using std::cout;
53 /// \cond CLASSIMP
54 ClassImp(AliMUONTrackReconstructorK) // Class implementation in ROOT context
55 /// \endcond
56
57   //__________________________________________________________________________
58 AliMUONTrackReconstructorK::AliMUONTrackReconstructorK(const AliMUONRecoParam* recoParam, AliMUONVClusterServer* clusterServer)
59   : AliMUONVTrackReconstructor(recoParam, clusterServer)
60 {
61   /// Constructor
62 }
63
64   //__________________________________________________________________________
65 AliMUONTrackReconstructorK::~AliMUONTrackReconstructorK()
66 {
67 /// Destructor
68
69
70   //__________________________________________________________________________
71 Bool_t AliMUONTrackReconstructorK::MakeTrackCandidates(AliMUONVClusterStore& clusterStore)
72 {
73   /// To make track candidates (assuming linear propagation if AliMUONRecoParam::MakeTrackCandidatesFast() return kTRUE):
74   /// Start with segments station(1..) 4 or 5 then follow track in station 5 or 4.
75   /// Good candidates are made of at least three clusters if both stations are requested (two otherwise).
76   /// Keep only best candidates or all of them according to the flag AliMUONRecoParam::TrackAllTracks().
77   
78   TClonesArray *segments;
79   AliMUONObjectPair *segment;
80   AliMUONTrack *track;
81   Int_t iCandidate = 0;
82   Bool_t clusterFound;
83
84   AliDebug(1,"Enter MakeTrackCandidates");
85
86   // Unless we're doing combined tracking, we'll clusterize all stations at once
87   Int_t firstChamber(0);
88   Int_t lastChamber(9);
89   
90   if (GetRecoParam()->CombineClusterTrackReco()) {
91     // ... Here's the exception : ask the clustering to reconstruct
92     // clusters *only* in station 4 and 5 for combined tracking
93     firstChamber = 6;
94   }
95   
96   for (Int_t i = firstChamber; i <= lastChamber; ++i ) 
97   {
98     if (fClusterServer && GetRecoParam()->UseChamber(i)) fClusterServer->Clusterize(i, clusterStore, AliMpArea(), GetRecoParam());
99   }
100   
101   // Loop over stations(1..) 5 and 4 and make track candidates
102   for (Int_t istat=4; istat>=3; istat--) {
103     
104     // Make segments in the station
105     segments = MakeSegmentsBetweenChambers(clusterStore, 2*istat, 2*istat+1);
106     
107     // Loop over segments
108     for (Int_t iSegment=0; iSegment<segments->GetEntriesFast(); iSegment++) {
109       AliDebug(1,Form("Making primary candidate(1..) %d",++iCandidate));
110       segment = (AliMUONObjectPair*) segments->UncheckedAt(iSegment);
111       
112       // Transform segments to tracks and put them at the end of fRecTracksPtr
113       track = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(segment,GetRecoParam()->GetBendingVertexDispersion());
114       fNRecTracks++;
115       
116       // Look for compatible cluster(s) in the other station
117       if (GetRecoParam()->MakeTrackCandidatesFast()) clusterFound = FollowLinearTrackInStation(*track, clusterStore, 7-istat);
118       else clusterFound = FollowTrackInStation(*track, clusterStore, 7-istat);
119       
120       // Remove track if no cluster found on a requested station
121       // or abort tracking if there are too many candidates
122       if (GetRecoParam()->RequestStation(7-istat)) {
123         if (!clusterFound) {
124           fRecTracksPtr->Remove(track);
125           fNRecTracks--;
126         } else if (fNRecTracks > GetRecoParam()->GetMaxTrackCandidates()) {
127           AliError(Form("Too many track candidates (%d tracks). Stop tracking.", fNRecTracks));
128           return kFALSE;
129         }
130       } else {
131         if ((fNRecTracks + segments->GetEntriesFast() - iSegment - 1) > GetRecoParam()->GetMaxTrackCandidates()) {
132           AliError(Form("Too many track candidates (%d tracks). Stop tracking.", fNRecTracks + segments->GetEntriesFast() - iSegment - 1));
133           return kFALSE;
134         }
135       }
136       
137     }
138     
139   }
140   
141   // Keep all different tracks if required
142   if (GetRecoParam()->TrackAllTracks()) RemoveIdenticalTracks();
143   
144   // Retrace tracks using Kalman filter and select them if needed
145   Int_t nCurrentTracks = fRecTracksPtr->GetLast()+1;
146   for (Int_t iRecTrack = 0; iRecTrack < nCurrentTracks; iRecTrack++) {
147     track = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iRecTrack);
148     
149     // skip empty slots
150     if(!track) continue;
151     
152     // retrace tracks using Kalman filter and remove the ones for which extrap failed or that are out of limits
153     if (!RetraceTrack(*track,kTRUE) || !IsAcceptable(*((AliMUONTrackParam*)track->GetTrackParamAtCluster()->First()))) {
154       fRecTracksPtr->Remove(track);
155       fNRecTracks--;
156     }
157     
158   }
159   
160   // Keep only the best tracks if required
161   if (!GetRecoParam()->TrackAllTracks()) RemoveDoubleTracks();
162   else fRecTracksPtr->Compress();
163   
164   AliDebug(1,Form("Number of good candidates = %d",fNRecTracks));
165   
166   return kTRUE;
167   
168 }
169
170   //__________________________________________________________________________
171 Bool_t AliMUONTrackReconstructorK::MakeMoreTrackCandidates(AliMUONVClusterStore& clusterStore)
172 {
173   /// To make extra track candidates assuming linear propagation:
174   /// clustering is supposed to be already done
175   /// Start with segments made of 1 cluster in each of the stations 4 and 5 then follow track in remaining chambers.
176   /// Good candidates are made of at least three clusters if both stations are requested (two otherwise).
177   /// Keep only best candidates or all of them according to the flag fgkTrackAllTracks.
178   
179   TClonesArray *segments;
180   AliMUONObjectPair *segment;
181   AliMUONTrack *track;
182   Int_t iCandidate = 0, iCurrentTrack, nCurrentTracks;
183   Int_t initialNRecTracks = fNRecTracks;
184   Bool_t clusterFound;
185   
186   AliDebug(1,"Enter MakeMoreTrackCandidates");
187   
188   // Double loop over chambers in stations(1..) 4 and 5 to make track candidates
189   for (Int_t ich1 = 6; ich1 <= 7; ich1++) {
190     for (Int_t ich2 = 8; ich2 <= 9; ich2++) {
191       
192       // Make segments between ch1 and ch2
193       segments = MakeSegmentsBetweenChambers(clusterStore, ich1, ich2);
194       
195       /// Remove segments already attached to a track
196       RemoveUsedSegments(*segments);
197       
198       // Loop over segments
199       for (Int_t iSegment=0; iSegment<segments->GetEntriesFast(); iSegment++) {
200         AliDebug(1,Form("Making primary candidate(1..) %d",++iCandidate));
201         segment = (AliMUONObjectPair*) segments->UncheckedAt(iSegment);
202         
203         // Transform segments to tracks and put them at the end of fRecTracksPtr
204         iCurrentTrack = fRecTracksPtr->GetLast()+1;
205         track = new ((*fRecTracksPtr)[iCurrentTrack]) AliMUONTrack(segment,GetRecoParam()->GetBendingVertexDispersion());
206         fNRecTracks++;
207         
208         // Look for compatible cluster(s) in the second chamber of station 5
209         clusterFound = FollowLinearTrackInChamber(*track, clusterStore, 17-ich2);
210         
211         // skip the original track in case it has been removed
212         if (GetRecoParam()->TrackAllTracks() && clusterFound) iCurrentTrack++;
213         
214         // loop over every new tracks
215         nCurrentTracks = fRecTracksPtr->GetLast()+1;
216         while (iCurrentTrack < nCurrentTracks) {
217           track = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iCurrentTrack);
218           
219           // Look for compatible cluster(s) in the second chamber of station 4
220           FollowLinearTrackInChamber(*track, clusterStore, 13-ich1);
221           
222           iCurrentTrack++;
223         }
224         
225         // abort tracking if there are too many candidates
226         if ((fNRecTracks + segments->GetEntriesFast() - iSegment - 1) > GetRecoParam()->GetMaxTrackCandidates()) {
227           AliError(Form("Too many track candidates (%d tracks). Stop tracking.", fNRecTracks + segments->GetEntriesFast() - iSegment - 1));
228           return kFALSE;
229         }
230         
231       }
232       
233     }
234   }
235   
236   // Retrace tracks using Kalman filter (also compute track chi2) and select them
237   nCurrentTracks = fRecTracksPtr->GetLast()+1;
238   for (Int_t iRecTrack = initialNRecTracks; iRecTrack < nCurrentTracks; iRecTrack++) {
239     track = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iRecTrack);
240     
241     // skip empty slots
242     if(!track) continue;
243     
244     // retrace tracks using Kalman filter and remove the ones for which extrap failed or that are out of limits
245     if (!RetraceTrack(*track,kTRUE) || !IsAcceptable(*((AliMUONTrackParam*)track->GetTrackParamAtCluster()->First()))) {
246       fRecTracksPtr->Remove(track);
247       fNRecTracks--;
248     }
249     
250   }
251   
252   // Keep only the best tracks if required
253   if (!GetRecoParam()->TrackAllTracks()) RemoveDoubleTracks();
254   else fRecTracksPtr->Compress();
255   
256   AliDebug(1,Form("Number of good candidates = %d",fNRecTracks));
257   
258   return kTRUE;
259   
260 }
261
262   //__________________________________________________________________________
263 Bool_t AliMUONTrackReconstructorK::RetraceTrack(AliMUONTrack &trackCandidate, Bool_t resetSeed)
264 {
265   /// Re-run the kalman filter from the most downstream cluster to the most uptream one
266   /// Return kFALSE in case of failure (i.e. extrapolation problem)
267   AliDebug(1,"Enter RetraceTrack");
268   
269   AliMUONTrackParam* lastTrackParam = (AliMUONTrackParam*) trackCandidate.GetTrackParamAtCluster()->Last();
270   
271   // Reset the "seed" (= track parameters and their covariances at last cluster) if required
272   if (resetSeed) {
273     
274     // parameters at last cluster
275     AliMUONVCluster* cluster2 = lastTrackParam->GetClusterPtr();
276     Double_t x2 = cluster2->GetX();
277     Double_t y2 = cluster2->GetY();
278     Double_t z2 = cluster2->GetZ();
279     
280     // parameters at last but one cluster
281     AliMUONTrackParam* previousTrackParam = (AliMUONTrackParam*) trackCandidate.GetTrackParamAtCluster()->Before(lastTrackParam);
282     AliMUONVCluster* cluster1 = previousTrackParam->GetClusterPtr();
283     // make sure it is on the previous chamber (can have 2 clusters in the same chamber after "ComplementTrack")
284     if (cluster2->GetChamberId() == cluster1->GetChamberId()) {
285       previousTrackParam = (AliMUONTrackParam*) trackCandidate.GetTrackParamAtCluster()->Before(previousTrackParam);
286       cluster1 = previousTrackParam->GetClusterPtr();
287     }
288     Double_t x1 = cluster1->GetX();
289     Double_t y1 = cluster1->GetY();
290     Double_t z1 = cluster1->GetZ();
291     
292     // reset track parameters
293     Double_t dZ = z1 - z2;
294     lastTrackParam->SetNonBendingCoor(x2);
295     lastTrackParam->SetBendingCoor(y2);
296     lastTrackParam->SetZ(z2);
297     lastTrackParam->SetNonBendingSlope((x1 - x2) / dZ);
298     lastTrackParam->SetBendingSlope((y1 - y2) / dZ);
299     Double_t bendingImpact = y2 - z2 * lastTrackParam->GetBendingSlope();
300     Double_t inverseBendingMomentum = 1. / AliMUONTrackExtrap::GetBendingMomentumFromImpactParam(bendingImpact);
301     lastTrackParam->SetInverseBendingMomentum(inverseBendingMomentum);
302     
303     // => Reset track parameter covariances at last cluster (as if the other clusters did not exist)
304     TMatrixD lastParamCov(5,5);
305     lastParamCov.Zero();
306     // Non bending plane
307     lastParamCov(0,0) = cluster2->GetErrX2();
308     lastParamCov(0,1) = - cluster2->GetErrX2() / dZ;
309     lastParamCov(1,0) = lastParamCov(0,1);
310     lastParamCov(1,1) = ( 1000. * cluster1->GetErrX2() + cluster2->GetErrX2() ) / dZ / dZ;
311     // Bending plane
312     lastParamCov(2,2) = cluster2->GetErrY2();
313     lastParamCov(2,3) = - cluster2->GetErrY2() / dZ;
314     lastParamCov(3,2) = lastParamCov(2,3);
315     lastParamCov(3,3) = ( 1000. * cluster1->GetErrY2() + cluster2->GetErrY2() ) / dZ / dZ;
316     // Inverse bending momentum (vertex resolution + bending slope resolution + 10% error on dipole parameters+field)
317     if (AliMUONTrackExtrap::IsFieldON()) {
318       lastParamCov(4,4) = ((GetRecoParam()->GetBendingVertexDispersion() *
319                             GetRecoParam()->GetBendingVertexDispersion() +
320                             (z1 * z1 * cluster2->GetErrY2() + z2 * z2 * 1000. * cluster1->GetErrY2()) / dZ / dZ) /
321                            bendingImpact / bendingImpact + 0.1 * 0.1) * inverseBendingMomentum * inverseBendingMomentum;
322       lastParamCov(2,4) = z1 * cluster2->GetErrY2() * inverseBendingMomentum / bendingImpact / dZ;
323       lastParamCov(4,2) = lastParamCov(2,4);
324       lastParamCov(3,4) = - (z1 * cluster2->GetErrY2() + z2 * 1000. * cluster1->GetErrY2()) *
325       inverseBendingMomentum / bendingImpact / dZ / dZ;
326       lastParamCov(4,3) = lastParamCov(3,4);
327     } else lastParamCov(4,4) = inverseBendingMomentum*inverseBendingMomentum;
328     lastTrackParam->SetCovariances(lastParamCov);
329     
330     // Reset the track chi2
331     lastTrackParam->SetTrackChi2(0.);
332   
333   }
334   
335   // Redo the tracking
336   return RetracePartialTrack(trackCandidate, lastTrackParam);
337   
338 }
339
340   //__________________________________________________________________________
341 Bool_t AliMUONTrackReconstructorK::RetracePartialTrack(AliMUONTrack &trackCandidate, const AliMUONTrackParam* startingTrackParam)
342 {
343   /// Re-run the kalman filter from the cluster attached to startingTrackParam to the most uptream cluster
344   /// Return kFALSE in case of failure (i.e. extrapolation problem)
345   AliDebug(1,"Enter RetracePartialTrack");
346   
347   // Printout for debuging
348   if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
349     cout << "RetracePartialTrack: track chi2 before re-tracking: " << trackCandidate.GetGlobalChi2() << endl;
350   }
351   
352   // Reset the track chi2
353   trackCandidate.SetGlobalChi2(startingTrackParam->GetTrackChi2());
354   
355   // loop over attached clusters until the first one and recompute track parameters and covariances using kalman filter
356   Bool_t extrapStatus = kTRUE;
357   Int_t expectedChamber = startingTrackParam->GetClusterPtr()->GetChamberId() - 1;
358   Int_t currentChamber;
359   Double_t addChi2TrackAtCluster;
360   AliMUONTrackParam* trackParamAtCluster = (AliMUONTrackParam*) trackCandidate.GetTrackParamAtCluster()->Before(startingTrackParam); 
361   while (trackParamAtCluster) {
362     
363     // reset track parameters and their covariances
364     trackParamAtCluster->SetParameters(startingTrackParam->GetParameters());
365     trackParamAtCluster->SetZ(startingTrackParam->GetZ());
366     trackParamAtCluster->SetCovariances(startingTrackParam->GetCovariances());
367     
368     // add MCS effect
369     AliMUONTrackExtrap::AddMCSEffect(trackParamAtCluster,AliMUONConstants::ChamberThicknessInX0(expectedChamber+1),-1.);
370     
371     // reset propagator for smoother
372     if (GetRecoParam()->UseSmoother()) trackParamAtCluster->ResetPropagator();
373     
374     // add MCS in missing chambers if any
375     currentChamber = trackParamAtCluster->GetClusterPtr()->GetChamberId();
376     while (currentChamber < expectedChamber) {
377       // extrapolation to the missing chamber (update the propagator)
378       if (!AliMUONTrackExtrap::ExtrapToZCov(trackParamAtCluster, AliMUONConstants::DefaultChamberZ(expectedChamber),
379                                             GetRecoParam()->UseSmoother())) extrapStatus = kFALSE;
380       // add MCS effect
381       AliMUONTrackExtrap::AddMCSEffect(trackParamAtCluster,AliMUONConstants::ChamberThicknessInX0(expectedChamber),-1.);
382       expectedChamber--;
383     }
384     
385     // extrapolation to the plane of the cluster attached to the current trackParamAtCluster (update the propagator)
386     if (!AliMUONTrackExtrap::ExtrapToZCov(trackParamAtCluster, trackParamAtCluster->GetClusterPtr()->GetZ(),
387                                           GetRecoParam()->UseSmoother())) extrapStatus = kFALSE;
388     
389     if (GetRecoParam()->UseSmoother()) {
390       // save extrapolated parameters for smoother
391       trackParamAtCluster->SetExtrapParameters(trackParamAtCluster->GetParameters());
392       
393       // save extrapolated covariance matrix for smoother
394       trackParamAtCluster->SetExtrapCovariances(trackParamAtCluster->GetCovariances());
395     }
396     
397     // Compute new track parameters using kalman filter
398     addChi2TrackAtCluster = RunKalmanFilter(*trackParamAtCluster);
399     
400     // Update the track chi2
401     trackCandidate.SetGlobalChi2(trackCandidate.GetGlobalChi2() + addChi2TrackAtCluster);
402     trackParamAtCluster->SetTrackChi2(trackCandidate.GetGlobalChi2());
403     
404     // prepare next step
405     expectedChamber = currentChamber - 1;
406     startingTrackParam = trackParamAtCluster;
407     trackParamAtCluster = (AliMUONTrackParam*) (trackCandidate.GetTrackParamAtCluster()->Before(startingTrackParam)); 
408   }
409   
410   // Printout for debuging
411   if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
412     cout << "RetracePartialTrack: track chi2 after re-tracking: " << trackCandidate.GetGlobalChi2() << endl;
413   }
414   
415   // set global chi2 to max value in case of problem during track extrapolation
416   if (!extrapStatus) trackCandidate.SetGlobalChi2(2.*AliMUONTrack::MaxChi2());
417   return extrapStatus;
418   
419 }
420
421   //__________________________________________________________________________
422 Bool_t AliMUONTrackReconstructorK::FollowTracks(AliMUONVClusterStore& clusterStore)
423 {
424   /// Follow tracks in stations(1..) 3, 2 and 1
425   AliDebug(1,"Enter FollowTracks");
426   
427   AliMUONTrack *track;
428   Int_t currentNRecTracks;
429   
430   for (Int_t station = 2; station >= 0; station--) {
431     
432     // Save the actual number of reconstructed track in case of
433     // tracks are added or suppressed during the tracking procedure
434     // !! Do not compress fRecTracksPtr until the end of the loop over tracks !!
435     currentNRecTracks = fNRecTracks;
436     
437     for (Int_t iRecTrack = 0; iRecTrack <currentNRecTracks; iRecTrack++) {
438       AliDebug(1,Form("FollowTracks: track candidate(1..) %d", iRecTrack+1));
439       
440       track = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iRecTrack);
441       
442       // Look for compatible cluster(s) in station(0..) "station"
443       if (!FollowTrackInStation(*track, clusterStore, station)) {
444         
445         // Try to recover track if required
446         if (GetRecoParam()->RecoverTracks()) {
447           
448           // work on a copy of the track if this station is not required
449           // to keep the case where no cluster is reconstructed as a possible candidate
450           if (!GetRecoParam()->RequestStation(station)) {
451             track = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(*track);
452             fNRecTracks++;
453           }
454           
455           // try to recover
456           if (!RecoverTrack(*track, clusterStore, station)) {
457             // remove track if no cluster found
458             fRecTracksPtr->Remove(track);
459             fNRecTracks--;
460           }
461           
462         } else if (GetRecoParam()->RequestStation(station)) {
463           // remove track if no cluster found
464           fRecTracksPtr->Remove(track);
465           fNRecTracks--;
466         } 
467         
468       }
469       
470       // abort tracking if there are too many candidates
471       if (GetRecoParam()->RequestStation(station)) {
472         if (fNRecTracks > GetRecoParam()->GetMaxTrackCandidates()) {
473           AliError(Form("Too many track candidates (%d tracks). Stop tracking.", fNRecTracks));
474           return kFALSE;
475         }
476       } else {
477         if ((fNRecTracks + currentNRecTracks - iRecTrack - 1) > GetRecoParam()->GetMaxTrackCandidates()) {
478           AliError(Form("Too many track candidates (%d tracks). Stop tracking.", fNRecTracks + currentNRecTracks - iRecTrack - 1));
479           return kFALSE;
480         }
481       }
482       
483     }
484     
485     fRecTracksPtr->Compress(); // this is essential before checking tracks
486     
487     // Keep only the best tracks if required
488     if (!GetRecoParam()->TrackAllTracks()) RemoveDoubleTracks();
489     
490   }
491   
492   return kTRUE;
493   
494 }
495
496   //__________________________________________________________________________
497 Bool_t AliMUONTrackReconstructorK::FollowTrackInChamber(AliMUONTrack &trackCandidate, AliMUONVClusterStore& clusterStore, Int_t nextChamber)
498 {
499   /// Follow trackCandidate in chamber(0..) nextChamber and search for compatible cluster(s)
500   /// Keep all possibilities or only the best one(s) according to the flag fgkTrackAllTracks:
501   /// kTRUE:  duplicate "trackCandidate" if there are several possibilities and add the new tracks at the end of
502   ///         fRecTracksPtr to avoid conficts with other track candidates at this current stage of the tracking procedure.
503   ///         Remove the obsolete "trackCandidate" at the end.
504   /// kFALSE: add only the best cluster(s) to the "trackCandidate". Try to add a couple of clusters in priority.
505   /// return kTRUE if new cluster(s) have been found (otherwise return kFALSE)
506   AliDebug(1,Form("Enter FollowTrackInChamber(1..) %d", nextChamber+1));
507   
508   Double_t chi2OfCluster;
509   Double_t maxChi2OfCluster = 2. * GetRecoParam()->GetSigmaCutForTracking() *
510                                    GetRecoParam()->GetSigmaCutForTracking(); // 2 because 2 quantities in chi2
511   Double_t addChi2TrackAtCluster;
512   Double_t bestAddChi2TrackAtCluster = AliMUONTrack::MaxChi2();
513   Bool_t foundOneCluster = kFALSE;
514   AliMUONTrack *newTrack = 0x0;
515   AliMUONVCluster *cluster;
516   AliMUONTrackParam extrapTrackParamAtCh;
517   AliMUONTrackParam extrapTrackParamAtCluster;
518   AliMUONTrackParam bestTrackParamAtCluster;
519   
520   // Get track parameters according to the propagation direction
521   if (nextChamber > 7) extrapTrackParamAtCh = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->Last();
522   else extrapTrackParamAtCh = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->First();
523   
524   // Printout for debuging
525   if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 2) || (AliLog::GetGlobalDebugLevel() >= 2)) {
526     cout<<endl<<"Track parameters and covariances at first cluster:"<<endl;
527     extrapTrackParamAtCh.GetParameters().Print();
528     extrapTrackParamAtCh.GetCovariances().Print();
529   }
530   
531   // Add MCS effect
532   Int_t currentChamber = extrapTrackParamAtCh.GetClusterPtr()->GetChamberId();
533   AliMUONTrackExtrap::AddMCSEffect(&extrapTrackParamAtCh,AliMUONConstants::ChamberThicknessInX0(currentChamber),-1.);
534   
535   // reset propagator for smoother
536   if (GetRecoParam()->UseSmoother()) extrapTrackParamAtCh.ResetPropagator();
537   
538   // Add MCS in the missing chamber(s) if any
539   while (currentChamber > nextChamber + 1) {
540     // extrapolation to the missing chamber
541     currentChamber--;
542     if (!AliMUONTrackExtrap::ExtrapToZCov(&extrapTrackParamAtCh, AliMUONConstants::DefaultChamberZ(currentChamber),
543                                           GetRecoParam()->UseSmoother())) return kFALSE;
544     // add MCS effect
545     AliMUONTrackExtrap::AddMCSEffect(&extrapTrackParamAtCh,AliMUONConstants::ChamberThicknessInX0(currentChamber),-1.);
546   }
547   
548   //Extrapolate trackCandidate to chamber
549   if (!AliMUONTrackExtrap::ExtrapToZCov(&extrapTrackParamAtCh, AliMUONConstants::DefaultChamberZ(nextChamber),
550                                         GetRecoParam()->UseSmoother())) return kFALSE;
551   
552   // Printout for debuging
553   if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 2) || (AliLog::GetGlobalDebugLevel() >= 2)) {
554     cout<<endl<<"Track parameters and covariances at first cluster extrapolated to z = "<<AliMUONConstants::DefaultChamberZ(nextChamber)<<":"<<endl;
555     extrapTrackParamAtCh.GetParameters().Print();
556     extrapTrackParamAtCh.GetCovariances().Print();
557   }
558   
559   // Printout for debuging
560   if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
561     cout << "FollowTrackInChamber: look for clusters in chamber(1..): " << nextChamber+1 << endl;
562   }
563   
564   // Ask the clustering to reconstruct new clusters around the track position in the current chamber
565   // except for station 4 and 5 that are already entirely clusterized
566   if (GetRecoParam()->CombineClusterTrackReco()) {
567     if (nextChamber < 6) AskForNewClustersInChamber(extrapTrackParamAtCh, clusterStore, nextChamber);
568   }
569   
570   // Create iterators to loop over clusters in both chambers
571   TIter next(clusterStore.CreateChamberIterator(nextChamber,nextChamber));
572   
573   // look for candidates in chamber
574   while ( ( cluster = static_cast<AliMUONVCluster*>(next()) ) ) {
575     
576     // try to add the current cluster fast
577     if (!TryOneClusterFast(extrapTrackParamAtCh, cluster)) continue;
578     
579     // try to add the current cluster accuratly
580     chi2OfCluster = TryOneCluster(extrapTrackParamAtCh, cluster, extrapTrackParamAtCluster,
581                                   GetRecoParam()->UseSmoother());
582     
583     // if good chi2 then consider to add cluster
584     if (chi2OfCluster < maxChi2OfCluster) {
585       
586       if (GetRecoParam()->UseSmoother()) {
587         // save extrapolated parameters for smoother
588         extrapTrackParamAtCluster.SetExtrapParameters(extrapTrackParamAtCluster.GetParameters());
589         
590         // save extrapolated covariance matrix for smoother
591         extrapTrackParamAtCluster.SetExtrapCovariances(extrapTrackParamAtCluster.GetCovariances());
592       }
593       
594       // Compute new track parameters including new cluster using kalman filter
595       addChi2TrackAtCluster = RunKalmanFilter(extrapTrackParamAtCluster);
596       
597       // skip track out of limits
598       if (!IsAcceptable(extrapTrackParamAtCluster)) continue;
599       
600       // remember a cluster was found
601       foundOneCluster = kTRUE;
602       
603       // Printout for debuging
604       if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
605         cout << "FollowTrackInChamber: found one cluster in chamber(1..): " << nextChamber+1
606         << " (Chi2 = " << chi2OfCluster << ")" << endl;
607         cluster->Print();
608       }
609       
610       if (GetRecoParam()->TrackAllTracks()) {
611         // copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new cluster
612         newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
613         UpdateTrack(*newTrack,extrapTrackParamAtCluster,addChi2TrackAtCluster);
614         fNRecTracks++;
615         
616         // Printout for debuging
617         if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
618           cout << "FollowTrackInChamber: added one cluster in chamber(1..): " << nextChamber+1 << endl;
619           if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
620         }
621         
622       } else if (addChi2TrackAtCluster < bestAddChi2TrackAtCluster) {
623         // keep track of the best cluster
624         bestAddChi2TrackAtCluster = addChi2TrackAtCluster;
625         bestTrackParamAtCluster = extrapTrackParamAtCluster;
626       }
627       
628     }
629     
630   }
631   
632   // fill out the best track if required else clean up the fRecTracksPtr array
633   if (!GetRecoParam()->TrackAllTracks()) {
634     if (foundOneCluster) {
635       UpdateTrack(trackCandidate,bestTrackParamAtCluster,bestAddChi2TrackAtCluster);
636       
637       // Printout for debuging
638       if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
639         cout << "FollowTrackInChamber: added the best cluster in chamber(1..): " << bestTrackParamAtCluster.GetClusterPtr()->GetChamberId()+1 << endl;
640         if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
641       }
642       
643     } else return kFALSE;
644     
645   } else if (foundOneCluster) {
646     
647     // remove obsolete track
648     fRecTracksPtr->Remove(&trackCandidate);
649     fNRecTracks--;
650     
651   } else return kFALSE;
652   
653   return kTRUE;
654   
655 }
656
657   //__________________________________________________________________________
658 Bool_t AliMUONTrackReconstructorK::FollowTrackInStation(AliMUONTrack &trackCandidate, AliMUONVClusterStore& clusterStore, Int_t nextStation)
659 {
660   /// Follow trackCandidate in station(0..) nextStation and search for compatible cluster(s)
661   /// Keep all possibilities or only the best one(s) according to the flag fgkTrackAllTracks:
662   /// kTRUE:  duplicate "trackCandidate" if there are several possibilities and add the new tracks at the end of
663   ///         fRecTracksPtr to avoid conficts with other track candidates at this current stage of the tracking procedure.
664   ///         Remove the obsolete "trackCandidate" at the end.
665   /// kFALSE: add only the best cluster(s) to the "trackCandidate". Try to add a couple of clusters in priority.
666   /// return kTRUE if new cluster(s) have been found (otherwise return kFALSE)
667   AliDebug(1,Form("Enter FollowTrackInStation(1..) %d", nextStation+1));
668   
669   // Order the chamber according to the propagation direction (tracking starts with chamber 2):
670   // - nextStation == station(1...) 5 => forward propagation
671   // - nextStation < station(1...) 5 => backward propagation
672   Int_t ch1, ch2;
673   if (nextStation==4) {
674     ch1 = 2*nextStation+1;
675     ch2 = 2*nextStation;
676   } else {
677     ch1 = 2*nextStation;
678     ch2 = 2*nextStation+1;
679   }
680   
681   Double_t chi2OfCluster;
682   Double_t maxChi2OfCluster = 2. * GetRecoParam()->GetSigmaCutForTracking() *
683                                    GetRecoParam()->GetSigmaCutForTracking(); // 2 because 2 quantities in chi2
684   Double_t addChi2TrackAtCluster1;
685   Double_t addChi2TrackAtCluster2;
686   Double_t bestAddChi2TrackAtCluster1 = AliMUONTrack::MaxChi2();
687   Double_t bestAddChi2TrackAtCluster2 = AliMUONTrack::MaxChi2();
688   Bool_t foundOneCluster = kFALSE;
689   Bool_t foundTwoClusters = kFALSE;
690   AliMUONTrack *newTrack = 0x0;
691   AliMUONVCluster *clusterCh1, *clusterCh2;
692   AliMUONTrackParam extrapTrackParam;
693   AliMUONTrackParam extrapTrackParamAtCh;
694   AliMUONTrackParam extrapTrackParamAtCluster1;
695   AliMUONTrackParam extrapTrackParamAtCluster2;
696   AliMUONTrackParam bestTrackParamAtCluster1;
697   AliMUONTrackParam bestTrackParamAtCluster2;
698   
699   // Get track parameters according to the propagation direction
700   if (nextStation==4) extrapTrackParamAtCh = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->Last();
701   else extrapTrackParamAtCh = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->First();
702   
703   // Printout for debuging
704   if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 2) || (AliLog::GetGlobalDebugLevel() >= 2)) {
705     cout<<endl<<"Track parameters and covariances at first cluster:"<<endl;
706     extrapTrackParamAtCh.GetParameters().Print();
707     extrapTrackParamAtCh.GetCovariances().Print();
708   }
709   
710   // Add MCS effect
711   Int_t currentChamber = extrapTrackParamAtCh.GetClusterPtr()->GetChamberId();
712   AliMUONTrackExtrap::AddMCSEffect(&extrapTrackParamAtCh,AliMUONConstants::ChamberThicknessInX0(currentChamber),-1.);
713   
714   // reset propagator for smoother
715   if (GetRecoParam()->UseSmoother()) extrapTrackParamAtCh.ResetPropagator();
716   
717   // Add MCS in the missing chamber(s) if any
718   while (ch1 < ch2 && currentChamber > ch2 + 1) {
719     // extrapolation to the missing chamber
720     currentChamber--;
721     if (!AliMUONTrackExtrap::ExtrapToZCov(&extrapTrackParamAtCh, AliMUONConstants::DefaultChamberZ(currentChamber),
722                                           GetRecoParam()->UseSmoother())) return kFALSE;
723     // add MCS effect
724     AliMUONTrackExtrap::AddMCSEffect(&extrapTrackParamAtCh,AliMUONConstants::ChamberThicknessInX0(currentChamber),-1.);
725   }
726   
727   //Extrapolate trackCandidate to chamber "ch2"
728   if (!AliMUONTrackExtrap::ExtrapToZCov(&extrapTrackParamAtCh, AliMUONConstants::DefaultChamberZ(ch2),
729                                         GetRecoParam()->UseSmoother())) return kFALSE;
730   
731   // Printout for debuging
732   if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 2) || (AliLog::GetGlobalDebugLevel() >= 2)) {
733     cout<<endl<<"Track parameters and covariances at first cluster extrapolated to z = "<<AliMUONConstants::DefaultChamberZ(ch2)<<":"<<endl;
734     extrapTrackParamAtCh.GetParameters().Print();
735     extrapTrackParamAtCh.GetCovariances().Print();
736   }
737   
738   // Printout for debuging
739   if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
740     cout << "FollowTrackInStation: look for clusters in chamber(1..): " << ch2+1 << endl;
741   }
742   
743   // Ask the clustering to reconstruct new clusters around the track position in the current station
744   // except for station 4 and 5 that are already entirely clusterized
745   if (GetRecoParam()->CombineClusterTrackReco()) {
746     if (nextStation < 3) AskForNewClustersInStation(extrapTrackParamAtCh, clusterStore, nextStation);
747   }
748   
749   Int_t nClusters = clusterStore.GetSize();
750   Bool_t *clusterCh1Used = new Bool_t[nClusters];
751   for (Int_t i = 0; i < nClusters; i++) clusterCh1Used[i] = kFALSE;
752   Int_t iCluster1;
753   
754   // Create iterators to loop over clusters in both chambers
755   TIter nextInCh1(clusterStore.CreateChamberIterator(ch1,ch1));
756   TIter nextInCh2(clusterStore.CreateChamberIterator(ch2,ch2));
757   
758   // look for candidates in chamber 2
759   while ( ( clusterCh2 = static_cast<AliMUONVCluster*>(nextInCh2()) ) ) {
760     
761     // try to add the current cluster fast
762     if (!TryOneClusterFast(extrapTrackParamAtCh, clusterCh2)) continue;
763     
764     // try to add the current cluster accuratly
765     chi2OfCluster = TryOneCluster(extrapTrackParamAtCh, clusterCh2, extrapTrackParamAtCluster2,
766                                   GetRecoParam()->UseSmoother());
767     
768     // if good chi2 then try to attach a cluster in the other chamber too
769     if (chi2OfCluster < maxChi2OfCluster) {
770       
771       if (GetRecoParam()->UseSmoother()) {
772         // save extrapolated parameters for smoother
773         extrapTrackParamAtCluster2.SetExtrapParameters(extrapTrackParamAtCluster2.GetParameters());
774         
775         // save extrapolated covariance matrix for smoother
776         extrapTrackParamAtCluster2.SetExtrapCovariances(extrapTrackParamAtCluster2.GetCovariances());
777       }
778       
779       // Compute new track parameters including "clusterCh2" using kalman filter
780       addChi2TrackAtCluster2 = RunKalmanFilter(extrapTrackParamAtCluster2);
781       
782       // skip track out of limits
783       if (!IsAcceptable(extrapTrackParamAtCluster2)) continue;
784       
785       // Printout for debuging
786       if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
787         cout << "FollowTrackInStation: found one cluster in chamber(1..): " << ch2+1
788         << " (Chi2 = " << chi2OfCluster << ")" << endl;
789         clusterCh2->Print();
790         cout << "                      look for second clusters in chamber(1..): " << ch1+1 << " ..." << endl;
791       }
792       
793       // copy new track parameters for next step
794       extrapTrackParam = extrapTrackParamAtCluster2;
795       
796       // add MCS effect
797       AliMUONTrackExtrap::AddMCSEffect(&extrapTrackParam,AliMUONConstants::ChamberThicknessInX0(ch2),-1.);
798       
799       // reset propagator for smoother
800       if (GetRecoParam()->UseSmoother()) extrapTrackParam.ResetPropagator();
801       
802       //Extrapolate track parameters to chamber "ch1"
803       Bool_t normalExtrap = AliMUONTrackExtrap::ExtrapToZCov(&extrapTrackParam, AliMUONConstants::DefaultChamberZ(ch1),
804                                                              GetRecoParam()->UseSmoother());
805       
806       // reset cluster iterator of chamber 1
807       nextInCh1.Reset();
808       iCluster1 = -1;
809       
810       // look for second candidates in chamber 1
811       Bool_t foundSecondCluster = kFALSE;
812       if (normalExtrap) while ( ( clusterCh1 = static_cast<AliMUONVCluster*>(nextInCh1()) ) ) {
813         iCluster1++;
814         
815         // try to add the current cluster fast
816         if (!TryOneClusterFast(extrapTrackParam, clusterCh1)) continue;
817         
818         // try to add the current cluster accuratly
819         chi2OfCluster = TryOneCluster(extrapTrackParam, clusterCh1, extrapTrackParamAtCluster1,
820                                       GetRecoParam()->UseSmoother());
821         
822         // if good chi2 then consider to add the 2 clusters to the "trackCandidate"
823         if (chi2OfCluster < maxChi2OfCluster) {
824           
825           if (GetRecoParam()->UseSmoother()) {
826             // save extrapolated parameters for smoother
827             extrapTrackParamAtCluster1.SetExtrapParameters(extrapTrackParamAtCluster1.GetParameters());
828             
829             // save extrapolated covariance matrix for smoother
830             extrapTrackParamAtCluster1.SetExtrapCovariances(extrapTrackParamAtCluster1.GetCovariances());
831           }
832           
833           // Compute new track parameters including "clusterCh1" using kalman filter
834           addChi2TrackAtCluster1 = RunKalmanFilter(extrapTrackParamAtCluster1);
835           
836           // skip track out of limits
837           if (!IsAcceptable(extrapTrackParamAtCluster1)) continue;
838           
839           // remember a second cluster was found
840           foundSecondCluster = kTRUE;
841           foundTwoClusters = kTRUE;
842           
843           // Printout for debuging
844           if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
845             cout << "FollowTrackInStation: found one cluster in chamber(1..): " << ch1+1
846             << " (Chi2 = " << chi2OfCluster << ")" << endl;
847             clusterCh1->Print();
848           }
849           
850           if (GetRecoParam()->TrackAllTracks()) {
851             // copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new clusters
852             newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
853             UpdateTrack(*newTrack,extrapTrackParamAtCluster1,extrapTrackParamAtCluster2,addChi2TrackAtCluster1,addChi2TrackAtCluster2);
854             fNRecTracks++;
855             
856             // Tag clusterCh1 as used
857             clusterCh1Used[iCluster1] = kTRUE;
858             
859             // Printout for debuging
860             if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
861               cout << "FollowTrackInStation: added two clusters in station(1..): " << nextStation+1 << endl;
862               if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
863             }
864             
865           } else if (addChi2TrackAtCluster1+addChi2TrackAtCluster2 < bestAddChi2TrackAtCluster1+bestAddChi2TrackAtCluster2) {
866             // keep track of the best couple of clusters
867             bestAddChi2TrackAtCluster1 = addChi2TrackAtCluster1;
868             bestAddChi2TrackAtCluster2 = addChi2TrackAtCluster2;
869             bestTrackParamAtCluster1 = extrapTrackParamAtCluster1;
870             bestTrackParamAtCluster2 = extrapTrackParamAtCluster2;
871           }
872           
873         }
874         
875       }
876       
877       // if no clusterCh1 found then consider to add clusterCh2 only
878       if (!foundSecondCluster) {
879         
880         // remember a cluster was found
881         foundOneCluster = kTRUE;
882         
883         if (GetRecoParam()->TrackAllTracks()) {
884           // copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new cluster
885           newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
886           UpdateTrack(*newTrack,extrapTrackParamAtCluster2,addChi2TrackAtCluster2);
887           fNRecTracks++;
888           
889           // Printout for debuging
890           if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
891             cout << "FollowTrackInStation: added one cluster in chamber(1..): " << ch2+1 << endl;
892             if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
893           }
894           
895         } else if (!foundTwoClusters && addChi2TrackAtCluster2 < bestAddChi2TrackAtCluster1) {
896           // keep track of the best single cluster except if a couple of clusters has already been found
897           bestAddChi2TrackAtCluster1 = addChi2TrackAtCluster2;
898           bestTrackParamAtCluster1 = extrapTrackParamAtCluster2;
899         }
900         
901       }
902       
903     }
904     
905   }
906   
907   // look for candidates in chamber 1 not already attached to a track
908   // if we want to keep all possible tracks or if no good couple of clusters has been found
909   if (GetRecoParam()->TrackAllTracks() || !foundTwoClusters) {
910     
911     // add MCS effect for next step
912     AliMUONTrackExtrap::AddMCSEffect(&extrapTrackParamAtCh,AliMUONConstants::ChamberThicknessInX0(ch2),-1.);
913     
914     //Extrapolate trackCandidate to chamber "ch1"
915     Bool_t normalExtrap = AliMUONTrackExtrap::ExtrapToZCov(&extrapTrackParamAtCh, AliMUONConstants::DefaultChamberZ(ch1),
916                                                            GetRecoParam()->UseSmoother());
917     
918     // Printout for debuging
919     if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 2) || (AliLog::GetGlobalDebugLevel() >= 2)) {
920       cout<<endl<<"Track parameters and covariances at first cluster extrapolated to z = "<<AliMUONConstants::DefaultChamberZ(ch1)<<":"<<endl;
921       extrapTrackParamAtCh.GetParameters().Print();
922       extrapTrackParamAtCh.GetCovariances().Print();
923     }
924     
925     // Printout for debuging
926     if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
927       cout << "FollowTrackInStation: look for single clusters in chamber(1..): " << ch1+1 << endl;
928     }
929     
930     // reset cluster iterator of chamber 1
931     nextInCh1.Reset();
932     iCluster1 = -1;
933     
934     // look for second candidates in chamber 1
935     if (normalExtrap) while ( ( clusterCh1 = static_cast<AliMUONVCluster*>(nextInCh1()) ) ) {
936       iCluster1++;
937       
938       if (clusterCh1Used[iCluster1]) continue; // Skip clusters already used
939       
940       // try to add the current cluster fast
941       if (!TryOneClusterFast(extrapTrackParamAtCh, clusterCh1)) continue;
942       
943       // try to add the current cluster accuratly
944       chi2OfCluster = TryOneCluster(extrapTrackParamAtCh, clusterCh1, extrapTrackParamAtCluster1,
945                                     GetRecoParam()->UseSmoother());
946       
947       // if good chi2 then consider to add clusterCh1
948       // We do not try to attach a cluster in the other chamber too since it has already been done above
949       if (chi2OfCluster < maxChi2OfCluster) {
950         
951         if (GetRecoParam()->UseSmoother()) {
952           // save extrapolated parameters for smoother
953           extrapTrackParamAtCluster1.SetExtrapParameters(extrapTrackParamAtCluster1.GetParameters());
954           
955           // save extrapolated covariance matrix for smoother
956           extrapTrackParamAtCluster1.SetExtrapCovariances(extrapTrackParamAtCluster1.GetCovariances());
957         }
958         
959         // Compute new track parameters including "clusterCh1" using kalman filter
960         addChi2TrackAtCluster1 = RunKalmanFilter(extrapTrackParamAtCluster1);
961         
962         // skip track out of limits
963         if (!IsAcceptable(extrapTrackParamAtCluster1)) continue;
964         
965         // remember a cluster was found
966         foundOneCluster = kTRUE;
967         
968         // Printout for debuging
969         if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
970           cout << "FollowTrackInStation: found one cluster in chamber(1..): " << ch1+1
971           << " (Chi2 = " << chi2OfCluster << ")" << endl;
972           clusterCh1->Print();
973         }
974         
975         if (GetRecoParam()->TrackAllTracks()) {
976           // copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new cluster
977           newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
978           UpdateTrack(*newTrack,extrapTrackParamAtCluster1,addChi2TrackAtCluster1);
979           fNRecTracks++;
980           
981           // Printout for debuging
982           if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
983             cout << "FollowTrackInStation: added one cluster in chamber(1..): " << ch1+1 << endl;
984             if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
985           }
986           
987         } else if (addChi2TrackAtCluster1 < bestAddChi2TrackAtCluster1) {
988           // keep track of the best single cluster except if a couple of clusters has already been found
989           bestAddChi2TrackAtCluster1 = addChi2TrackAtCluster1;
990           bestTrackParamAtCluster1 = extrapTrackParamAtCluster1;
991         }
992         
993       }
994       
995     }
996     
997   }
998   
999   // fill out the best track if required else clean up the fRecTracksPtr array
1000   if (!GetRecoParam()->TrackAllTracks()) {
1001     if (foundTwoClusters) {
1002       UpdateTrack(trackCandidate,bestTrackParamAtCluster1,bestTrackParamAtCluster2,bestAddChi2TrackAtCluster1,bestAddChi2TrackAtCluster2);
1003       
1004       // Printout for debuging
1005       if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
1006         cout << "FollowTrackInStation: added the two best clusters in station(1..): " << nextStation+1 << endl;
1007         if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
1008       }
1009       
1010     } else if (foundOneCluster) {
1011       UpdateTrack(trackCandidate,bestTrackParamAtCluster1,bestAddChi2TrackAtCluster1);
1012       
1013       // Printout for debuging
1014       if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
1015         cout << "FollowTrackInStation: added the best cluster in chamber(1..): " << bestTrackParamAtCluster1.GetClusterPtr()->GetChamberId()+1 << endl;
1016         if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
1017       }
1018       
1019     } else {
1020       delete [] clusterCh1Used;
1021       return kFALSE;
1022     }
1023     
1024   } else if (foundOneCluster || foundTwoClusters) {
1025     
1026     // remove obsolete track
1027     fRecTracksPtr->Remove(&trackCandidate);
1028     fNRecTracks--;
1029     
1030   } else {
1031     delete [] clusterCh1Used;
1032     return kFALSE;
1033   }
1034   
1035   delete [] clusterCh1Used;
1036   return kTRUE;
1037   
1038 }
1039
1040   //__________________________________________________________________________
1041 Double_t AliMUONTrackReconstructorK::RunKalmanFilter(AliMUONTrackParam &trackParamAtCluster)
1042 {
1043   /// Compute new track parameters and their covariances including new cluster using kalman filter
1044   /// return the additional track chi2
1045   AliDebug(1,"Enter RunKalmanFilter");
1046   
1047   // Get actual track parameters (p)
1048   TMatrixD param(trackParamAtCluster.GetParameters());
1049   
1050   // Get new cluster parameters (m)
1051   AliMUONVCluster *cluster = trackParamAtCluster.GetClusterPtr();
1052   TMatrixD clusterParam(5,1);
1053   clusterParam.Zero();
1054   clusterParam(0,0) = cluster->GetX();
1055   clusterParam(2,0) = cluster->GetY();
1056   
1057   // Compute the actual parameter weight (W)
1058   TMatrixD paramWeight(trackParamAtCluster.GetCovariances());
1059   if (paramWeight.Determinant() != 0) {
1060     paramWeight.Invert();
1061   } else {
1062     AliWarning(" Determinant = 0");
1063     return 2.*AliMUONTrack::MaxChi2();
1064   }
1065   
1066   // Compute the new cluster weight (U)
1067   TMatrixD clusterWeight(5,5);
1068   clusterWeight.Zero();
1069   clusterWeight(0,0) = 1. / cluster->GetErrX2();
1070   clusterWeight(2,2) = 1. / cluster->GetErrY2();
1071
1072   // Compute the new parameters covariance matrix ( (W+U)^-1 )
1073   TMatrixD newParamCov(paramWeight,TMatrixD::kPlus,clusterWeight);
1074   if (newParamCov.Determinant() != 0) {
1075     newParamCov.Invert();
1076   } else {
1077     AliWarning(" Determinant = 0");
1078     return 2.*AliMUONTrack::MaxChi2();
1079   }
1080   
1081   // Save the new parameters covariance matrix
1082   trackParamAtCluster.SetCovariances(newParamCov);
1083   
1084   // Compute the new parameters (p' = ((W+U)^-1)U(m-p) + p)
1085   TMatrixD tmp(clusterParam,TMatrixD::kMinus,param);
1086   TMatrixD tmp2(clusterWeight,TMatrixD::kMult,tmp); // U(m-p)
1087   TMatrixD newParam(newParamCov,TMatrixD::kMult,tmp2); // ((W+U)^-1)U(m-p)
1088   newParam += param; // ((W+U)^-1)U(m-p) + p
1089   
1090   // Save the new parameters
1091   trackParamAtCluster.SetParameters(newParam);
1092   
1093   // Compute the additional chi2 (= ((p'-p)^-1)W(p'-p) + ((p'-m)^-1)U(p'-m))
1094   tmp = newParam; // p'
1095   tmp -= param; // (p'-p)
1096   TMatrixD tmp3(paramWeight,TMatrixD::kMult,tmp); // W(p'-p)
1097   TMatrixD addChi2Track(tmp,TMatrixD::kTransposeMult,tmp3); // ((p'-p)^-1)W(p'-p)
1098   tmp = newParam; // p'
1099   tmp -= clusterParam; // (p'-m)
1100   TMatrixD tmp4(clusterWeight,TMatrixD::kMult,tmp); // U(p'-m)
1101   addChi2Track += TMatrixD(tmp,TMatrixD::kTransposeMult,tmp4); // ((p'-p)^-1)W(p'-p) + ((p'-m)^-1)U(p'-m)
1102   
1103   return addChi2Track(0,0);
1104   
1105 }
1106
1107   //__________________________________________________________________________
1108 void AliMUONTrackReconstructorK::UpdateTrack(AliMUONTrack &track, AliMUONTrackParam &trackParamAtCluster, Double_t addChi2)
1109 {
1110   /// Add 1 cluster to the track candidate
1111   /// Update chi2 of the track 
1112   
1113   // Flag cluster as being (not) removable
1114   if (GetRecoParam()->RequestStation(trackParamAtCluster.GetClusterPtr()->GetChamberId()/2))
1115     trackParamAtCluster.SetRemovable(kFALSE);
1116   else trackParamAtCluster.SetRemovable(kTRUE);
1117   trackParamAtCluster.SetLocalChi2(0.); // --> Local chi2 not used
1118   
1119   // Update the track chi2 into trackParamAtCluster
1120   trackParamAtCluster.SetTrackChi2(track.GetGlobalChi2() + addChi2);
1121   
1122   // Update the chi2 of the new track
1123   track.SetGlobalChi2(trackParamAtCluster.GetTrackChi2());
1124   
1125   // Update array of TrackParamAtCluster
1126   track.AddTrackParamAtCluster(trackParamAtCluster,*(trackParamAtCluster.GetClusterPtr()));
1127   
1128 }
1129
1130   //__________________________________________________________________________
1131 void AliMUONTrackReconstructorK::UpdateTrack(AliMUONTrack &track, AliMUONTrackParam &trackParamAtCluster1, AliMUONTrackParam &trackParamAtCluster2,
1132                                              Double_t addChi2AtCluster1, Double_t addChi2AtCluster2)
1133 {
1134   /// Add 2 clusters to the track candidate (order is important)
1135   /// Update track and local chi2
1136   
1137   // Update local chi2 at first cluster
1138   AliMUONVCluster* cluster1 = trackParamAtCluster1.GetClusterPtr();
1139   Double_t deltaX = trackParamAtCluster1.GetNonBendingCoor() - cluster1->GetX();
1140   Double_t deltaY = trackParamAtCluster1.GetBendingCoor() - cluster1->GetY();
1141   Double_t localChi2AtCluster1 = deltaX*deltaX / cluster1->GetErrX2() +
1142                                  deltaY*deltaY / cluster1->GetErrY2();
1143   trackParamAtCluster1.SetLocalChi2(localChi2AtCluster1);
1144   
1145   // Flag first cluster as being removable
1146   trackParamAtCluster1.SetRemovable(kTRUE);
1147   
1148   // Update local chi2 at second cluster
1149   AliMUONVCluster* cluster2 = trackParamAtCluster2.GetClusterPtr();
1150   AliMUONTrackParam extrapTrackParamAtCluster2(trackParamAtCluster1);
1151   AliMUONTrackExtrap::ExtrapToZ(&extrapTrackParamAtCluster2, trackParamAtCluster2.GetZ());
1152   deltaX = extrapTrackParamAtCluster2.GetNonBendingCoor() - cluster2->GetX();
1153   deltaY = extrapTrackParamAtCluster2.GetBendingCoor() - cluster2->GetY();
1154   Double_t localChi2AtCluster2 = deltaX*deltaX / cluster2->GetErrX2() +
1155                                  deltaY*deltaY / cluster2->GetErrY2();
1156   trackParamAtCluster2.SetLocalChi2(localChi2AtCluster2);
1157   
1158   // Flag second cluster as being removable
1159   trackParamAtCluster2.SetRemovable(kTRUE);
1160   
1161   // Update the track chi2 into trackParamAtCluster2
1162   trackParamAtCluster2.SetTrackChi2(track.GetGlobalChi2() + addChi2AtCluster2);
1163   
1164   // Update the track chi2 into trackParamAtCluster1
1165   trackParamAtCluster1.SetTrackChi2(trackParamAtCluster2.GetTrackChi2() + addChi2AtCluster1);
1166   
1167   // Update the chi2 of the new track
1168   track.SetGlobalChi2(trackParamAtCluster1.GetTrackChi2());
1169   
1170   // Update array of trackParamAtCluster
1171   track.AddTrackParamAtCluster(trackParamAtCluster1,*cluster1);
1172   track.AddTrackParamAtCluster(trackParamAtCluster2,*cluster2);
1173   
1174 }
1175
1176   //__________________________________________________________________________
1177 Bool_t AliMUONTrackReconstructorK::RecoverTrack(AliMUONTrack &trackCandidate, AliMUONVClusterStore& clusterStore, Int_t nextStation)
1178 {
1179   /// Try to recover the track candidate in the next station
1180   /// by removing the worst of the two clusters attached in the current station
1181   /// Return kTRUE if recovering succeeds
1182   AliDebug(1,"Enter RecoverTrack");
1183   
1184   // Do not try to recover track until we have attached cluster(s) on station(1..) 3
1185   if (nextStation > 1) return kFALSE;
1186   
1187   Int_t worstClusterNumber = -1;
1188   Double_t localChi2, worstLocalChi2 = -1.;
1189   
1190   // Look for the cluster to remove
1191   for (Int_t clusterNumber = 0; clusterNumber < 2; clusterNumber++) {
1192     AliMUONTrackParam *trackParamAtCluster = (AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->UncheckedAt(clusterNumber);
1193     
1194     // check if current cluster is in the previous station
1195     if (trackParamAtCluster->GetClusterPtr()->GetChamberId()/2 != nextStation+1) break;
1196     
1197     // check if current cluster is removable
1198     if (!trackParamAtCluster->IsRemovable()) return kFALSE;
1199     
1200     // reset the current cluster as being not removable if it is on a required station
1201     if (GetRecoParam()->RequestStation(nextStation+1)) trackParamAtCluster->SetRemovable(kFALSE);
1202     
1203     // Pick up cluster with the worst chi2
1204     localChi2 = trackParamAtCluster->GetLocalChi2();
1205     if (localChi2 > worstLocalChi2) {
1206       worstLocalChi2 = localChi2;
1207       worstClusterNumber = clusterNumber;
1208     }
1209   }
1210   
1211   // check if worst cluster found
1212   if (worstClusterNumber < 0) return kFALSE;
1213   
1214   // Remove the worst cluster
1215   trackCandidate.RemoveTrackParamAtCluster((AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->UncheckedAt(worstClusterNumber));
1216   
1217   // Re-calculate track parameters at the (new) first cluster
1218   if (!RetracePartialTrack(trackCandidate,(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->UncheckedAt(1))) return kFALSE;
1219   
1220   // skip track out of limits
1221   if (!IsAcceptable(*((AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->First()))) return kFALSE;
1222   
1223   // Look for new cluster(s) in next station
1224   return FollowTrackInStation(trackCandidate, clusterStore, nextStation);
1225   
1226 }
1227
1228   //__________________________________________________________________________
1229 Bool_t AliMUONTrackReconstructorK::RunSmoother(AliMUONTrack &track)
1230 {
1231   /// Compute new track parameters and their covariances using smoother
1232   AliDebug(1,"Enter UseSmoother");
1233   
1234   AliMUONTrackParam *previousTrackParam = (AliMUONTrackParam*) track.GetTrackParamAtCluster()->First();
1235   
1236   // Smoothed parameters and covariances at first cluster = filtered parameters and covariances
1237   previousTrackParam->SetSmoothParameters(previousTrackParam->GetParameters());
1238   previousTrackParam->SetSmoothCovariances(previousTrackParam->GetCovariances());
1239
1240   AliMUONTrackParam *currentTrackParam = (AliMUONTrackParam*) track.GetTrackParamAtCluster()->After(previousTrackParam);
1241   
1242   // Save local chi2 at first cluster = last additional chi2 provided by Kalman
1243   previousTrackParam->SetLocalChi2(previousTrackParam->GetTrackChi2() - currentTrackParam->GetTrackChi2());
1244   
1245   // if the track contains only 2 clusters simply copy the filtered parameters
1246   if (track.GetNClusters() == 2) {
1247     currentTrackParam->SetSmoothParameters(currentTrackParam->GetParameters());
1248     currentTrackParam->SetSmoothCovariances(currentTrackParam->GetCovariances());
1249     currentTrackParam->SetLocalChi2(currentTrackParam->GetTrackChi2());
1250     return kTRUE;
1251   }
1252   
1253   while (currentTrackParam) {
1254     
1255     // Get variables
1256     const TMatrixD &extrapParameters          = previousTrackParam->GetExtrapParameters();  // X(k+1 k)
1257     const TMatrixD &filteredParameters        = currentTrackParam->GetParameters();         // X(k k)
1258     const TMatrixD &previousSmoothParameters  = previousTrackParam->GetSmoothParameters();  // X(k+1 n)
1259     const TMatrixD &propagator                = previousTrackParam->GetPropagator();        // F(k)
1260     const TMatrixD &extrapCovariances         = previousTrackParam->GetExtrapCovariances(); // C(k+1 k)
1261     const TMatrixD &filteredCovariances       = currentTrackParam->GetCovariances();        // C(k k)
1262     const TMatrixD &previousSmoothCovariances = previousTrackParam->GetSmoothCovariances(); // C(k+1 n)
1263     
1264     // Compute smoother gain: A(k) = C(kk) * F(k)^t * (C(k+1 k))^-1
1265     TMatrixD extrapWeight(extrapCovariances);
1266     if (extrapWeight.Determinant() != 0) {
1267       extrapWeight.Invert(); // (C(k+1 k))^-1
1268     } else {
1269       AliWarning(" Determinant = 0");
1270       return kFALSE;
1271     }
1272     TMatrixD smootherGain(filteredCovariances,TMatrixD::kMultTranspose,propagator); // C(kk) * F(k)^t
1273     smootherGain *= extrapWeight; // C(kk) * F(k)^t * (C(k+1 k))^-1
1274     
1275     // Compute smoothed parameters: X(k n) = X(k k) + A(k) * (X(k+1 n) - X(k+1 k))
1276     TMatrixD tmpParam(previousSmoothParameters,TMatrixD::kMinus,extrapParameters); // X(k+1 n) - X(k+1 k)
1277     TMatrixD smoothParameters(smootherGain,TMatrixD::kMult,tmpParam); // A(k) * (X(k+1 n) - X(k+1 k))
1278     smoothParameters += filteredParameters; // X(k k) + A(k) * (X(k+1 n) - X(k+1 k))
1279     
1280     // Save smoothed parameters
1281     currentTrackParam->SetSmoothParameters(smoothParameters);
1282     
1283     // Compute smoothed covariances: C(k n) = C(k k) + A(k) * (C(k+1 n) - C(k+1 k)) * (A(k))^t
1284     TMatrixD tmpCov(previousSmoothCovariances,TMatrixD::kMinus,extrapCovariances); // C(k+1 n) - C(k+1 k)
1285     TMatrixD tmpCov2(tmpCov,TMatrixD::kMultTranspose,smootherGain); // (C(k+1 n) - C(k+1 k)) * (A(k))^t
1286     TMatrixD smoothCovariances(smootherGain,TMatrixD::kMult,tmpCov2); // A(k) * (C(k+1 n) - C(k+1 k)) * (A(k))^t
1287     smoothCovariances += filteredCovariances; // C(k k) + A(k) * (C(k+1 n) - C(k+1 k)) * (A(k))^t
1288     
1289     // Save smoothed covariances
1290     currentTrackParam->SetSmoothCovariances(smoothCovariances);
1291     
1292     // Compute smoothed residual: r(k n) = cluster - X(k n)
1293     AliMUONVCluster* cluster = currentTrackParam->GetClusterPtr();
1294     TMatrixD smoothResidual(2,1);
1295     smoothResidual.Zero();
1296     smoothResidual(0,0) = cluster->GetX() - smoothParameters(0,0);
1297     smoothResidual(1,0) = cluster->GetY() - smoothParameters(2,0);
1298     
1299     // Compute weight of smoothed residual: W(k n) = (clusterCov - C(k n))^-1
1300     TMatrixD smoothResidualWeight(2,2);
1301     smoothResidualWeight(0,0) = cluster->GetErrX2() - smoothCovariances(0,0);
1302     smoothResidualWeight(0,1) = - smoothCovariances(0,2);
1303     smoothResidualWeight(1,0) = - smoothCovariances(2,0);
1304     smoothResidualWeight(1,1) = cluster->GetErrY2() - smoothCovariances(2,2);
1305     if (smoothResidualWeight.Determinant() != 0) {
1306       smoothResidualWeight.Invert();
1307     } else {
1308       AliWarning(" Determinant = 0");
1309       return kFALSE;
1310     }
1311     
1312     // Compute local chi2 = (r(k n))^t * W(k n) * r(k n)
1313     TMatrixD tmpChi2(smoothResidual,TMatrixD::kTransposeMult,smoothResidualWeight); // (r(k n))^t * W(k n)
1314     TMatrixD localChi2(tmpChi2,TMatrixD::kMult,smoothResidual); // (r(k n))^t * W(k n) * r(k n)
1315     
1316     // Save local chi2
1317     currentTrackParam->SetLocalChi2(localChi2(0,0));
1318     
1319     previousTrackParam = currentTrackParam;
1320     currentTrackParam = (AliMUONTrackParam*) track.GetTrackParamAtCluster()->After(previousTrackParam);
1321   }
1322   
1323   return kTRUE;
1324   
1325 }
1326
1327   //__________________________________________________________________________
1328 Bool_t AliMUONTrackReconstructorK::ComplementTracks(const AliMUONVClusterStore& clusterStore)
1329 {
1330   /// Complete tracks by adding missing clusters (if there is an overlap between
1331   /// two detection elements, the track may have two clusters in the same chamber).
1332   /// Recompute track parameters and covariances at each clusters.
1333   /// Remove tracks getting abnormal (i.e. extrapolation failed...) after being complemented.
1334   /// Return kTRUE if one or more tracks have been complemented or removed.
1335   AliDebug(1,"Enter ComplementTracks");
1336   
1337   Int_t chamberId, detElemId;
1338   Double_t chi2OfCluster, addChi2TrackAtCluster, bestAddChi2TrackAtCluster;
1339   Double_t maxChi2OfCluster = 2. * GetRecoParam()->GetSigmaCutForTracking() *
1340                                    GetRecoParam()->GetSigmaCutForTracking(); // 2 because 2 quantities in chi2
1341   Bool_t foundOneCluster, trackModified, hasChanged = kFALSE;
1342   AliMUONVCluster *cluster;
1343   AliMUONTrackParam *trackParam, *previousTrackParam, *nextTrackParam, trackParamAtCluster, bestTrackParamAtCluster;
1344   AliMUONTrack *nextTrack;
1345   
1346   AliMUONTrack *track = (AliMUONTrack*) fRecTracksPtr->First();
1347   while (track) {
1348     trackModified = kFALSE;
1349     
1350     trackParam = (AliMUONTrackParam*)track->GetTrackParamAtCluster()->First();
1351     previousTrackParam = trackParam;
1352     while (trackParam) {
1353       foundOneCluster = kFALSE;
1354       bestAddChi2TrackAtCluster = AliMUONTrack::MaxChi2();
1355       chamberId = trackParam->GetClusterPtr()->GetChamberId();
1356       detElemId = trackParam->GetClusterPtr()->GetDetElemId();
1357       
1358       // prepare nextTrackParam before adding new cluster because of the sorting
1359       nextTrackParam = (AliMUONTrackParam*)track->GetTrackParamAtCluster()->After(trackParam);
1360       
1361       // Create iterators to loop over clusters in current chamber
1362       TIter nextInCh(clusterStore.CreateChamberIterator(chamberId,chamberId));
1363       
1364       // look for one second candidate in the same chamber
1365       while ( ( cluster = static_cast<AliMUONVCluster*>(nextInCh()) ) ) {
1366         
1367         // look for a cluster in another detection element
1368         if (cluster->GetDetElemId() == detElemId) continue;
1369         
1370         // try to add the current cluster fast
1371         if (!TryOneClusterFast(*trackParam, cluster)) continue;
1372         
1373         // try to add the current cluster accurately
1374         // never use track parameters at last cluster because the covariance matrix is meaningless
1375         if (nextTrackParam) chi2OfCluster = TryOneCluster(*trackParam, cluster, trackParamAtCluster);
1376         else chi2OfCluster = TryOneCluster(*previousTrackParam, cluster, trackParamAtCluster);
1377         
1378         // if good chi2 then consider to add this cluster to the track
1379         if (chi2OfCluster < maxChi2OfCluster) {
1380           
1381           // Compute local track parameters including current cluster using kalman filter
1382           addChi2TrackAtCluster = RunKalmanFilter(trackParamAtCluster);
1383           
1384           // keep track of the best cluster
1385           if (addChi2TrackAtCluster < bestAddChi2TrackAtCluster) {
1386             bestAddChi2TrackAtCluster = addChi2TrackAtCluster;
1387             bestTrackParamAtCluster = trackParamAtCluster;
1388             foundOneCluster = kTRUE;
1389           }
1390           
1391         }
1392         
1393       }
1394       
1395       // add new cluster if any
1396       if (foundOneCluster) {
1397         
1398         // Printout for debuging
1399         if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
1400           cout << "ComplementTracks: found one cluster in chamber(1..): " << chamberId+1 << endl;
1401           bestTrackParamAtCluster.GetClusterPtr()->Print();
1402           if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 2) || (AliLog::GetGlobalDebugLevel() >= 2)) {
1403             cout<<endl<<"Track parameters and covariances at cluster:"<<endl;
1404             bestTrackParamAtCluster.GetParameters().Print();
1405             bestTrackParamAtCluster.GetCovariances().Print();
1406           }
1407         }
1408         
1409         trackParam->SetRemovable(kTRUE);
1410         bestTrackParamAtCluster.SetRemovable(kTRUE);
1411         track->AddTrackParamAtCluster(bestTrackParamAtCluster,*(bestTrackParamAtCluster.GetClusterPtr()));
1412         trackModified = kTRUE;
1413         hasChanged = kTRUE;
1414       }
1415       
1416       previousTrackParam = trackParam;
1417       trackParam = nextTrackParam;
1418     }
1419     
1420     // prepare next track
1421     nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track);
1422     
1423     // re-compute track parameters using kalman filter if needed
1424     if (trackModified && !RetraceTrack(*track,kTRUE)) {
1425       AliWarning("track modified but problem occur during refitting --> remove track");
1426       fRecTracksPtr->Remove(track);
1427       fNRecTracks--;
1428     }
1429     
1430     track = nextTrack;
1431   }
1432   
1433   return hasChanged;
1434   
1435 }
1436
1437 //__________________________________________________________________________
1438 void AliMUONTrackReconstructorK::ImproveTrack(AliMUONTrack &track)
1439 {
1440   /// Improve the given track by removing removable clusters with local chi2 highter than the defined cut
1441   /// Removable clusters are identified by the method AliMUONTrack::TagRemovableClusters()
1442   /// Recompute track parameters and covariances at the remaining clusters
1443   /// and if something goes wrong (i.e. extrapolation failed...) set track chi2 to max value
1444   AliDebug(1,"Enter ImproveTrack");
1445   
1446   Double_t localChi2, worstLocalChi2;
1447   AliMUONTrackParam *trackParamAtCluster, *worstTrackParamAtCluster, *nextTrackParam, *next2nextTrackParam;
1448   Int_t nextChamber, next2nextChamber;
1449   Bool_t smoothed;
1450   Double_t sigmaCut2 = GetRecoParam()->GetSigmaCutForImprovement() *
1451                        GetRecoParam()->GetSigmaCutForImprovement();
1452   
1453   while (!track.IsImproved()) {
1454     
1455     // identify removable clusters
1456     track.TagRemovableClusters(GetRecoParam()->RequestedStationMask());
1457     
1458     // Run smoother if required
1459     smoothed = kFALSE;
1460     if (GetRecoParam()->UseSmoother()) smoothed = RunSmoother(track);
1461     
1462     // Use standard procedure to compute local chi2 if smoother not required or not working
1463     if (!smoothed) {
1464       
1465       // Update track parameters and covariances
1466       if (!track.UpdateCovTrackParamAtCluster()) {
1467         AliWarning("unable to update track parameters and covariances --> stop improvement");
1468         // restore the kalman parameters
1469         RetraceTrack(track,kTRUE);
1470         break;
1471       }
1472       
1473       // Compute local chi2 of each clusters
1474       track.ComputeLocalChi2(kTRUE);
1475     }
1476     
1477     // Look for the cluster to remove
1478     worstTrackParamAtCluster = 0x0;
1479     worstLocalChi2 = -1.;
1480     trackParamAtCluster = (AliMUONTrackParam*)track.GetTrackParamAtCluster()->First();
1481     while (trackParamAtCluster) {
1482       
1483       // save parameters into smooth parameters in case of smoother did not work properly
1484       if (GetRecoParam()->UseSmoother() && !smoothed) {
1485         trackParamAtCluster->SetSmoothParameters(trackParamAtCluster->GetParameters());
1486         trackParamAtCluster->SetSmoothCovariances(trackParamAtCluster->GetCovariances());
1487       }
1488       
1489       // Pick up cluster with the worst chi2
1490       localChi2 = trackParamAtCluster->GetLocalChi2();
1491       if (localChi2 > worstLocalChi2) {
1492         worstLocalChi2 = localChi2;
1493         worstTrackParamAtCluster = trackParamAtCluster;
1494       }
1495       
1496       trackParamAtCluster = (AliMUONTrackParam*)track.GetTrackParamAtCluster()->After(trackParamAtCluster);
1497     }
1498     
1499     // Check whether the worst chi2 is under requirement or not
1500     if (worstLocalChi2 < 2. * sigmaCut2) { // 2 because 2 quantities in chi2
1501       track.SetImproved(kTRUE);
1502       break;
1503     }
1504     
1505     // if the worst cluster is not removable then stop improvement
1506     if (!worstTrackParamAtCluster->IsRemovable()) {
1507       // restore the kalman parameters in case they have been lost
1508       if (!smoothed) RetraceTrack(track,kTRUE);
1509       break;
1510     }
1511     
1512     // get track parameters at cluster next to the one to be removed
1513     nextTrackParam = (AliMUONTrackParam*) track.GetTrackParamAtCluster()->After(worstTrackParamAtCluster);
1514     
1515     // Remove the worst cluster
1516     track.RemoveTrackParamAtCluster(worstTrackParamAtCluster);
1517     
1518     // Re-calculate track parameters
1519     // - from the cluster immediately downstream the one suppressed
1520     // - or from the begining - if parameters have been re-computed using the standard method (kalman parameters have been lost)
1521     //                        - or if the removed cluster was used to compute the tracking seed
1522     Bool_t normalExtrap;
1523     if (smoothed && nextTrackParam) {
1524       
1525       nextChamber = nextTrackParam->GetClusterPtr()->GetChamberId();
1526       next2nextTrackParam = nextTrackParam;
1527       do {
1528         
1529         next2nextChamber = next2nextTrackParam->GetClusterPtr()->GetChamberId();
1530         next2nextTrackParam = (AliMUONTrackParam*) track.GetTrackParamAtCluster()->After(next2nextTrackParam);
1531         
1532       } while (next2nextTrackParam && (next2nextChamber == nextChamber));
1533       
1534       if (next2nextChamber == nextChamber) normalExtrap = RetraceTrack(track,kTRUE);
1535       else normalExtrap = RetracePartialTrack(track,nextTrackParam);
1536       
1537     } else normalExtrap = RetraceTrack(track,kTRUE);
1538     
1539     // stop in case of extrapolation problem
1540     if (!normalExtrap) {
1541       AliWarning("track partially improved but problem occur during refitting --> stop improvement");
1542       break;
1543     }
1544     
1545     // Printout for debuging
1546     if ((AliLog::GetDebugLevel("MUON","AliMUONTrackReconstructorK") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
1547       cout << "ImproveTracks: track " << fRecTracksPtr->IndexOf(&track)+1 << " improved " << endl;
1548     }
1549     
1550   }
1551   
1552 }
1553
1554 //__________________________________________________________________________
1555 Bool_t AliMUONTrackReconstructorK::FinalizeTrack(AliMUONTrack &track)
1556 {
1557   /// Update track parameters and covariances at each attached cluster
1558   /// using smoother if required, if not already done
1559   /// return kFALSE if the track cannot be extrapolated uo to the last chamber
1560   
1561   AliMUONTrackParam *trackParamAtCluster;
1562   Bool_t smoothed = kFALSE;
1563   
1564   // update track parameters (using smoother if required) if not already done
1565   if (track.IsImproved()) smoothed = GetRecoParam()->UseSmoother();
1566   else {
1567     if (GetRecoParam()->UseSmoother()) smoothed = RunSmoother(track);
1568     if (!smoothed) {
1569       if (track.UpdateCovTrackParamAtCluster()) track.ComputeLocalChi2(kTRUE);
1570       else {
1571         AliWarning("finalization failed due to extrapolation problem");
1572         return kFALSE;
1573       }
1574     }
1575   }
1576   
1577   // copy smoothed parameters and covariances if any
1578   if (smoothed) {
1579     
1580     trackParamAtCluster = (AliMUONTrackParam*) (track.GetTrackParamAtCluster()->First());
1581     while (trackParamAtCluster) {
1582       
1583       trackParamAtCluster->SetParameters(trackParamAtCluster->GetSmoothParameters());
1584       trackParamAtCluster->SetCovariances(trackParamAtCluster->GetSmoothCovariances());
1585       
1586       trackParamAtCluster = (AliMUONTrackParam*) (track.GetTrackParamAtCluster()->After(trackParamAtCluster));
1587     }
1588     
1589   }
1590   
1591   return kTRUE;
1592   
1593 }
1594
1595   //__________________________________________________________________________
1596 Bool_t AliMUONTrackReconstructorK::RefitTrack(AliMUONTrack &track, Bool_t enableImprovement)
1597 {
1598   /// re-fit the given track
1599   AliDebug(1,"Enter RefitTrack");
1600   
1601   // check validity of the track (i.e. at least 2 chambers hit on stations 4 and 5)
1602   if (!track.IsValid(0)) {
1603     AliWarning("the track is not valid --> unable to refit");
1604     return kFALSE;
1605   }
1606   
1607   // re-compute track parameters and covariances using Kalman filter
1608   if (!RetraceTrack(track,kTRUE)) {
1609     AliWarning("bad track refitting due to extrapolation failure");
1610     return kFALSE;
1611   }
1612   
1613   // Improve the reconstructed tracks if required
1614   track.SetImproved(kFALSE);
1615   if (enableImprovement && GetRecoParam()->ImproveTracks()) ImproveTrack(track);
1616   
1617   // Fill AliMUONTrack data members
1618   if (track.GetGlobalChi2() < AliMUONTrack::MaxChi2()) return FinalizeTrack(track);
1619   else {
1620     AliWarning("track not finalized due to extrapolation failure");
1621     return kFALSE;
1622   }
1623   
1624 }
1625