]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MFT/AliMFTTrackerMU.cxx
MFT classes modified to include MUON+MFT tracks in ESD&AOD. Clusterization algorithm...
[u/mrichter/AliRoot.git] / MFT / AliMFTTrackerMU.cxx
1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
7 * Permission to use, copy, modify and distribute this software and its   *
8 * documentation strictly for non-commercial purposes is hereby granted   *
9 * without fee, provided that the above copyright notice appears in all   *
10 * copies and that both the copyright notice and this permission notice   *
11 * appear in the supporting documentation. The authors make no claims     *
12 * about the suitability of this software for any purpose. It is          *
13 * provided "as is" without express or implied warranty.                  *
14 **************************************************************************/
15
16 //====================================================================================================================================================
17 //
18 //                          MFT tracker
19 //
20 // Class for the creation of the "global muon tracks" built from the tracks reconstructed in the 
21 // muon spectrometer and the clusters of the Muon Forward Tracker
22 //
23 // Contact author: antonio.uras@cern.ch
24 //
25 //====================================================================================================================================================
26
27 #include "TTree.h"
28 #include "AliLog.h"
29 #include "AliGeomManager.h"
30 #include "AliESDEvent.h"
31 #include "AliESDMuonTrack.h"
32 // #include "AliESDMuonGlobalTrack.h"
33 #include "AliMFTTrackerMU.h"
34 #include "TMath.h"
35 #include "AliRun.h"
36 #include "AliMFT.h"
37 #include "AliMUONTrackExtrap.h"
38 #include "AliMUONTrack.h"
39 #include "AliMUONESDInterface.h"
40 #include "AliMuonForwardTrack.h"
41
42 ClassImp(AliMFTTrackerMU)
43
44 const Double_t AliMFTTrackerMU::fRadLengthSi = AliMFTConstants::fRadLengthSi;
45
46 //====================================================================================================================================================
47
48 AliMFTTrackerMU::AliMFTTrackerMU() : 
49   AliTracker(),
50   fESD(0),
51   fMFT(0),
52   fSegmentation(0),
53   fNPlanesMFT(0),
54   fNPlanesMFTAnalyzed(0),
55   fSigmaClusterCut(2),
56   fScaleSigmaClusterCut(1.),
57   fNMaxMissingMFTClusters(0),
58   fGlobalTrackingDiverged(kFALSE),
59   fCandidateTracks(0),
60   fMUONTrack(0),
61   fCurrentTrack(0),
62   fFinalBestCandidate(0),
63   fVertexErrorX(0.015),
64   fVertexErrorY(0.015),
65   fVertexErrorZ(0.010),
66   fBransonCorrection(kFALSE)
67 {
68
69   //--------------------------------------------------------------------
70   // This is the AliMFTTrackerMU constructor
71   //--------------------------------------------------------------------
72
73   fMFT = (AliMFT*) gAlice->GetDetector("MFT");
74   fSegmentation = fMFT->GetSegmentation();
75   SetNPlanesMFT(fSegmentation->GetNPlanes());
76   AliMUONTrackExtrap::SetField();                 // set the magnetic field for track extrapolations
77
78   AliInfo(Form("fMFT = %p, fSegmentation = %p", fMFT, fSegmentation));
79
80   for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
81     fMFTClusterArray[iPlane]      = new TClonesArray("AliMFTCluster");
82     fMFTClusterArrayFront[iPlane] = new TClonesArray("AliMFTCluster");
83     fMFTClusterArrayBack[iPlane]  = new TClonesArray("AliMFTCluster");
84     fMFTClusterArray[iPlane]      -> SetOwner(kTRUE);
85     fMFTClusterArrayFront[iPlane] -> SetOwner(kTRUE);
86     fMFTClusterArrayBack[iPlane]  -> SetOwner(kTRUE);
87     fMinResearchRadiusAtPlane[iPlane] = 0.;
88   }
89
90   fCandidateTracks = new TClonesArray("AliMuonForwardTrack",50000);
91
92 }
93
94 //====================================================================================================================================================
95
96 AliMFTTrackerMU::~AliMFTTrackerMU() {
97
98   // destructor
99
100   for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
101     fMFTClusterArray[iPlane] -> Delete();
102     delete fMFTClusterArray[iPlane];
103     delete fMFTClusterArrayFront[iPlane];
104     delete fMFTClusterArrayBack[iPlane];
105   }
106
107   delete fCandidateTracks;
108
109 }
110
111 //====================================================================================================================================================
112
113 Int_t AliMFTTrackerMU::LoadClusters(TTree *cTree) {
114
115   //--------------------------------------------------------------------
116   // This function loads the MFT clusters
117   //--------------------------------------------------------------------
118  
119   for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
120     AliDebug(1, Form("Setting Address for Branch Plane_%02d", iPlane)); 
121     cTree->SetBranchAddress(Form("Plane_%02d",iPlane), &fMFTClusterArray[iPlane]);
122   }
123  
124   if (!cTree->GetEvent()) return kFALSE;
125   for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
126     AliInfo(Form("plane %02d: nClusters = %d\n", iPlane, fMFTClusterArray[iPlane]->GetEntries()));
127   }
128   SeparateFrontBackClusters();
129
130   return 0;
131
132 }
133
134 //====================================================================================================================================================
135
136 void AliMFTTrackerMU::UnloadClusters() {
137   
138 //   //--------------------------------------------------------------------
139 //   // This function unloads MFT clusters
140 //   //--------------------------------------------------------------------
141
142   for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
143     fMFTClusterArray[iPlane]      -> Clear("C");
144     fMFTClusterArrayFront[iPlane] -> Clear("C");
145     fMFTClusterArrayBack[iPlane]  -> Clear("C");
146   }
147
148 }
149
150 //====================================================================================================================================================
151
152 Int_t AliMFTTrackerMU::Clusters2Tracks(AliESDEvent *event) {
153
154   //--------------------------------------------------------------------
155   // This functions reconstructs the Muon Forward Tracks
156   // The clusters must be already loaded !
157   //--------------------------------------------------------------------
158
159   // Tree of AliMuonForwardTrack objects. Created outside the ESD framework for cross-check purposes
160
161   TFile *outputFileMuonGlobalTracks = new TFile("MuonGlobalTracks.root", "update");
162   TTree *outputTreeMuonGlobalTracks = new TTree("AliMuonForwardTracks", "Tree of AliMuonForwardTracks");
163   TClonesArray *muonForwardTracks = new TClonesArray("AliMuonForwardTrack");
164   outputTreeMuonGlobalTracks -> Branch("tracks", &muonForwardTracks);
165  
166   //--------------------------------------------------------------------
167
168   fESD = event;
169
170   // get the ITS primary vertex
171
172   Double_t vertex[3] = {0., 0., 0.};
173   const AliESDVertex* esdVert = fESD->GetVertex(); 
174   if (esdVert->GetNContributors() > 0 || !strcmp(esdVert->GetTitle(),"vertexer: smearMC")) {
175     esdVert->GetXYZ(vertex);
176     AliDebug(1,Form("found vertex (%e,%e,%e)",vertex[0],vertex[1],vertex[2]));
177   }
178   
179   //----------- Read ESD MUON tracks -------------------
180
181   Int_t nTracksMUON = event->GetNumberOfMuonTracks();
182
183   AliInfo(Form("Number of ESD MUON tracks: %d\n", nTracksMUON));
184
185   Int_t iTrack=0;
186   while (iTrack<nTracksMUON) {
187
188     fNPlanesMFTAnalyzed = 0;
189
190     AliInfo("****************************************************************************************");
191     AliInfo(Form("***************************   MUON TRACK %3d/%d   ***************************************", iTrack, nTracksMUON));
192     AliInfo("****************************************************************************************");
193     
194     fCandidateTracks -> Delete();
195     
196     fNPlanesMFTAnalyzed = 0;
197     
198     const AliESDMuonTrack *esdTrack = event->GetMuonTrack(iTrack);
199     if (fMUONTrack) delete fMUONTrack;
200     fMUONTrack = new AliMUONTrack();
201
202     AliMUONESDInterface::ESDToMUON(*esdTrack, *fMUONTrack, kTRUE);     // last arguments should be kTRUE or kFALSE??
203
204     // the track we are going to build, starting from fMUONTrack and adding the MFT clusters
205     AliMuonForwardTrack *track = new ((*fCandidateTracks)[0]) AliMuonForwardTrack();
206     track -> SetMUONTrack(new AliMUONTrack(*fMUONTrack));
207     track -> SetMCLabel(fMUONTrack->GetMCLabel());
208     track -> SetMatchTrigger(fMUONTrack->GetMatchTrigger());
209
210     //------------------------- NOW THE CYCLE OVER THE MFT PLANES STARTS ---------------------------------------
211
212     for (Int_t iPlane=fNPlanesMFT-1; iPlane>=0; iPlane--) {   /* *** do not reverse the order of this cycle!!! 
213                                                                  *** this reflects the fact that the extrapolation is performed 
214                                                                  *** starting from the last MFT plane back to the origin */
215       
216       // --------- updating the array of candidates according to the clusters available in the i-th plane ---------
217       
218       fNPlanesMFTAnalyzed++;
219       
220       Int_t nCandidates = fCandidateTracks->GetEntriesFast();
221       for (Int_t iCandidate=0; iCandidate<nCandidates; iCandidate++) {
222         fCurrentTrack = (AliMuonForwardTrack*) fCandidateTracks->UncheckedAt(iCandidate);
223
224         // if the old track is compatible with the new cluster, the track is updated and inserted as new track in the array 
225         // (several new tracks can be created for one old track)
226         if (FindClusterInPlane(iPlane) == kDiverged) {
227           fGlobalTrackingDiverged = kTRUE;
228           break;
229         }
230
231         if ((fNPlanesMFTAnalyzed-fCurrentTrack->GetNMFTClusters())>fNMaxMissingMFTClusters || fIsPlaneMandatory[iPlane]) {
232           fCandidateTracks->Remove(fCurrentTrack);     // the old track is removed after the check;
233         }
234       }
235       if (fGlobalTrackingDiverged) {
236         if (fScaleSigmaClusterCut>0) fScaleSigmaClusterCut -= 0.1;
237         continue;
238       }
239
240       fCandidateTracks->Compress();
241       
242     }      
243
244     // -------------------------- END OF THE CYCLE OVER THE MFT PLANES --------------------------------------------
245     
246     fGlobalTrackingDiverged = kFALSE;
247     fScaleSigmaClusterCut = 1.0;
248     
249     AliDebug(1, "Finished cycle over planes");
250     
251     iTrack++;
252
253     // If we have several final tracks, we must find the best candidate:
254     
255     Int_t nFinalTracks = fCandidateTracks->GetEntriesFast();
256     AliInfo(Form("nFinalTracks = %d", nFinalTracks));
257
258     Int_t nGoodClustersBestCandidate =  0;
259     Int_t idBestCandidate            =  0;
260     Double_t bestChi2                = -1.;  // variable defining the best candidate
261
262     for (Int_t iFinalCandidate=0; iFinalCandidate<nFinalTracks; iFinalCandidate++) {
263       
264       AliMuonForwardTrack *finalTrack = (AliMuonForwardTrack*) fCandidateTracks->UncheckedAt(iFinalCandidate);
265       Int_t nMFTClusters  = finalTrack->GetNMFTClusters();
266
267       Double_t chi2 = 0;
268       for (Int_t iCluster=0; iCluster<nMFTClusters; iCluster++) {
269         AliMFTCluster *localCluster = finalTrack->GetMFTCluster(iCluster);
270         chi2 += localCluster->GetLocalChi2();
271       }
272       chi2 /= nMFTClusters;
273
274       // now comparing the tracks in order to find the best one
275       
276       if (chi2<bestChi2 || bestChi2<0) {
277         bestChi2 = chi2;
278         idBestCandidate = iFinalCandidate;
279       }
280       
281     }
282     
283     if (nFinalTracks) {
284
285       AliMuonForwardTrack *newTrack = (AliMuonForwardTrack*) fCandidateTracks->UncheckedAt(idBestCandidate);
286       newTrack -> SetNWrongClustersMC(newTrack->GetNMFTClusters() - nGoodClustersBestCandidate);
287
288       new ((*muonForwardTracks)[muonForwardTracks->GetEntries()]) AliMuonForwardTrack(*newTrack);
289
290       //----------------------- Save the information to the AliESDMuonGlobalTrack object
291
292       newTrack -> EvalKinem(vertex[2]);     // we evaluate the kinematics at the primary vertex
293
294       AliDebug(2,"Creating a new Muon Global Track");
295       AliESDMuonGlobalTrack *myESDTrack = event->NewMuonGlobalTrack();
296       myESDTrack -> SetPxPyPz(newTrack->Px(), newTrack->Py(), newTrack->Pz());
297       myESDTrack -> SetChi2OverNdf(newTrack->GetChi2OverNdf());
298       myESDTrack -> SetCharge(newTrack->GetCharge());
299       myESDTrack -> SetMatchTrigger(newTrack->GetMatchTrigger());
300       myESDTrack -> SetFirstTrackingPoint(newTrack->GetMFTCluster(0)->GetX(), newTrack->GetMFTCluster(0)->GetY(), newTrack->GetMFTCluster(0)->GetZ());
301       myESDTrack -> SetXYAtVertex(newTrack->GetOffsetX(vertex[0], vertex[2]), newTrack->GetOffsetX(vertex[1], vertex[2]));
302       myESDTrack -> SetRAtAbsorberEnd(newTrack->GetRAtAbsorberEnd());
303       myESDTrack -> SetChi2MatchTrigger(esdTrack->GetChi2MatchTrigger());
304       myESDTrack -> SetMuonClusterMap(esdTrack->GetMuonClusterMap());
305       myESDTrack -> SetHitsPatternInTrigCh(esdTrack->GetHitsPatternInTrigCh());
306       myESDTrack -> SetHitsPatternInTrigChTrk(esdTrack->GetHitsPatternInTrigChTrk());
307       myESDTrack -> Connected(esdTrack->IsConnected());
308
309       //---------------------------------------------------------------------------------
310
311     }
312
313     fFinalBestCandidate = NULL;
314    
315   }
316
317   outputTreeMuonGlobalTracks -> Fill();
318
319   Int_t myEventID = 0;
320   while (outputFileMuonGlobalTracks->cd(Form("Event%d",myEventID))) myEventID++;
321   outputFileMuonGlobalTracks -> mkdir(Form("Event%d",myEventID));
322   outputFileMuonGlobalTracks -> cd(Form("Event%d",myEventID));
323   outputTreeMuonGlobalTracks -> Write();
324   outputFileMuonGlobalTracks -> Close();
325
326   muonForwardTracks -> Delete();
327   delete muonForwardTracks;
328
329   return 0;
330
331 }
332
333 //=========================================================================================================================================
334
335 void AliMFTTrackerMU::SeparateFrontBackClusters() {
336
337   for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
338     fMFTClusterArrayFront[iPlane]->Delete();
339     fMFTClusterArrayBack[iPlane] ->Delete();
340     for (Int_t iCluster=0; iCluster<fMFTClusterArray[iPlane]->GetEntries(); iCluster++) {
341       AliMFTCluster *cluster = (AliMFTCluster*) fMFTClusterArray[iPlane]->At(iCluster);
342       if (TMath::Abs(cluster->GetZ())<TMath::Abs(fSegmentation->GetPlane(iPlane)->GetZCenter())) {
343         new ((*fMFTClusterArrayFront[iPlane])[fMFTClusterArrayFront[iPlane]->GetEntries()]) AliMFTCluster(*cluster);
344       }
345       else {
346         new ((*fMFTClusterArrayBack[iPlane])[fMFTClusterArrayBack[iPlane]->GetEntries()]) AliMFTCluster(*cluster);
347       }
348     }
349   }
350
351 }
352
353 //==========================================================================================================================================
354
355 Int_t AliMFTTrackerMU::FindClusterInPlane(Int_t planeId) { 
356   
357   // !!!!!!!!! coordinates and errors on the interaction vertex should be taken from the event itself (ITS) if available
358
359   // propagate track to plane #planeId (both to front and back active sensors)
360   // look for compatible clusters
361   // update TrackParam at found cluster (if any) using Kalman Filter
362
363   AliMUONTrackParam currentParamFront, currentParamBack, currentParamForResearchFront, currentParamForResearchBack;
364
365   if (planeId == fNPlanesMFT-1) {      // last plane of the telecope
366     currentParamFront = (*((AliMUONTrackParam*)(fMUONTrack->GetTrackParamAtCluster()->First())));
367     currentParamBack  = (*((AliMUONTrackParam*)(fMUONTrack->GetTrackParamAtCluster()->First())));
368     currentParamForResearchFront = currentParamFront;
369     currentParamForResearchBack  = currentParamBack;
370     Double_t xExtrap = gRandom->Gaus(0,fVertexErrorX);
371     Double_t yExtrap = gRandom->Gaus(0,fVertexErrorY);
372     Double_t zExtrap = gRandom->Gaus(0,fVertexErrorZ);
373     if (fBransonCorrection) {
374       AliMUONTrackExtrap::ExtrapToVertex(&currentParamFront, xExtrap, yExtrap, zExtrap, fVertexErrorX, fVertexErrorY); 
375       AliMUONTrackExtrap::ExtrapToVertex(&currentParamBack,  xExtrap, yExtrap, zExtrap, fVertexErrorX, fVertexErrorY); 
376     }
377     else {
378       AliMUONTrackExtrap::ExtrapToVertexWithoutBranson(&currentParamFront, zExtrap);
379       AliMUONTrackExtrap::ExtrapToVertexWithoutBranson(&currentParamBack,  zExtrap);
380     }
381     AliMUONTrackExtrap::ExtrapToVertex(&currentParamForResearchFront, xExtrap, yExtrap, zExtrap, fVertexErrorX, fVertexErrorY); 
382     AliMUONTrackExtrap::ExtrapToVertex(&currentParamForResearchBack,  xExtrap, yExtrap, zExtrap, fVertexErrorX, fVertexErrorY); 
383   }
384   else {          // MFT planes others than the last one: mult. scattering correction because of the upstream MFT planes is performed
385     currentParamFront = (*((AliMUONTrackParam*)(fCurrentTrack->GetTrackParamAtCluster()->First())));
386     currentParamBack  = (*((AliMUONTrackParam*)(fCurrentTrack->GetTrackParamAtCluster()->First())));
387     currentParamForResearchFront = currentParamFront;
388     currentParamForResearchBack  = currentParamBack;
389     AliMUONTrackExtrap::AddMCSEffect(&currentParamFront,           (fSegmentation->GetPlane(planeId+1)->GetEquivalentSilicon()+
390                                                                     fSegmentation->GetPlane(planeId)->GetEquivalentSiliconBeforeFront())/fRadLengthSi,-1.);
391     AliMUONTrackExtrap::AddMCSEffect(&currentParamForResearchFront,(fSegmentation->GetPlane(planeId+1)->GetEquivalentSilicon()+
392                                                                     fSegmentation->GetPlane(planeId)->GetEquivalentSiliconBeforeFront())/fRadLengthSi,-1.);
393     AliMUONTrackExtrap::AddMCSEffect(&currentParamBack,            (fSegmentation->GetPlane(planeId+1)->GetEquivalentSilicon()+
394                                                                     fSegmentation->GetPlane(planeId)->GetEquivalentSiliconBeforeBack())/fRadLengthSi,-1.);
395     AliMUONTrackExtrap::AddMCSEffect(&currentParamForResearchBack, (fSegmentation->GetPlane(planeId+1)->GetEquivalentSilicon()+
396                                                                     fSegmentation->GetPlane(planeId)->GetEquivalentSiliconBeforeBack())/fRadLengthSi,-1.);
397   }
398   // for all planes: extrapolation to the Z of the plane
399   AliMUONTrackExtrap::ExtrapToZCov(&currentParamFront,            -1.*fSegmentation->GetPlane(planeId)->GetZCenterActiveFront());   
400   AliMUONTrackExtrap::ExtrapToZCov(&currentParamForResearchFront, -1.*fSegmentation->GetPlane(planeId)->GetZCenterActiveFront());
401   AliMUONTrackExtrap::ExtrapToZCov(&currentParamBack,             -1.*fSegmentation->GetPlane(planeId)->GetZCenterActiveBack());   
402   AliMUONTrackExtrap::ExtrapToZCov(&currentParamForResearchBack,  -1.*fSegmentation->GetPlane(planeId)->GetZCenterActiveBack());
403
404   //---------------------------------------------------------------------------------------
405
406   TMatrixD covFront(5,5); covFront = currentParamForResearchFront.GetCovariances();
407   TMatrixD covBack(5,5);  covBack  = currentParamForResearchBack.GetCovariances();
408   
409   Double_t squaredError_X_Front = covFront(0,0);
410   Double_t squaredError_Y_Front = covFront(2,2);
411   Double_t squaredError_X_Back  = covBack(0,0);
412   Double_t squaredError_Y_Back  = covBack(2,2);
413
414   Double_t corrFact = 1.0;
415
416   Double_t researchRadiusFront = TMath::Sqrt(squaredError_X_Front + squaredError_Y_Front);
417   Double_t researchRadiusBack  = TMath::Sqrt(squaredError_X_Back  + squaredError_Y_Back);
418   if (0.5*(researchRadiusFront+researchRadiusBack)<fMinResearchRadiusAtPlane[planeId]) {
419     corrFact = fMinResearchRadiusAtPlane[planeId]/(0.5*(researchRadiusFront+researchRadiusBack));
420   }
421
422   //---------------------------------------------------------------------------------------
423
424   Double_t chi2cut = 2.*fScaleSigmaClusterCut*fScaleSigmaClusterCut*fSigmaClusterCut*fSigmaClusterCut;     // depends on the number of variables (here, 2)
425   
426   // Analyizing the clusters: FRONT ACTIVE ELEMENTS
427   
428   Int_t nClustersFront = fMFTClusterArrayFront[planeId]->GetEntries();
429   AliDebug(2, Form("There are %3d clusters in plane %02d FRONT\n", nClustersFront, planeId));
430   
431   for (Int_t iCluster=0; iCluster<nClustersFront; iCluster++) {
432
433     Bool_t isGoodChi2 = kFALSE;
434
435     AliMFTCluster *cluster = (AliMFTCluster*) fMFTClusterArrayFront[planeId]->At(iCluster); 
436     Double_t chi2 = (1./(corrFact*corrFact)) * TryOneCluster(currentParamForResearchFront, cluster);     // describes the compatibility between the track and the cluster
437     if (chi2<chi2cut) isGoodChi2 = kTRUE;
438
439     if (isGoodChi2) {
440       AliDebug(3, Form("accepting cluster: chi2=%f (cut = %f)\n", chi2, chi2cut));
441       AliMuonForwardTrack *newTrack = new ((*fCandidateTracks)[fCandidateTracks->GetEntriesFast()]) AliMuonForwardTrack(*fCurrentTrack);
442       if (fCandidateTracks->GetEntriesFast() > fMaxNCandidates) return kDiverged;
443       newTrack->AddTrackParamAtMFTCluster(currentParamFront, *cluster);    // creating new track param and attaching the cluster
444       AliDebug(2, Form("After plane %02d: newTrack->GetNMFTClusters() = %d (fCurrentTrack->GetNMFTClusters() = %d)", 
445                        planeId, newTrack->GetNMFTClusters(), fCurrentTrack->GetNMFTClusters()));
446       newTrack->SetPlaneExists(planeId);
447     }
448     else AliDebug(3, Form("discarding cluster: chi2=%f (cut = %f)\n", chi2, chi2cut));
449
450   }
451
452   // Analyizing the clusters: BACK ACTIVE ELEMENTS
453   
454   Int_t nClustersBack = fMFTClusterArrayBack[planeId]->GetEntries();
455   AliDebug(2, Form("There are %3d clusters in plane %02d BACK\n", nClustersBack, planeId));
456   
457   for (Int_t iCluster=0; iCluster<nClustersBack; iCluster++) {
458
459     Bool_t isGoodChi2 = kFALSE;
460
461     AliMFTCluster *cluster = (AliMFTCluster*) fMFTClusterArrayBack[planeId]->At(iCluster); 
462     Double_t chi2 = (1./(corrFact*corrFact)) * TryOneCluster(currentParamForResearchBack, cluster);     // describes the compatibility between the track and the cluster
463     if (chi2<chi2cut) isGoodChi2 = kTRUE;
464
465     if (isGoodChi2) {
466       AliDebug(3,Form("accepting cluster: chi2=%f (cut = %f)\n", chi2, chi2cut));
467       AliMuonForwardTrack *newTrack = new ((*fCandidateTracks)[fCandidateTracks->GetEntriesFast()]) AliMuonForwardTrack(*fCurrentTrack);
468       if (fCandidateTracks->GetEntriesFast() > fMaxNCandidates) return kDiverged;
469       newTrack->AddTrackParamAtMFTCluster(currentParamBack, *cluster);    // creating new track param and attaching the cluster
470       AliDebug(2, Form("After plane %02d: newTrack->GetNMFTClusters() = %d (fCurrentTrack->GetNMFTClusters() = %d)", 
471                        planeId, newTrack->GetNMFTClusters(), fCurrentTrack->GetNMFTClusters()));
472       newTrack->SetPlaneExists(planeId);
473     }
474     else AliDebug(3,Form("discarding cluster: chi2=%f (cut = %f)\n", chi2, chi2cut));
475
476   }
477
478   //---------------------------------------------------------------------------------------------
479
480   return kConverged;
481   
482 }
483
484 //==========================================================================================================================================
485
486 Double_t AliMFTTrackerMU::TryOneCluster(const AliMUONTrackParam &trackParam, AliMFTCluster *cluster) {
487
488   // Test the compatibility between the track and the cluster (using trackParam's covariance matrix):
489   // return the corresponding Chi2
490   // assume the track parameters are given at the Z of the cluster
491   
492   // Set differences between trackParam and cluster in the bending and non bending directions
493   Double_t dX = cluster->GetX() - trackParam.GetNonBendingCoor();
494   Double_t dY = cluster->GetY() - trackParam.GetBendingCoor();
495   AliDebug(3,Form("dX = %f, dY = %f\n", dX, dY));
496   
497   // Calculate errors and covariances
498   const TMatrixD& kParamCov = trackParam.GetCovariances();
499   Double_t sigmaX2 = kParamCov(0,0) + cluster->GetErrX2();
500   Double_t sigmaY2 = kParamCov(2,2) + cluster->GetErrY2();
501   AliDebug(3, Form("dX2 = %f, dY2 = %f\n", sigmaX2, sigmaY2));
502   Double_t covXY   = kParamCov(0,2);
503   Double_t det     = sigmaX2 * sigmaY2 - covXY * covXY;
504   
505   // Compute chi2
506   if (det==0.) return 1.e10;
507   return (dX*dX*sigmaY2 + dY*dY*sigmaX2 - 2.*dX*dY*covXY) / det;
508   
509 }
510
511 //=========================================================================================================================================
512