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