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