Line data Source code
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 : #include "AliMUONGeometryTransformer.h"
78 : #include "AliMUONVDigit.h"
79 :
80 : #include "AliMpDEManager.h"
81 : #include "AliMpArea.h"
82 :
83 : #include "AliMpDDLStore.h"
84 : #include "AliMpVSegmentation.h"
85 : #include "AliMpSegmentation.h"
86 : #include "AliMpPad.h"
87 : #include "AliMpDetElement.h"
88 : #include "AliMpCathodType.h"
89 :
90 : #include "AliLog.h"
91 : #include "AliCodeTimer.h"
92 : #include "AliTracker.h"
93 :
94 : #include <TClonesArray.h>
95 : #include <TMath.h>
96 : #include <TMatrixD.h>
97 : #include <TVector2.h>
98 :
99 : #include <Riostream.h>
100 :
101 : using std::cout;
102 : using std::endl;
103 : /// \cond CLASSIMP
104 18 : ClassImp(AliMUONVTrackReconstructor) // Class implementation in ROOT context
105 : /// \endcond
106 :
107 : //__________________________________________________________________________
108 : AliMUONVTrackReconstructor::AliMUONVTrackReconstructor(const AliMUONRecoParam* recoParam,
109 : AliMUONVClusterServer* clusterServer,
110 : const AliMUONGeometryTransformer* transformer)
111 2 : : TObject(),
112 2 : fRecTracksPtr(0x0),
113 2 : fNRecTracks(0),
114 2 : fClusterServer(clusterServer),
115 2 : fkRecoParam(recoParam),
116 2 : fkTransformer(transformer),
117 2 : fMaxMCSAngle2(0x0)
118 6 : {
119 : /// Constructor for class AliMUONVTrackReconstructor
120 : /// WARNING: if clusterServer=0x0, no clusterization will be possible at this level
121 :
122 : // Memory allocation for the TClonesArray of reconstructed tracks
123 6 : fRecTracksPtr = new TClonesArray("AliMUONTrack", 100);
124 :
125 : // set the magnetic field for track extrapolations
126 2 : AliMUONTrackExtrap::SetField();
127 :
128 : // set the maximum MCS angle in chamber from the minimum acceptable momentum
129 2 : AliMUONTrackParam param;
130 6 : Double_t inverseBendingP = (GetRecoParam()->GetMinBendingMomentum() > 0.) ? 1./GetRecoParam()->GetMinBendingMomentum() : 1.;
131 2 : param.SetInverseBendingMomentum(inverseBendingP);
132 6 : fMaxMCSAngle2 = new Double_t [AliMUONConstants::NTrackingCh()];
133 66 : for (Int_t iCh=0; iCh<AliMUONConstants::NTrackingCh(); iCh++)
134 40 : fMaxMCSAngle2[iCh] = AliMUONTrackExtrap::GetMCSAngle2(param, AliMUONConstants::ChamberThicknessInX0(iCh), 1.);
135 :
136 2 : }
137 :
138 : //__________________________________________________________________________
139 : AliMUONVTrackReconstructor::~AliMUONVTrackReconstructor()
140 4 : {
141 : /// Destructor for class AliMUONVTrackReconstructor
142 4 : delete fRecTracksPtr;
143 4 : delete[] fMaxMCSAngle2;
144 2 : }
145 :
146 : //__________________________________________________________________________
147 : void AliMUONVTrackReconstructor::ResetTracks()
148 : {
149 : /// To reset the TClonesArray of reconstructed tracks
150 24 : if (fRecTracksPtr) fRecTracksPtr->Clear("C");
151 8 : fNRecTracks = 0;
152 8 : return;
153 : }
154 :
155 : //__________________________________________________________________________
156 : void AliMUONVTrackReconstructor::EventReconstruct(AliMUONVClusterStore& clusterStore, AliMUONVTrackStore& trackStore)
157 : {
158 : /// To reconstruct one event
159 32 : AliDebug(1,"");
160 8 : AliCodeTimerAuto("",0);
161 :
162 : // Reset array of tracks
163 8 : ResetTracks();
164 :
165 : // Look for candidates from clusters in stations(1..) 4 and 5 (abort in case of failure)
166 16 : if (!MakeTrackCandidates(clusterStore)) return;
167 :
168 : // Look for extra candidates from clusters in stations(1..) 4 and 5 (abort in case of failure)
169 8 : if (GetRecoParam()->MakeMoreTrackCandidates()) {
170 0 : if (!MakeMoreTrackCandidates(clusterStore)) return;
171 : }
172 :
173 : // Stop tracking if no candidate found
174 16 : if (fRecTracksPtr->GetEntriesFast() == 0) return;
175 :
176 : // Follow tracks in stations(1..) 3, 2 and 1 (abort in case of failure)
177 16 : if (!FollowTracks(clusterStore)) return;
178 :
179 : // Complement the reconstructed tracks
180 8 : if (GetRecoParam()->ComplementTracks()) {
181 18 : if (ComplementTracks(clusterStore)) RemoveIdenticalTracks();
182 : }
183 :
184 : // Improve the reconstructed tracks
185 16 : if (GetRecoParam()->ImproveTracks()) ImproveTracks();
186 :
187 : // Remove connected tracks
188 8 : RemoveConnectedTracks(3, 4, kFALSE);
189 8 : RemoveConnectedTracks(2, 2, kFALSE);
190 8 : if (GetRecoParam()->RemoveConnectedTracksInSt12()) RemoveConnectedTracks(0, 1, kFALSE);
191 :
192 : // Fill AliMUONTrack data members
193 8 : Finalize();
194 16 : if (!GetRecoParam()->RemoveConnectedTracksInSt12()) TagConnectedTracks(0, 1, kTRUE);
195 :
196 : // Make sure there is no bad track left
197 8 : RemoveBadTracks();
198 :
199 : // Refit the reconstructed tracks with a different resolution for mono-cathod clusters
200 8 : if (GetRecoParam()->DiscardMonoCathodClusters()) DiscardMonoCathodClusters();
201 :
202 : // Add tracks to MUON data container
203 48 : for (Int_t i=0; i<fNRecTracks; ++i)
204 : {
205 32 : AliMUONTrack * track = (AliMUONTrack*) fRecTracksPtr->At(i);
206 16 : track->SetUniqueID(i+1);
207 16 : trackStore.Add(*track);
208 : }
209 :
210 16 : }
211 :
212 : //__________________________________________________________________________
213 : Bool_t AliMUONVTrackReconstructor::IsAcceptable(AliMUONTrackParam &trackParam)
214 : {
215 : /// Return kTRUE if the track is within given limits on momentum/angle/origin
216 :
217 412 : const TMatrixD& kParamCov = trackParam.GetCovariances();
218 206 : Int_t chamber = trackParam.GetClusterPtr()->GetChamberId();
219 206 : Double_t z = trackParam.GetZ();
220 206 : Double_t sigmaCut = GetRecoParam()->GetSigmaCutForTracking();
221 :
222 : // MCS dipersion
223 : Double_t angleMCS2 = 0.;
224 : Double_t impactMCS2 = 0.;
225 206 : if (AliMUONTrackExtrap::IsFieldON() && chamber < 6) {
226 :
227 : // track momentum is known
228 1056 : for (Int_t iCh=0; iCh<=chamber; iCh++) {
229 410 : Double_t localMCS2 = AliMUONTrackExtrap::GetMCSAngle2(trackParam, AliMUONConstants::ChamberThicknessInX0(iCh), 1.);
230 410 : angleMCS2 += localMCS2;
231 410 : impactMCS2 += AliMUONConstants::DefaultChamberZ(chamber) * AliMUONConstants::DefaultChamberZ(chamber) * localMCS2;
232 : }
233 :
234 118 : } else {
235 :
236 : // track momentum is unknown
237 1620 : for (Int_t iCh=0; iCh<=chamber; iCh++) {
238 722 : angleMCS2 += fMaxMCSAngle2[iCh];
239 722 : impactMCS2 += AliMUONConstants::DefaultChamberZ(chamber) * AliMUONConstants::DefaultChamberZ(chamber) * fMaxMCSAngle2[iCh];
240 : }
241 :
242 : }
243 :
244 : // ------ track selection in non bending direction ------
245 206 : if (GetRecoParam()->SelectOnTrackSlope()) {
246 :
247 : // check if non bending slope is within tolerances
248 0 : Double_t nonBendingSlopeErr = TMath::Sqrt(kParamCov(1,1) + angleMCS2);
249 0 : if ((TMath::Abs(trackParam.GetNonBendingSlope()) - sigmaCut * nonBendingSlopeErr) > GetRecoParam()->GetMaxNonBendingSlope()) return kFALSE;
250 :
251 0 : } else {
252 :
253 : // or check if non bending impact parameter is within tolerances
254 206 : Double_t nonBendingImpactParam = TMath::Abs(trackParam.GetNonBendingCoor() - z * trackParam.GetNonBendingSlope());
255 206 : Double_t nonBendingImpactParamErr = TMath::Sqrt(kParamCov(0,0) + z * z * kParamCov(1,1) - 2. * z * kParamCov(0,1) + impactMCS2);
256 206 : if ((nonBendingImpactParam - sigmaCut * nonBendingImpactParamErr) > (3. * GetRecoParam()->GetNonBendingVertexDispersion())) return kFALSE;
257 :
258 206 : }
259 :
260 : // ------ track selection in bending direction ------
261 206 : if (AliMUONTrackExtrap::IsFieldON()) { // depending whether the field is ON or OFF
262 :
263 : // check if bending momentum is within tolerances
264 206 : Double_t bendingMomentum = TMath::Abs(1. / trackParam.GetInverseBendingMomentum());
265 206 : Double_t bendingMomentumErr = TMath::Sqrt(kParamCov(4,4)) * bendingMomentum * bendingMomentum;
266 324 : if (chamber < 6 && (bendingMomentum + sigmaCut * bendingMomentumErr) < GetRecoParam()->GetMinBendingMomentum()) return kFALSE;
267 206 : else if ((bendingMomentum + 3. * bendingMomentumErr) < GetRecoParam()->GetMinBendingMomentum()) return kFALSE;
268 :
269 206 : } else {
270 :
271 0 : if (GetRecoParam()->SelectOnTrackSlope()) {
272 :
273 : // check if bending slope is within tolerances
274 0 : Double_t bendingSlopeErr = TMath::Sqrt(kParamCov(3,3) + angleMCS2);
275 0 : if ((TMath::Abs(trackParam.GetBendingSlope()) - sigmaCut * bendingSlopeErr) > GetRecoParam()->GetMaxBendingSlope()) return kFALSE;
276 :
277 0 : } else {
278 :
279 : // or check if bending impact parameter is within tolerances
280 0 : Double_t bendingImpactParam = TMath::Abs(trackParam.GetBendingCoor() - z * trackParam.GetBendingSlope());
281 0 : Double_t bendingImpactParamErr = TMath::Sqrt(kParamCov(2,2) + z * z * kParamCov(3,3) - 2. * z * kParamCov(2,3) + impactMCS2);
282 0 : if ((bendingImpactParam - sigmaCut * bendingImpactParamErr) > (3. * GetRecoParam()->GetBendingVertexDispersion())) return kFALSE;
283 :
284 0 : }
285 :
286 : }
287 :
288 206 : return kTRUE;
289 :
290 206 : }
291 :
292 : //__________________________________________________________________________
293 : TClonesArray* AliMUONVTrackReconstructor::MakeSegmentsBetweenChambers(const AliMUONVClusterStore& clusterStore, Int_t ch1, Int_t ch2)
294 : {
295 : /// To make the list of segments from the list of clusters in the 2 given chambers.
296 : /// Return a TClonesArray of new segments (segments made in a previous call of this function are removed).
297 64 : AliDebug(1,Form("Enter MakeSegmentsBetweenChambers (1..) %d-%d", ch1+1, ch2+1));
298 16 : AliCodeTimerAuto("",0);
299 :
300 : AliMUONVCluster *cluster1, *cluster2;
301 : AliMUONObjectPair *segment;
302 : Double_t z1 = 0., z2 = 0., dZ = 0.;
303 : Double_t nonBendingSlope = 0., nonBendingSlopeErr = 0., nonBendingImpactParam = 0., nonBendingImpactParamErr = 0.;
304 : Double_t bendingSlope = 0., bendingSlopeErr = 0., bendingImpactParam = 0., bendingImpactParamErr = 0., bendingImpactParamErr2 = 0.;
305 : Double_t bendingMomentum = 0., bendingMomentumErr = 0.;
306 16 : Double_t bendingVertexDispersion2 = GetRecoParam()->GetBendingVertexDispersion() * GetRecoParam()->GetBendingVertexDispersion();
307 : Double_t angleMCS2 = 0.; // maximum angular dispersion**2 due to MCS in chamber
308 : Double_t impactMCS2 = 0.; // maximum impact parameter dispersion**2 due to MCS in chamber
309 288 : for (Int_t iCh=0; iCh<=ch1; iCh++) {
310 128 : angleMCS2 += fMaxMCSAngle2[iCh];
311 128 : impactMCS2 += AliMUONConstants::DefaultChamberZ(iCh) * AliMUONConstants::DefaultChamberZ(iCh) * fMaxMCSAngle2[iCh];
312 : }
313 16 : Double_t sigmaCut = GetRecoParam()->GetSigmaCutForTracking();
314 :
315 : // Create iterators to loop over clusters in both chambers
316 32 : TIter nextInCh1(clusterStore.CreateChamberIterator(ch1,ch1));
317 32 : TIter nextInCh2(clusterStore.CreateChamberIterator(ch2,ch2));
318 :
319 : // list of segments
320 24 : static TClonesArray *segments = new TClonesArray("AliMUONObjectPair", 100);
321 16 : segments->Clear("C");
322 :
323 : // Loop over clusters in the first chamber of the station
324 100 : while ( ( cluster1 = static_cast<AliMUONVCluster*>(nextInCh1()) ) ) {
325 34 : z1 = cluster1->GetZ();
326 :
327 : // reset cluster iterator of chamber 2
328 34 : nextInCh2.Reset();
329 :
330 : // Loop over clusters in the second chamber of the station
331 204 : while ( ( cluster2 = static_cast<AliMUONVCluster*>(nextInCh2()) ) ) {
332 68 : z2 = cluster2->GetZ();
333 68 : dZ = z1 - z2;
334 :
335 : // ------ track selection in non bending direction ------
336 204 : nonBendingSlope = (cluster1->GetX() - cluster2->GetX()) / dZ;
337 136 : if (GetRecoParam()->SelectOnTrackSlope()) {
338 :
339 : // check if non bending slope is within tolerances
340 68 : nonBendingSlopeErr = TMath::Sqrt((cluster1->GetErrX2() + cluster2->GetErrX2()) / dZ / dZ + angleMCS2);
341 0 : if ((TMath::Abs(nonBendingSlope) - sigmaCut * nonBendingSlopeErr) > GetRecoParam()->GetMaxNonBendingSlope()) continue;
342 :
343 : } else {
344 :
345 : // or check if non bending impact parameter is within tolerances
346 204 : nonBendingImpactParam = TMath::Abs(cluster1->GetX() - cluster1->GetZ() * nonBendingSlope);
347 204 : nonBendingImpactParamErr = TMath::Sqrt((z1 * z1 * cluster2->GetErrX2() + z2 * z2 * cluster1->GetErrX2()) / dZ / dZ + impactMCS2);
348 68 : if ((nonBendingImpactParam - sigmaCut * nonBendingImpactParamErr) > (3. * GetRecoParam()->GetNonBendingVertexDispersion())) continue;
349 :
350 : }
351 :
352 : // ------ track selection in bending direction ------
353 102 : bendingSlope = (cluster1->GetY() - cluster2->GetY()) / dZ;
354 34 : if (AliMUONTrackExtrap::IsFieldON()) { // depending whether the field is ON or OFF
355 :
356 : // check if bending momentum is within tolerances
357 102 : bendingImpactParam = cluster1->GetY() - cluster1->GetZ() * bendingSlope;
358 102 : bendingImpactParamErr2 = (z1 * z1 * cluster2->GetErrY2() + z2 * z2 * cluster1->GetErrY2()) / dZ / dZ + impactMCS2;
359 68 : bendingMomentum = TMath::Abs(AliMUONTrackExtrap::GetBendingMomentumFromImpactParam(bendingImpactParam));
360 68 : bendingMomentumErr = TMath::Sqrt((bendingVertexDispersion2 + bendingImpactParamErr2) /
361 68 : bendingImpactParam / bendingImpactParam + 0.01) * bendingMomentum;
362 34 : if ((bendingMomentum + 3. * bendingMomentumErr) < GetRecoParam()->GetMinBendingMomentum()) continue;
363 :
364 : } else {
365 :
366 0 : if (GetRecoParam()->SelectOnTrackSlope()) {
367 :
368 : // check if bending slope is within tolerances
369 0 : bendingSlopeErr = TMath::Sqrt((cluster1->GetErrY2() + cluster2->GetErrY2()) / dZ / dZ + angleMCS2);
370 0 : if ((TMath::Abs(bendingSlope) - sigmaCut * bendingSlopeErr) > GetRecoParam()->GetMaxBendingSlope()) continue;
371 :
372 : } else {
373 :
374 : // or check if bending impact parameter is within tolerances
375 0 : bendingImpactParam = TMath::Abs(cluster1->GetY() - cluster1->GetZ() * bendingSlope);
376 0 : bendingImpactParamErr = TMath::Sqrt((z1 * z1 * cluster2->GetErrY2() + z2 * z2 * cluster1->GetErrY2()) / dZ / dZ + impactMCS2);
377 0 : if ((bendingImpactParam - sigmaCut * bendingImpactParamErr) > (3. * GetRecoParam()->GetBendingVertexDispersion())) continue;
378 :
379 : }
380 :
381 : }
382 :
383 : // make new segment
384 136 : segment = new ((*segments)[segments->GetLast()+1]) AliMUONObjectPair(cluster1, cluster2, kFALSE, kFALSE);
385 :
386 : // Printout for debug
387 68 : if (AliLog::GetGlobalDebugLevel() > 1) {
388 0 : cout << "segmentIndex(0...): " << segments->GetLast() << endl;
389 0 : segment->Dump();
390 0 : cout << "Cluster in first chamber" << endl;
391 0 : cluster1->Print();
392 0 : cout << "Cluster in second chamber" << endl;
393 0 : cluster2->Print();
394 : }
395 :
396 : }
397 :
398 : }
399 :
400 : // Printout for debug
401 80 : AliDebug(1,Form("chambers%d-%d: NSegments = %d ", ch1+1, ch2+1, segments->GetEntriesFast()));
402 :
403 16 : return segments;
404 16 : }
405 :
406 : //__________________________________________________________________________
407 : void AliMUONVTrackReconstructor::RemoveUsedSegments(TClonesArray& segments)
408 : {
409 : /// To remove pairs of clusters already attached to a track
410 0 : AliDebug(1,"Enter RemoveUsedSegments");
411 0 : Int_t nSegments = segments.GetEntriesFast();
412 0 : Int_t nTracks = fRecTracksPtr->GetEntriesFast();
413 : AliMUONObjectPair *segment;
414 : AliMUONTrack *track;
415 : AliMUONVCluster *cluster, *cluster1, *cluster2;
416 : Bool_t foundCluster1, foundCluster2, removeSegment;
417 :
418 : // Loop over segments
419 0 : for (Int_t iSegment=0; iSegment<nSegments; iSegment++) {
420 0 : segment = (AliMUONObjectPair*) segments.UncheckedAt(iSegment);
421 :
422 0 : cluster1 = (AliMUONVCluster*) segment->First();
423 0 : cluster2 = (AliMUONVCluster*) segment->Second();
424 : removeSegment = kFALSE;
425 :
426 : // Loop over tracks
427 0 : for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
428 0 : track = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iTrack);
429 :
430 : // skip empty slot
431 0 : if (!track) continue;
432 :
433 : foundCluster1 = kFALSE;
434 : foundCluster2 = kFALSE;
435 :
436 : // Loop over clusters
437 0 : Int_t nClusters = track->GetNClusters();
438 0 : for (Int_t iCluster = 0; iCluster < nClusters; iCluster++) {
439 0 : cluster = ((AliMUONTrackParam*) track->GetTrackParamAtCluster()->UncheckedAt(iCluster))->GetClusterPtr();
440 :
441 : // check if both clusters are in that track
442 0 : if (cluster == cluster1) foundCluster1 = kTRUE;
443 0 : else if (cluster == cluster2) foundCluster2 = kTRUE;
444 :
445 0 : if (foundCluster1 && foundCluster2) {
446 : removeSegment = kTRUE;
447 0 : break;
448 : }
449 :
450 : }
451 :
452 0 : if (removeSegment) break;
453 :
454 0 : }
455 :
456 0 : if (removeSegment) segments.RemoveAt(iSegment);
457 :
458 : }
459 :
460 0 : segments.Compress();
461 :
462 : // Printout for debug
463 0 : AliDebug(1,Form("NSegments = %d ", segments.GetEntriesFast()));
464 0 : }
465 :
466 : //__________________________________________________________________________
467 : void AliMUONVTrackReconstructor::RemoveIdenticalTracks()
468 : {
469 : /// To remove identical tracks:
470 : /// Tracks are considered identical if they have all their clusters in common.
471 : /// One keeps the track with the larger number of clusters if need be
472 : AliMUONTrack *track1, *track2;
473 20 : Int_t nTracks = fRecTracksPtr->GetEntriesFast();
474 : Int_t clustersInCommon, nClusters1, nClusters2;
475 : // Loop over first track of the pair
476 172 : for (Int_t iTrack1 = 0; iTrack1 < nTracks; iTrack1++) {
477 76 : track1 = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iTrack1);
478 : // skip empty slot
479 76 : if (!track1) continue;
480 22 : nClusters1 = track1->GetNClusters();
481 : // Loop over second track of the pair
482 286 : for (Int_t iTrack2 = iTrack1+1; iTrack2 < nTracks; iTrack2++) {
483 110 : track2 = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iTrack2);
484 : // skip empty slot
485 110 : if (!track2) continue;
486 46 : nClusters2 = track2->GetNClusters();
487 : // number of clusters in common between two tracks
488 46 : clustersInCommon = track1->ClustersInCommon(track2);
489 : // check for identical tracks
490 72 : if ((clustersInCommon == nClusters1) || (clustersInCommon == nClusters2)) {
491 : // decide which track to remove
492 40 : if (nClusters2 > nClusters1) {
493 : // remove track1 and continue the first loop with the track next to track1
494 20 : fRecTracksPtr->RemoveAt(iTrack1);
495 0 : fNRecTracks--;
496 0 : break;
497 : } else {
498 : // remove track2 and continue the second loop with the track next to track2
499 20 : fRecTracksPtr->RemoveAt(iTrack2);
500 20 : fNRecTracks--;
501 : }
502 20 : }
503 : } // track2
504 22 : } // track1
505 10 : fRecTracksPtr->Compress(); // this is essential to retrieve the TClonesArray afterwards
506 10 : }
507 :
508 : //__________________________________________________________________________
509 : void AliMUONVTrackReconstructor::RemoveDoubleTracks()
510 : {
511 : /// To remove double tracks:
512 : /// Tracks are considered identical if more than half of the clusters of the track
513 : /// which has the smaller number of clusters are in common with the other track.
514 : /// Among two identical tracks, one keeps the track with the larger number of clusters
515 : /// or, if these numbers are equal, the track with the minimum chi2.
516 : AliMUONTrack *track1, *track2;
517 0 : Int_t nTracks = fRecTracksPtr->GetEntriesFast();
518 : Int_t clustersInCommon2, nClusters1, nClusters2;
519 : // Loop over first track of the pair
520 0 : for (Int_t iTrack1 = 0; iTrack1 < nTracks; iTrack1++) {
521 0 : track1 = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iTrack1);
522 : // skip empty slot
523 0 : if (!track1) continue;
524 0 : nClusters1 = track1->GetNClusters();
525 : // Loop over second track of the pair
526 0 : for (Int_t iTrack2 = iTrack1+1; iTrack2 < nTracks; iTrack2++) {
527 0 : track2 = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iTrack2);
528 : // skip empty slot
529 0 : if (!track2) continue;
530 0 : nClusters2 = track2->GetNClusters();
531 : // number of clusters in common between two tracks
532 0 : clustersInCommon2 = 2 * track1->ClustersInCommon(track2);
533 : // check for identical tracks
534 0 : if (clustersInCommon2 > nClusters1 || clustersInCommon2 > nClusters2) {
535 : // decide which track to remove
536 0 : if ((nClusters1 > nClusters2) || ((nClusters1 == nClusters2) && (track1->GetGlobalChi2() <= track2->GetGlobalChi2()))) {
537 : // remove track2 and continue the second loop with the track next to track2
538 0 : fRecTracksPtr->RemoveAt(iTrack2);
539 0 : fNRecTracks--;
540 : } else {
541 : // else remove track1 and continue the first loop with the track next to track1
542 0 : fRecTracksPtr->RemoveAt(iTrack1);
543 0 : fNRecTracks--;
544 0 : break;
545 : }
546 0 : }
547 : } // track2
548 0 : } // track1
549 0 : fRecTracksPtr->Compress(); // this is essential to retrieve the TClonesArray afterwards
550 0 : }
551 :
552 : //__________________________________________________________________________
553 : void AliMUONVTrackReconstructor::RemoveBadTracks()
554 : {
555 : /// Remove tracks for which a problem occured somewhere during the tracking
556 :
557 : AliMUONTrack *track, *nextTrack;
558 : Bool_t trackRemoved = kFALSE;
559 :
560 16 : track = (AliMUONTrack*) fRecTracksPtr->First();
561 48 : while (track) {
562 :
563 16 : nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track);
564 :
565 16 : if (track->GetGlobalChi2() >= AliMUONTrack::MaxChi2()) {
566 0 : AliWarning("problem occured somewhere during the tracking --> discard track");
567 0 : fRecTracksPtr->Remove(track);
568 0 : fNRecTracks--;
569 : trackRemoved = kTRUE;
570 0 : }
571 :
572 : track = nextTrack;
573 :
574 : }
575 :
576 : // compress array of tracks if needed
577 8 : if (trackRemoved) fRecTracksPtr->Compress();
578 :
579 8 : }
580 :
581 : //__________________________________________________________________________
582 : void AliMUONVTrackReconstructor::RemoveConnectedTracks(Int_t stMin, Int_t stMax, Bool_t all)
583 : {
584 : /// Find and remove tracks sharing 1 cluster or more in station(s) [stMin, stMax].
585 : /// For each couple of connected tracks, one removes the one with the smallest
586 : /// number of clusters or with the highest chi2 value in case of equality.
587 : /// If all=kTRUE: both tracks are removed.
588 :
589 : // tag the tracks to be removed
590 16 : TagConnectedTracks(stMin, stMax, all);
591 :
592 : // remove them
593 16 : Int_t nTracks = fRecTracksPtr->GetEntriesFast();
594 100 : for (Int_t i = 0; i < nTracks; i++) {
595 34 : if (((AliMUONTrack*) fRecTracksPtr->UncheckedAt(i))->IsConnected()) {
596 2 : fRecTracksPtr->RemoveAt(i);
597 2 : fNRecTracks--;
598 2 : }
599 : }
600 :
601 : // remove holes in the array if any
602 16 : fRecTracksPtr->Compress();
603 16 : }
604 :
605 : //__________________________________________________________________________
606 : void AliMUONVTrackReconstructor::TagConnectedTracks(Int_t stMin, Int_t stMax, Bool_t all)
607 : {
608 : /// Find and tag tracks sharing 1 cluster or more in station(s) [stMin, stMax].
609 : /// For each couple of connected tracks, one tags the one with the smallest
610 : /// number of clusters or with the highest chi2 value in case of equality.
611 : /// If all=kTRUE: both tracks are tagged.
612 :
613 : AliMUONTrack *track1, *track2;
614 : Int_t nClusters1, nClusters2;
615 24 : Int_t nTracks = fRecTracksPtr->GetEntriesFast();
616 :
617 : // reset the tags
618 148 : for (Int_t i = 0; i < nTracks; i++) ((AliMUONTrack*) fRecTracksPtr->UncheckedAt(i))->Connected(kFALSE);
619 :
620 : // Loop over first track of the pair
621 148 : for (Int_t iTrack1 = 0; iTrack1 < nTracks; iTrack1++) {
622 50 : track1 = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iTrack1);
623 :
624 : // Loop over second track of the pair
625 156 : for (Int_t iTrack2 = iTrack1+1; iTrack2 < nTracks; iTrack2++) {
626 28 : track2 = (AliMUONTrack*) fRecTracksPtr->UncheckedAt(iTrack2);
627 :
628 : // check for connected tracks and tag them
629 28 : if (track1->ClustersInCommon(track2, stMin, stMax) > 0) {
630 :
631 2 : if (all) {
632 :
633 : // tag both tracks
634 0 : track1->Connected();
635 0 : track2->Connected();
636 :
637 0 : } else {
638 :
639 : // tag only the worst track
640 2 : nClusters1 = track1->GetNClusters();
641 2 : nClusters2 = track2->GetNClusters();
642 4 : if ((nClusters1 > nClusters2) || ((nClusters1 == nClusters2) && (track1->GetGlobalChi2() <= track2->GetGlobalChi2()))) {
643 0 : track2->Connected();
644 0 : } else {
645 2 : track1->Connected();
646 : }
647 :
648 : }
649 :
650 : }
651 :
652 : }
653 :
654 : }
655 :
656 24 : }
657 :
658 : //__________________________________________________________________________
659 : void AliMUONVTrackReconstructor::AskForNewClustersInChamber(const AliMUONTrackParam &trackParam,
660 : AliMUONVClusterStore& clusterStore, Int_t chamber)
661 : {
662 : /// Ask the clustering to reconstruct new clusters around the track candidate position
663 :
664 : // check if the current chamber is useable
665 0 : if (!fClusterServer || !GetRecoParam()->UseChamber(chamber)) return;
666 :
667 : // maximum distance between the center of the chamber and a detection element
668 : // (accounting for the inclination of the chamber)
669 : static const Double_t kMaxDZ = 15.; // 15 cm
670 :
671 : // extrapolate track parameters to the chamber
672 0 : AliMUONTrackParam extrapTrackParam(trackParam);
673 0 : if (!AliMUONTrackExtrap::ExtrapToZCov(&extrapTrackParam, AliMUONConstants::DefaultChamberZ(chamber))) return;
674 :
675 : // build the searching area using the track and chamber resolutions and the maximum-distance-to-track value
676 0 : const TMatrixD& kParamCov = extrapTrackParam.GetCovariances();
677 0 : Double_t errX2 = kParamCov(0,0) + kMaxDZ * kMaxDZ * kParamCov(1,1) + 2. * kMaxDZ * TMath::Abs(kParamCov(0,1)) +
678 0 : GetRecoParam()->GetDefaultNonBendingReso(chamber) * GetRecoParam()->GetDefaultNonBendingReso(chamber);
679 0 : Double_t errY2 = kParamCov(2,2) + kMaxDZ * kMaxDZ * kParamCov(3,3) + 2. * kMaxDZ * TMath::Abs(kParamCov(2,3)) +
680 0 : GetRecoParam()->GetDefaultBendingReso(chamber) * GetRecoParam()->GetDefaultBendingReso(chamber);
681 0 : Double_t dX = TMath::Abs(trackParam.GetNonBendingSlope()) * kMaxDZ +
682 0 : GetRecoParam()->GetMaxNonBendingDistanceToTrack() +
683 0 : GetRecoParam()->GetSigmaCutForTracking() * TMath::Sqrt(2. * errX2);
684 0 : Double_t dY = TMath::Abs(trackParam.GetBendingSlope()) * kMaxDZ +
685 0 : GetRecoParam()->GetMaxBendingDistanceToTrack() +
686 0 : GetRecoParam()->GetSigmaCutForTracking() * TMath::Sqrt(2. * errY2);
687 0 : AliMpArea area(extrapTrackParam.GetNonBendingCoor(),
688 0 : extrapTrackParam.GetBendingCoor(),
689 : dX, dY);
690 :
691 : // ask to cluterize in the given area of the given chamber
692 0 : fClusterServer->Clusterize(chamber, clusterStore, area, GetRecoParam());
693 :
694 0 : }
695 :
696 : //__________________________________________________________________________
697 : void AliMUONVTrackReconstructor::AskForNewClustersInStation(const AliMUONTrackParam &trackParam,
698 : AliMUONVClusterStore& clusterStore, Int_t station)
699 : {
700 : /// Ask the clustering to reconstruct new clusters around the track candidate position
701 : /// in the 2 chambers of the given station
702 0 : AskForNewClustersInChamber(trackParam, clusterStore, 2*station+1);
703 0 : AskForNewClustersInChamber(trackParam, clusterStore, 2*station);
704 0 : }
705 :
706 : //__________________________________________________________________________
707 : Double_t AliMUONVTrackReconstructor::TryOneCluster(const AliMUONTrackParam &trackParam, AliMUONVCluster* cluster,
708 : AliMUONTrackParam &trackParamAtCluster, Bool_t updatePropagator)
709 : {
710 : /// Test the compatibility between the track and the cluster (using trackParam's covariance matrix):
711 : /// return the corresponding Chi2
712 : /// return trackParamAtCluster
713 :
714 : // extrapolate track parameters and covariances at the z position of the tested cluster
715 : // and set pointer to cluster into trackParamAtCluster
716 198 : trackParamAtCluster = trackParam;
717 198 : trackParamAtCluster.SetClusterPtr(cluster);
718 198 : if (!AliMUONTrackExtrap::ExtrapToZCov(&trackParamAtCluster, cluster->GetZ(), updatePropagator))
719 0 : return 2.*AliMUONTrack::MaxChi2();
720 :
721 : // Set differences between trackParam and cluster in the bending and non bending directions
722 198 : Double_t dX = cluster->GetX() - trackParamAtCluster.GetNonBendingCoor();
723 198 : Double_t dY = cluster->GetY() - trackParamAtCluster.GetBendingCoor();
724 :
725 : // Calculate errors and covariances
726 198 : const TMatrixD& kParamCov = trackParamAtCluster.GetCovariances();
727 198 : Double_t sigmaX2 = kParamCov(0,0) + cluster->GetErrX2();
728 198 : Double_t sigmaY2 = kParamCov(2,2) + cluster->GetErrY2();
729 198 : Double_t covXY = kParamCov(0,2);
730 198 : Double_t det = sigmaX2 * sigmaY2 - covXY * covXY;
731 :
732 : // Compute chi2
733 198 : if (det == 0.) return 2.*AliMUONTrack::MaxChi2();
734 198 : return (dX * dX * sigmaY2 + dY * dY * sigmaX2 - 2. * dX * dY * covXY) / det;
735 :
736 198 : }
737 :
738 : //__________________________________________________________________________
739 : Bool_t AliMUONVTrackReconstructor::TryOneClusterFast(const AliMUONTrackParam &trackParam, const AliMUONVCluster* cluster)
740 : {
741 : /// Test the compatibility between the track and the cluster
742 : /// given the track and cluster resolutions + the maximum-distance-to-track value
743 : /// and assuming linear propagation of the track:
744 : /// return kTRUE if they are compatibles
745 :
746 1288 : Double_t dZ = cluster->GetZ() - trackParam.GetZ();
747 644 : Double_t dX = cluster->GetX() - (trackParam.GetNonBendingCoor() + trackParam.GetNonBendingSlope() * dZ);
748 644 : Double_t dY = cluster->GetY() - (trackParam.GetBendingCoor() + trackParam.GetBendingSlope() * dZ);
749 644 : const TMatrixD& kParamCov = trackParam.GetCovariances();
750 644 : Double_t errX2 = kParamCov(0,0) + dZ * dZ * kParamCov(1,1) + 2. * dZ * kParamCov(0,1) + cluster->GetErrX2();
751 644 : Double_t errY2 = kParamCov(2,2) + dZ * dZ * kParamCov(3,3) + 2. * dZ * kParamCov(2,3) + cluster->GetErrY2();
752 :
753 1288 : Double_t dXmax = GetRecoParam()->GetSigmaCutForTracking() * TMath::Sqrt(2. * errX2) +
754 644 : GetRecoParam()->GetMaxNonBendingDistanceToTrack();
755 1288 : Double_t dYmax = GetRecoParam()->GetSigmaCutForTracking() * TMath::Sqrt(2. * errY2) +
756 644 : GetRecoParam()->GetMaxBendingDistanceToTrack();
757 :
758 1288 : if (TMath::Abs(dX) > dXmax || TMath::Abs(dY) > dYmax) return kFALSE;
759 :
760 198 : return kTRUE;
761 :
762 644 : }
763 :
764 : //__________________________________________________________________________
765 : Double_t AliMUONVTrackReconstructor::TryTwoClustersFast(const AliMUONTrackParam &trackParamAtCluster1, AliMUONVCluster* cluster2,
766 : AliMUONTrackParam &trackParamAtCluster2)
767 : {
768 : /// Test the compatibility between the track and the 2 clusters together (using trackParam's covariance matrix)
769 : /// assuming linear propagation between the two clusters:
770 : /// return the corresponding Chi2 accounting for covariances between the 2 clusters
771 : /// return trackParamAtCluster2
772 :
773 : // extrapolate linearly track parameters and covariances at the z position of the second cluster
774 0 : trackParamAtCluster2 = trackParamAtCluster1;
775 0 : AliMUONTrackExtrap::LinearExtrapToZCov(&trackParamAtCluster2, cluster2->GetZ());
776 :
777 : // set pointer to cluster2 into trackParamAtCluster2
778 0 : trackParamAtCluster2.SetClusterPtr(cluster2);
779 :
780 : // Set differences between track and clusters in the bending and non bending directions
781 0 : AliMUONVCluster* cluster1 = trackParamAtCluster1.GetClusterPtr();
782 0 : Double_t dX1 = cluster1->GetX() - trackParamAtCluster1.GetNonBendingCoor();
783 0 : Double_t dX2 = cluster2->GetX() - trackParamAtCluster2.GetNonBendingCoor();
784 0 : Double_t dY1 = cluster1->GetY() - trackParamAtCluster1.GetBendingCoor();
785 0 : Double_t dY2 = cluster2->GetY() - trackParamAtCluster2.GetBendingCoor();
786 :
787 : // Calculate errors and covariances
788 0 : const TMatrixD& kParamCov1 = trackParamAtCluster1.GetCovariances();
789 0 : const TMatrixD& kParamCov2 = trackParamAtCluster2.GetCovariances();
790 0 : Double_t dZ = trackParamAtCluster2.GetZ() - trackParamAtCluster1.GetZ();
791 0 : Double_t sigma2X1 = kParamCov1(0,0) + cluster1->GetErrX2();
792 0 : Double_t sigma2X2 = kParamCov2(0,0) + cluster2->GetErrX2();
793 0 : Double_t covX1X2 = kParamCov1(0,0) + dZ * kParamCov1(0,1);
794 0 : Double_t sigma2Y1 = kParamCov1(2,2) + cluster1->GetErrY2();
795 0 : Double_t sigma2Y2 = kParamCov2(2,2) + cluster2->GetErrY2();
796 0 : Double_t covY1Y2 = kParamCov1(2,2) + dZ * kParamCov1(2,3);
797 :
798 : // Compute chi2
799 0 : Double_t detX = sigma2X1 * sigma2X2 - covX1X2 * covX1X2;
800 0 : Double_t detY = sigma2Y1 * sigma2Y2 - covY1Y2 * covY1Y2;
801 0 : if (detX == 0. || detY == 0.) return 2.*AliMUONTrack::MaxChi2();
802 0 : return (dX1 * dX1 * sigma2X2 + dX2 * dX2 * sigma2X1 - 2. * dX1 * dX2 * covX1X2) / detX
803 0 : + (dY1 * dY1 * sigma2Y2 + dY2 * dY2 * sigma2Y1 - 2. * dY1 * dY2 * covY1Y2) / detY;
804 :
805 0 : }
806 :
807 : //__________________________________________________________________________
808 : Bool_t AliMUONVTrackReconstructor::FollowLinearTrackInChamber(AliMUONTrack &trackCandidate, const AliMUONVClusterStore& clusterStore,
809 : Int_t nextChamber)
810 : {
811 : /// Follow trackCandidate in chamber(0..) nextChamber assuming linear propagation, and search for compatible cluster(s)
812 : /// Keep all possibilities or only the best one(s) according to the flag fgkTrackAllTracks:
813 : /// kTRUE: duplicate "trackCandidate" if there are several possibilities and add the new tracks at the end of
814 : /// fRecTracksPtr to avoid conficts with other track candidates at this current stage of the tracking procedure.
815 : /// Remove the obsolete "trackCandidate" at the end.
816 : /// kFALSE: add only the best cluster(s) to the "trackCandidate". Try to add a couple of clusters in priority.
817 : /// return kTRUE if new cluster(s) have been found (otherwise return kFALSE)
818 0 : AliDebug(1,Form("Enter FollowLinearTrackInChamber(1..) %d", nextChamber+1));
819 :
820 0 : Double_t chi2WithOneCluster = AliMUONTrack::MaxChi2();
821 0 : Double_t maxChi2WithOneCluster = 2. * GetRecoParam()->GetSigmaCutForTracking() *
822 0 : GetRecoParam()->GetSigmaCutForTracking(); // 2 because 2 quantities in chi2
823 : Double_t bestChi2WithOneCluster = maxChi2WithOneCluster;
824 : Bool_t foundOneCluster = kFALSE;
825 : AliMUONTrack *newTrack = 0x0;
826 : AliMUONVCluster *cluster;
827 0 : AliMUONTrackParam trackParam;
828 0 : AliMUONTrackParam extrapTrackParamAtCluster;
829 0 : AliMUONTrackParam bestTrackParamAtCluster;
830 :
831 : // Get track parameters according to the propagation direction
832 0 : if (nextChamber > 7) trackParam = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->Last();
833 0 : else trackParam = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->First();
834 :
835 : // Printout for debuging
836 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 2) || (AliLog::GetGlobalDebugLevel() >= 2)) {
837 0 : cout<<endl<<"Track parameters and covariances at first cluster:"<<endl;
838 0 : trackParam.GetParameters().Print();
839 0 : trackParam.GetCovariances().Print();
840 : }
841 :
842 : // Add MCS effect
843 0 : AliMUONTrackExtrap::AddMCSEffect(&trackParam,AliMUONConstants::ChamberThicknessInX0(trackParam.GetClusterPtr()->GetChamberId()),-1.);
844 :
845 : // Printout for debuging
846 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
847 0 : cout << "FollowLinearTrackInChamber: look for cluster in chamber(1..): " << nextChamber+1 << endl;
848 : }
849 :
850 : // Create iterators to loop over clusters in chamber
851 0 : TIter next(clusterStore.CreateChamberIterator(nextChamber,nextChamber));
852 :
853 : // look for candidates in chamber
854 0 : while ( ( cluster = static_cast<AliMUONVCluster*>(next()) ) ) {
855 :
856 : // try to add the current cluster fast
857 0 : if (!TryOneClusterFast(trackParam, cluster)) continue;
858 :
859 : // try to add the current cluster accuratly
860 0 : extrapTrackParamAtCluster = trackParam;
861 0 : AliMUONTrackExtrap::LinearExtrapToZCov(&extrapTrackParamAtCluster, cluster->GetZ());
862 0 : chi2WithOneCluster = TryOneCluster(extrapTrackParamAtCluster, cluster, extrapTrackParamAtCluster);
863 :
864 : // if good chi2 then consider to add cluster
865 0 : if (chi2WithOneCluster < maxChi2WithOneCluster) {
866 : foundOneCluster = kTRUE;
867 :
868 : // Printout for debuging
869 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
870 0 : cout << "FollowLinearTrackInChamber: found one cluster in chamber(1..): " << nextChamber+1
871 0 : << " (Chi2 = " << chi2WithOneCluster << ")" << endl;
872 0 : cluster->Print();
873 : }
874 :
875 0 : if (GetRecoParam()->TrackAllTracks()) {
876 : // copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new cluster
877 0 : newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
878 0 : if (GetRecoParam()->RequestStation(nextChamber/2))
879 0 : extrapTrackParamAtCluster.SetRemovable(kFALSE);
880 0 : else extrapTrackParamAtCluster.SetRemovable(kTRUE);
881 0 : newTrack->AddTrackParamAtCluster(extrapTrackParamAtCluster,*cluster);
882 0 : newTrack->SetGlobalChi2(trackCandidate.GetGlobalChi2()+chi2WithOneCluster);
883 0 : fNRecTracks++;
884 :
885 : // Printout for debuging
886 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
887 0 : cout << "FollowLinearTrackInChamber: added one cluster in chamber(1..): " << nextChamber+1 << endl;
888 0 : if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
889 : }
890 :
891 0 : } else if (chi2WithOneCluster < bestChi2WithOneCluster) {
892 : // keep track of the best cluster
893 : bestChi2WithOneCluster = chi2WithOneCluster;
894 0 : bestTrackParamAtCluster = extrapTrackParamAtCluster;
895 : }
896 :
897 : }
898 :
899 : }
900 :
901 : // fill out the best track if required else clean up the fRecTracksPtr array
902 0 : if (!GetRecoParam()->TrackAllTracks()) {
903 0 : if (foundOneCluster) {
904 0 : if (GetRecoParam()->RequestStation(nextChamber/2))
905 0 : bestTrackParamAtCluster.SetRemovable(kFALSE);
906 0 : else bestTrackParamAtCluster.SetRemovable(kTRUE);
907 0 : trackCandidate.AddTrackParamAtCluster(bestTrackParamAtCluster,*(bestTrackParamAtCluster.GetClusterPtr()));
908 0 : trackCandidate.SetGlobalChi2(trackCandidate.GetGlobalChi2()+bestChi2WithOneCluster);
909 :
910 : // Printout for debuging
911 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
912 0 : cout << "FollowLinearTrackInChamber: added the best cluster in chamber(1..): " << bestTrackParamAtCluster.GetClusterPtr()->GetChamberId()+1 << endl;
913 0 : if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
914 : }
915 :
916 0 : } else return kFALSE;
917 :
918 0 : } else if (foundOneCluster) {
919 :
920 : // remove obsolete track
921 0 : fRecTracksPtr->Remove(&trackCandidate);
922 0 : fNRecTracks--;
923 :
924 0 : } else return kFALSE;
925 :
926 0 : return kTRUE;
927 :
928 0 : }
929 :
930 : //__________________________________________________________________________
931 : Bool_t AliMUONVTrackReconstructor::FollowLinearTrackInStation(AliMUONTrack &trackCandidate, const AliMUONVClusterStore& clusterStore,
932 : Int_t nextStation)
933 : {
934 : /// Follow trackCandidate in station(0..) nextStation assuming linear propagation, and search for compatible cluster(s)
935 : /// Keep all possibilities or only the best one(s) according to the flag fgkTrackAllTracks:
936 : /// kTRUE: duplicate "trackCandidate" if there are several possibilities and add the new tracks at the end of
937 : /// fRecTracksPtr to avoid conficts with other track candidates at this current stage of the tracking procedure.
938 : /// Remove the obsolete "trackCandidate" at the end.
939 : /// kFALSE: add only the best cluster(s) to the "trackCandidate". Try to add a couple of clusters in priority.
940 : /// return kTRUE if new cluster(s) have been found (otherwise return kFALSE)
941 0 : AliDebug(1,Form("Enter FollowLinearTrackInStation(1..) %d", nextStation+1));
942 :
943 : // Order the chamber according to the propagation direction (tracking starts with chamber 2):
944 : // - nextStation == station(1...) 5 => forward propagation
945 : // - nextStation < station(1...) 5 => backward propagation
946 : Int_t ch1, ch2;
947 0 : if (nextStation==4) {
948 0 : ch1 = 2*nextStation+1;
949 : ch2 = 2*nextStation;
950 0 : } else {
951 : ch1 = 2*nextStation;
952 0 : ch2 = 2*nextStation+1;
953 : }
954 :
955 0 : Double_t chi2WithOneCluster = AliMUONTrack::MaxChi2();
956 0 : Double_t chi2WithTwoClusters = AliMUONTrack::MaxChi2();
957 0 : Double_t maxChi2WithOneCluster = 2. * GetRecoParam()->GetSigmaCutForTracking() *
958 0 : GetRecoParam()->GetSigmaCutForTracking(); // 2 because 2 quantities in chi2
959 0 : Double_t maxChi2WithTwoClusters = 4. * GetRecoParam()->GetSigmaCutForTracking() *
960 0 : GetRecoParam()->GetSigmaCutForTracking(); // 4 because 4 quantities in chi2
961 : Double_t bestChi2WithOneCluster = maxChi2WithOneCluster;
962 : Double_t bestChi2WithTwoClusters = maxChi2WithTwoClusters;
963 : Bool_t foundOneCluster = kFALSE;
964 : Bool_t foundTwoClusters = kFALSE;
965 : AliMUONTrack *newTrack = 0x0;
966 : AliMUONVCluster *clusterCh1, *clusterCh2;
967 0 : AliMUONTrackParam trackParam;
968 0 : AliMUONTrackParam extrapTrackParamAtCluster1;
969 0 : AliMUONTrackParam extrapTrackParamAtCluster2;
970 0 : AliMUONTrackParam bestTrackParamAtCluster1;
971 0 : AliMUONTrackParam bestTrackParamAtCluster2;
972 :
973 0 : Int_t nClusters = clusterStore.GetSize();
974 0 : Bool_t *clusterCh1Used = new Bool_t[nClusters];
975 0 : for (Int_t i = 0; i < nClusters; i++) clusterCh1Used[i] = kFALSE;
976 : Int_t iCluster1;
977 :
978 : // Get track parameters according to the propagation direction
979 0 : if (nextStation==4) trackParam = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->Last();
980 0 : else trackParam = *(AliMUONTrackParam*)trackCandidate.GetTrackParamAtCluster()->First();
981 :
982 : // Printout for debuging
983 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 2) || (AliLog::GetGlobalDebugLevel() >= 2)) {
984 0 : cout<<endl<<"Track parameters and covariances at first cluster:"<<endl;
985 0 : trackParam.GetParameters().Print();
986 0 : trackParam.GetCovariances().Print();
987 : }
988 :
989 : // Add MCS effect
990 0 : AliMUONTrackExtrap::AddMCSEffect(&trackParam,AliMUONConstants::ChamberThicknessInX0(trackParam.GetClusterPtr()->GetChamberId()),-1.);
991 :
992 : // Printout for debuging
993 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
994 0 : cout << "FollowLinearTrackInStation: look for clusters in chamber(1..): " << ch2+1 << endl;
995 : }
996 :
997 : // Create iterators to loop over clusters in both chambers
998 0 : TIter nextInCh1(clusterStore.CreateChamberIterator(ch1,ch1));
999 0 : TIter nextInCh2(clusterStore.CreateChamberIterator(ch2,ch2));
1000 :
1001 : // look for candidates in chamber 2
1002 0 : while ( ( clusterCh2 = static_cast<AliMUONVCluster*>(nextInCh2()) ) ) {
1003 :
1004 : // try to add the current cluster fast
1005 0 : if (!TryOneClusterFast(trackParam, clusterCh2)) continue;
1006 :
1007 : // try to add the current cluster accuratly
1008 0 : extrapTrackParamAtCluster2 = trackParam;
1009 0 : AliMUONTrackExtrap::LinearExtrapToZCov(&extrapTrackParamAtCluster2, clusterCh2->GetZ());
1010 0 : chi2WithOneCluster = TryOneCluster(extrapTrackParamAtCluster2, clusterCh2, extrapTrackParamAtCluster2);
1011 :
1012 : // if good chi2 then try to attach a cluster in the other chamber too
1013 0 : if (chi2WithOneCluster < maxChi2WithOneCluster) {
1014 : Bool_t foundSecondCluster = kFALSE;
1015 :
1016 : // Printout for debuging
1017 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
1018 0 : cout << "FollowLinearTrackInStation: found one cluster in chamber(1..): " << ch2+1
1019 0 : << " (Chi2 = " << chi2WithOneCluster << ")" << endl;
1020 0 : clusterCh2->Print();
1021 0 : cout << " look for second clusters in chamber(1..): " << ch1+1 << " ..." << endl;
1022 : }
1023 :
1024 : // add MCS effect
1025 0 : AliMUONTrackExtrap::AddMCSEffect(&extrapTrackParamAtCluster2,AliMUONConstants::ChamberThicknessInX0(ch2),-1.);
1026 :
1027 : // reset cluster iterator of chamber 1
1028 0 : nextInCh1.Reset();
1029 : iCluster1 = -1;
1030 :
1031 : // look for second candidates in chamber 1
1032 0 : while ( ( clusterCh1 = static_cast<AliMUONVCluster*>(nextInCh1()) ) ) {
1033 0 : iCluster1++;
1034 :
1035 : // try to add the current cluster fast
1036 0 : if (!TryOneClusterFast(extrapTrackParamAtCluster2, clusterCh1)) continue;
1037 :
1038 : // try to add the current cluster in addition to the one found in the previous chamber
1039 0 : chi2WithTwoClusters = TryTwoClustersFast(extrapTrackParamAtCluster2, clusterCh1, extrapTrackParamAtCluster1);
1040 :
1041 : // if good chi2 then consider to add the 2 clusters to the "trackCandidate"
1042 0 : if (chi2WithTwoClusters < maxChi2WithTwoClusters) {
1043 : foundSecondCluster = kTRUE;
1044 : foundTwoClusters = kTRUE;
1045 :
1046 : // Printout for debuging
1047 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
1048 0 : cout << "FollowLinearTrackInStation: found one cluster in chamber(1..): " << ch1+1
1049 0 : << " (Chi2 = " << chi2WithTwoClusters << ")" << endl;
1050 0 : clusterCh1->Print();
1051 : }
1052 :
1053 0 : if (GetRecoParam()->TrackAllTracks()) {
1054 : // copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new clusters
1055 0 : newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
1056 0 : extrapTrackParamAtCluster1.SetRemovable(kTRUE);
1057 0 : newTrack->AddTrackParamAtCluster(extrapTrackParamAtCluster1,*clusterCh1);
1058 0 : extrapTrackParamAtCluster2.SetRemovable(kTRUE);
1059 0 : newTrack->AddTrackParamAtCluster(extrapTrackParamAtCluster2,*clusterCh2);
1060 0 : newTrack->SetGlobalChi2(newTrack->GetGlobalChi2()+chi2WithTwoClusters);
1061 0 : fNRecTracks++;
1062 :
1063 : // Tag clusterCh1 as used
1064 0 : clusterCh1Used[iCluster1] = kTRUE;
1065 :
1066 : // Printout for debuging
1067 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
1068 0 : cout << "FollowLinearTrackInStation: added two clusters in station(1..): " << nextStation+1 << endl;
1069 0 : if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
1070 : }
1071 :
1072 0 : } else if (chi2WithTwoClusters < bestChi2WithTwoClusters) {
1073 : // keep track of the best couple of clusters
1074 : bestChi2WithTwoClusters = chi2WithTwoClusters;
1075 0 : bestTrackParamAtCluster1 = extrapTrackParamAtCluster1;
1076 0 : bestTrackParamAtCluster2 = extrapTrackParamAtCluster2;
1077 : }
1078 :
1079 : }
1080 :
1081 : }
1082 :
1083 : // if no cluster found in chamber1 then consider to add clusterCh2 only
1084 0 : if (!foundSecondCluster) {
1085 : foundOneCluster = kTRUE;
1086 :
1087 0 : if (GetRecoParam()->TrackAllTracks()) {
1088 : // copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new cluster
1089 0 : newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
1090 0 : if (GetRecoParam()->RequestStation(nextStation))
1091 0 : extrapTrackParamAtCluster2.SetRemovable(kFALSE);
1092 0 : else extrapTrackParamAtCluster2.SetRemovable(kTRUE);
1093 0 : newTrack->AddTrackParamAtCluster(extrapTrackParamAtCluster2,*clusterCh2);
1094 0 : newTrack->SetGlobalChi2(newTrack->GetGlobalChi2()+chi2WithOneCluster);
1095 0 : fNRecTracks++;
1096 :
1097 : // Printout for debuging
1098 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
1099 0 : cout << "FollowLinearTrackInStation: added one cluster in chamber(1..): " << ch2+1 << endl;
1100 0 : if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
1101 : }
1102 :
1103 0 : } else if (!foundTwoClusters && chi2WithOneCluster < bestChi2WithOneCluster) {
1104 : // keep track of the best cluster except if a couple of clusters has already been found
1105 : bestChi2WithOneCluster = chi2WithOneCluster;
1106 0 : bestTrackParamAtCluster1 = extrapTrackParamAtCluster2;
1107 : }
1108 :
1109 : }
1110 :
1111 0 : }
1112 :
1113 : }
1114 :
1115 : // look for candidates in chamber 1 not already attached to a track
1116 : // if we want to keep all possible tracks or if no good couple of clusters has been found
1117 0 : if (GetRecoParam()->TrackAllTracks() || !foundTwoClusters) {
1118 :
1119 : // Printout for debuging
1120 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
1121 0 : cout << "FollowLinearTrackInStation: look for single cluster in chamber(1..): " << ch1+1 << endl;
1122 : }
1123 :
1124 : //Extrapolate trackCandidate to chamber "ch2"
1125 0 : AliMUONTrackExtrap::LinearExtrapToZCov(&trackParam, AliMUONConstants::DefaultChamberZ(ch2));
1126 :
1127 : // add MCS effect for next step
1128 0 : AliMUONTrackExtrap::AddMCSEffect(&trackParam,AliMUONConstants::ChamberThicknessInX0(ch2),-1.);
1129 :
1130 : // reset cluster iterator of chamber 1
1131 0 : nextInCh1.Reset();
1132 : iCluster1 = -1;
1133 :
1134 : // look for second candidates in chamber 1
1135 0 : while ( ( clusterCh1 = static_cast<AliMUONVCluster*>(nextInCh1()) ) ) {
1136 0 : iCluster1++;
1137 :
1138 0 : if (clusterCh1Used[iCluster1]) continue; // Skip clusters already used
1139 :
1140 : // try to add the current cluster fast
1141 0 : if (!TryOneClusterFast(trackParam, clusterCh1)) continue;
1142 :
1143 : // try to add the current cluster accuratly
1144 0 : extrapTrackParamAtCluster1 = trackParam;
1145 0 : AliMUONTrackExtrap::LinearExtrapToZCov(&extrapTrackParamAtCluster1, clusterCh1->GetZ());
1146 0 : chi2WithOneCluster = TryOneCluster(extrapTrackParamAtCluster1, clusterCh1, extrapTrackParamAtCluster1);
1147 :
1148 : // if good chi2 then consider to add clusterCh1
1149 : // We do not try to attach a cluster in the other chamber too since it has already been done above
1150 0 : if (chi2WithOneCluster < maxChi2WithOneCluster) {
1151 : foundOneCluster = kTRUE;
1152 :
1153 : // Printout for debuging
1154 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
1155 0 : cout << "FollowLinearTrackInStation: found one cluster in chamber(1..): " << ch1+1
1156 0 : << " (Chi2 = " << chi2WithOneCluster << ")" << endl;
1157 0 : clusterCh1->Print();
1158 : }
1159 :
1160 0 : if (GetRecoParam()->TrackAllTracks()) {
1161 : // copy trackCandidate into a new track put at the end of fRecTracksPtr and add the new cluster
1162 0 : newTrack = new ((*fRecTracksPtr)[fRecTracksPtr->GetLast()+1]) AliMUONTrack(trackCandidate);
1163 0 : if (GetRecoParam()->RequestStation(nextStation))
1164 0 : extrapTrackParamAtCluster1.SetRemovable(kFALSE);
1165 0 : else extrapTrackParamAtCluster1.SetRemovable(kTRUE);
1166 0 : newTrack->AddTrackParamAtCluster(extrapTrackParamAtCluster1,*clusterCh1);
1167 0 : newTrack->SetGlobalChi2(newTrack->GetGlobalChi2()+chi2WithOneCluster);
1168 0 : fNRecTracks++;
1169 :
1170 : // Printout for debuging
1171 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
1172 0 : cout << "FollowLinearTrackInStation: added one cluster in chamber(1..): " << ch1+1 << endl;
1173 0 : if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
1174 : }
1175 :
1176 0 : } else if (chi2WithOneCluster < bestChi2WithOneCluster) {
1177 : // keep track of the best cluster except if a couple of clusters has already been found
1178 : bestChi2WithOneCluster = chi2WithOneCluster;
1179 0 : bestTrackParamAtCluster1 = extrapTrackParamAtCluster1;
1180 : }
1181 :
1182 : }
1183 :
1184 : }
1185 :
1186 : }
1187 :
1188 : // fill out the best track if required else clean up the fRecTracksPtr array
1189 0 : if (!GetRecoParam()->TrackAllTracks()) {
1190 0 : if (foundTwoClusters) {
1191 0 : bestTrackParamAtCluster1.SetRemovable(kTRUE);
1192 0 : trackCandidate.AddTrackParamAtCluster(bestTrackParamAtCluster1,*(bestTrackParamAtCluster1.GetClusterPtr()));
1193 0 : bestTrackParamAtCluster2.SetRemovable(kTRUE);
1194 0 : trackCandidate.AddTrackParamAtCluster(bestTrackParamAtCluster2,*(bestTrackParamAtCluster2.GetClusterPtr()));
1195 0 : trackCandidate.SetGlobalChi2(trackCandidate.GetGlobalChi2()+bestChi2WithTwoClusters);
1196 :
1197 : // Printout for debuging
1198 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
1199 0 : cout << "FollowLinearTrackInStation: added the two best clusters in station(1..): " << nextStation+1 << endl;
1200 0 : if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
1201 : }
1202 :
1203 0 : } else if (foundOneCluster) {
1204 0 : if (GetRecoParam()->RequestStation(nextStation))
1205 0 : bestTrackParamAtCluster1.SetRemovable(kFALSE);
1206 0 : else bestTrackParamAtCluster1.SetRemovable(kTRUE);
1207 0 : trackCandidate.AddTrackParamAtCluster(bestTrackParamAtCluster1,*(bestTrackParamAtCluster1.GetClusterPtr()));
1208 0 : trackCandidate.SetGlobalChi2(trackCandidate.GetGlobalChi2()+bestChi2WithOneCluster);
1209 :
1210 : // Printout for debuging
1211 0 : if ((AliLog::GetDebugLevel("MUON","AliMUONVTrackReconstructor") >= 1) || (AliLog::GetGlobalDebugLevel() >= 1)) {
1212 0 : cout << "FollowLinearTrackInStation: added the best cluster in chamber(1..): " << bestTrackParamAtCluster1.GetClusterPtr()->GetChamberId()+1 << endl;
1213 0 : if (AliLog::GetGlobalDebugLevel() >= 3) newTrack->RecursiveDump();
1214 : }
1215 :
1216 : } else {
1217 0 : delete [] clusterCh1Used;
1218 0 : return kFALSE;
1219 : }
1220 :
1221 0 : } else if (foundOneCluster || foundTwoClusters) {
1222 :
1223 : // remove obsolete track
1224 0 : fRecTracksPtr->Remove(&trackCandidate);
1225 0 : fNRecTracks--;
1226 :
1227 : } else {
1228 0 : delete [] clusterCh1Used;
1229 0 : return kFALSE;
1230 : }
1231 :
1232 0 : delete [] clusterCh1Used;
1233 0 : return kTRUE;
1234 :
1235 0 : }
1236 :
1237 : //__________________________________________________________________________
1238 : void AliMUONVTrackReconstructor::ImproveTracks()
1239 : {
1240 : /// Improve tracks by removing clusters with local chi2 highter than the defined cut
1241 : /// Recompute track parameters and covariances at the remaining clusters
1242 32 : AliDebug(1,"Enter ImproveTracks");
1243 :
1244 : AliMUONTrack *track, *nextTrack;
1245 :
1246 8 : track = (AliMUONTrack*) fRecTracksPtr->First();
1247 78 : while (track) {
1248 :
1249 : // prepare next track in case the actual track is suppressed
1250 44 : nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track);
1251 :
1252 18 : ImproveTrack(*track);
1253 :
1254 : // remove track if improvement failed
1255 18 : if (!track->IsImproved()) {
1256 0 : fRecTracksPtr->Remove(track);
1257 0 : fNRecTracks--;
1258 0 : }
1259 :
1260 : track = nextTrack;
1261 : }
1262 :
1263 : // compress the array in case of some tracks have been removed
1264 8 : fRecTracksPtr->Compress();
1265 :
1266 8 : }
1267 :
1268 : //__________________________________________________________________________
1269 : void AliMUONVTrackReconstructor::Finalize()
1270 : {
1271 : /// Recompute track parameters and covariances at each attached cluster
1272 : /// Set the label pointing to the corresponding MC track
1273 : /// Remove the track if finalization failed
1274 :
1275 : AliMUONTrack *track, *nextTrack;
1276 : Bool_t trackRemoved = kFALSE;
1277 :
1278 16 : track = (AliMUONTrack*) fRecTracksPtr->First();
1279 48 : while (track) {
1280 :
1281 16 : nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track);
1282 :
1283 32 : if (FinalizeTrack(*track)) track->FindMCLabel();
1284 : else {
1285 0 : fRecTracksPtr->Remove(track);
1286 0 : fNRecTracks--;
1287 : trackRemoved = kTRUE;
1288 : }
1289 :
1290 : track = nextTrack;
1291 :
1292 : }
1293 :
1294 : // compress array of tracks if needed
1295 8 : if (trackRemoved) fRecTracksPtr->Compress();
1296 :
1297 8 : }
1298 :
1299 : //__________________________________________________________________________
1300 : void AliMUONVTrackReconstructor::DiscardMonoCathodClusters()
1301 : {
1302 : /// Assign a different resolution to the mono-cathod clusters
1303 : /// in the direction of the missing plane and refit the track
1304 : /// Remove the track in case of failure
1305 :
1306 0 : if (!fkTransformer) AliFatal("missing geometry transformer");
1307 :
1308 : AliMUONTrack *track, *nextTrack;
1309 : Bool_t trackRemoved = kFALSE;
1310 :
1311 0 : track = (AliMUONTrack*) fRecTracksPtr->First();
1312 0 : while (track) {
1313 :
1314 0 : nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track);
1315 :
1316 0 : ChangeMonoCathodClusterRes(*track);
1317 :
1318 0 : if (!RefitTrack(*track) || (GetRecoParam()->ImproveTracks() && !track->IsImproved())) {
1319 0 : fRecTracksPtr->Remove(track);
1320 0 : fNRecTracks--;
1321 : trackRemoved = kTRUE;
1322 0 : }
1323 :
1324 : track = nextTrack;
1325 :
1326 : }
1327 :
1328 : // compress array of tracks if needed
1329 0 : if (trackRemoved) fRecTracksPtr->Compress();
1330 :
1331 0 : }
1332 :
1333 : //__________________________________________________________________________
1334 : void AliMUONVTrackReconstructor::ChangeMonoCathodClusterRes(AliMUONTrack &track)
1335 : {
1336 : /// Assign a different resolution to the mono-cathod clusters
1337 : /// in the direction of the missing plane and refit the track
1338 :
1339 : // Loop over clusters
1340 : AliMUONVCluster *cluster;
1341 0 : Int_t nClusters = track.GetNClusters();
1342 0 : for (Int_t iCluster = 0; iCluster < nClusters; iCluster++) {
1343 0 : cluster = ((AliMUONTrackParam*) track.GetTrackParamAtCluster()->UncheckedAt(iCluster))->GetClusterPtr();
1344 :
1345 : // do it only for stations 3, 4 & 5
1346 0 : if (cluster->GetChamberId() < 4) continue;
1347 :
1348 : // get the cathod corresponding to the bending/non-bending plane
1349 0 : Int_t deId = cluster->GetDetElemId();
1350 0 : AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(deId, kFALSE);
1351 0 : if (!de) continue;
1352 0 : AliMp::CathodType cath1 = de->GetCathodType(AliMp::kBendingPlane);
1353 0 : AliMp::CathodType cath2 = de->GetCathodType(AliMp::kNonBendingPlane);
1354 :
1355 : // get the corresponding segmentation
1356 0 : const AliMpVSegmentation* seg1 = AliMpSegmentation::Instance()->GetMpSegmentation(deId, cath1);
1357 0 : const AliMpVSegmentation* seg2 = AliMpSegmentation::Instance()->GetMpSegmentation(deId, cath2);
1358 0 : if (!seg1 || !seg2) continue;
1359 :
1360 : // get local coordinate of the cluster
1361 0 : Double_t lX,lY,lZ;
1362 0 : Double_t gX = cluster->GetX();
1363 0 : Double_t gY = cluster->GetY();
1364 0 : Double_t gZ = cluster->GetZ();
1365 0 : fkTransformer->Global2Local(deId,gX,gY,gZ,lX,lY,lZ);
1366 :
1367 : // find pads below the cluster
1368 0 : AliMpPad pad1 = seg1->PadByPosition(lX, lY, kFALSE);
1369 0 : AliMpPad pad2 = seg2->PadByPosition(lX, lY, kFALSE);
1370 :
1371 : // build their ID if pads are valid
1372 0 : UInt_t padId1 = (pad1.IsValid()) ? AliMUONVDigit::BuildUniqueID(deId, pad1.GetManuId(), pad1.GetManuChannel(), cath1) : 0;
1373 0 : UInt_t padId2 = (pad2.IsValid()) ? AliMUONVDigit::BuildUniqueID(deId, pad2.GetManuId(), pad2.GetManuChannel(), cath2) : 0;
1374 :
1375 : // check if the cluster contains these pads
1376 : Bool_t hasNonBending = kFALSE;
1377 : Bool_t hasBending = kFALSE;
1378 0 : for (Int_t iDigit = 0; iDigit < cluster->GetNDigits(); iDigit++) {
1379 :
1380 0 : UInt_t digitId = cluster->GetDigitId(iDigit);
1381 :
1382 0 : if (digitId == padId1) {
1383 :
1384 : hasBending = kTRUE;
1385 0 : if (hasNonBending) break;
1386 :
1387 0 : } else if (digitId == padId2) {
1388 :
1389 : hasNonBending = kTRUE;
1390 0 : if (hasBending) break;
1391 :
1392 : }
1393 :
1394 0 : }
1395 :
1396 : // modify the cluster resolution if needed
1397 0 : if (!hasNonBending) cluster->SetErrXY(GetRecoParam()->GetMonoCathodClNonBendingRes(), cluster->GetErrY());
1398 0 : if (!hasBending) cluster->SetErrXY(cluster->GetErrX(), GetRecoParam()->GetMonoCathodClBendingRes());
1399 :
1400 0 : }
1401 :
1402 0 : }
1403 :
1404 : //__________________________________________________________________________
1405 : void AliMUONVTrackReconstructor::ValidateTracksWithTrigger(AliMUONVTrackStore& trackStore,
1406 : const AliMUONVTriggerTrackStore& triggerTrackStore,
1407 : const AliMUONVTriggerStore& triggerStore,
1408 : const AliMUONTrackHitPattern& trackHitPattern)
1409 : {
1410 : /// Try to match track from tracking system with trigger track
1411 16 : AliCodeTimerAuto("",0);
1412 :
1413 8 : trackHitPattern.ExecuteValidation(trackStore, triggerTrackStore, triggerStore);
1414 8 : }
1415 :
1416 :
1417 : //__________________________________________________________________________
1418 : void AliMUONVTrackReconstructor::EventReconstructTrigger(const AliMUONTriggerCircuit& circuit,
1419 : const AliMUONVTriggerStore& triggerStore,
1420 : AliMUONVTriggerTrackStore& triggerTrackStore)
1421 : {
1422 : /// Fill trigger track store from local trigger
1423 32 : AliDebug(1, "");
1424 8 : AliCodeTimerAuto("",0);
1425 :
1426 8 : AliMUONGlobalTrigger* globalTrigger = triggerStore.Global();
1427 :
1428 : UChar_t gloTrigPat = 0;
1429 :
1430 8 : if (globalTrigger)
1431 : {
1432 8 : gloTrigPat = globalTrigger->GetGlobalResponse();
1433 8 : }
1434 :
1435 8 : AliMUONTriggerTrack triggerTrack;
1436 :
1437 16 : TIter next(triggerStore.CreateIterator());
1438 : AliMUONLocalTrigger* locTrg(0x0);
1439 :
1440 118 : while ( ( locTrg = static_cast<AliMUONLocalTrigger*>(next()) ) )
1441 : {
1442 124 : if ( locTrg->IsTrigX() && locTrg->IsTrigY() )
1443 : { // make Trigger Track if trigger in X and Y
1444 :
1445 28 : if (TriggerToTrack(circuit, *locTrg, triggerTrack, gloTrigPat))
1446 14 : triggerTrackStore.Add(triggerTrack);
1447 :
1448 : } // board is fired
1449 : } // end of loop on Local Trigger
1450 8 : }
1451 :
1452 : //__________________________________________________________________________
1453 : Bool_t AliMUONVTrackReconstructor::TriggerToTrack(const AliMUONTriggerCircuit& circuit,
1454 : const AliMUONLocalTrigger& locTrg,
1455 : AliMUONTriggerTrack& triggerTrack,
1456 : UChar_t globalTriggerPattern)
1457 : {
1458 : /// To make the trigger tracks from Local Trigger
1459 28 : const Double_t kTrigNonBendReso = AliMUONConstants::TriggerNonBendingReso();
1460 14 : const Double_t kTrigBendReso = AliMUONConstants::TriggerBendingReso();
1461 14 : const Double_t kSqrt12 = TMath::Sqrt(12.);
1462 :
1463 14 : TMatrixD trigCov(3,3);
1464 :
1465 14 : Int_t localBoardId = locTrg.LoCircuit();
1466 :
1467 14 : Float_t y11 = circuit.GetY11Pos(localBoardId, locTrg.LoStripX());
1468 14 : Float_t z11 = circuit.GetZ11Pos(localBoardId, locTrg.LoStripX());
1469 : // need first to convert deviation to [0-30]
1470 : // (see AliMUONLocalTriggerBoard::LocalTrigger)
1471 14 : Int_t deviation = locTrg.GetDeviation();
1472 14 : Int_t stripX21 = locTrg.LoStripX()+deviation+1;
1473 14 : Float_t y21 = circuit.GetY21Pos(localBoardId, stripX21);
1474 14 : Float_t z21 = circuit.GetZ21Pos(localBoardId, stripX21);
1475 14 : Float_t x11 = circuit.GetX11Pos(localBoardId, locTrg.LoStripY());
1476 :
1477 70 : AliDebug(1, Form(" MakeTriggerTrack %3d %2d %2d %2d (%f %f %f) (%f %f)\n",locTrg.LoCircuit(),
1478 : locTrg.LoStripX(),locTrg.LoStripX()+deviation+1,locTrg.LoStripY(),x11, y11, z11, y21, z21));
1479 :
1480 14 : if (TMath::Abs(z11) < 0.00001) return kFALSE;
1481 :
1482 14 : Double_t deltaZ = z11 - z21;
1483 :
1484 14 : Float_t slopeX = x11/z11;
1485 14 : Float_t slopeY = (y11-y21) / deltaZ;
1486 :
1487 28 : Float_t sigmaX = circuit.GetX11Width(localBoardId, locTrg.LoStripY()) / kSqrt12;
1488 28 : Float_t sigmaY = circuit.GetY11Width(localBoardId, locTrg.LoStripX()) / kSqrt12;
1489 28 : Float_t sigmaY21 = circuit.GetY21Width(localBoardId, locTrg.LoStripX()) / kSqrt12;
1490 :
1491 14 : trigCov.Zero();
1492 28 : trigCov(0,0) = kTrigNonBendReso * kTrigNonBendReso + sigmaX * sigmaX;
1493 28 : trigCov(1,1) = kTrigBendReso * kTrigBendReso + sigmaY * sigmaY;
1494 28 : trigCov(2,2) =
1495 14 : (2. * kTrigBendReso * kTrigBendReso + sigmaY * sigmaY + sigmaY21 * sigmaY21 ) / deltaZ / deltaZ;
1496 56 : trigCov(1,2) = trigCov(2,1) = trigCov(1,1) / deltaZ;
1497 :
1498 14 : triggerTrack.SetX11(x11);
1499 14 : triggerTrack.SetY11(y11);
1500 14 : triggerTrack.SetZ11(z11);
1501 14 : triggerTrack.SetZ21(z21);
1502 14 : triggerTrack.SetSlopeX(slopeX);
1503 14 : triggerTrack.SetSlopeY(slopeY);
1504 14 : triggerTrack.SetGTPattern(globalTriggerPattern);
1505 14 : triggerTrack.SetLoTrgNum(localBoardId);
1506 14 : triggerTrack.SetCovariances(trigCov);
1507 28 : triggerTrack.SetUniqueID(locTrg.GetUniqueID());
1508 :
1509 : return kTRUE;
1510 :
1511 14 : }
|