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