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