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