]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONTrackReconstructor.cxx
1cb6010d4c9353302992868a4a4edfd2f2172ef0
[u/mrichter/AliRoot.git] / MUON / AliMUONTrackReconstructor.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 /// \class AliMUONTrackReconstructor
19 /// MUON track reconstructor using the original method
20 ///
21 /// This class contains as data:
22 /// - the parameters for the track reconstruction
23 ///
24 /// It contains as methods, among others:
25 /// - MakeTracks to build the tracks
26 ///
27
28 #include <stdlib.h>
29 #include <Riostream.h>
30 #include <TMatrixD.h>
31
32 #include "AliMUONVTrackReconstructor.h"
33 #include "AliMUONTrackReconstructor.h"
34 #include "AliMUONData.h"
35 #include "AliMUONConstants.h"
36 #include "AliMUONRawCluster.h"
37 #include "AliMUONHitForRec.h"
38 #include "AliMUONSegment.h"
39 #include "AliMUONTrack.h"
40 #include "AliLog.h"
41 #include <TVirtualFitter.h>
42
43 // Functions to be minimized with Minuit
44 void TrackChi2(Int_t &NParam, Double_t *Gradient, Double_t &Chi2, Double_t *Param, Int_t Flag);
45 void TrackChi2MCS(Int_t &NParam, Double_t *Gradient, Double_t &Chi2, Double_t *Param, Int_t Flag);
46
47 void mnvertLocal(Double_t* a, Int_t l, Int_t m, Int_t n, Int_t& ifail);
48
49 Double_t MultipleScatteringAngle2(AliMUONTrackParam *param);
50
51 /// \cond CLASSIMP
52 ClassImp(AliMUONTrackReconstructor) // Class implementation in ROOT context
53 /// \endcond
54
55 //************* Defaults parameters for reconstruction
56 const Double_t AliMUONTrackReconstructor::fgkDefaultMaxChi2 = 100.0;
57
58 TVirtualFitter* AliMUONTrackReconstructor::fgFitter = 0x0; 
59
60 //__________________________________________________________________________
61 AliMUONTrackReconstructor::AliMUONTrackReconstructor(AliMUONData* data)
62   : AliMUONVTrackReconstructor(data),
63     fMaxChi2(fgkDefaultMaxChi2)
64 {
65   /// Constructor for class AliMUONTrackReconstructor
66   
67   // Memory allocation for the TClonesArray of reconstructed tracks
68   fRecTracksPtr = new TClonesArray("AliMUONTrack", 10);
69 }
70
71   //__________________________________________________________________________
72 AliMUONTrackReconstructor::~AliMUONTrackReconstructor(void)
73 {
74   /// Destructor for class AliMUONTrackReconstructor
75   delete fRecTracksPtr;
76 }
77
78   //__________________________________________________________________________
79 void AliMUONTrackReconstructor::AddHitsForRecFromRawClusters()
80 {
81   /// To add to the list of hits for reconstruction all the raw clusters
82   TTree *TR = fMUONData->TreeR();
83   AliMUONHitForRec *hitForRec;
84   AliMUONRawCluster *clus;
85   Int_t iclus, nclus, nTRentries;
86   TClonesArray *rawclusters;
87   AliDebug(1,"Enter AddHitsForRecFromRawClusters");
88   
89   fMUONData->SetTreeAddress("RC");
90   nTRentries = Int_t(TR->GetEntries());
91   if (nTRentries != 1) {
92     AliError(Form("nTRentries = %d not equal to 1 ",nTRentries));
93     exit(0);
94   }
95   fMUONData->GetRawClusters(); // only one entry  
96   
97   // Loop over tracking chambers
98   for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
99     // number of HitsForRec to 0 for the chamber
100     fNHitsForRecPerChamber[ch] = 0;
101     // index of first HitForRec for the chamber
102     if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0;
103     else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec;
104     rawclusters =fMUONData->RawClusters(ch);
105     nclus = (Int_t) (rawclusters->GetEntries());
106     // Loop over (cathode correlated) raw clusters
107     for (iclus = 0; iclus < nclus; iclus++) {
108       clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(iclus);
109       // new AliMUONHitForRec from raw cluster
110       // and increment number of AliMUONHitForRec's (total and in chamber)
111       hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(clus);
112       fNHitsForRec++;
113       (fNHitsForRecPerChamber[ch])++;
114       // more information into HitForRec
115       hitForRec->SetBendingReso2(clus->GetErrY() * clus->GetErrY());
116       hitForRec->SetNonBendingReso2(clus->GetErrX() * clus->GetErrX());
117       //  original raw cluster
118       hitForRec->SetChamberNumber(ch);
119       hitForRec->SetHitNumber(iclus);
120       // Z coordinate of the raw cluster (cm)
121       hitForRec->SetZ(clus->GetZ(0));
122       
123       StdoutToAliDebug(3,
124                        cout << "Chamber " << ch <<
125                        " raw cluster  " << iclus << " : " << endl;
126                        clus->Print("full");
127                        cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
128                        hitForRec->Print("full");
129                        );
130     } // end of cluster loop
131   } // end of chamber loop
132   SortHitsForRecWithIncreasingChamber(); 
133   
134   AliDebug(1,"End of AddHitsForRecFromRawClusters");
135     if (AliLog::GetGlobalDebugLevel() > 0) {
136       AliDebug(1, Form("NHitsForRec: %d",fNHitsForRec));
137       for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
138         AliDebug(1, Form("Chamber(0...): %d",ch));
139         AliDebug(1, Form("NHitsForRec: %d", fNHitsForRecPerChamber[ch]));
140         AliDebug(1, Form("Index(first HitForRec): %d", fIndexOfFirstHitForRecPerChamber[ch]));
141         for (Int_t hit = fIndexOfFirstHitForRecPerChamber[ch];
142              hit < fIndexOfFirstHitForRecPerChamber[ch] + fNHitsForRecPerChamber[ch];
143              hit++) {
144           AliDebug(1, Form("HitForRec index(0...): %d",hit));
145           ((*fHitsForRecPtr)[hit])->Dump();
146       }
147     }
148   }
149   
150   return;
151 }
152
153   //__________________________________________________________________________
154 void AliMUONTrackReconstructor::MakeSegments(void)
155 {
156   /// To make the list of segments in all stations,
157   /// from the list of hits to be reconstructed
158   AliDebug(1,"Enter MakeSegments");
159   // Loop over stations
160   for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) MakeSegmentsPerStation(st); 
161   
162   StdoutToAliDebug(3,
163     cout << "end of MakeSegments" << endl;
164     for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) 
165     {
166       cout << "station " << st
167             << "  has " << fNSegments[st] << " segments:"
168             << endl;
169       for (Int_t seg = 0; seg < fNSegments[st]; seg++) 
170       {
171               ((*fSegmentsPtr[st])[seg])->Print();
172       }
173     }
174                    );
175   return;
176 }
177
178   //__________________________________________________________________________
179 void AliMUONTrackReconstructor::MakeTracks(void)
180 {
181   /// To make the tracks from the list of segments and points in all stations
182   AliDebug(1,"Enter MakeTracks");
183   // Look for candidates from at least 3 aligned points in stations(1..) 4 and 5
184   MakeTrackCandidates();
185   // Follow tracks in stations(1..) 3, 2 and 1
186   FollowTracks();
187   // Remove double tracks
188   RemoveDoubleTracks();
189   UpdateHitForRecAtHit();
190 }
191
192   //__________________________________________________________________________
193 void AliMUONTrackReconstructor::MakeTrackCandidates(void)
194 {
195   /// To make track candidates
196   /// with at least 3 aligned points in stations(1..) 4 and 5
197   /// (two Segment's or one Segment and one HitForRec)
198   Int_t begStation, iBegSegment, nbCan1Seg1Hit, nbCan2Seg;
199   AliMUONSegment *begSegment;
200   AliDebug(1,"Enter MakeTrackCandidates");
201   // Loop over stations(1..) 5 and 4 for the beginning segment
202   for (begStation = 4; begStation > 2; begStation--) {
203     // Loop over segments in the beginning station
204     for (iBegSegment = 0; iBegSegment < fNSegments[begStation]; iBegSegment++) {
205       // pointer to segment
206       begSegment = (AliMUONSegment*) ((*fSegmentsPtr[begStation])[iBegSegment]);
207       AliDebug(2,Form("Look for TrackCandidate's with Segment %d  in Station(0..) %d", iBegSegment, begStation));
208       // Look for track candidates with two segments,
209       // "begSegment" and all compatible segments in other station.
210       // Only for beginning station(1..) 5
211       // because candidates with 2 segments have to looked for only once.
212       if (begStation == 4)
213         nbCan2Seg = MakeTrackCandidatesWithTwoSegments(begSegment);
214       // Look for track candidates with one segment and one point,
215       // "begSegment" and all compatible HitForRec's in other station.
216       // Only if "begSegment" does not belong already to a track candidate.
217       // Is that a too strong condition ????
218       if (!(begSegment->GetInTrack()))
219         nbCan1Seg1Hit = MakeTrackCandidatesWithOneSegmentAndOnePoint(begSegment);
220     } // for (iBegSegment = 0;...
221   } // for (begStation = 4;...
222   return;
223 }
224
225   //__________________________________________________________________________
226 Int_t AliMUONTrackReconstructor::MakeTrackCandidatesWithTwoSegments(AliMUONSegment *BegSegment)
227 {
228   /// To make track candidates with two segments in stations(1..) 4 and 5,
229   /// the first segment being pointed to by "BegSegment".
230   /// Returns the number of such track candidates.
231   Int_t endStation, iEndSegment, nbCan2Seg;
232   AliMUONSegment *endSegment;
233   AliMUONSegment *extrapSegment = NULL;
234   AliMUONTrack *recTrack;
235   Double_t mcsFactor;
236   AliDebug(1,"Enter MakeTrackCandidatesWithTwoSegments");
237   // Station for the end segment
238   endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
239   // multiple scattering factor corresponding to one chamber
240   mcsFactor = 0.0136 /
241     GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
242   mcsFactor     = fChamberThicknessInX0 * mcsFactor * mcsFactor;
243   // linear extrapolation to end station
244   // number of candidates with 2 segments to 0
245   nbCan2Seg = 0;
246   // Loop over segments in the end station
247   for (iEndSegment = 0; iEndSegment < fNSegments[endStation]; iEndSegment++) {
248     // pointer to segment
249     endSegment = (AliMUONSegment*) ((*fSegmentsPtr[endStation])[iEndSegment]);
250     // test compatibility between current segment and "extrapSegment"
251     // 4 because 4 quantities in chi2
252     extrapSegment =
253       BegSegment->CreateSegmentFromLinearExtrapToStation(endSegment->GetZ(), mcsFactor);
254     if ((endSegment->
255          NormalizedChi2WithSegment(extrapSegment,
256                                    fMaxSigma2Distance)) <= 4.0) {
257       // both segments compatible:
258       // make track candidate from "begSegment" and "endSegment"
259       AliDebug(2,Form("TrackCandidate with Segment %d in Station(0..) %d", iEndSegment, endStation));
260       // flag for both segments in one track:
261       // to be done in track constructor ????
262       BegSegment->SetInTrack(kTRUE);
263       endSegment->SetInTrack(kTRUE);
264       recTrack = new ((*fRecTracksPtr)[fNRecTracks]) AliMUONTrack(BegSegment, endSegment);
265       // Set track parameters at vertex from last stations 4 & 5
266       CalcTrackParamAtVertex(recTrack);
267       fNRecTracks++;
268       if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump();
269       // increment number of track candidates with 2 segments
270       nbCan2Seg++;
271     }
272     delete extrapSegment; // should not delete HitForRec's it points to !!!!
273   } // for (iEndSegment = 0;...
274   return nbCan2Seg;
275 }
276
277   //__________________________________________________________________________
278 Int_t AliMUONTrackReconstructor::MakeTrackCandidatesWithOneSegmentAndOnePoint(AliMUONSegment *BegSegment)
279 {
280   /// To make track candidates with one segment and one point
281   /// in stations(1..) 4 and 5,
282   /// the segment being pointed to by "BegSegment".
283   Int_t ch, ch1, ch2, endStation, iHit, iHitMax, iHitMin, nbCan1Seg1Hit;
284   AliMUONHitForRec *extrapHitForRec= NULL;
285   AliMUONHitForRec *hit;
286   AliMUONTrack *recTrack;
287   Double_t mcsFactor;
288   AliDebug(1,"Enter MakeTrackCandidatesWithOneSegmentAndOnePoint");
289   // station for the end point
290   endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
291   // multiple scattering factor corresponding to one chamber
292   mcsFactor = 0.0136 /
293     GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
294   mcsFactor     = fChamberThicknessInX0 * mcsFactor * mcsFactor;
295   // first and second chambers(0..) in the end station
296   ch1 = 2 * endStation;
297   ch2 = ch1 + 1;
298   // number of candidates to 0
299   nbCan1Seg1Hit = 0;
300   // Loop over chambers of the end station
301   for (ch = ch2; ch >= ch1; ch--) {
302     // limits for the hit index in the loop
303     iHitMin = fIndexOfFirstHitForRecPerChamber[ch];
304     iHitMax = iHitMin + fNHitsForRecPerChamber[ch];
305     // Loop over HitForRec's in the chamber
306     for (iHit = iHitMin; iHit < iHitMax; iHit++) {
307       // pointer to HitForRec
308       hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
309       // test compatibility between current HitForRec and "extrapHitForRec"
310       // 2 because 2 quantities in chi2
311       // linear extrapolation to chamber
312       extrapHitForRec =
313         BegSegment->CreateHitForRecFromLinearExtrapToChamber( hit->GetZ(), mcsFactor);
314       if ((hit->
315            NormalizedChi2WithHitForRec(extrapHitForRec,
316                                        fMaxSigma2Distance)) <= 2.0) {
317         // both HitForRec's compatible:
318         // make track candidate from begSegment and current HitForRec
319         AliDebug(2, Form("TrackCandidate with HitForRec  %d in Chamber(0..) %d", iHit, ch));
320         // flag for beginning segments in one track:
321         // to be done in track constructor ????
322         BegSegment->SetInTrack(kTRUE);
323         recTrack = new ((*fRecTracksPtr)[fNRecTracks]) AliMUONTrack(BegSegment, hit);
324         // Set track parameters at vertex from last stations 4 & 5
325         CalcTrackParamAtVertex(recTrack);
326         // the right place to eliminate "double counting" ???? how ????
327         fNRecTracks++;
328         if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump();
329         // increment number of track candidates
330         nbCan1Seg1Hit++;
331       }
332       delete extrapHitForRec;
333     } // for (iHit = iHitMin;...
334   } // for (ch = ch2;...
335   return nbCan1Seg1Hit;
336 }
337
338   //__________________________________________________________________________
339 void AliMUONTrackReconstructor::CalcTrackParamAtVertex(AliMUONTrack *Track) const
340 {
341   /// Set track parameters at vertex.
342   /// TrackHit's are assumed to be only in stations(1..) 4 and 5,
343   /// and sorted according to increasing Z.
344   /// Parameters are calculated from information in HitForRec's
345   /// of first and last TrackHit's.
346   AliMUONTrackParam *trackParamVertex = Track->GetTrackParamAtVertex(); // pointer to track parameters at vertex
347   // Pointer to HitForRec attached to first TrackParamAtHit
348   AliMUONHitForRec *firstHit = ((AliMUONTrackParam*) (Track->GetTrackParamAtHit()->First()))->GetHitForRecPtr();
349   // Pointer to HitForRec attached to last TrackParamAtHit
350   AliMUONHitForRec *lastHit = ((AliMUONTrackParam*) (Track->GetTrackParamAtHit()->Last()))->GetHitForRecPtr();
351   // Z difference between first and last hits
352   Double_t deltaZ = firstHit->GetZ() - lastHit->GetZ();
353   // bending slope in stations(1..) 4 and 5
354   Double_t bendingSlope = (firstHit->GetBendingCoor() - lastHit->GetBendingCoor()) / deltaZ;
355   trackParamVertex->SetBendingSlope(bendingSlope);
356   // impact parameter
357   Double_t impactParam = firstHit->GetBendingCoor() - bendingSlope * firstHit->GetZ();
358   // signed bending momentum
359   trackParamVertex->SetInverseBendingMomentum(1.0 / GetBendingMomentumFromImpactParam(impactParam));
360   // bending slope at vertex
361   trackParamVertex->SetBendingSlope(bendingSlope + impactParam / fSimpleBPosition);
362   // non bending slope
363   trackParamVertex->SetNonBendingSlope((firstHit->GetNonBendingCoor() - lastHit->GetNonBendingCoor()) / deltaZ);
364   // vertex coordinates at (0,0,0)
365   trackParamVertex->SetZ(0.0);
366   trackParamVertex->SetBendingCoor(0.0);
367   trackParamVertex->SetNonBendingCoor(0.0);
368 }
369
370   //__________________________________________________________________________
371 void AliMUONTrackReconstructor::FollowTracks(void)
372 {
373   /// Follow tracks in stations(1..) 3, 2 and 1
374   // too long: should be made more modular !!!!
375   AliMUONHitForRec *bestHit, *extrapHit, *hit;
376   AliMUONSegment *bestSegment, *extrapSegment, *segment;
377   AliMUONTrack *track, *nextTrack;
378   AliMUONTrackParam *trackParam1, trackParam[2], trackParamVertex;
379   // -1 to avoid compilation warnings
380   Int_t ch = -1, chInStation, chBestHit = -1, iHit, iSegment, station, trackIndex; 
381   Double_t bestChi2, chi2, dZ1, dZ2, dZ3, maxSigma2Distance, mcsFactor;
382   Double_t bendingMomentum, chi2Norm = 0.;
383
384
385   // local maxSigma2Distance, for easy increase in testing
386   maxSigma2Distance = fMaxSigma2Distance;
387   AliDebug(2,"Enter FollowTracks");
388   // Loop over track candidates
389   track = (AliMUONTrack*) fRecTracksPtr->First();
390   trackIndex = -1;
391   while (track) {
392     trackIndex++;
393     nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track); // prepare next track
394     AliDebug(2,Form("FollowTracks: track candidate(0..): %d", trackIndex));
395     // Fit track candidate from parameters at vertex
396     // -> with 3 parameters (X_vertex and Y_vertex are fixed)
397     // without multiple Coulomb scattering
398     Fit(track,0,0);
399     if (AliLog::GetGlobalDebugLevel()> 2) {
400       cout << "FollowTracks: track candidate(0..): " << trackIndex
401            << " after fit in stations(0..) 3 and 4" << endl;
402       track->RecursiveDump();
403     }
404     // Loop over stations(1..) 3, 2 and 1
405     // something SPECIAL for stations 2 and 1 for majority 3 coincidence ????
406     // otherwise: majority coincidence 2 !!!!
407     for (station = 2; station >= 0; station--) {
408       // Track parameters at first track hit (smallest Z)
409       trackParam1 = (AliMUONTrackParam*) (track->GetTrackParamAtHit()->First());
410       // extrapolation to station
411       trackParam1->ExtrapToStation(station, trackParam);
412       extrapSegment = new AliMUONSegment(); //  empty segment
413       // multiple scattering factor corresponding to one chamber
414       // and momentum in bending plane (not total)
415       mcsFactor = 0.0136 * trackParam1->GetInverseBendingMomentum();
416       mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
417       // Z difference from previous station
418       dZ1 = AliMUONConstants::DefaultChamberZ(2 * station) -
419             AliMUONConstants::DefaultChamberZ(2 * station + 2);
420       // Z difference between the two previous stations
421       dZ2 = AliMUONConstants::DefaultChamberZ(2 * station + 2) -
422             AliMUONConstants::DefaultChamberZ(2 * station + 4);
423       // Z difference between the two chambers in the previous station
424       dZ3 = AliMUONConstants::DefaultChamberZ(2 * station) -
425             AliMUONConstants::DefaultChamberZ(2 * station + 1);
426       extrapSegment->SetBendingCoorReso2(fBendingResolution * fBendingResolution);
427       extrapSegment->SetNonBendingCoorReso2(fNonBendingResolution * fNonBendingResolution);
428       extrapSegment->UpdateFromStationTrackParam(trackParam, mcsFactor, dZ1, dZ2, dZ3, station,
429                                                  trackParam1->GetInverseBendingMomentum());
430       bestChi2 = 5.0;
431       bestSegment = NULL;
432       if (AliLog::GetGlobalDebugLevel() > 2) {
433         cout << "FollowTracks: track candidate(0..): " << trackIndex
434              << " Look for segment in station(0..): " << station << endl;
435       }
436
437       // Loop over segments in station
438       for (iSegment = 0; iSegment < fNSegments[station]; iSegment++) {
439         // Look for best compatible Segment in station
440         // should consider all possibilities ????
441         // multiple scattering ????
442         // separation in 2 functions: Segment and HitForRec ????
443         segment = (AliMUONSegment*) ((*fSegmentsPtr[station])[iSegment]);
444         // correction of corrected segment (fBendingCoor and fNonBendingCoor)
445         // according to real Z value of "segment" and slopes of "extrapSegment"
446         trackParam[0].ExtrapToZ(segment->GetZ());
447         trackParam[1].ExtrapToZ(segment->GetZ()); // now same as trackParam[0] !?!?!?!?!?!
448         extrapSegment->SetBendingCoor((&(trackParam[0]))->GetBendingCoor());
449         extrapSegment->SetNonBendingCoor((&(trackParam[0]))->GetNonBendingCoor());
450         extrapSegment->SetBendingSlope((&(trackParam[0]))->GetBendingSlope());
451         extrapSegment->SetNonBendingSlope((&(trackParam[0]))->GetNonBendingSlope());
452         chi2 = segment->NormalizedChi2WithSegment(extrapSegment, maxSigma2Distance);
453         if (chi2 < bestChi2) {
454           // update best Chi2 and Segment if better found
455           bestSegment = segment;
456           bestChi2 = chi2;
457         }
458       }
459       if (bestSegment) {
460         // best segment found: add it to track candidate
461         trackParam[0].ExtrapToZ(bestSegment->GetZ());
462         track->AddTrackParamAtHit(&(trackParam[0]),bestSegment->GetHitForRec1());
463         trackParam[1].ExtrapToZ(bestSegment->GetZ()); // now same as trackParam[0] !?!?!?!?!?!
464         track->AddTrackParamAtHit(&(trackParam[1]),bestSegment->GetHitForRec2());
465         AliDebug(3, Form("FollowTracks: track candidate(0..): %d  Added segment in station(0..): %d", trackIndex, station));
466         if (AliLog::GetGlobalDebugLevel()>2) track->RecursiveDump();
467       } else {
468         // No best segment found:
469         // Look for best compatible HitForRec in station:
470         // should consider all possibilities ????
471         // multiple scattering ???? do about like for extrapSegment !!!!
472         extrapHit = new AliMUONHitForRec(); //  empty hit
473         bestChi2 = 3.0;
474         bestHit = NULL;
475         AliDebug(3, Form("FollowTracks: track candidate(0..): %d Segment not found, look for hit in station(0..): %d ", 
476                          trackIndex, station));
477         
478         // Loop over chambers of the station
479         for (chInStation = 0; chInStation < 2; chInStation++) {
480           ch = 2 * station + chInStation;
481           for (iHit = fIndexOfFirstHitForRecPerChamber[ch]; iHit < fIndexOfFirstHitForRecPerChamber[ch]+fNHitsForRecPerChamber[ch]; iHit++) {
482             hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
483             // coordinates of extrapolated hit
484             trackParam[chInStation].ExtrapToZ(hit->GetZ());
485             extrapHit->SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor());
486             extrapHit->SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor());
487             // resolutions from "extrapSegment"
488             extrapHit->SetBendingReso2(extrapSegment->GetBendingCoorReso2());
489             extrapHit->SetNonBendingReso2(extrapSegment->GetNonBendingCoorReso2());
490             // Loop over hits in the chamber
491             // condition for hit not already in segment ????
492             chi2 = hit->NormalizedChi2WithHitForRec(extrapHit, maxSigma2Distance);
493             if (chi2 < bestChi2) {
494               // update best Chi2 and HitForRec if better found
495               bestHit = hit;
496               bestChi2 = chi2;
497               chBestHit = chInStation;
498             }
499           }
500         }
501         if (bestHit) {
502           // best hit found: add it to track candidate
503           trackParam[chBestHit].ExtrapToZ(bestHit->GetZ());
504           track->AddTrackParamAtHit(&(trackParam[chBestHit]),bestHit);
505           if (AliLog::GetGlobalDebugLevel() > 2) {
506             cout << "FollowTracks: track candidate(0..): " << trackIndex
507                  << " Added hit in station(0..): " << station << endl;
508             track->RecursiveDump();
509           }
510         } else {
511           // Remove current track candidate
512           // and corresponding TrackHit's, ...
513           fRecTracksPtr->Remove(track);
514           fNRecTracks--;
515           delete extrapSegment;
516           delete extrapHit;
517           break; // stop the search for this candidate:
518           // exit from the loop over station
519         }
520         delete extrapHit;
521       }
522       delete extrapSegment;
523       // Sort TrackParamAtHit according to increasing Z
524       track->GetTrackParamAtHit()->Sort();
525       // Update track parameters at first track hit (smallest Z)
526       trackParam1 = (AliMUONTrackParam*) (track->GetTrackParamAtHit()->First());
527       bendingMomentum = 0.;
528       if (TMath::Abs(trackParam1->GetInverseBendingMomentum()) > 0.)
529         bendingMomentum = TMath::Abs(1/(trackParam1->GetInverseBendingMomentum()));
530       // Track removed if bendingMomentum not in window [min, max]
531       if ((bendingMomentum < fMinBendingMomentum) || (bendingMomentum > fMaxBendingMomentum)) {
532         fRecTracksPtr->Remove(track);
533         fNRecTracks--;
534         break; // stop the search for this candidate:
535         // exit from the loop over station 
536       }
537       // Track fit from parameters at first hit
538       // -> with 5 parameters (momentum and position)
539       // with multiple Coulomb scattering if all stations
540       if (station == 0) Fit(track,1,1);
541       // without multiple Coulomb scattering if not all stations
542       else Fit(track,1,0);
543       Double_t numberOfDegFree = (2.0 * track->GetNTrackHits() - 5);
544       if (numberOfDegFree > 0) {
545         chi2Norm =  track->GetFitFMin() / numberOfDegFree;
546       } else {
547         chi2Norm = 1.e10;
548       }
549       // Track removed if normalized chi2 too high
550       if (chi2Norm > fMaxChi2) {
551         fRecTracksPtr->Remove(track);
552         fNRecTracks--;
553         break; // stop the search for this candidate:
554         // exit from the loop over station 
555       }
556       if (AliLog::GetGlobalDebugLevel() > 2) {
557         cout << "FollowTracks: track candidate(0..): " << trackIndex
558              << " after fit from station(0..): " << station << " to 4" << endl;
559         track->RecursiveDump();
560       }
561       // Track extrapolation to the vertex through the absorber (Branson)
562       // after going through the first station
563       if (station == 0) {
564         trackParamVertex = *((AliMUONTrackParam*) (track->GetTrackParamAtHit()->First()));
565         (&trackParamVertex)->ExtrapToVertex(0.,0.,0.);
566         track->SetTrackParamAtVertex(&trackParamVertex);
567         if (AliLog::GetGlobalDebugLevel() > 0) {
568           cout << "FollowTracks: track candidate(0..): " << trackIndex
569                << " after extrapolation to vertex" << endl;
570           track->RecursiveDump();
571         }
572       }
573     } // for (station = 2;...
574     // go really to next track
575     track = nextTrack;
576   } // while (track)
577   // Compression of track array (necessary after Remove)
578   fRecTracksPtr->Compress();
579   return;
580 }
581
582   //__________________________________________________________________________
583 void AliMUONTrackReconstructor::Fit(AliMUONTrack *Track, Int_t FitStart, Int_t FitMCS)
584 {
585   /// Fit the track "Track",
586   /// with or without multiple Coulomb scattering according to "FitMCS",
587   /// starting, according to "FitStart",
588   /// with track parameters at vertex or at the first TrackHit.
589   
590   if ((FitStart != 0) && (FitStart != 1)) {
591     cout << "ERROR in AliMUONTrackReconstructor::Fit(...)" << endl;
592     cout << "FitStart = " << FitStart << " is neither 0 nor 1" << endl;
593     exit(0);
594   }
595   if ((FitMCS != 0) && (FitMCS != 1)) {
596     cout << "ERROR in AliMUONTrackReconstructor::Fit(...)" << endl;
597     cout << "FitMCS = " << FitMCS << " is neither 0 nor 1" << endl;
598     exit(0);
599   }
600   
601   Double_t arg[1], benC, errorParam, invBenP, lower, nonBenC, upper, x, y;
602   char parName[50];
603   AliMUONTrackParam *trackParam;
604   // Check if Minuit is initialized...
605   fgFitter = TVirtualFitter::Fitter(Track,5);
606   fgFitter->Clear(); // necessary ???? probably yes
607   // how to reset the printout number at every fit ????
608   // is there any risk to leave it like that ????
609   // how to go faster ???? choice of Minuit parameters like EDM ????
610   // choice of function to be minimized according to fFitMCS
611   if (FitMCS == 0) fgFitter->SetFCN(TrackChi2);
612   else fgFitter->SetFCN(TrackChi2MCS);
613   // Switch off printout
614   arg[0] = -1;
615   fgFitter->ExecuteCommand("SET PRINT", arg, 1); // More printing !!!!
616   // No warnings
617   fgFitter->ExecuteCommand("SET NOW", arg, 0);
618   // Parameters according to "fFitStart"
619   // (should be a function to be used at every place where needed ????)
620   if (FitStart == 0) trackParam = Track->GetTrackParamAtVertex();
621   else trackParam = (AliMUONTrackParam*) (Track->GetTrackParamAtHit()->First());
622   // set first 3 Minuit parameters
623   // could be tried with no limits for the search (min=max=0) ????
624   fgFitter->SetParameter(0, "InvBenP",
625                          trackParam->GetInverseBendingMomentum(),
626                          0.003, -0.4, 0.4);
627   fgFitter->SetParameter(1, "BenS",
628                          trackParam->GetBendingSlope(),
629                          0.001, -0.5, 0.5);
630   fgFitter->SetParameter(2, "NonBenS",
631                          trackParam->GetNonBendingSlope(),
632                          0.001, -0.5, 0.5);
633   if (FitStart == 1) {
634     // set last 2 Minuit parameters when we start from first track hit
635     // mandatory limits in Bending to avoid NaN values of parameters
636     fgFitter->SetParameter(3, "X",
637                            trackParam->GetNonBendingCoor(),
638                            0.03, -500.0, 500.0);
639     // mandatory limits in non Bending to avoid NaN values of parameters
640     fgFitter->SetParameter(4, "Y",
641                            trackParam->GetBendingCoor(),
642                            0.10, -500.0, 500.0);
643   }
644   // search without gradient calculation in the function
645   fgFitter->ExecuteCommand("SET NOGRADIENT", arg, 0);
646   // minimization
647   fgFitter->ExecuteCommand("MINIMIZE", arg, 0);
648   // exit from Minuit
649   //  fgFitter->ExecuteCommand("EXIT", arg, 0); // necessary ????
650   // get results into "invBenP", "benC", "nonBenC" ("x", "y")
651   fgFitter->GetParameter(0, parName, invBenP, errorParam, lower, upper);
652   trackParam->SetInverseBendingMomentum(invBenP);
653   fgFitter->GetParameter(1, parName, benC, errorParam, lower, upper);
654   trackParam->SetBendingSlope(benC);
655   fgFitter->GetParameter(2, parName, nonBenC, errorParam, lower, upper);
656   trackParam->SetNonBendingSlope(nonBenC);
657   if (FitStart == 1) {
658     fgFitter->GetParameter(3, parName, x, errorParam, lower, upper);
659     trackParam->SetNonBendingCoor(x);
660     fgFitter->GetParameter(4, parName, y, errorParam, lower, upper);
661     trackParam->SetBendingCoor(y);
662   }
663   // global result of the fit
664   Double_t fedm, errdef, fitFMin;
665   Int_t npari, nparx;
666   fgFitter->GetStats(fitFMin, fedm, errdef, npari, nparx);
667   Track->SetFitFMin(fitFMin);
668 }
669
670   //__________________________________________________________________________
671 void TrackChi2(Int_t &NParam, Double_t * /*Gradient*/, Double_t &Chi2, Double_t *Param, Int_t /*Flag*/)
672 {
673   /// Return the "Chi2" to be minimized with Minuit for track fitting,
674   /// with "NParam" parameters
675   /// and their current values in array pointed to by "Param".
676   /// Assumes that the track hits are sorted according to increasing Z.
677   /// Track parameters at each TrackHit are updated accordingly.
678   /// Multiple Coulomb scattering is not taken into account
679   AliMUONTrack *trackBeingFitted;
680   AliMUONTrackParam param1;
681   AliMUONTrackParam* trackParamAtHit;
682   AliMUONHitForRec* hitForRec;
683   Chi2 = 0.0; // initialize Chi2
684   // copy of track parameters to be fitted
685   trackBeingFitted = (AliMUONTrack*) AliMUONTrackReconstructor::Fitter()->GetObjectFit();
686   // 3 parameters means fit track candidate from parameters at vertex (X_vertex and Y_vertex are fixed)
687   if (NParam == 3) param1 = *(trackBeingFitted->GetTrackParamAtVertex());
688   else param1 = *((AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->First()));
689   // Minuit parameters to be fitted into this copy
690   param1.SetInverseBendingMomentum(Param[0]);
691   param1.SetBendingSlope(Param[1]);
692   param1.SetNonBendingSlope(Param[2]);
693   if (NParam == 5) {
694     param1.SetNonBendingCoor(Param[3]);
695     param1.SetBendingCoor(Param[4]);
696   }
697   // Follow track through all planes of track hits
698   trackParamAtHit = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->First());
699   while (trackParamAtHit) {
700     hitForRec = trackParamAtHit->GetHitForRecPtr();
701     // extrapolation to the plane of the hitForRec attached to the current trackParamAtHit
702     param1.ExtrapToZ(hitForRec->GetZ());
703     // update track parameters of the current hit
704     trackParamAtHit->SetTrackParam(param1);
705     // Increment Chi2
706     // done hit per hit, with hit resolution,
707     // and not with point and angle like in "reco_muon.F" !!!!
708     // Needs to add multiple scattering contribution ????
709     Double_t dX = hitForRec->GetNonBendingCoor() - param1.GetNonBendingCoor();
710     Double_t dY = hitForRec->GetBendingCoor() - param1.GetBendingCoor();
711     Chi2 = Chi2 + dX * dX / hitForRec->GetNonBendingReso2() + dY * dY / hitForRec->GetBendingReso2();
712     trackParamAtHit = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->After(trackParamAtHit));
713   }
714 }
715
716   //__________________________________________________________________________
717 void TrackChi2MCS(Int_t &NParam, Double_t * /*Gradient*/, Double_t &Chi2, Double_t *Param, Int_t /*Flag*/)
718 {
719   /// Return the "Chi2" to be minimized with Minuit for track fitting,
720   /// with "NParam" parameters
721   /// and their current values in array pointed to by "Param".
722   /// Assumes that the track hits are sorted according to increasing Z.
723   /// Track parameters at each TrackHit are updated accordingly.
724   /// Multiple Coulomb scattering is taken into account with covariance matrix.
725   AliMUONTrack *trackBeingFitted;
726   AliMUONTrackParam param1;
727   AliMUONTrackParam* trackParamAtHit;
728   AliMUONHitForRec* hitForRec;
729   Chi2 = 0.0; // initialize Chi2
730   // copy of track parameters to be fitted
731   trackBeingFitted = (AliMUONTrack*) AliMUONTrackReconstructor::Fitter()->GetObjectFit();
732   // 3 parameters means fit track candidate from parameters at vertex (X_vertex and Y_vertex are fixed)
733   if (NParam == 3) param1 = *(trackBeingFitted->GetTrackParamAtVertex());
734   else param1 = *((AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->First()));
735   // Minuit parameters to be fitted into this copy
736   param1.SetInverseBendingMomentum(Param[0]);
737   param1.SetBendingSlope(Param[1]);
738   param1.SetNonBendingSlope(Param[2]);
739   if (NParam == 5) {
740     param1.SetNonBendingCoor(Param[3]);
741     param1.SetBendingCoor(Param[4]);
742   }
743
744   Int_t chCurrent, chPrev = 0, hitNumber, hitNumber1, hitNumber2, hitNumber3;
745   Double_t z1, z2, z3;
746   AliMUONTrackParam *trackParamAtHit1, *trackParamAtHit2, *trackParamAtHit3;
747   AliMUONHitForRec *hitForRec1, *hitForRec2;
748   Double_t hbc1, hbc2, pbc1, pbc2;
749   Double_t hnbc1, hnbc2, pnbc1, pnbc2;
750   Int_t numberOfHit = trackBeingFitted->GetNTrackHits();
751   TMatrixD *covBending = new TMatrixD(numberOfHit, numberOfHit);
752   TMatrixD *covNonBending = new TMatrixD(numberOfHit, numberOfHit);
753   Double_t *msa2 = new Double_t[numberOfHit];
754
755   // Predicted coordinates and  multiple scattering angles are first calculated
756   for (hitNumber = 0; hitNumber < numberOfHit; hitNumber++) {
757     trackParamAtHit = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->UncheckedAt(hitNumber));
758     hitForRec = trackParamAtHit->GetHitForRecPtr();
759     // extrapolation to the plane of the hitForRec attached to the current trackParamAtHit
760     param1.ExtrapToZ(hitForRec->GetZ());
761     // update track parameters of the current hit
762     trackParamAtHit->SetTrackParam(param1);
763     // square of multiple scattering angle at current hit, with one chamber
764     msa2[hitNumber] = MultipleScatteringAngle2(&param1);
765     // correction for eventual missing hits or multiple hits in a chamber,
766     // according to the number of chambers
767     // between the current hit and the previous one
768     chCurrent = hitForRec->GetChamberNumber();
769     if (hitNumber > 0) msa2[hitNumber] = msa2[hitNumber] * (chCurrent - chPrev);
770     chPrev = chCurrent;
771   }
772
773   // Calculates the covariance matrix
774   for (hitNumber1 = 0; hitNumber1 < numberOfHit; hitNumber1++) { 
775     trackParamAtHit1 = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->UncheckedAt(hitNumber1));
776     hitForRec1 = trackParamAtHit1->GetHitForRecPtr();
777     z1 = hitForRec1->GetZ();
778     for (hitNumber2 = hitNumber1; hitNumber2 < numberOfHit; hitNumber2++) {
779       trackParamAtHit2 = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->UncheckedAt(hitNumber2));
780       z2 = trackParamAtHit2->GetHitForRecPtr()->GetZ();
781       // initialization to 0 (diagonal plus upper triangular part)
782       (*covBending)(hitNumber2, hitNumber1) = 0.0;
783       // contribution from multiple scattering in bending plane:
784       // loop over upstream hits
785       for (hitNumber3 = 0; hitNumber3 < hitNumber1; hitNumber3++) {     
786         trackParamAtHit3 = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->UncheckedAt(hitNumber3));
787         z3 = trackParamAtHit3->GetHitForRecPtr()->GetZ();
788         (*covBending)(hitNumber2, hitNumber1) = (*covBending)(hitNumber2, hitNumber1) + ((z1 - z3) * (z2 - z3) * msa2[hitNumber3]); 
789       }
790       // equal contribution from multiple scattering in non bending plane
791       (*covNonBending)(hitNumber2, hitNumber1) = (*covBending)(hitNumber2, hitNumber1);
792       if (hitNumber1 == hitNumber2) {
793         // Diagonal elements: add contribution from position measurements
794         // in bending plane
795         (*covBending)(hitNumber2, hitNumber1) = (*covBending)(hitNumber2, hitNumber1) + hitForRec1->GetBendingReso2();
796         // and in non bending plane
797         (*covNonBending)(hitNumber2, hitNumber1) = (*covNonBending)(hitNumber2, hitNumber1) + hitForRec1->GetNonBendingReso2();
798       } else {
799         // Non diagonal elements: symmetrization
800         // for bending plane
801         (*covBending)(hitNumber1, hitNumber2) = (*covBending)(hitNumber2, hitNumber1);
802         // and non bending plane
803         (*covNonBending)(hitNumber1, hitNumber2) = (*covNonBending)(hitNumber2, hitNumber1);
804       }
805     } // for (hitNumber2 = hitNumber1;...
806   } // for (hitNumber1 = 0;...
807     
808   // Inversion of covariance matrices
809   // with "mnvertLocal", local "mnvert" function of Minuit.
810   // One cannot use directly "mnvert" since "TVirtualFitter" does not know it.
811   // One will have to replace this local function by the right inversion function
812   // from a specialized Root package for symmetric positive definite matrices,
813   // when available!!!!
814   Int_t ifailBending;
815   mnvertLocal(&((*covBending)(0,0)), numberOfHit, numberOfHit, numberOfHit, ifailBending);
816   Int_t ifailNonBending;
817   mnvertLocal(&((*covNonBending)(0,0)), numberOfHit, numberOfHit, numberOfHit, ifailNonBending);
818
819   // It would be worth trying to calculate the inverse of the covariance matrix
820   // only once per fit, since it cannot change much in principle,
821   // and it would save a lot of computing time !!!!
822   
823   // Calculates Chi2
824   if ((ifailBending == 0) && (ifailNonBending == 0)) {
825     // with Multiple Scattering if inversion correct
826     for (hitNumber1 = 0; hitNumber1 < numberOfHit ; hitNumber1++) { 
827       trackParamAtHit1 = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->UncheckedAt(hitNumber1));
828       hitForRec1 = trackParamAtHit1->GetHitForRecPtr();
829       hbc1 = hitForRec1->GetBendingCoor();
830       pbc1 = trackParamAtHit1->GetBendingCoor();
831       hnbc1 = hitForRec1->GetNonBendingCoor();
832       pnbc1 = trackParamAtHit1->GetNonBendingCoor();
833       for (hitNumber2 = 0; hitNumber2 < numberOfHit; hitNumber2++) {
834         trackParamAtHit2 = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->UncheckedAt(hitNumber2));
835         hitForRec2 = trackParamAtHit2->GetHitForRecPtr();
836         hbc2 = hitForRec2->GetBendingCoor();
837         pbc2 = trackParamAtHit2->GetBendingCoor();
838         hnbc2 = hitForRec2->GetNonBendingCoor();
839         pnbc2 = trackParamAtHit2->GetNonBendingCoor();
840         Chi2 += ((*covBending)(hitNumber2, hitNumber1) * (hbc1 - pbc1) * (hbc2 - pbc2)) +
841                 ((*covNonBending)(hitNumber2, hitNumber1) * (hnbc1 - pnbc1) * (hnbc2 - pnbc2));
842       }
843     }
844   } else {
845     // without Multiple Scattering if inversion impossible
846     for (hitNumber1 = 0; hitNumber1 < numberOfHit ; hitNumber1++) { 
847       trackParamAtHit1 = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->UncheckedAt(hitNumber1));
848       hitForRec1 = trackParamAtHit1->GetHitForRecPtr();
849       hbc1 = hitForRec1->GetBendingCoor();
850       pbc1 = trackParamAtHit1->GetBendingCoor();
851       hnbc1 = hitForRec1->GetNonBendingCoor();
852       pnbc1 = trackParamAtHit1->GetNonBendingCoor();
853       Chi2 += ((hbc1 - pbc1) * (hbc1 - pbc1) / hitForRec1->GetBendingReso2()) +
854               ((hnbc1 - pnbc1) * (hnbc1 - pnbc1) / hitForRec1->GetNonBendingReso2());
855     }
856   }
857   
858   delete covBending;
859   delete covNonBending;
860   delete [] msa2;
861 }
862
863 Double_t MultipleScatteringAngle2(AliMUONTrackParam *param)
864 {
865   /// Returns square of multiple Coulomb scattering angle
866   /// from TrackParamAtHit pointed to by "param"
867   Double_t slopeBending, slopeNonBending, radiationLength, inverseBendingMomentum2, inverseTotalMomentum2;
868   Double_t varMultipleScatteringAngle;
869   // Better implementation in AliMUONTrack class ????
870   slopeBending = param->GetBendingSlope();
871   slopeNonBending = param->GetNonBendingSlope();
872   // thickness in radiation length for the current track,
873   // taking local angle into account
874   radiationLength = AliMUONConstants::DefaultChamberThicknessInX0() *
875                     TMath::Sqrt(1.0 + slopeBending*slopeBending + slopeNonBending*slopeNonBending);
876   inverseBendingMomentum2 =  param->GetInverseBendingMomentum() * param->GetInverseBendingMomentum();
877   inverseTotalMomentum2 = inverseBendingMomentum2 * (1.0 + slopeBending * slopeBending) /
878                           (1.0 + slopeBending *slopeBending + slopeNonBending * slopeNonBending); 
879   varMultipleScatteringAngle = 0.0136 * (1.0 + 0.038 * TMath::Log(radiationLength));
880   // The velocity is assumed to be 1 !!!!
881   varMultipleScatteringAngle = inverseTotalMomentum2 * radiationLength * varMultipleScatteringAngle * varMultipleScatteringAngle;
882   return varMultipleScatteringAngle;
883 }
884
885 //______________________________________________________________________________
886  void mnvertLocal(Double_t *a, Int_t l, Int_t, Int_t n, Int_t &ifail)
887 {
888 ///*-*-*-*-*-*-*-*-*-*-*-*Inverts a symmetric matrix*-*-*-*-*-*-*-*-*-*-*-*-*
889 ///*-*                    ==========================
890 ///*-*        inverts a symmetric matrix.   matrix is first scaled to
891 ///*-*        have all ones on the diagonal (equivalent to change of units)
892 ///*-*        but no pivoting is done since matrix is positive-definite.
893 ///*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
894
895   // taken from TMinuit package of Root (l>=n)
896   // fVERTs, fVERTq and fVERTpp changed to localVERTs, localVERTq and localVERTpp
897   //  Double_t localVERTs[n], localVERTq[n], localVERTpp[n];
898   Double_t * localVERTs = new Double_t[n];
899   Double_t * localVERTq = new Double_t[n];
900   Double_t * localVERTpp = new Double_t[n];
901   // fMaxint changed to localMaxint
902   Int_t localMaxint = n;
903
904     /* System generated locals */
905     Int_t aOffset;
906
907     /* Local variables */
908     Double_t si;
909     Int_t i, j, k, kp1, km1;
910
911     /* Parameter adjustments */
912     aOffset = l + 1;
913     a -= aOffset;
914
915     /* Function Body */
916     ifail = 0;
917     if (n < 1) goto L100;
918     if (n > localMaxint) goto L100;
919 //*-*-                  scale matrix by sqrt of diag elements
920     for (i = 1; i <= n; ++i) {
921         si = a[i + i*l];
922         if (si <= 0) goto L100;
923         localVERTs[i-1] = 1 / TMath::Sqrt(si);
924     }
925     for (i = 1; i <= n; ++i) {
926         for (j = 1; j <= n; ++j) {
927             a[i + j*l] = a[i + j*l]*localVERTs[i-1]*localVERTs[j-1];
928         }
929     }
930 //*-*-                                       . . . start main loop . . . .
931     for (i = 1; i <= n; ++i) {
932         k = i;
933 //*-*-                  preparation for elimination step1
934         if (a[k + k*l] != 0) localVERTq[k-1] = 1 / a[k + k*l];
935         else goto L100;
936         localVERTpp[k-1] = 1;
937         a[k + k*l] = 0;
938         kp1 = k + 1;
939         km1 = k - 1;
940         if (km1 < 0) goto L100;
941         else if (km1 == 0) goto L50;
942         else               goto L40;
943 L40:
944         for (j = 1; j <= km1; ++j) {
945             localVERTpp[j-1] = a[j + k*l];
946             localVERTq[j-1]  = a[j + k*l]*localVERTq[k-1];
947             a[j + k*l]   = 0;
948         }
949 L50:
950         if (k - n < 0) goto L51;
951         else if (k - n == 0) goto L60;
952         else                goto L100;
953 L51:
954         for (j = kp1; j <= n; ++j) {
955             localVERTpp[j-1] = a[k + j*l];
956             localVERTq[j-1]  = -a[k + j*l]*localVERTq[k-1];
957             a[k + j*l]   = 0;
958         }
959 //*-*-                  elimination proper
960 L60:
961         for (j = 1; j <= n; ++j) {
962             for (k = j; k <= n; ++k) { a[j + k*l] += localVERTpp[j-1]*localVERTq[k-1]; }
963         }
964     }
965 //*-*-                  elements of left diagonal and unscaling
966     for (j = 1; j <= n; ++j) {
967         for (k = 1; k <= j; ++k) {
968             a[k + j*l] = a[k + j*l]*localVERTs[k-1]*localVERTs[j-1];
969             a[j + k*l] = a[k + j*l];
970         }
971     }
972     delete [] localVERTs;
973     delete [] localVERTq;
974     delete [] localVERTpp;
975     return;
976 //*-*-                  failure return
977 L100:
978     delete [] localVERTs;
979     delete [] localVERTq;
980     delete [] localVERTpp;
981     ifail = 1;
982 } /* mnvertLocal */
983
984   //__________________________________________________________________________
985 void AliMUONTrackReconstructor::RemoveDoubleTracks(void)
986 {
987   /// To remove double tracks.
988   /// Tracks are considered identical
989   /// if they have at least half of their hits in common.
990   /// Among two identical tracks, one keeps the track with the larger number of hits
991   /// or, if these numbers are equal, the track with the minimum Chi2.
992   AliMUONTrack *track1, *track2, *trackToRemove;
993   Int_t hitsInCommon, nHits1, nHits2;
994   Bool_t removedTrack1;
995   // Loop over first track of the pair
996   track1 = (AliMUONTrack*) fRecTracksPtr->First();
997   while (track1) {
998     removedTrack1 = kFALSE;
999     nHits1 = track1->GetNTrackHits();
1000     // Loop over second track of the pair
1001     track2 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1002     while (track2) {
1003       nHits2 = track2->GetNTrackHits();
1004       // number of hits in common between two tracks
1005       hitsInCommon = track1->HitsInCommon(track2);
1006       // check for identical tracks
1007       if ((4 * hitsInCommon) >= (nHits1 + nHits2)) {
1008         // decide which track to remove
1009         if ((nHits1 > nHits2) || ((nHits1 == nHits2) && (track1->GetFitFMin() < track2->GetFitFMin()))) {
1010           // remove track2 and continue the second loop with the track next to track2
1011           trackToRemove = track2;
1012           track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
1013           fRecTracksPtr->Remove(trackToRemove);
1014           fNRecTracks--;
1015           fRecTracksPtr->Compress(); // this is essential to retrieve the TClonesArray afterwards
1016         } else {
1017           // else remove track1 and continue the first loop with the track next to track1
1018           trackToRemove = track1;
1019           track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1020           fRecTracksPtr->Remove(trackToRemove);
1021           fNRecTracks--;
1022           fRecTracksPtr->Compress(); // this is essential to retrieve the TClonesArray afterwards
1023           removedTrack1 = kTRUE;
1024           break;
1025         }
1026       } else track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
1027     } // track2
1028     if (removedTrack1) continue;
1029     track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1030   } // track1
1031   return;
1032 }
1033
1034   //__________________________________________________________________________
1035 void AliMUONTrackReconstructor::UpdateHitForRecAtHit()
1036 {
1037   /// Set cluster parameters after track fitting. Fill fHitForRecAtHit of AliMUONTrack's
1038   AliMUONTrack *track;
1039   AliMUONTrackParam *trackParamAtHit;
1040   track = (AliMUONTrack*) fRecTracksPtr->First();
1041   while (track) {
1042     trackParamAtHit = (AliMUONTrackParam*) (track->GetTrackParamAtHit()->First());
1043     while (trackParamAtHit) {
1044       track->AddHitForRecAtHit(trackParamAtHit->GetHitForRecPtr());
1045       trackParamAtHit = (AliMUONTrackParam*) (track->GetTrackParamAtHit()->After(trackParamAtHit)); 
1046     }
1047     track = (AliMUONTrack*) fRecTracksPtr->After(track);
1048   }
1049   return;
1050 }
1051
1052   //__________________________________________________________________________
1053 void AliMUONTrackReconstructor::EventDump(void)
1054 {
1055   /// Dump reconstructed event (track parameters at vertex and at first hit),
1056   /// and the particle parameters
1057   AliMUONTrack *track;
1058   AliMUONTrackParam *trackParam, *trackParam1;
1059   Double_t bendingSlope, nonBendingSlope, pYZ;
1060   Double_t pX, pY, pZ, x, y, z, c;
1061   Int_t trackIndex, nTrackHits;
1062  
1063   AliDebug(1,"****** enter EventDump ******");
1064   AliDebug(1, Form("Number of Reconstructed tracks : %d", fNRecTracks)); 
1065   
1066   fRecTracksPtr->Compress(); // for simple loop without "Next" since no hole
1067   // Loop over reconstructed tracks
1068   for (trackIndex = 0; trackIndex < fNRecTracks; trackIndex++) {
1069     AliDebug(1, Form("track number: %d", trackIndex));
1070     // function for each track for modularity ????
1071     track = (AliMUONTrack*) ((*fRecTracksPtr)[trackIndex]);
1072     nTrackHits = track->GetNTrackHits();
1073     AliDebug(1, Form("Number of track hits: %d ", nTrackHits));
1074     // track parameters at Vertex
1075     trackParam = track->GetTrackParamAtVertex();
1076     x = trackParam->GetNonBendingCoor();
1077     y = trackParam->GetBendingCoor();
1078     z = trackParam->GetZ();
1079     bendingSlope = trackParam->GetBendingSlope();
1080     nonBendingSlope = trackParam->GetNonBendingSlope();
1081     pYZ = 1/TMath::Abs(trackParam->GetInverseBendingMomentum());
1082     pZ = pYZ/TMath::Sqrt(1+bendingSlope*bendingSlope);
1083     pX = pZ * nonBendingSlope;
1084     pY = pZ * bendingSlope;
1085     c = TMath::Sign(1.0, trackParam->GetInverseBendingMomentum());
1086     AliDebug(1, Form("Track parameters at Vertex z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1087                      z, x, y, pX, pY, pZ, c));
1088
1089     // track parameters at first hit
1090     trackParam1 = (AliMUONTrackParam*) track->GetTrackParamAtHit()->First();
1091     x = trackParam1->GetNonBendingCoor();
1092     y = trackParam1->GetBendingCoor();
1093     z = trackParam1->GetZ();
1094     bendingSlope = trackParam1->GetBendingSlope();
1095     nonBendingSlope = trackParam1->GetNonBendingSlope();
1096     pYZ = 1/TMath::Abs(trackParam1->GetInverseBendingMomentum());
1097     pZ = pYZ/TMath::Sqrt(1.0 + bendingSlope * bendingSlope);
1098     pX = pZ * nonBendingSlope;
1099     pY = pZ * bendingSlope;
1100     c = TMath::Sign(1.0, trackParam1->GetInverseBendingMomentum());
1101     AliDebug(1, Form("track parameters at z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1102                      z, x, y, pX, pY, pZ, c));
1103   }
1104   // informations about generated particles NO !!!!!!!!
1105   
1106 //    for (Int_t iPart = 0; iPart < np; iPart++) {
1107 //      p = gAlice->Particle(iPart);
1108 //      printf(" particle %d: type= %d px= %f py= %f pz= %f pdg= %d\n",
1109 //         iPart, p->GetPdgCode(), p->Px(), p->Py(), p->Pz(), p->GetPdgCode());    
1110 //    }
1111   return;
1112 }
1113
1114