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