]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONTrackReconstructor.cxx
Adding AliMUONGeometryMisAligner
[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 ////////////////////////////////////
19 //
20 // MUON track reconstructor in ALICE (class renamed from AliMUONEventReconstructor)
21 //
22 // This class contains as data:
23 // * the parameters for the track reconstruction
24 // * a pointer to the array of hits to be reconstructed (the event)
25 // * a pointer to the array of segments made with these hits inside each station
26 // * a pointer to the array of reconstructed tracks
27 //
28 // It contains as methods, among others:
29 // * MakeEventToBeReconstructed to build the array of hits to be reconstructed
30 // * MakeSegments to build the segments
31 // * MakeTracks to build the tracks
32 //
33 ////////////////////////////////////
34
35 #include <stdlib.h>
36 #include <Riostream.h>
37 #include <TDirectory.h>
38 #include <TFile.h>
39 #include <TMatrixD.h>
40
41 #include "AliMUONTrackReconstructor.h"
42 #include "AliMUON.h"
43 #include "AliMUONConstants.h"
44 #include "AliMUONHitForRec.h"
45 #include "AliMUONTriggerTrack.h"
46 #include "AliMUONTriggerCircuit.h"
47 #include "AliMUONRawCluster.h"
48 #include "AliMUONLocalTrigger.h"
49 #include "AliMUONGlobalTrigger.h"
50 #include "AliMUONSegment.h"
51 #include "AliMUONTrack.h"
52 #include "AliMUONTrackHit.h"
53 #include "AliMagF.h"
54 #include "AliRun.h" // for gAlice
55 #include "AliRunLoader.h"
56 #include "AliLoader.h"
57 #include "AliMUONTrackK.h" 
58 #include "AliMC.h"
59 #include "AliLog.h"
60 #include "AliTrackReference.h"
61
62 //************* Defaults parameters for reconstruction
63 const Double_t AliMUONTrackReconstructor::fgkDefaultMinBendingMomentum = 3.0;
64 const Double_t AliMUONTrackReconstructor::fgkDefaultMaxBendingMomentum = 3000.0;
65 const Double_t AliMUONTrackReconstructor::fgkDefaultMaxChi2 = 100.0;
66 const Double_t AliMUONTrackReconstructor::fgkDefaultMaxSigma2Distance = 16.0;
67 const Double_t AliMUONTrackReconstructor::fgkDefaultBendingResolution = 0.01;
68 const Double_t AliMUONTrackReconstructor::fgkDefaultNonBendingResolution = 0.144;
69 const Double_t AliMUONTrackReconstructor::fgkDefaultChamberThicknessInX0 = 0.03;
70 // Simple magnetic field:
71 // Value taken from macro MUONtracking.C: 0.7 T, hence 7 kG
72 // Length and Position from reco_muon.F, with opposite sign:
73 // Length = ZMAGEND-ZCOIL
74 // Position = (ZMAGEND+ZCOIL)/2
75 // to be ajusted differently from real magnetic field ????
76 const Double_t AliMUONTrackReconstructor::fgkDefaultSimpleBValue = 7.0;
77 const Double_t AliMUONTrackReconstructor::fgkDefaultSimpleBLength = 428.0;
78 const Double_t AliMUONTrackReconstructor::fgkDefaultSimpleBPosition = 1019.0;
79 const Int_t    AliMUONTrackReconstructor::fgkDefaultRecTrackRefHits = 0;
80 const Double_t AliMUONTrackReconstructor::fgkDefaultEfficiency = 0.95;
81
82 ClassImp(AliMUONTrackReconstructor) // Class implementation in ROOT context
83
84   //__________________________________________________________________________
85 AliMUONTrackReconstructor::AliMUONTrackReconstructor(AliLoader* loader)
86   : TObject()
87 {
88   // Constructor for class AliMUONTrackReconstructor
89   SetReconstructionParametersToDefaults();
90   fTrackMethod = 1; //AZ - tracking method (1-default, 2-Kalman)
91   // Memory allocation for the TClonesArray of hits for reconstruction
92   // Is 10000 the right size ????
93   fHitsForRecPtr = new TClonesArray("AliMUONHitForRec", 10000);
94   fNHitsForRec = 0; // really needed or GetEntriesFast sufficient ????
95   // Memory allocation for the TClonesArray's of segments in stations
96   // Is 2000 the right size ????
97   for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) {
98     fSegmentsPtr[st] = new TClonesArray("AliMUONSegment", 2000);
99     fNSegments[st] = 0; // really needed or GetEntriesFast sufficient ????
100   }
101   // Memory allocation for the TClonesArray of reconstructed tracks
102   // Is 10 the right size ????
103   fRecTracksPtr = new TClonesArray("AliMUONTrack", 10);
104   fNRecTracks = 0; // really needed or GetEntriesFast sufficient ????
105   // Memory allocation for the TClonesArray of hits on reconstructed tracks
106   // Is 100 the right size ????
107   fRecTrackHitsPtr = new TClonesArray("AliMUONTrackHit", 100);
108   fNRecTrackHits = 0; // really needed or GetEntriesFast sufficient ????
109
110   // Sign of fSimpleBValue according to sign of Bx value at (50,50,-950).
111   Float_t b[3], x[3];
112   x[0] = 50.; x[1] = 50.; x[2] = -950.;
113   gAlice->Field()->Field(x, b);
114   fSimpleBValue = TMath::Sign(fSimpleBValue,(Double_t) b[0]);
115   fSimpleBPosition = TMath::Sign(fSimpleBPosition,(Double_t) x[2]);
116   // See how to get fSimple(BValue, BLength, BPosition)
117   // automatically calculated from the actual magnetic field ????
118
119   AliDebug(1,"AliMUONTrackReconstructor constructed with defaults"); 
120   if ( AliLog::GetGlobalDebugLevel()>0) Dump();
121   AliDebug(1,"Magnetic field from root file:");
122   if ( AliLog::GetGlobalDebugLevel()>0) gAlice->Field()->Dump();
123
124   
125   // Initializions for track ref. background events
126   fBkgTrackRefFile = 0;
127   fBkgTrackRefTK = 0;
128   fBkgTrackRefParticles = 0;
129   fBkgTrackRefTTR = 0;
130   fBkgTrackRefEventNumber = -1;
131    
132   // initialize loader's
133   fLoader = loader;
134
135   // initialize container
136   fMUONData  = new AliMUONData(fLoader,"MUON","MUON");
137
138   return;
139 }
140   //__________________________________________________________________________
141 AliMUONTrackReconstructor::AliMUONTrackReconstructor (const AliMUONTrackReconstructor& rhs)
142   : TObject(rhs)
143 {
144 // Protected copy constructor
145
146   AliFatal("Not implemented.");
147 }
148
149 AliMUONTrackReconstructor & 
150 AliMUONTrackReconstructor::operator=(const AliMUONTrackReconstructor& rhs)
151 {
152 // Protected assignement operator
153
154   if (this == &rhs) return *this;
155
156   AliFatal("Not implemented.");
157     
158   return *this;  
159 }
160
161   //__________________________________________________________________________
162 AliMUONTrackReconstructor::~AliMUONTrackReconstructor(void)
163 {
164   // Destructor for class AliMUONTrackReconstructor
165   delete fHitsForRecPtr; // Correct destruction of everything ???? or delete [] ????
166   for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++)
167     delete fSegmentsPtr[st]; // Correct destruction of everything ????
168   return;
169 }
170
171   //__________________________________________________________________________
172 void AliMUONTrackReconstructor::SetReconstructionParametersToDefaults(void)
173 {
174   // Set reconstruction parameters to default values
175   // Would be much more convenient with a structure (or class) ????
176   fMinBendingMomentum = fgkDefaultMinBendingMomentum;
177   fMaxBendingMomentum = fgkDefaultMaxBendingMomentum;
178   fMaxChi2 = fgkDefaultMaxChi2;
179   fMaxSigma2Distance = fgkDefaultMaxSigma2Distance;
180
181   // ******** Parameters for making HitsForRec
182   // minimum radius,
183   // like in TRACKF_STAT:
184   // 2 degrees for stations 1 and 2, or ch(0...) from 0 to 3;
185   // 30 cm for stations 3 to 5, or ch(0...) from 4 to 9
186   for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
187     if (ch < 4) fRMin[ch] = TMath::Abs(AliMUONConstants::DefaultChamberZ(ch)) *
188                   2.0 * TMath::Pi() / 180.0;
189     else fRMin[ch] = 30.0;
190     // maximum radius at 10 degrees and Z of chamber
191     fRMax[ch] = TMath::Abs(AliMUONConstants::DefaultChamberZ(ch)) *
192                   10.0 * TMath::Pi() / 180.0;
193   }
194
195   // ******** Parameters for making segments
196   // should be parametrized ????
197   // according to interval between chambers in a station ????
198   // Maximum distance in non bending plane
199   // 5 * 0.22 just to remember the way it was made in TRACKF_STAT
200   // SIGCUT*DYMAX(IZ)
201   for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++)
202     fSegmentMaxDistNonBending[st] = 5. * 0.22;
203   // Maximum distance in bending plane:
204   // values from TRACKF_STAT, corresponding to (J psi 20cm),
205   // scaled to the real distance between chambers in a station
206   fSegmentMaxDistBending[0] = TMath::Abs( 1.5 *
207                                           (AliMUONConstants::DefaultChamberZ(1) - AliMUONConstants::DefaultChamberZ(0)) / 20.0);
208   fSegmentMaxDistBending[1] = TMath::Abs( 1.5 *
209                                           (AliMUONConstants::DefaultChamberZ(3) - AliMUONConstants::DefaultChamberZ(2)) / 20.0);
210   fSegmentMaxDistBending[2] = TMath::Abs( 3.0 *
211                                           (AliMUONConstants::DefaultChamberZ(5) - AliMUONConstants::DefaultChamberZ(4)) / 20.0);
212   fSegmentMaxDistBending[3] = TMath::Abs( 6.0 *
213                                           (AliMUONConstants::DefaultChamberZ(7) - AliMUONConstants::DefaultChamberZ(6)) / 20.0);
214   fSegmentMaxDistBending[4] = TMath::Abs( 6.0 *
215                                           (AliMUONConstants::DefaultChamberZ(9) - AliMUONConstants::DefaultChamberZ(8)) / 20.0);
216
217   
218   fBendingResolution = fgkDefaultBendingResolution;
219   fNonBendingResolution = fgkDefaultNonBendingResolution;
220   fChamberThicknessInX0 = fgkDefaultChamberThicknessInX0;
221   fSimpleBValue = fgkDefaultSimpleBValue;
222   fSimpleBLength = fgkDefaultSimpleBLength;
223   fSimpleBPosition = fgkDefaultSimpleBPosition;
224   fRecTrackRefHits = fgkDefaultRecTrackRefHits;
225   fEfficiency = fgkDefaultEfficiency;
226   return;
227 }
228
229 //__________________________________________________________________________
230 Double_t AliMUONTrackReconstructor::GetImpactParamFromBendingMomentum(Double_t BendingMomentum) const
231 {
232   // Returns impact parameter at vertex in bending plane (cm),
233   // from the signed bending momentum "BendingMomentum" in bending plane (GeV/c),
234   // using simple values for dipole magnetic field.
235   // The sign of "BendingMomentum" is the sign of the charge.
236   return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
237           BendingMomentum);
238 }
239
240 //__________________________________________________________________________
241 Double_t AliMUONTrackReconstructor::GetBendingMomentumFromImpactParam(Double_t ImpactParam) const
242 {
243   // Returns signed bending momentum in bending plane (GeV/c),
244   // the sign being the sign of the charge for particles moving forward in Z,
245   // from the impact parameter "ImpactParam" at vertex in bending plane (cm),
246   // using simple values for dipole magnetic field.
247   return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
248           ImpactParam);
249 }
250
251 //__________________________________________________________________________
252 void AliMUONTrackReconstructor::SetBkgTrackRefFile(Text_t *BkgTrackRefFileName)
253 {
254   // Set background file ... for track ref. hits
255   // Must be called after having loaded the firts signal event
256   AliDebug(1,Form("Enter SetBkgTrackRefFile with BkgTrackRefFileName %s",BkgTrackRefFileName));
257   
258   if (strlen(BkgTrackRefFileName)) {
259     // BkgTrackRefFileName not empty: try to open the file
260
261     if(AliLog::GetGlobalDebugLevel()>1) {
262       cout << "Before File(Bkg)" << endl; gDirectory->Dump();
263     }
264     fBkgTrackRefFile = new TFile(BkgTrackRefFileName);
265     if(AliLog::GetGlobalDebugLevel()>1) {
266       cout << "After File(Bkg)" << endl; gDirectory->Dump();
267     }
268     if (fBkgTrackRefFile-> IsOpen()) {
269       if(AliLog::GetGlobalDebugLevel()>0) {
270         cout << "Background for Track ref. hits in file: ``" << BkgTrackRefFileName
271              << "'' successfully opened" << endl;}
272     }
273     else {
274       cout << "Background for Track Ref. hits in file: " << BkgTrackRefFileName << endl;
275       cout << "NOT FOUND: EXIT" << endl;
276       exit(0); // right instruction for exit ????
277     }
278     // Arrays for "particles" and "hits"
279     fBkgTrackRefParticles = new TClonesArray("TParticle", 200);
280     // Event number to -1 for initialization
281     fBkgTrackRefEventNumber = -1;
282     // Back to the signal file:
283     // first signal event must have been loaded previously,
284     // otherwise, Segmentation violation at the next instruction
285     // How is it possible to do smething better ????
286     ((gAlice->TreeK())->GetCurrentFile())->cd();
287     if(AliLog::GetGlobalDebugLevel()>1) cout << "After cd(gAlice)" << endl; gDirectory->Dump();
288   }
289   return;
290 }
291
292 //__________________________________________________________________________
293 void AliMUONTrackReconstructor::NextBkgTrackRefEvent(void)
294 {
295   // Get next event in background file for track ref. hits
296   // Goes back to event number 0 when end of file is reached
297   char treeName[20];
298   
299   AliDebug(1,"Enter NextBkgTrackRefEvent");
300   // Clean previous event
301   if(fBkgTrackRefTK) delete fBkgTrackRefTK;
302   fBkgTrackRefTK = NULL;
303   if(fBkgTrackRefParticles) fBkgTrackRefParticles->Clear();
304   if(fBkgTrackRefTTR) delete fBkgTrackRefTTR;
305   fBkgTrackRefTTR = NULL;
306   // Increment event number
307   fBkgTrackRefEventNumber++;
308   // Get access to Particles and Hits for event from background file
309   if (AliLog::GetGlobalDebugLevel()>1) cout << "Before cd(Bkg)" << endl; gDirectory->Dump();
310   fBkgTrackRefFile->cd();
311   if (AliLog::GetGlobalDebugLevel()>1) cout << "After cd(Bkg)" << endl; gDirectory->Dump();
312   // Particles: TreeK for event and branch "Particles"
313   sprintf(treeName, "TreeK%d", fBkgTrackRefEventNumber);
314   fBkgTrackRefTK = (TTree*)gDirectory->Get(treeName);
315   if (!fBkgTrackRefTK) {
316    
317     AliDebug(1,Form("Cannot find Kine Tree for background event: %d",fBkgTrackRefEventNumber));
318     AliDebug(1,"Goes back to event 0");
319     
320     fBkgTrackRefEventNumber = 0;
321     sprintf(treeName, "TreeK%d", fBkgTrackRefEventNumber);
322     fBkgTrackRefTK = (TTree*)gDirectory->Get(treeName);
323     if (!fBkgTrackRefTK) {
324                 AliError(Form("cannot find Kine Tree for background event: %d",fBkgTrackRefEventNumber));
325         exit(0);
326     }
327   }
328   if (fBkgTrackRefTK) 
329     fBkgTrackRefTK->SetBranchAddress("Particles", &fBkgTrackRefParticles);
330   fBkgTrackRefTK->GetEvent(0); // why event 0 ???? necessary ????
331   // Hits: TreeH for event and branch "MUON"
332   sprintf(treeName, "TreeTR%d", fBkgTrackRefEventNumber);
333   fBkgTrackRefTTR = (TTree*)gDirectory->Get(treeName);
334   if (!fBkgTrackRefTTR) {
335           AliError(Form("cannot find Hits Tree for background event: %d",fBkgTrackRefEventNumber));
336       exit(0);
337   }
338   fBkgTrackRefTTR->GetEntries(); // necessary ????
339   // Back to the signal file
340   ((gAlice->TreeK())->GetCurrentFile())->cd();
341   if (AliLog::GetGlobalDebugLevel()>1) 
342     cout << "After cd(gAlice)" << endl; gDirectory->Dump();
343   return;
344 }
345
346 //__________________________________________________________________________
347 void AliMUONTrackReconstructor::EventReconstruct(void)
348 {
349   // To reconstruct one event
350   AliDebug(1,"Enter EventReconstruct");
351   MakeEventToBeReconstructed();
352   MakeSegments();
353   MakeTracks();
354   if (fMUONData->IsTriggerTrackBranchesInTree()) 
355     ValidateTracksWithTrigger(); 
356
357   // Add tracks to MUON data container 
358   for(Int_t i=0; i<GetNRecTracks(); i++) {
359     AliMUONTrack * track = (AliMUONTrack*) GetRecTracksPtr()->At(i);
360     fMUONData->AddRecTrack(*track);
361   }
362
363   return;
364 }
365
366 //__________________________________________________________________________
367 void AliMUONTrackReconstructor::EventReconstructTrigger(void)
368 {
369   // To reconstruct one event
370   AliDebug(1,"Enter EventReconstructTrigger");
371   MakeTriggerTracks();  
372   return;
373 }
374
375   //__________________________________________________________________________
376 void AliMUONTrackReconstructor::ResetHitsForRec(void)
377 {
378   // To reset the array and the number of HitsForRec,
379   // and also the number of HitsForRec
380   // and the index of the first HitForRec per chamber
381   if (fHitsForRecPtr) fHitsForRecPtr->Clear();
382   fNHitsForRec = 0;
383   for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++)
384     fNHitsForRecPerChamber[ch] = fIndexOfFirstHitForRecPerChamber[ch] = 0;
385   return;
386 }
387
388   //__________________________________________________________________________
389 void AliMUONTrackReconstructor::ResetSegments(void)
390 {
391   // To reset the TClonesArray of segments and the number of Segments
392   // for all stations
393   for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) {
394     if (fSegmentsPtr[st]) fSegmentsPtr[st]->Clear();
395     fNSegments[st] = 0;
396   }
397   return;
398 }
399
400   //__________________________________________________________________________
401 void AliMUONTrackReconstructor::ResetTracks(void)
402 {
403   // To reset the TClonesArray of reconstructed tracks
404   if (fRecTracksPtr) fRecTracksPtr->Delete();
405   // Delete in order that the Track destructors are called,
406   // hence the space for the TClonesArray of pointers to TrackHit's is freed
407   fNRecTracks = 0;
408   return;
409 }
410
411
412   //__________________________________________________________________________
413 void AliMUONTrackReconstructor::ResetTrackHits(void)
414 {
415   // To reset the TClonesArray of hits on reconstructed tracks
416   if (fRecTrackHitsPtr) fRecTrackHitsPtr->Clear();
417   fNRecTrackHits = 0;
418   return;
419 }
420
421   //__________________________________________________________________________
422 void AliMUONTrackReconstructor::MakeEventToBeReconstructed(void)
423 {
424   // To make the list of hits to be reconstructed,
425   // either from the track ref. hits or from the raw clusters
426   // according to the parameter set for the reconstructor
427 //   TString evfoldname = AliConfig::GetDefaultEventFolderName();//to be interfaced properly
428   
429 //   AliRunLoader* rl = AliRunLoader::GetRunLoader(evfoldname);
430 //   if (rl == 0x0)
431 //    {
432 //      Error("MakeEventToBeReconstructed",
433 //            "Can not find Run Loader in Event Folder named %s.",
434 //            evfoldname.Data());
435 //      return;
436 //    }
437 //   AliLoader* gime = rl->GetLoader("MUONLoader");
438 //   if (gime == 0x0)
439 //    {
440 //      Error("MakeEventToBeReconstructed","Can not get MUON Loader from Run Loader.");
441 //      return;
442 //    }
443   AliRunLoader *runLoader = fLoader->GetRunLoader();
444
445   AliDebug(1,"Enter MakeEventToBeReconstructed");
446   ResetHitsForRec();
447   if (fRecTrackRefHits == 1) {
448     // Reconstruction from track ref. hits
449     // Back to the signal file
450     TTree* treeTR = runLoader->TreeTR();
451     if (treeTR == 0x0)
452       {
453         Int_t retval = runLoader->LoadTrackRefs("READ");
454         if ( retval)
455           {
456             AliError("Error occured while loading hits.");
457             return;
458           }
459         treeTR = runLoader->TreeTR();
460         if (treeTR == 0x0)
461           {
462             AliError("Can not get TreeTR");
463             return;
464           }
465       }
466     
467     AddHitsForRecFromTrackRef(treeTR,1);
468     
469     // Background hits
470     AddHitsForRecFromTrackRef(fBkgTrackRefTTR,0);
471     // Sort HitsForRec in increasing order with respect to chamber number
472     SortHitsForRecWithIncreasingChamber();
473   }
474   else {
475     // Reconstruction from raw clusters
476     // AliMUON *MUON  = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
477     // Security on MUON ????
478     // TreeR assumed to be be "prepared" in calling function
479     // by "MUON->GetTreeR(nev)" ????
480     TTree *treeR = fLoader->TreeR();
481     //AZ? fMUONData->SetTreeAddress("RC");
482     AddHitsForRecFromRawClusters(treeR);
483     // No sorting: it is done automatically in the previous function
484   }
485  
486   AliDebug(1,"End of MakeEventToBeReconstructed");
487     if (AliLog::GetGlobalDebugLevel() > 0) {
488       AliDebug(1, Form("NHitsForRec: %d",fNHitsForRec));
489       for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
490         AliDebug(1, Form("Chamber(0...): %d",ch));
491         AliDebug(1, Form("NHitsForRec: %d", fNHitsForRecPerChamber[ch]));
492         AliDebug(1, Form("Index(first HitForRec): %d", fIndexOfFirstHitForRecPerChamber[ch]));
493         for (Int_t hit = fIndexOfFirstHitForRecPerChamber[ch];
494              hit < fIndexOfFirstHitForRecPerChamber[ch] + fNHitsForRecPerChamber[ch];
495              hit++) {
496           AliDebug(1, Form("HitForRec index(0...): %d",hit));
497           ((*fHitsForRecPtr)[hit])->Dump();
498       }
499     }
500   }
501   return;
502 }
503
504   //__________________________________________________________________________
505 void AliMUONTrackReconstructor::AddHitsForRecFromTrackRef(TTree *TTR, Int_t signal)
506 {
507   // To add to the list of hits for reconstruction
508   // the signal hits from a track reference tree TreeTR.
509   TClonesArray *listOfTrackRefs = NULL;
510   AliTrackReference *trackRef;
511   
512   Int_t track = 0;
513   AliDebug(2,Form("Enter AddHitsForRecFromTrackRef with TTR: %d",  TTR));
514   if (TTR == NULL) return;
515   
516   listOfTrackRefs = CleanTrackRefs(TTR);
517
518   Int_t ntracks = listOfTrackRefs->GetEntriesFast();
519
520   AliDebug(2,Form("ntracks: %d", ntracks));
521
522   for (Int_t index = 0; index < ntracks; index++) {
523     trackRef = (AliTrackReference*) listOfTrackRefs->At(index);
524     track =  trackRef->GetTrack();
525     
526     NewHitForRecFromTrackRef(trackRef,track,signal);
527   } // end of track ref.
528   
529   listOfTrackRefs->Delete();
530   delete listOfTrackRefs;
531   return;
532 }
533
534
535   //__________________________________________________________________________
536 AliMUONHitForRec* AliMUONTrackReconstructor::NewHitForRecFromTrackRef(AliTrackReference* Hit, Int_t TrackNumber, Int_t Signal)
537 {
538   // To make a new hit for reconstruction from a track ref. hit pointed to by "Hit",
539   // with  the track numbered "TrackNumber",
540   // either from signal ("Signal" = 1) or background ("Signal" = 0) event.
541   // Selects hits in tracking (not trigger) chambers.
542   // Takes into account the efficiency (fEfficiency)
543   // and the smearing from resolution (fBendingResolution and fNonBendingResolution).
544   // Adds a condition on the radius between RMin and RMax
545   // to better simulate the real chambers.
546   // Returns the pointer to the new hit for reconstruction,
547   // or NULL in case of inefficiency or non tracking chamber or bad radius.
548   // No condition on at most 20 cm from a muon from a resonance
549   // like in Fortran TRACKF_STAT.
550   AliMUONHitForRec* hitForRec;
551   Double_t bendCoor, nonBendCoor, radius;
552   Int_t chamber = AliMUONConstants::ChamberNumber(Hit->Z()); // chamber(0...)
553   if (chamber < 0) return NULL;
554   // only in tracking chambers (fChamber starts at 1)
555   if (chamber >= AliMUONConstants::NTrackingCh()) return NULL;
556   // only if hit is efficient (keep track for checking ????)
557   if (gRandom->Rndm() > fEfficiency) return NULL;
558   // only if radius between RMin and RMax
559   bendCoor = Hit->Y();
560   nonBendCoor = Hit->X();
561   radius = TMath::Sqrt((bendCoor * bendCoor) + (nonBendCoor * nonBendCoor));
562   // This cut is not needed with a realistic chamber geometry !!!!
563 //   if ((radius < fRMin[chamber]) || (radius > fRMax[chamber])) return NULL;
564   // new AliMUONHitForRec from track ref. hit and increment number of AliMUONHitForRec's
565   hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(Hit);
566   fNHitsForRec++;
567   // add smearing from resolution
568   hitForRec->SetBendingCoor(bendCoor + gRandom->Gaus(0., fBendingResolution));
569   hitForRec->SetNonBendingCoor(nonBendCoor
570                                + gRandom->Gaus(0., fNonBendingResolution));
571 //   // !!!! without smearing
572 //   hitForRec->SetBendingCoor(bendCoor);
573 //   hitForRec->SetNonBendingCoor(nonBendCoor);
574   // more information into HitForRec
575   //  resolution: angular effect to be added here ????
576   hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
577   hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
578   //   track ref. info
579   hitForRec->SetTTRTrack(TrackNumber);
580   hitForRec->SetTrackRefSignal(Signal);
581   if (AliLog::GetGlobalDebugLevel()> 1) {
582     AliDebug(2,Form("Track: %d", TrackNumber));
583     Hit->Dump();
584     cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
585     hitForRec->Dump();
586   }
587   return hitForRec;
588 }
589   //__________________________________________________________________________
590 TClonesArray* AliMUONTrackReconstructor::CleanTrackRefs(TTree *treeTR)
591 {
592   // Make hits from track ref..  
593   // Re-calculate hits parameters because two AliTrackReferences are recorded for
594   // each chamber (one when particle is entering + one when particle is leaving 
595   // the sensitive volume) 
596  
597   AliTrackReference *trackReference;
598   Float_t x1, y1, z1, pX1, pY1, pZ1;
599   Float_t x2, y2, z2, pX2, pY2, pZ2;
600   Int_t track1, track2;
601   Int_t nRec = 0;
602   Float_t maxGasGap = 1.; // cm 
603   Int_t iHit1 = 0;
604   Int_t iHitMin = 0;
605  
606   AliTrackReference *trackReferenceNew = new AliTrackReference();
607   
608   TClonesArray* trackRefs = new TClonesArray("AliTrackReference", 10);
609   TClonesArray* cleanTrackRefs = new TClonesArray("AliTrackReference", 10);
610
611   if (treeTR == NULL) return NULL;
612   TBranch* branch = treeTR->GetBranch("MUON");
613   if (branch == NULL) return NULL;
614   branch->SetAddress(&trackRefs);
615
616   Int_t nTrack = (Int_t)treeTR->GetEntries();
617   for (Int_t iTrack  = 0; iTrack < nTrack; iTrack++) {
618     treeTR->GetEntry(iTrack);
619     iHitMin = 0;
620     iHit1 = 0;
621     while (iHit1 < trackRefs->GetEntries()) {
622       trackReference = (AliTrackReference*)trackRefs->At(iHit1);
623       x1 = trackReference->X();
624       y1 = trackReference->Y();
625       z1 = trackReference->Z();
626       pX1 = trackReference->Px();
627       pY1 = trackReference->Py();
628       pZ1 = trackReference->Pz();
629       track1 = trackReference->GetTrack();
630       nRec = 1;
631       iHitMin = iHit1+1;
632       for (Int_t iHit2 = iHit1+1; iHit2 < trackRefs->GetEntries(); iHit2++) {
633         trackReference = (AliTrackReference*)trackRefs->At(iHit2);
634         x2 = trackReference->X();
635         y2 = trackReference->Y();
636         z2 = trackReference->Z();
637         pX2 = trackReference->Px();
638         pY2 = trackReference->Py();
639         pZ2 = trackReference->Pz();
640         track2 = trackReference->GetTrack();
641         if (track2 == track1 && TMath::Abs(z2-z1) < maxGasGap ) {
642           nRec++;
643           x1 += x2;
644           y1 += y2;
645           z1 += z2;                             
646           pX1 += pX2;
647           pY1 += pY2;
648           pZ1 += pZ2;
649           iHitMin = iHit2+1;
650         }
651         
652       } // for iHit2
653       x1 /= (Float_t)nRec;
654       y1 /= (Float_t)nRec;
655       z1 /= (Float_t)nRec;
656       pX1 /= (Float_t)nRec;
657       pY1 /= (Float_t)nRec;
658       pZ1 /= (Float_t)nRec;
659       trackReferenceNew->SetPosition(x1,y1,z1);
660       trackReferenceNew->SetMomentum(pX1,pY1,pZ1);
661       trackReferenceNew->SetTrack(track1);
662       {new ((*cleanTrackRefs)[cleanTrackRefs->GetEntriesFast()]) AliTrackReference(*trackReferenceNew);}       
663       iHit1 = iHitMin;
664     } // while iHit1
665   } // for track
666
667   trackRefs->Delete();
668   delete trackRefs;
669   delete trackReferenceNew;
670   return cleanTrackRefs;
671  
672 }
673   //__________________________________________________________________________
674 void AliMUONTrackReconstructor::SortHitsForRecWithIncreasingChamber()
675 {
676   // Sort HitsForRec's in increasing order with respect to chamber number.
677   // Uses the function "Compare".
678   // Update the information for HitsForRec per chamber too.
679   Int_t ch, nhits, prevch;
680   fHitsForRecPtr->Sort();
681   for (ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
682     fNHitsForRecPerChamber[ch] = 0;
683     fIndexOfFirstHitForRecPerChamber[ch] = 0;
684   }
685   prevch = 0; // previous chamber
686   nhits = 0; // number of hits in current chamber
687   // Loop over HitsForRec
688   for (Int_t hit = 0; hit < fNHitsForRec; hit++) {
689     // chamber number (0...)
690     ch = ((AliMUONHitForRec*)  ((*fHitsForRecPtr)[hit]))->GetChamberNumber();
691     // increment number of hits in current chamber
692     (fNHitsForRecPerChamber[ch])++;
693     // update index of first HitForRec in current chamber
694     // if chamber number different from previous one
695     if (ch != prevch) {
696       fIndexOfFirstHitForRecPerChamber[ch] = hit;
697       prevch = ch;
698     }
699   }
700   return;
701 }
702
703   //__________________________________________________________________________
704 void AliMUONTrackReconstructor::AddHitsForRecFromRawClusters(TTree* TR)
705 {
706   // To add to the list of hits for reconstruction all the raw clusters
707   // No condition added, like in Fortran TRACKF_STAT,
708   // on the radius between RMin and RMax.
709   AliMUONHitForRec *hitForRec;
710   AliMUONRawCluster *clus;
711   Int_t iclus, nclus, nTRentries;
712   TClonesArray *rawclusters;
713   AliDebug(1,"Enter AddHitsForRecFromRawClusters");
714
715   if (fTrackMethod != 3) { //AZ
716     fMUONData->SetTreeAddress("RC"); //AZ
717     nTRentries = Int_t(TR->GetEntries());
718     if (nTRentries != 1) {
719       AliError(Form("nTRentries = %d not equal to 1 ",nTRentries));
720       exit(0);
721     }
722     fLoader->TreeR()->GetEvent(0); // only one entry  
723   }
724
725   // Loop over tracking chambers
726   for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
727     // number of HitsForRec to 0 for the chamber
728     fNHitsForRecPerChamber[ch] = 0;
729     // index of first HitForRec for the chamber
730     if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0;
731     else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec;
732     rawclusters =fMUONData->RawClusters(ch);
733 //     pMUON->ResetRawClusters();
734 //     TR->GetEvent((Int_t) (TR->GetEntries()) - 1); // to be checked ????
735     nclus = (Int_t) (rawclusters->GetEntries());
736     // Loop over (cathode correlated) raw clusters
737     for (iclus = 0; iclus < nclus; iclus++) {
738       clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(iclus);
739       // new AliMUONHitForRec from raw cluster
740       // and increment number of AliMUONHitForRec's (total and in chamber)
741       hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(clus);
742       fNHitsForRec++;
743       (fNHitsForRecPerChamber[ch])++;
744       // more information into HitForRec
745       //  resolution: info should be already in raw cluster and taken from it ????
746       hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
747       hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
748       //hitForRec->SetBendingReso2(clus->GetErrY() * clus->GetErrY());
749       //hitForRec->SetNonBendingReso2(clus->GetErrX() * clus->GetErrX());
750       //  original raw cluster
751       hitForRec->SetChamberNumber(ch);
752       hitForRec->SetHitNumber(iclus);
753       // Z coordinate of the raw cluster (cm)
754       hitForRec->SetZ(clus->GetZ(0));
755       if (AliLog::GetGlobalDebugLevel() > 1) {
756         cout << "chamber (0...): " << ch <<
757           " raw cluster (0...): " << iclus << endl;
758         clus->Dump();
759         cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
760         hitForRec->Dump();}
761     } // end of cluster loop
762   } // end of chamber loop
763   SortHitsForRecWithIncreasingChamber(); 
764   return;
765 }
766
767   //__________________________________________________________________________
768 void AliMUONTrackReconstructor::MakeSegments(void)
769 {
770   // To make the list of segments in all stations,
771   // from the list of hits to be reconstructed
772   AliDebug(1,"Enter MakeSegments");
773   ResetSegments();
774   // Loop over stations
775   Int_t nb = (fTrackMethod != 1) ? 3 : 0; //AZ
776   for (Int_t st = nb; st < AliMUONConstants::NTrackingCh()/2; st++) 
777     MakeSegmentsPerStation(st); 
778   if (AliLog::GetGlobalDebugLevel() > 1) {
779     cout << "end of MakeSegments" << endl;
780     for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) {
781       cout << "station(0...): " << st
782            << "  Segments: " << fNSegments[st]
783            << endl;
784       for (Int_t seg = 0;
785            seg < fNSegments[st];
786            seg++) {
787         cout << "Segment index(0...): " << seg << endl;
788         ((*fSegmentsPtr[st])[seg])->Dump();
789       }
790     }
791   }
792   return;
793 }
794
795   //__________________________________________________________________________
796 void AliMUONTrackReconstructor::MakeSegmentsPerStation(Int_t Station)
797 {
798   // To make the list of segments in station number "Station" (0...)
799   // from the list of hits to be reconstructed.
800   // Updates "fNSegments"[Station].
801   // Segments in stations 4 and 5 are sorted
802   // according to increasing absolute value of "impact parameter"
803   AliMUONHitForRec *hit1Ptr, *hit2Ptr;
804   AliMUONSegment *segment;
805   Bool_t last2st;
806   Double_t bendingSlope, distBend, distNonBend, extBendCoor, extNonBendCoor,
807       impactParam = 0., maxImpactParam = 0., minImpactParam = 0.; // =0 to avoid compilation warnings.
808   AliDebug(1,Form("Enter MakeSegmentsPerStation (0...) %d",Station));
809   // first and second chambers (0...) in the station
810   Int_t ch1 = 2 * Station;
811   Int_t ch2 = ch1 + 1;
812   // variable true for stations downstream of the dipole:
813   // Station(0..4) equal to 3 or 4
814   if ((Station == 3) || (Station == 4)) {
815     last2st = kTRUE;
816     // maximum impact parameter (cm) according to fMinBendingMomentum...
817     maxImpactParam =
818       TMath::Abs(GetImpactParamFromBendingMomentum(fMinBendingMomentum));
819     // minimum impact parameter (cm) according to fMaxBendingMomentum...
820     minImpactParam =
821       TMath::Abs(GetImpactParamFromBendingMomentum(fMaxBendingMomentum));
822   }
823   else last2st = kFALSE;
824   // extrapolation factor from Z of first chamber to Z of second chamber
825   // dZ to be changed to take into account fine structure of chambers ????
826   Double_t extrapFact;
827   // index for current segment
828   Int_t segmentIndex = 0;
829   // Loop over HitsForRec in the first chamber of the station
830   for (Int_t hit1 = fIndexOfFirstHitForRecPerChamber[ch1];
831        hit1 < fIndexOfFirstHitForRecPerChamber[ch1] + fNHitsForRecPerChamber[ch1];
832        hit1++) {
833     // pointer to the HitForRec
834     hit1Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit1]);
835     // extrapolation,
836     // on the straight line joining the HitForRec to the vertex (0,0,0),
837     // to the Z of the second chamber of the station
838     // Loop over HitsForRec in the second chamber of the station
839     for (Int_t hit2 = fIndexOfFirstHitForRecPerChamber[ch2];
840          hit2 < fIndexOfFirstHitForRecPerChamber[ch2] + fNHitsForRecPerChamber[ch2];
841          hit2++) {
842       // pointer to the HitForRec
843       hit2Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit2]);
844       // absolute values of distances, in bending and non bending planes,
845       // between the HitForRec in the second chamber
846       // and the previous extrapolation
847       extrapFact = hit2Ptr->GetZ()/ hit1Ptr->GetZ();
848       extBendCoor = extrapFact * hit1Ptr->GetBendingCoor();
849       extNonBendCoor = extrapFact * hit1Ptr->GetNonBendingCoor();
850       distBend = TMath::Abs(hit2Ptr->GetBendingCoor() - extBendCoor);
851       distNonBend = TMath::Abs(hit2Ptr->GetNonBendingCoor() - extNonBendCoor);
852       if (last2st) {
853         // bending slope
854         bendingSlope = (hit1Ptr->GetBendingCoor() - hit2Ptr->GetBendingCoor()) /
855           (hit1Ptr->GetZ() - hit2Ptr->GetZ());
856         // absolute value of impact parameter
857         impactParam =
858           TMath::Abs(hit1Ptr->GetBendingCoor() - hit1Ptr->GetZ() * bendingSlope);
859       }
860       // check for distances not too large,
861       // and impact parameter not too big if stations downstream of the dipole.
862       // Conditions "distBend" and "impactParam" correlated for these stations ????
863       if ((distBend < fSegmentMaxDistBending[Station]) &&
864           (distNonBend < fSegmentMaxDistNonBending[Station]) &&
865           (!last2st || (impactParam < maxImpactParam)) &&
866           (!last2st || (impactParam > minImpactParam))) {
867         // make new segment
868         segment = new ((*fSegmentsPtr[Station])[segmentIndex])
869           AliMUONSegment(hit1Ptr, hit2Ptr);
870         // update "link" to this segment from the hit in the first chamber
871         if (hit1Ptr->GetNSegments() == 0)
872           hit1Ptr->SetIndexOfFirstSegment(segmentIndex);
873         hit1Ptr->SetNSegments(hit1Ptr->GetNSegments() + 1);
874         if (AliLog::GetGlobalDebugLevel() > 1) {
875           cout << "segmentIndex(0...): " << segmentIndex
876                << "  distBend: " << distBend
877                << "  distNonBend: " << distNonBend
878                << endl;
879           segment->Dump();
880           cout << "HitForRec in first chamber" << endl;
881           hit1Ptr->Dump();
882           cout << "HitForRec in second chamber" << endl;
883           hit2Ptr->Dump();
884         };
885         // increment index for current segment
886         segmentIndex++;
887       }
888     } //for (Int_t hit2
889   } // for (Int_t hit1...
890   fNSegments[Station] = segmentIndex;
891   // Sorting according to "impact parameter" if station(1..5) 4 or 5,
892   // i.e. Station(0..4) 3 or 4, using the function "Compare".
893   // After this sorting, it is impossible to use
894   // the "fNSegments" and "fIndexOfFirstSegment"
895   // of the HitForRec in the first chamber to explore all segments formed with it.
896   // Is this sorting really needed ????
897   if ((Station == 3) || (Station == 4)) (fSegmentsPtr[Station])->Sort();
898   AliDebug(1,Form("Station: %d  NSegments:  %d ", Station, fNSegments[Station]));
899   return;
900 }
901
902   //__________________________________________________________________________
903 void AliMUONTrackReconstructor::MakeTracks(void)
904 {
905   // To make the tracks,
906   // from the list of segments and points in all stations
907   AliDebug(1,"Enter MakeTracks");
908   // The order may be important for the following Reset's
909   ResetTracks();
910   ResetTrackHits();
911   if (fTrackMethod != 1) { //AZ - Kalman filter
912     MakeTrackCandidatesK();
913     if (fRecTracksPtr->GetEntriesFast() == 0) return;
914     // Follow tracks in stations(1..) 3, 2 and 1
915     FollowTracksK();
916     // Remove double tracks
917     RemoveDoubleTracksK();
918     // Propagate tracks to the vertex thru absorber
919     GoToVertex();
920     // Fill AliMUONTrack data members
921     FillMUONTrack();
922   } else { 
923     // Look for candidates from at least 3 aligned points in stations(1..) 4 and 5
924     MakeTrackCandidates();
925     // Follow tracks in stations(1..) 3, 2 and 1
926     FollowTracks();
927     // Remove double tracks
928     RemoveDoubleTracks();
929     UpdateTrackParamAtHit();
930     UpdateHitForRecAtHit();
931   }
932   return;
933 }
934
935   //__________________________________________________________________________
936 void AliMUONTrackReconstructor::ValidateTracksWithTrigger(void)
937 {
938   // Try to match track from tracking system with trigger track
939   AliMUONTrack *track;
940   TClonesArray *recTriggerTracks;
941   
942   fMUONData->SetTreeAddress("RL");
943   fMUONData->GetRecTriggerTracks();
944   recTriggerTracks = fMUONData->RecTriggerTracks();
945
946   track = (AliMUONTrack*) fRecTracksPtr->First();
947   while (track) {
948     track->MatchTriggerTrack(recTriggerTracks);
949     track = (AliMUONTrack*) fRecTracksPtr->After(track);
950   }
951
952 }
953
954   //__________________________________________________________________________
955 Bool_t AliMUONTrackReconstructor::MakeTriggerTracks(void)
956 {
957     // To make the trigger tracks from Local Trigger
958   AliDebug(1, "Enter MakeTriggerTracks");
959     
960     Int_t nTRentries;
961     Long_t gloTrigPat;
962     TClonesArray *localTrigger;
963     TClonesArray *globalTrigger;
964     AliMUONLocalTrigger *locTrg;
965     AliMUONGlobalTrigger *gloTrg;
966     AliMUONTriggerCircuit *circuit;
967     AliMUONTriggerTrack *recTriggerTrack = 0;
968
969     TTree* treeR = fLoader->TreeR();
970     
971     // Loading MUON subsystem
972     AliMUON * pMUON = (AliMUON *) gAlice->GetDetector("MUON");
973     
974     nTRentries = Int_t(treeR->GetEntries());
975      
976     treeR->GetEvent(0); // only one entry  
977
978     if (!(fMUONData->IsTriggerBranchesInTree())) {
979       AliWarning(Form("Trigger information is not avalaible, nTRentries = %d not equal to 1",nTRentries));
980       return kFALSE;
981     }
982
983     fMUONData->SetTreeAddress("TC");
984     fMUONData->GetTrigger();
985
986     // global trigger for trigger pattern
987     gloTrigPat = 0;
988     globalTrigger = fMUONData->GlobalTrigger(); 
989     gloTrg = (AliMUONGlobalTrigger*)globalTrigger->UncheckedAt(0);      
990     if (gloTrg) {
991       if (gloTrg->SinglePlusLpt())  gloTrigPat|= 0x1;
992       if (gloTrg->SinglePlusHpt())  gloTrigPat|= 0x2;
993       if (gloTrg->SinglePlusApt())  gloTrigPat|= 0x4;
994       
995       if (gloTrg->SingleMinusLpt()) gloTrigPat|= 0x8;
996       if (gloTrg->SingleMinusHpt()) gloTrigPat|= 0x10;
997       if (gloTrg->SingleMinusApt()) gloTrigPat|= 0x20;
998       
999       if (gloTrg->SingleUndefLpt()) gloTrigPat|= 0x40;
1000       if (gloTrg->SingleUndefHpt()) gloTrigPat|= 0x80;
1001       if (gloTrg->SingleUndefApt()) gloTrigPat|= 0x100;
1002       
1003       if (gloTrg->PairUnlikeLpt())  gloTrigPat|= 0x200;
1004       if (gloTrg->PairUnlikeHpt())  gloTrigPat|= 0x400;
1005       if (gloTrg->PairUnlikeApt())  gloTrigPat|= 0x800;
1006       
1007       if (gloTrg->PairLikeLpt())    gloTrigPat|= 0x1000;
1008       if (gloTrg->PairLikeHpt())    gloTrigPat|= 0x2000;
1009       if (gloTrg->PairLikeApt())    gloTrigPat|= 0x4000;
1010     }
1011
1012     // local trigger for tracking 
1013     localTrigger = fMUONData->LocalTrigger();    
1014     Int_t nlocals = (Int_t) (localTrigger->GetEntries());
1015
1016     Float_t z11 = AliMUONConstants::DefaultChamberZ(10);
1017     Float_t z21 = AliMUONConstants::DefaultChamberZ(12);
1018
1019     for (Int_t i=0; i<nlocals; i++) { // loop on Local Trigger
1020       locTrg = (AliMUONLocalTrigger*)localTrigger->UncheckedAt(i);      
1021       circuit = &(pMUON->TriggerCircuit(locTrg->LoCircuit()));
1022       Float_t y11 = circuit->GetY11Pos(locTrg->LoStripX()); 
1023       Int_t stripX21 = locTrg->LoStripX()+locTrg->LoDev()+1;
1024       Float_t y21 = circuit->GetY21Pos(stripX21);       
1025       Float_t x11 = circuit->GetX11Pos(locTrg->LoStripY());
1026       Float_t thetax = TMath::ATan2( x11 , z11 );
1027       Float_t thetay = TMath::ATan2( (y21-y11) , (z21-z11) );
1028
1029       recTriggerTrack = new AliMUONTriggerTrack(x11,y11,thetax,thetay,gloTrigPat);
1030       
1031       // since static statement does not work, set gloTrigPat for each track
1032
1033       fMUONData->AddRecTriggerTrack(*recTriggerTrack);
1034       delete recTriggerTrack;
1035     } // end of loop on Local Trigger
1036     return kTRUE;    
1037 }
1038
1039   //__________________________________________________________________________
1040 Int_t AliMUONTrackReconstructor::MakeTrackCandidatesWithTwoSegments(AliMUONSegment *BegSegment)
1041 {
1042   // To make track candidates with two segments in stations(1..) 4 and 5,
1043   // the first segment being pointed to by "BegSegment".
1044   // Returns the number of such track candidates.
1045   Int_t endStation, iEndSegment, nbCan2Seg;
1046   AliMUONSegment *endSegment;
1047   AliMUONSegment *extrapSegment = NULL;
1048   AliMUONTrack *recTrack;
1049   Double_t mcsFactor;
1050   AliDebug(1,"Enter MakeTrackCandidatesWithTwoSegments");
1051   // Station for the end segment
1052   endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
1053   // multiple scattering factor corresponding to one chamber
1054   mcsFactor = 0.0136 /
1055     GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
1056   mcsFactor     = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1057   // linear extrapolation to end station
1058   // number of candidates with 2 segments to 0
1059   nbCan2Seg = 0;
1060   // Loop over segments in the end station
1061   for (iEndSegment = 0; iEndSegment < fNSegments[endStation]; iEndSegment++) {
1062     // pointer to segment
1063     endSegment = (AliMUONSegment*) ((*fSegmentsPtr[endStation])[iEndSegment]);
1064     // test compatibility between current segment and "extrapSegment"
1065     // 4 because 4 quantities in chi2
1066     extrapSegment =
1067       BegSegment->CreateSegmentFromLinearExtrapToStation(endSegment->GetZ(), mcsFactor);
1068     if ((endSegment->
1069          NormalizedChi2WithSegment(extrapSegment,
1070                                    fMaxSigma2Distance)) <= 4.0) {
1071       // both segments compatible:
1072       // make track candidate from "begSegment" and "endSegment"
1073       AliDebug(2,Form("TrackCandidate with Segment %d in Station(0..) %d", iEndSegment, endStation));
1074       // flag for both segments in one track:
1075       // to be done in track constructor ????
1076       BegSegment->SetInTrack(kTRUE);
1077       endSegment->SetInTrack(kTRUE);
1078       recTrack = new ((*fRecTracksPtr)[fNRecTracks])
1079         AliMUONTrack(BegSegment, endSegment, this);
1080       fNRecTracks++;
1081       if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump();
1082       // increment number of track candidates with 2 segments
1083       nbCan2Seg++;
1084     }
1085     delete extrapSegment; // should not delete HitForRec's it points to !!!!
1086   } // for (iEndSegment = 0;...
1087   return nbCan2Seg;
1088 }
1089
1090   //__________________________________________________________________________
1091 Int_t AliMUONTrackReconstructor::MakeTrackCandidatesWithOneSegmentAndOnePoint(AliMUONSegment *BegSegment)
1092 {
1093   // To make track candidates with one segment and one point
1094   // in stations(1..) 4 and 5,
1095   // the segment being pointed to by "BegSegment".
1096   Int_t ch, ch1, ch2, endStation, iHit, iHitMax, iHitMin, nbCan1Seg1Hit;
1097   AliMUONHitForRec *extrapHitForRec= NULL;
1098   AliMUONHitForRec *hit;
1099   AliMUONTrack *recTrack;
1100   Double_t mcsFactor;
1101   AliDebug(1,"Enter MakeTrackCandidatesWithOneSegmentAndOnePoint");
1102   // station for the end point
1103   endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
1104   // multiple scattering factor corresponding to one chamber
1105   mcsFactor = 0.0136 /
1106     GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
1107   mcsFactor     = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1108   // first and second chambers(0..) in the end station
1109   ch1 = 2 * endStation;
1110   ch2 = ch1 + 1;
1111   // number of candidates to 0
1112   nbCan1Seg1Hit = 0;
1113   // Loop over chambers of the end station
1114   for (ch = ch2; ch >= ch1; ch--) {
1115     // limits for the hit index in the loop
1116     iHitMin = fIndexOfFirstHitForRecPerChamber[ch];
1117     iHitMax = iHitMin + fNHitsForRecPerChamber[ch];
1118     // Loop over HitForRec's in the chamber
1119     for (iHit = iHitMin; iHit < iHitMax; iHit++) {
1120       // pointer to HitForRec
1121       hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1122       // test compatibility between current HitForRec and "extrapHitForRec"
1123       // 2 because 2 quantities in chi2
1124       // linear extrapolation to chamber
1125       extrapHitForRec =
1126         BegSegment->CreateHitForRecFromLinearExtrapToChamber( hit->GetZ(), mcsFactor);
1127       if ((hit->
1128            NormalizedChi2WithHitForRec(extrapHitForRec,
1129                                        fMaxSigma2Distance)) <= 2.0) {
1130         // both HitForRec's compatible:
1131         // make track candidate from begSegment and current HitForRec
1132         AliDebug(2, Form("TrackCandidate with HitForRec  %d in Chamber(0..) %d", iHit, ch));
1133         // flag for beginning segments in one track:
1134         // to be done in track constructor ????
1135         BegSegment->SetInTrack(kTRUE);
1136         recTrack = new ((*fRecTracksPtr)[fNRecTracks])
1137           AliMUONTrack(BegSegment, hit, this);
1138         // the right place to eliminate "double counting" ???? how ????
1139         fNRecTracks++;
1140         if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump();
1141         // increment number of track candidates
1142         nbCan1Seg1Hit++;
1143       }
1144       delete extrapHitForRec;
1145     } // for (iHit = iHitMin;...
1146   } // for (ch = ch2;...
1147   return nbCan1Seg1Hit;
1148 }
1149
1150   //__________________________________________________________________________
1151 void AliMUONTrackReconstructor::MakeTrackCandidates(void)
1152 {
1153   // To make track candidates
1154   // with at least 3 aligned points in stations(1..) 4 and 5
1155   // (two Segment's or one Segment and one HitForRec)
1156   Int_t begStation, iBegSegment, nbCan1Seg1Hit, nbCan2Seg;
1157   AliMUONSegment *begSegment;
1158   AliDebug(1,"Enter MakeTrackCandidates");
1159   // Loop over stations(1..) 5 and 4 for the beginning segment
1160   for (begStation = 4; begStation > 2; begStation--) {
1161     // Loop over segments in the beginning station
1162     for (iBegSegment = 0; iBegSegment < fNSegments[begStation]; iBegSegment++) {
1163       // pointer to segment
1164       begSegment = (AliMUONSegment*) ((*fSegmentsPtr[begStation])[iBegSegment]);
1165       AliDebug(2,Form("Look for TrackCandidate's with Segment %d  in Station(0..) %d", iBegSegment, begStation));
1166       // Look for track candidates with two segments,
1167       // "begSegment" and all compatible segments in other station.
1168       // Only for beginning station(1..) 5
1169       // because candidates with 2 segments have to looked for only once.
1170       if (begStation == 4)
1171         nbCan2Seg = MakeTrackCandidatesWithTwoSegments(begSegment);
1172       // Look for track candidates with one segment and one point,
1173       // "begSegment" and all compatible HitForRec's in other station.
1174       // Only if "begSegment" does not belong already to a track candidate.
1175       // Is that a too strong condition ????
1176       if (!(begSegment->GetInTrack()))
1177         nbCan1Seg1Hit = MakeTrackCandidatesWithOneSegmentAndOnePoint(begSegment);
1178     } // for (iBegSegment = 0;...
1179   } // for (begStation = 4;...
1180   return;
1181 }
1182
1183   //__________________________________________________________________________
1184 void AliMUONTrackReconstructor::FollowTracks(void)
1185 {
1186   // Follow tracks in stations(1..) 3, 2 and 1
1187   // too long: should be made more modular !!!!
1188   AliMUONHitForRec *bestHit, *extrapHit, *hit;
1189   AliMUONSegment *bestSegment, *extrapSegment, *segment;
1190   AliMUONTrack *track, *nextTrack;
1191   AliMUONTrackParam *trackParam1, trackParam[2], trackParamVertex;
1192   // -1 to avoid compilation warnings
1193   Int_t ch = -1, chInStation, chBestHit = -1, iHit, iSegment, station, trackIndex; 
1194   Double_t bestChi2, chi2, dZ1, dZ2, dZ3, maxSigma2Distance, mcsFactor;
1195   Double_t bendingMomentum, chi2Norm = 0.;
1196
1197   // local maxSigma2Distance, for easy increase in testing
1198   maxSigma2Distance = fMaxSigma2Distance;
1199   AliDebug(2,"Enter FollowTracks");
1200   // Loop over track candidates
1201   track = (AliMUONTrack*) fRecTracksPtr->First();
1202   trackIndex = -1;
1203   while (track) {
1204     // Follow function for each track candidate ????
1205     trackIndex++;
1206     nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track); // prepare next track
1207     AliDebug(2,Form("FollowTracks: track candidate(0..): %d", trackIndex));
1208     // Fit track candidate
1209     track->SetFitMCS(0); // without multiple Coulomb scattering
1210     track->SetFitNParam(3); // with 3 parameters (X = Y = 0)
1211     track->SetFitStart(0); // from parameters at vertex
1212     track->Fit();
1213     if (AliLog::GetGlobalDebugLevel()> 2) {
1214       cout << "FollowTracks: track candidate(0..): " << trackIndex
1215            << " after fit in stations(0..) 3 and 4" << endl;
1216       track->RecursiveDump();
1217     }
1218     // Loop over stations(1..) 3, 2 and 1
1219     // something SPECIAL for stations 2 and 1 for majority 3 coincidence ????
1220     // otherwise: majority coincidence 2 !!!!
1221     for (station = 2; station >= 0; station--) {
1222       // Track parameters at first track hit (smallest Z)
1223       trackParam1 = ((AliMUONTrackHit*)
1224                      (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1225       // extrapolation to station
1226       trackParam1->ExtrapToStation(station, trackParam);
1227       extrapSegment = new AliMUONSegment(); //  empty segment
1228       // multiple scattering factor corresponding to one chamber
1229       // and momentum in bending plane (not total)
1230       mcsFactor = 0.0136 * trackParam1->GetInverseBendingMomentum();
1231       mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1232       // Z difference from previous station
1233       dZ1 = AliMUONConstants::DefaultChamberZ(2 * station) -
1234             AliMUONConstants::DefaultChamberZ(2 * station + 2);
1235       // Z difference between the two previous stations
1236       dZ2 = AliMUONConstants::DefaultChamberZ(2 * station + 2) -
1237             AliMUONConstants::DefaultChamberZ(2 * station + 4);
1238       // Z difference between the two chambers in the previous station
1239       dZ3 = AliMUONConstants::DefaultChamberZ(2 * station) -
1240             AliMUONConstants::DefaultChamberZ(2 * station + 1);
1241       extrapSegment->SetBendingCoorReso2(fBendingResolution * fBendingResolution);
1242       extrapSegment->
1243         SetNonBendingCoorReso2(fNonBendingResolution * fNonBendingResolution);
1244       extrapSegment->UpdateFromStationTrackParam
1245         (trackParam, mcsFactor, dZ1, dZ2, dZ3, station,
1246          trackParam1->GetInverseBendingMomentum());
1247       bestChi2 = 5.0;
1248       bestSegment = NULL;
1249       if (AliLog::GetGlobalDebugLevel() > 2) {
1250         cout << "FollowTracks: track candidate(0..): " << trackIndex
1251              << " Look for segment in station(0..): " << station << endl;
1252       }
1253       // Loop over segments in station
1254       for (iSegment = 0; iSegment < fNSegments[station]; iSegment++) {
1255         // Look for best compatible Segment in station
1256         // should consider all possibilities ????
1257         // multiple scattering ????
1258         // separation in 2 functions: Segment and HitForRec ????
1259         segment = (AliMUONSegment*) ((*fSegmentsPtr[station])[iSegment]);
1260         // correction of corrected segment (fBendingCoor and fNonBendingCoor)
1261         // according to real Z value of "segment" and slopes of "extrapSegment"
1262         (&(trackParam[0]))->ExtrapToZ(segment->GetZ());
1263         (&(trackParam[1]))->ExtrapToZ(segment->GetZ());
1264         extrapSegment->SetBendingCoor((&(trackParam[0]))->GetBendingCoor());
1265         extrapSegment->SetNonBendingCoor((&(trackParam[0]))->GetNonBendingCoor());
1266         extrapSegment->SetBendingSlope((&(trackParam[0]))->GetBendingSlope());
1267         extrapSegment->SetNonBendingSlope((&(trackParam[0]))->GetNonBendingSlope());
1268         chi2 = segment->
1269           NormalizedChi2WithSegment(extrapSegment, maxSigma2Distance);
1270         if (chi2 < bestChi2) {
1271           // update best Chi2 and Segment if better found
1272           bestSegment = segment;
1273           bestChi2 = chi2;
1274         }
1275       }
1276       if (bestSegment) {
1277         // best segment found: add it to track candidate
1278         (&(trackParam[0]))->ExtrapToZ(bestSegment->GetZ());
1279         (&(trackParam[1]))->ExtrapToZ(bestSegment->GetZ());
1280         track->AddSegment(bestSegment);
1281         // set track parameters at these two TrakHit's
1282         track->SetTrackParamAtHit(track->GetNTrackHits() - 2, &(trackParam[0]));
1283         track->SetTrackParamAtHit(track->GetNTrackHits() - 1, &(trackParam[1]));
1284         AliDebug(3, Form("FollowTracks: track candidate(0..): %d  Added segment in station(0..): %d", trackIndex, station));
1285         if (AliLog::GetGlobalDebugLevel()>2) track->RecursiveDump();
1286       }
1287       else {
1288         // No best segment found:
1289         // Look for best compatible HitForRec in station:
1290         // should consider all possibilities ????
1291         // multiple scattering ???? do about like for extrapSegment !!!!
1292         extrapHit = new AliMUONHitForRec(); //  empty hit
1293         bestChi2 = 3.0;
1294         bestHit = NULL;
1295         AliDebug(3, Form("FollowTracks: track candidate(0..): %d Segment not found, look for hit in station(0..): %d ", 
1296                          trackIndex, station));
1297         
1298         // Loop over chambers of the station
1299         for (chInStation = 0; chInStation < 2; chInStation++) {
1300           ch = 2 * station + chInStation;
1301           for (iHit = fIndexOfFirstHitForRecPerChamber[ch];
1302                iHit < fIndexOfFirstHitForRecPerChamber[ch] +
1303                  fNHitsForRecPerChamber[ch];
1304                iHit++) {
1305             hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1306             // coordinates of extrapolated hit
1307             (&(trackParam[chInStation]))->ExtrapToZ(hit->GetZ());
1308             extrapHit->
1309               SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor());
1310             extrapHit->
1311               SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor());
1312             // resolutions from "extrapSegment"
1313             extrapHit->SetBendingReso2(extrapSegment->GetBendingCoorReso2());
1314             extrapHit->SetNonBendingReso2(extrapSegment->GetNonBendingCoorReso2());
1315             // Loop over hits in the chamber
1316             // condition for hit not already in segment ????
1317             chi2 = hit->NormalizedChi2WithHitForRec(extrapHit, maxSigma2Distance);
1318             if (chi2 < bestChi2) {
1319               // update best Chi2 and HitForRec if better found
1320               bestHit = hit;
1321               bestChi2 = chi2;
1322               chBestHit = chInStation;
1323             }
1324           }
1325         }
1326         if (bestHit) {
1327           // best hit found: add it to track candidate
1328           (&(trackParam[chBestHit]))->ExtrapToZ(bestHit->GetZ());
1329           track->AddHitForRec(bestHit);
1330           // set track parameters at this TrackHit
1331           track->SetTrackParamAtHit(track->GetNTrackHits() - 1,
1332                                     &(trackParam[chBestHit]));
1333           if (AliLog::GetGlobalDebugLevel() > 2) {
1334             cout << "FollowTracks: track candidate(0..): " << trackIndex
1335                  << " Added hit in station(0..): " << station << endl;
1336             track->RecursiveDump();
1337           }
1338         }
1339         else {
1340           // Remove current track candidate
1341           // and corresponding TrackHit's, ...
1342           track->Remove();
1343           delete extrapSegment;
1344           delete extrapHit;
1345           break; // stop the search for this candidate:
1346           // exit from the loop over station
1347         }
1348         delete extrapHit;
1349       }
1350       delete extrapSegment;
1351       // Sort track hits according to increasing Z
1352       track->GetTrackHitsPtr()->Sort();
1353       // Update track parameters at first track hit (smallest Z)
1354       trackParam1 = ((AliMUONTrackHit*)
1355                      (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1356       bendingMomentum = 0.;
1357       if (TMath::Abs(trackParam1->GetInverseBendingMomentum()) > 0.)
1358         bendingMomentum = TMath::Abs(1/(trackParam1->GetInverseBendingMomentum()));
1359       // Track removed if bendingMomentum not in window [min, max]
1360       if ((bendingMomentum < fMinBendingMomentum) || (bendingMomentum > fMaxBendingMomentum)) {
1361         track->Remove();
1362         break; // stop the search for this candidate:
1363         // exit from the loop over station 
1364       }
1365       // Track fit
1366       // with multiple Coulomb scattering if all stations
1367       if (station == 0) track->SetFitMCS(1);
1368       // without multiple Coulomb scattering if not all stations
1369       else track->SetFitMCS(0);
1370       track->SetFitNParam(5);  // with 5 parameters (momentum and position)
1371       track->SetFitStart(1);  // from parameters at first hit
1372       track->Fit();
1373       Double_t numberOfDegFree = (2.0 * track->GetNTrackHits() - 5);
1374       if (numberOfDegFree > 0) {
1375         chi2Norm =  track->GetFitFMin() / numberOfDegFree;
1376       } else {
1377         chi2Norm = 1.e10;
1378       }
1379       // Track removed if normalized chi2 too high
1380       if (chi2Norm > fMaxChi2) {
1381         track->Remove();
1382         break; // stop the search for this candidate:
1383         // exit from the loop over station 
1384       }
1385       if (AliLog::GetGlobalDebugLevel() > 2) {
1386         cout << "FollowTracks: track candidate(0..): " << trackIndex
1387              << " after fit from station(0..): " << station << " to 4" << endl;
1388         track->RecursiveDump();
1389       }
1390       // Track extrapolation to the vertex through the absorber (Branson)
1391       // after going through the first station
1392       if (station == 0) {
1393         trackParamVertex = *trackParam1;
1394         (&trackParamVertex)->ExtrapToVertex(0.,0.,0.);
1395         track->SetTrackParamAtVertex(&trackParamVertex);
1396         if (AliLog::GetGlobalDebugLevel() > 0) {
1397           cout << "FollowTracks: track candidate(0..): " << trackIndex
1398                << " after extrapolation to vertex" << endl;
1399           track->RecursiveDump();
1400         }
1401       }
1402     } // for (station = 2;...
1403     // go really to next track
1404     track = nextTrack;
1405   } // while (track)
1406   // Compression of track array (necessary after Remove ????)
1407   fRecTracksPtr->Compress();
1408   return;
1409 }
1410
1411   //__________________________________________________________________________
1412 void AliMUONTrackReconstructor::RemoveDoubleTracks(void)
1413 {
1414   // To remove double tracks.
1415   // Tracks are considered identical
1416   // if they have at least half of their hits in common.
1417   // Among two identical tracks, one keeps the track with the larger number of hits
1418   // or, if these numbers are equal, the track with the minimum Chi2.
1419   AliMUONTrack *track1, *track2, *trackToRemove;
1420   Bool_t identicalTracks;
1421   Int_t hitsInCommon, nHits1, nHits2;
1422   identicalTracks = kTRUE;
1423   while (identicalTracks) {
1424     identicalTracks = kFALSE;
1425     // Loop over first track of the pair
1426     track1 = (AliMUONTrack*) fRecTracksPtr->First();
1427     while (track1 && (!identicalTracks)) {
1428       nHits1 = track1->GetNTrackHits();
1429       // Loop over second track of the pair
1430       track2 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1431       while (track2 && (!identicalTracks)) {
1432         nHits2 = track2->GetNTrackHits();
1433         // number of hits in common between two tracks
1434         hitsInCommon = track1->HitsInCommon(track2);
1435         // check for identical tracks
1436         if ((4 * hitsInCommon) >= (nHits1 + nHits2)) {
1437           identicalTracks = kTRUE;
1438           // decide which track to remove
1439           if (nHits1 > nHits2) trackToRemove = track2;
1440           else if (nHits1 < nHits2) trackToRemove = track1;
1441           else if ((track1->GetFitFMin()) < (track2->GetFitFMin()))
1442             trackToRemove = track2;
1443           else trackToRemove = track1;
1444           // remove it
1445           trackToRemove->Remove();
1446         }
1447         track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
1448       } // track2
1449       track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1450     } // track1
1451   }
1452   return;
1453 }
1454
1455   //__________________________________________________________________________
1456 void AliMUONTrackReconstructor::UpdateTrackParamAtHit()
1457 {
1458   // Set track parameters after track fitting. Fill fTrackParamAtHit of AliMUONTrack's
1459   AliMUONTrack *track;
1460   AliMUONTrackHit *trackHit;
1461   AliMUONTrackParam *trackParam;
1462   track = (AliMUONTrack*) fRecTracksPtr->First();
1463   while (track) {
1464     trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->First();
1465     while (trackHit) {
1466       trackParam = trackHit->GetTrackParam();
1467       track->AddTrackParamAtHit(trackParam);
1468       trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->After(trackHit); 
1469     } // trackHit    
1470     track = (AliMUONTrack*) fRecTracksPtr->After(track);
1471   } // track
1472   return;
1473 }
1474
1475   //__________________________________________________________________________
1476 void AliMUONTrackReconstructor::UpdateHitForRecAtHit()
1477 {
1478   // Set cluster parameterss after track fitting. Fill fHitForRecAtHit of AliMUONTrack's
1479   AliMUONTrack *track;
1480   AliMUONTrackHit *trackHit;
1481   AliMUONHitForRec *hitForRec;
1482   track = (AliMUONTrack*) fRecTracksPtr->First();
1483   while (track) {
1484     trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->First();
1485     while (trackHit) {
1486       hitForRec = trackHit->GetHitForRecPtr();
1487       track->AddHitForRecAtHit(hitForRec);
1488       trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->After(trackHit); 
1489     } // trackHit    
1490     track = (AliMUONTrack*) fRecTracksPtr->After(track);
1491   } // track
1492   return;
1493 }
1494
1495   //__________________________________________________________________________
1496 void AliMUONTrackReconstructor::FillMUONTrack()
1497 {
1498   // Set track parameters at hits for Kalman track. Fill fTrackParamAtHit of AliMUONTrack's
1499   AliMUONTrackK *track;
1500   track = (AliMUONTrackK*) fRecTracksPtr->First();
1501   while (track) {
1502     track->FillMUONTrack();
1503     track = (AliMUONTrackK*) fRecTracksPtr->After(track);
1504   } 
1505   return;
1506 }
1507
1508   //__________________________________________________________________________
1509 void AliMUONTrackReconstructor::EventDump(void)
1510 {
1511   // Dump reconstructed event (track parameters at vertex and at first hit),
1512   // and the particle parameters
1513
1514   AliMUONTrack *track;
1515   AliMUONTrackParam *trackParam, *trackParam1;
1516   Double_t bendingSlope, nonBendingSlope, pYZ;
1517   Double_t pX, pY, pZ, x, y, z, c;
1518   Int_t np, trackIndex, nTrackHits;
1519  
1520   AliDebug(1,"****** enter EventDump ******");
1521   AliDebug(1, Form("Number of Reconstructed tracks : %d", fNRecTracks)); 
1522   
1523   fRecTracksPtr->Compress(); // for simple loop without "Next" since no hole
1524   // Loop over reconstructed tracks
1525   for (trackIndex = 0; trackIndex < fNRecTracks; trackIndex++) {
1526     if (fTrackMethod != 1) continue; //AZ - skip the rest for now
1527     AliDebug(1, Form("track number: %d", trackIndex));
1528     // function for each track for modularity ????
1529     track = (AliMUONTrack*) ((*fRecTracksPtr)[trackIndex]);
1530     nTrackHits = track->GetNTrackHits();
1531     AliDebug(1, Form("Number of track hits: %d ", nTrackHits));
1532     // track parameters at Vertex
1533     trackParam = track->GetTrackParamAtVertex();
1534     x = trackParam->GetNonBendingCoor();
1535     y = trackParam->GetBendingCoor();
1536     z = trackParam->GetZ();
1537     bendingSlope = trackParam->GetBendingSlope();
1538     nonBendingSlope = trackParam->GetNonBendingSlope();
1539     pYZ = 1/TMath::Abs(trackParam->GetInverseBendingMomentum());
1540     pZ = pYZ/TMath::Sqrt(1+bendingSlope*bendingSlope);
1541     pX = pZ * nonBendingSlope;
1542     pY = pZ * bendingSlope;
1543     c = TMath::Sign(1.0, trackParam->GetInverseBendingMomentum());
1544     AliDebug(1, Form("Track parameters at Vertex z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1545                      z, x, y, pX, pY, pZ, c));
1546
1547     // track parameters at first hit
1548     trackParam1 = ((AliMUONTrackHit*)
1549                    (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1550     x = trackParam1->GetNonBendingCoor();
1551     y = trackParam1->GetBendingCoor();
1552     z = trackParam1->GetZ();
1553     bendingSlope = trackParam1->GetBendingSlope();
1554     nonBendingSlope = trackParam1->GetNonBendingSlope();
1555     pYZ = 1/TMath::Abs(trackParam1->GetInverseBendingMomentum());
1556     pZ = pYZ/TMath::Sqrt(1.0 + bendingSlope * bendingSlope);
1557     pX = pZ * nonBendingSlope;
1558     pY = pZ * bendingSlope;
1559     c = TMath::Sign(1.0, trackParam1->GetInverseBendingMomentum());
1560     AliDebug(1, Form("track parameters at z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1561                      z, x, y, pX, pY, pZ, c));
1562   }
1563   // informations about generated particles
1564   np = gAlice->GetMCApp()->GetNtrack();
1565   printf(" **** number of generated particles: %d  \n", np);
1566   
1567 //    for (Int_t iPart = 0; iPart < np; iPart++) {
1568 //      p = gAlice->Particle(iPart);
1569 //      printf(" particle %d: type= %d px= %f py= %f pz= %f pdg= %d\n",
1570 //         iPart, p->GetPdgCode(), p->Px(), p->Py(), p->Pz(), p->GetPdgCode());    
1571 //    }
1572   return;
1573 }
1574
1575
1576 //__________________________________________________________________________
1577 void AliMUONTrackReconstructor::EventDumpTrigger(void)
1578 {
1579   // Dump reconstructed trigger event 
1580   // and the particle parameters
1581     
1582   AliMUONTriggerTrack *triggertrack ;
1583   Int_t nTriggerTracks = fMUONData->RecTriggerTracks()->GetEntriesFast();
1584  
1585   AliDebug(1, "****** enter EventDumpTrigger ******");
1586   AliDebug(1, Form("Number of Reconstructed tracks : %d ",  nTriggerTracks));
1587   
1588   // Loop over reconstructed tracks
1589   for (Int_t trackIndex = 0; trackIndex < nTriggerTracks; trackIndex++) {
1590     triggertrack = (AliMUONTriggerTrack*)fMUONData->RecTriggerTracks()->At(trackIndex);
1591       printf(" trigger track number %i x11=%f y11=%f thetax=%f thetay=%f \n",
1592              trackIndex,
1593              triggertrack->GetX11(),triggertrack->GetY11(),
1594              triggertrack->GetThetax(),triggertrack->GetThetay());      
1595   } 
1596 }
1597
1598 //__________________________________________________________________________
1599 void AliMUONTrackReconstructor::MakeTrackCandidatesK(void)
1600 {
1601   // To make initial tracks for Kalman filter from the list of segments
1602   Int_t istat, iseg;
1603   AliMUONSegment *segment;
1604   AliMUONTrackK *trackK;
1605
1606   AliDebug(1,"Enter MakeTrackCandidatesK");
1607   // Reset the TClonesArray of reconstructed tracks
1608   if (fRecTracksPtr) fRecTracksPtr->Delete();
1609   // Delete in order that the Track destructors are called,
1610   // hence the space for the TClonesArray of pointers to TrackHit's is freed
1611   fNRecTracks = 0;
1612
1613   AliMUONTrackK a(this, fHitsForRecPtr); // bad idea ???
1614   // Loop over stations(1...) 5 and 4
1615   for (istat=4; istat>=3; istat--) {
1616     // Loop over segments in the station
1617     for (iseg=0; iseg<fNSegments[istat]; iseg++) {
1618       // Transform segments to tracks and evaluate covariance matrix
1619       segment = (AliMUONSegment*) ((*fSegmentsPtr[istat])[iseg]);
1620       trackK = new ((*fRecTracksPtr)[fNRecTracks++]) AliMUONTrackK(segment);
1621     } // for (iseg=0;...)
1622   } // for (istat=4;...)
1623   return;
1624 }
1625
1626 //__________________________________________________________________________
1627 void AliMUONTrackReconstructor::FollowTracksK(void)
1628 {
1629   // Follow tracks using Kalman filter
1630   Bool_t ok;
1631   Int_t icand, ichamBeg = 0, ichamEnd, chamBits;
1632   Double_t zDipole1, zDipole2;
1633   AliMUONTrackK *trackK;
1634   AliMUONHitForRec *hit;
1635   AliMUONRawCluster *clus;
1636   TClonesArray *rawclusters;
1637   clus = 0; rawclusters = 0;
1638
1639   zDipole1 = GetSimpleBPosition() + GetSimpleBLength()/2;
1640   zDipole2 = zDipole1 - GetSimpleBLength();
1641
1642   // Print hits
1643   trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[0]);
1644   if (trackK->DebugLevel() > 0) {
1645     for (Int_t i1=0; i1<fNHitsForRec; i1++) {
1646       hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[i1]);
1647       //if (hit->GetTTRTrack() > 1 || hit->GetTrackRefSignal() == 0) continue;
1648       printf(" Hit # %d %10.4f %10.4f %10.4f",
1649              hit->GetChamberNumber(), hit->GetBendingCoor(),
1650              hit->GetNonBendingCoor(), hit->GetZ());
1651       if (fRecTrackRefHits) {
1652         // from track ref hits
1653         printf(" %3d %3d \n", hit->GetTrackRefSignal(), hit->GetTTRTrack());
1654       } else {
1655         // from raw clusters
1656         rawclusters = fMUONData->RawClusters(hit->GetChamberNumber());
1657         clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(hit->
1658                                                  GetHitNumber());
1659         printf(" %d", clus->GetTrack(1));
1660         if (clus->GetTrack(2) != -1) printf(" %d \n", clus->GetTrack(2));
1661         else printf("\n");
1662       } // if (fRecTrackRefHits)
1663     }
1664   } // if (trackK->DebugLevel() > 0)
1665
1666   icand = -1;
1667   Int_t nSeeds;
1668   nSeeds = fNRecTracks; // starting number of seeds
1669   // Loop over track candidates
1670   while (icand < fNRecTracks-1) {
1671     icand ++;
1672     if (trackK->DebugLevel()>0) cout << " *** Kalman track candidate No. " << icand << endl;
1673     trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]);
1674     if (trackK->GetRecover() < 0) continue; // failed track
1675
1676     // Discard candidate which will produce the double track
1677     /*
1678     if (icand > 0) {
1679       ok = CheckCandidateK(icand,nSeeds);
1680       if (!ok) {
1681         trackK->SetRecover(-1); // mark candidate to be removed
1682         continue;
1683       }
1684     }
1685     */
1686
1687     ok = kTRUE;
1688     if (trackK->GetRecover() == 0) hit = (AliMUONHitForRec*) 
1689                                    trackK->GetHitOnTrack()->Last(); // last hit
1690     //else hit = (AliMUONHitForRec*) (*trackK->GetHitOnTrack())[1]; // 2'nd hit
1691     else hit = trackK->GetHitLastOk(); // hit where track stopped
1692     if (hit) ichamBeg = hit->GetChamberNumber();
1693     ichamEnd = 0;
1694     // Check propagation direction
1695     if (!hit) { ichamBeg = ichamEnd; AliFatal(" ??? "); }
1696     else if (trackK->GetTrackDir() < 0) {
1697       ichamEnd = 9; // forward propagation
1698       ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kFALSE,zDipole1,zDipole2);
1699       if (ok) {
1700         ichamBeg = ichamEnd;
1701         ichamEnd = 6; // backward propagation
1702         // Change weight matrix and zero fChi2 for backpropagation
1703         trackK->StartBack();
1704         trackK->SetTrackDir(1);
1705         ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kTRUE,zDipole1,zDipole2);
1706         ichamBeg = ichamEnd;
1707         ichamEnd = 0;
1708       }
1709     } else {
1710       if (trackK->GetBPFlag()) {
1711         // backpropagation
1712         ichamEnd = 6; // backward propagation
1713         // Change weight matrix and zero fChi2 for backpropagation
1714         trackK->StartBack();
1715         ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kTRUE,zDipole1,zDipole2);
1716         ichamBeg = ichamEnd;
1717         ichamEnd = 0;
1718       }
1719     }
1720
1721     if (ok) {
1722       trackK->SetTrackDir(1);
1723       trackK->SetBPFlag(kFALSE);
1724       ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kFALSE,zDipole1,zDipole2);
1725     }
1726     if (!ok) { trackK->SetRecover(-1); continue; } // mark candidate to be removed
1727
1728     // Apply smoother
1729     if (trackK->GetRecover() >= 0) {
1730       ok = trackK->Smooth();
1731       if (!ok) trackK->SetRecover(-1); // mark candidate to be removed
1732     }
1733
1734     // Majority 3 of 4 in first 2 stations
1735     if (!ok) continue;
1736     chamBits = 0;
1737     Double_t chi2max = 0;
1738     for (Int_t i=0; i<trackK->GetNTrackHits(); i++) {
1739       hit = (AliMUONHitForRec*) (*trackK->GetHitOnTrack())[i];
1740       chamBits |= BIT(hit->GetChamberNumber());
1741       if (trackK->GetChi2PerPoint(i) > chi2max) chi2max = trackK->GetChi2PerPoint(i);
1742     }
1743     if (!((chamBits&3)==3 || (chamBits>>2&3)==3) && chi2max > 25) {
1744       //trackK->Recover();
1745       trackK->SetRecover(-1); //mark candidate to be removed
1746       continue;
1747     }
1748     if (ok) trackK->SetTrackQuality(0); // compute "track quality"
1749   } // while
1750
1751   for (Int_t i=0; i<fNRecTracks; i++) {
1752     trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[i]);
1753     if (trackK->GetRecover() < 0) fRecTracksPtr->RemoveAt(i);
1754   }
1755
1756   // Compress TClonesArray
1757   fRecTracksPtr->Compress();
1758   fNRecTracks = fRecTracksPtr->GetEntriesFast();
1759   return;
1760 }
1761
1762 //__________________________________________________________________________
1763 Bool_t AliMUONTrackReconstructor::CheckCandidateK(Int_t icand, Int_t nSeeds) const
1764 {
1765   // Discards track candidate if it will produce the double track (having
1766   // the same seed segment hits as hits of a good track found before)
1767   AliMUONTrackK *track1, *track2;
1768   AliMUONHitForRec *hit1, *hit2, *hit;
1769
1770   track1 = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]);
1771   hit1 = (AliMUONHitForRec*) (*track1->GetHitOnTrack())[0]; // 1'st hit
1772   hit2 = (AliMUONHitForRec*) (*track1->GetHitOnTrack())[1]; // 2'nd hit
1773
1774   for (Int_t i=0; i<icand; i++) {
1775     track2 = (AliMUONTrackK*) ((*fRecTracksPtr)[i]);
1776     //if (track2->GetRecover() < 0) continue;
1777     if (track2->GetRecover() < 0 && icand >= nSeeds) continue;
1778
1779     if (track1->GetStartSegment() == track2->GetStartSegment()) {
1780       return kFALSE;
1781     } else {
1782       Int_t nSame = 0;
1783       for (Int_t j=0; j<track2->GetNTrackHits(); j++) {
1784         hit = (AliMUONHitForRec*) (*track2->GetHitOnTrack())[j];
1785         if (hit == hit1 || hit == hit2) {
1786           nSame++;
1787           if (nSame == 2) return kFALSE;
1788         }
1789       } // for (Int_t j=0;
1790     }
1791   } // for (Int_t i=0;
1792   return kTRUE;
1793 }
1794
1795 //__________________________________________________________________________
1796 void AliMUONTrackReconstructor::RemoveDoubleTracksK(void)
1797 {
1798   // Removes double tracks (sharing more than half of their hits). Keeps
1799   // the track with higher quality
1800   AliMUONTrackK *track1, *track2, *trackToKill;
1801
1802   // Sort tracks according to their quality
1803   fRecTracksPtr->Sort();
1804
1805   // Loop over first track of the pair
1806   track1 = (AliMUONTrackK*) fRecTracksPtr->First();
1807   Int_t debug = track1->DebugLevel();
1808   while (track1) {
1809     // Loop over second track of the pair
1810     track2 = (AliMUONTrackK*) fRecTracksPtr->After(track1);
1811     while (track2) {
1812       // Check whether or not to keep track2
1813       if (!track2->KeepTrack(track1)) {
1814         if (debug >= 0) cout << " Killed track: " << 1/(*track2->GetTrackParameters())(4,0) <<
1815           " " << track2->GetTrackQuality() << endl;
1816         trackToKill = track2;
1817         track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2);
1818         trackToKill->Kill();
1819         fRecTracksPtr->Compress();
1820       } else track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2);
1821     } // track2
1822     track1 = (AliMUONTrackK*) fRecTracksPtr->After(track1);
1823   } // track1
1824
1825   fNRecTracks = fRecTracksPtr->GetEntriesFast();
1826   if (debug >= 0) cout << " Number of Kalman tracks: " << fNRecTracks << endl;
1827 }
1828
1829 //__________________________________________________________________________
1830 void AliMUONTrackReconstructor::GoToVertex(void)
1831 {
1832   // Propagates track to the vertex thru absorber
1833   // (using Branson correction for now)
1834
1835   Double_t zVertex;
1836   zVertex = 0;
1837   for (Int_t i=0; i<fNRecTracks; i++) {
1838     //((AliMUONTrackK*)(*fRecTracksPtr)[i])->Branson();
1839     ((AliMUONTrackK*)(*fRecTracksPtr)[i])->SetTrackQuality(1); // compute Chi2
1840     //((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToZ(zVertex); // w/out absorber
1841     ((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToVertex(1); // with absorber
1842   }
1843 }
1844
1845 //__________________________________________________________________________
1846 void AliMUONTrackReconstructor::SetTrackMethod(Int_t iTrackMethod)
1847 {
1848   // Set track method and recreate track container if necessary
1849   
1850   fTrackMethod = TMath::Min (iTrackMethod, 3);
1851   fTrackMethod = TMath::Max (fTrackMethod, 1);
1852   if (fTrackMethod != 1) {
1853     if (fRecTracksPtr) delete fRecTracksPtr;
1854     fRecTracksPtr = new TClonesArray("AliMUONTrackK", 10);
1855     if (fTrackMethod == 2) cout << " *** Tracking with the Kalman filter *** " << endl;
1856     else cout << " *** Combined cluster / track finder ***" << endl;
1857   } else cout << " *** Traditional tracking *** " << endl;
1858
1859 }