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