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