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