]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONEventReconstructor.cxx
memory leak fixed in track reconstruction and obsolete code removed
[u/mrichter/AliRoot.git] / MUON / AliMUONEventReconstructor.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 event reconstructor in ALICE
21 //
22 // This class contains as data:
23 // * the parameters for the event 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> //AZ
40
41 #include "AliMUONEventReconstructor.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" //AZ
58 #include "AliMC.h"
59 #include "AliLog.h"
60 #include "AliTrackReference.h"
61
62 //************* Defaults parameters for reconstruction
63 const Double_t AliMUONEventReconstructor::fgkDefaultMinBendingMomentum = 3.0;
64 const Double_t AliMUONEventReconstructor::fgkDefaultMaxBendingMomentum = 2000.0;
65 const Double_t AliMUONEventReconstructor::fgkDefaultMaxChi2 = 100.0;
66 const Double_t AliMUONEventReconstructor::fgkDefaultMaxSigma2Distance = 16.0;
67 const Double_t AliMUONEventReconstructor::fgkDefaultBendingResolution = 0.01;
68 const Double_t AliMUONEventReconstructor::fgkDefaultNonBendingResolution = 0.144;
69 const Double_t AliMUONEventReconstructor::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 AliMUONEventReconstructor::fgkDefaultSimpleBValue = 7.0;
77 const Double_t AliMUONEventReconstructor::fgkDefaultSimpleBLength = 428.0;
78 const Double_t AliMUONEventReconstructor::fgkDefaultSimpleBPosition = 1019.0;
79 const Int_t    AliMUONEventReconstructor::fgkDefaultRecTrackRefHits = 0;
80 const Double_t AliMUONEventReconstructor::fgkDefaultEfficiency = 0.95;
81
82 ClassImp(AliMUONEventReconstructor) // Class implementation in ROOT context
83
84   //__________________________________________________________________________
85 AliMUONEventReconstructor::AliMUONEventReconstructor(AliLoader* loader)
86   : TObject()
87 {
88   // Constructor for class AliMUONEventReconstructor
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   //AZ fRecTrackHitsPtr = new TClonesArray("AliMUONTrack", 100);
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,"AliMUONEventReconstructor 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
139   return;
140 }
141   //__________________________________________________________________________
142 AliMUONEventReconstructor::AliMUONEventReconstructor (const AliMUONEventReconstructor& rhs)
143   : TObject(rhs)
144 {
145 // Protected copy constructor
146
147   AliFatal("Not implemented.");
148 }
149
150 AliMUONEventReconstructor & 
151 AliMUONEventReconstructor::operator=(const AliMUONEventReconstructor& 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 AliMUONEventReconstructor::~AliMUONEventReconstructor(void)
164 {
165   // Destructor for class AliMUONEventReconstructor
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 AliMUONEventReconstructor::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 AliMUONEventReconstructor::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 AliMUONEventReconstructor::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 AliMUONEventReconstructor::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 AliMUONEventReconstructor::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 AliMUONEventReconstructor::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 AliMUONEventReconstructor::EventReconstructTrigger(void)
369 {
370   // To reconstruct one event
371   AliDebug(1,"Enter EventReconstructTrigger");
372   MakeTriggerTracks();  
373   return;
374 }
375
376   //__________________________________________________________________________
377 void AliMUONEventReconstructor::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 AliMUONEventReconstructor::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 AliMUONEventReconstructor::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 AliMUONEventReconstructor::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 AliMUONEventReconstructor::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     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 AliMUONEventReconstructor::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* AliMUONEventReconstructor::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* AliMUONEventReconstructor::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 AliMUONEventReconstructor::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 AliMUONEventReconstructor::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   nTRentries = Int_t(TR->GetEntries());
717   if (nTRentries != 1) {
718     AliError(Form("nTRentries = %d not equal to 1 ",nTRentries));
719     exit(0);
720   }
721   fLoader->TreeR()->GetEvent(0); // only one entry  
722
723   // Loop over tracking chambers
724   for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
725     // number of HitsForRec to 0 for the chamber
726     fNHitsForRecPerChamber[ch] = 0;
727     // index of first HitForRec for the chamber
728     if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0;
729     else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec;
730     rawclusters =fMUONData->RawClusters(ch);
731 //     pMUON->ResetRawClusters();
732 //     TR->GetEvent((Int_t) (TR->GetEntries()) - 1); // to be checked ????
733     nclus = (Int_t) (rawclusters->GetEntries());
734     // Loop over (cathode correlated) raw clusters
735     for (iclus = 0; iclus < nclus; iclus++) {
736       clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(iclus);
737       // new AliMUONHitForRec from raw cluster
738       // and increment number of AliMUONHitForRec's (total and in chamber)
739       hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(clus);
740       fNHitsForRec++;
741       (fNHitsForRecPerChamber[ch])++;
742       // more information into HitForRec
743       //  resolution: info should be already in raw cluster and taken from it ????
744       hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
745       hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
746       //  original raw cluster
747       hitForRec->SetChamberNumber(ch);
748       hitForRec->SetHitNumber(iclus);
749       // Z coordinate of the raw cluster (cm)
750       hitForRec->SetZ(clus->GetZ(0));
751       if (AliLog::GetGlobalDebugLevel() > 1) {
752         cout << "chamber (0...): " << ch <<
753           " raw cluster (0...): " << iclus << endl;
754         clus->Dump();
755         cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
756         hitForRec->Dump();}
757     } // end of cluster loop
758   } // end of chamber loop
759   SortHitsForRecWithIncreasingChamber(); //AZ 
760   return;
761 }
762
763   //__________________________________________________________________________
764 void AliMUONEventReconstructor::MakeSegments(void)
765 {
766   // To make the list of segments in all stations,
767   // from the list of hits to be reconstructed
768   AliDebug(1,"Enter MakeSegments");
769   ResetSegments();
770   // Loop over stations
771   Int_t nb = (fTrackMethod == 2) ? 3 : 0; //AZ
772   //AZ for (Int_t st = 0; st < fgkMaxMuonTrackingStations; st++)
773   for (Int_t st = nb; st < AliMUONConstants::NTrackingCh()/2; st++) //AZ
774     MakeSegmentsPerStation(st); 
775   if (AliLog::GetGlobalDebugLevel() > 1) {
776     cout << "end of MakeSegments" << endl;
777     for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) {
778       cout << "station(0...): " << st
779            << "  Segments: " << fNSegments[st]
780            << endl;
781       for (Int_t seg = 0;
782            seg < fNSegments[st];
783            seg++) {
784         cout << "Segment index(0...): " << seg << endl;
785         ((*fSegmentsPtr[st])[seg])->Dump();
786       }
787     }
788   }
789   return;
790 }
791
792   //__________________________________________________________________________
793 void AliMUONEventReconstructor::MakeSegmentsPerStation(Int_t Station)
794 {
795   // To make the list of segments in station number "Station" (0...)
796   // from the list of hits to be reconstructed.
797   // Updates "fNSegments"[Station].
798   // Segments in stations 4 and 5 are sorted
799   // according to increasing absolute value of "impact parameter"
800   AliMUONHitForRec *hit1Ptr, *hit2Ptr;
801   AliMUONSegment *segment;
802   Bool_t last2st;
803   Double_t bendingSlope, distBend, distNonBend, extBendCoor, extNonBendCoor,
804       impactParam = 0., maxImpactParam = 0., minImpactParam = 0.; // =0 to avoid compilation warnings.
805   AliDebug(1,Form("Enter MakeSegmentsPerStation (0...) %d",Station));
806   // first and second chambers (0...) in the station
807   Int_t ch1 = 2 * Station;
808   Int_t ch2 = ch1 + 1;
809   // variable true for stations downstream of the dipole:
810   // Station(0..4) equal to 3 or 4
811   if ((Station == 3) || (Station == 4)) {
812     last2st = kTRUE;
813     // maximum impact parameter (cm) according to fMinBendingMomentum...
814     maxImpactParam =
815       TMath::Abs(GetImpactParamFromBendingMomentum(fMinBendingMomentum));
816     // minimum impact parameter (cm) according to fMaxBendingMomentum...
817     minImpactParam =
818       TMath::Abs(GetImpactParamFromBendingMomentum(fMaxBendingMomentum));
819   }
820   else last2st = kFALSE;
821   // extrapolation factor from Z of first chamber to Z of second chamber
822   // dZ to be changed to take into account fine structure of chambers ????
823   Double_t extrapFact;
824   // index for current segment
825   Int_t segmentIndex = 0;
826   // Loop over HitsForRec in the first chamber of the station
827   for (Int_t hit1 = fIndexOfFirstHitForRecPerChamber[ch1];
828        hit1 < fIndexOfFirstHitForRecPerChamber[ch1] + fNHitsForRecPerChamber[ch1];
829        hit1++) {
830     // pointer to the HitForRec
831     hit1Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit1]);
832     // extrapolation,
833     // on the straight line joining the HitForRec to the vertex (0,0,0),
834     // to the Z of the second chamber of the station
835     // Loop over HitsForRec in the second chamber of the station
836     for (Int_t hit2 = fIndexOfFirstHitForRecPerChamber[ch2];
837          hit2 < fIndexOfFirstHitForRecPerChamber[ch2] + fNHitsForRecPerChamber[ch2];
838          hit2++) {
839       // pointer to the HitForRec
840       hit2Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit2]);
841       // absolute values of distances, in bending and non bending planes,
842       // between the HitForRec in the second chamber
843       // and the previous extrapolation
844       extrapFact = hit2Ptr->GetZ()/ hit1Ptr->GetZ();
845       extBendCoor = extrapFact * hit1Ptr->GetBendingCoor();
846       extNonBendCoor = extrapFact * hit1Ptr->GetNonBendingCoor();
847       distBend = TMath::Abs(hit2Ptr->GetBendingCoor() - extBendCoor);
848       distNonBend = TMath::Abs(hit2Ptr->GetNonBendingCoor() - extNonBendCoor);
849       if (last2st) {
850         // bending slope
851         bendingSlope = (hit1Ptr->GetBendingCoor() - hit2Ptr->GetBendingCoor()) /
852           (hit1Ptr->GetZ() - hit2Ptr->GetZ());
853         // absolute value of impact parameter
854         impactParam =
855           TMath::Abs(hit1Ptr->GetBendingCoor() - hit1Ptr->GetZ() * bendingSlope);
856       }
857       // check for distances not too large,
858       // and impact parameter not too big if stations downstream of the dipole.
859       // Conditions "distBend" and "impactParam" correlated for these stations ????
860       if ((distBend < fSegmentMaxDistBending[Station]) &&
861           (distNonBend < fSegmentMaxDistNonBending[Station]) &&
862           (!last2st || (impactParam < maxImpactParam)) &&
863           (!last2st || (impactParam > minImpactParam))) {
864         // make new segment
865         segment = new ((*fSegmentsPtr[Station])[segmentIndex])
866           AliMUONSegment(hit1Ptr, hit2Ptr);
867         // update "link" to this segment from the hit in the first chamber
868         if (hit1Ptr->GetNSegments() == 0)
869           hit1Ptr->SetIndexOfFirstSegment(segmentIndex);
870         hit1Ptr->SetNSegments(hit1Ptr->GetNSegments() + 1);
871         if (AliLog::GetGlobalDebugLevel() > 1) {
872           cout << "segmentIndex(0...): " << segmentIndex
873                << "  distBend: " << distBend
874                << "  distNonBend: " << distNonBend
875                << endl;
876           segment->Dump();
877           cout << "HitForRec in first chamber" << endl;
878           hit1Ptr->Dump();
879           cout << "HitForRec in second chamber" << endl;
880           hit2Ptr->Dump();
881         };
882         // increment index for current segment
883         segmentIndex++;
884       }
885     } //for (Int_t hit2
886   } // for (Int_t hit1...
887   fNSegments[Station] = segmentIndex;
888   // Sorting according to "impact parameter" if station(1..5) 4 or 5,
889   // i.e. Station(0..4) 3 or 4, using the function "Compare".
890   // After this sorting, it is impossible to use
891   // the "fNSegments" and "fIndexOfFirstSegment"
892   // of the HitForRec in the first chamber to explore all segments formed with it.
893   // Is this sorting really needed ????
894   if ((Station == 3) || (Station == 4)) (fSegmentsPtr[Station])->Sort();
895   AliDebug(1,Form("Station: %d  NSegments:  %d ", Station, fNSegments[Station]));
896   return;
897 }
898
899   //__________________________________________________________________________
900 void AliMUONEventReconstructor::MakeTracks(void)
901 {
902   // To make the tracks,
903   // from the list of segments and points in all stations
904   AliDebug(1,"Enter MakeTracks");
905   // The order may be important for the following Reset's
906   ResetTracks();
907   ResetTrackHits();
908   if (fTrackMethod == 2) { //AZ - Kalman filter
909     MakeTrackCandidatesK();
910     if (fRecTracksPtr->GetEntriesFast() == 0) return;
911     // Follow tracks in stations(1..) 3, 2 and 1
912     FollowTracksK();
913     // Remove double tracks
914     RemoveDoubleTracksK();
915     // Propagate tracks to the vertex thru absorber
916     GoToVertex();
917     // Fill AliMUONTrack data members
918     FillMUONTrack();
919   } else { 
920     // Look for candidates from at least 3 aligned points in stations(1..) 4 and 5
921     MakeTrackCandidates();
922     // Follow tracks in stations(1..) 3, 2 and 1
923     FollowTracks();
924     // Remove double tracks
925     RemoveDoubleTracks();
926     UpdateTrackParamAtHit();
927     UpdateHitForRecAtHit();
928   }
929   return;
930 }
931
932   //__________________________________________________________________________
933 void AliMUONEventReconstructor::ValidateTracksWithTrigger(void)
934 {
935   // Try to match track from tracking system with trigger track
936   AliMUONTrack *track;
937   TClonesArray *recTriggerTracks;
938   
939   fMUONData->SetTreeAddress("RL");
940   fMUONData->GetRecTriggerTracks();
941   recTriggerTracks = fMUONData->RecTriggerTracks();
942
943   track = (AliMUONTrack*) fRecTracksPtr->First();
944   while (track) {
945     track->MatchTriggerTrack(recTriggerTracks);
946     track = (AliMUONTrack*) fRecTracksPtr->After(track);
947   }
948
949 }
950
951   //__________________________________________________________________________
952 Bool_t AliMUONEventReconstructor::MakeTriggerTracks(void)
953 {
954     // To make the trigger tracks from Local Trigger
955   AliDebug(1, "Enter MakeTriggerTracks");
956     
957     Int_t nTRentries;
958     Long_t gloTrigPat;
959     TClonesArray *localTrigger;
960     TClonesArray *globalTrigger;
961     AliMUONLocalTrigger *locTrg;
962     AliMUONGlobalTrigger *gloTrg;
963     AliMUONTriggerCircuit *circuit;
964     AliMUONTriggerTrack *recTriggerTrack = 0;
965
966     TTree* treeR = fLoader->TreeR();
967     
968     // Loading MUON subsystem
969     AliMUON * pMUON = (AliMUON *) gAlice->GetDetector("MUON");
970     
971     nTRentries = Int_t(treeR->GetEntries());
972      
973     treeR->GetEvent(0); // only one entry  
974
975     if (!(fMUONData->IsTriggerBranchesInTree())) {
976       AliWarning(Form("Trigger information is not avalaible, nTRentries = %d not equal to 1",nTRentries));
977       return kFALSE;
978     }
979
980     fMUONData->SetTreeAddress("TC");
981     fMUONData->GetTrigger();
982
983     // global trigger for trigger pattern
984     gloTrigPat = 0;
985     globalTrigger = fMUONData->GlobalTrigger(); 
986     gloTrg = (AliMUONGlobalTrigger*)globalTrigger->UncheckedAt(0);      
987     if (gloTrg) {
988       if (gloTrg->SinglePlusLpt())  gloTrigPat|= 0x1;
989       if (gloTrg->SinglePlusHpt())  gloTrigPat|= 0x2;
990       if (gloTrg->SinglePlusApt())  gloTrigPat|= 0x4;
991       
992       if (gloTrg->SingleMinusLpt()) gloTrigPat|= 0x8;
993       if (gloTrg->SingleMinusHpt()) gloTrigPat|= 0x10;
994       if (gloTrg->SingleMinusApt()) gloTrigPat|= 0x20;
995       
996       if (gloTrg->SingleUndefLpt()) gloTrigPat|= 0x40;
997       if (gloTrg->SingleUndefHpt()) gloTrigPat|= 0x80;
998       if (gloTrg->SingleUndefApt()) gloTrigPat|= 0x100;
999       
1000       if (gloTrg->PairUnlikeLpt())  gloTrigPat|= 0x200;
1001       if (gloTrg->PairUnlikeHpt())  gloTrigPat|= 0x400;
1002       if (gloTrg->PairUnlikeApt())  gloTrigPat|= 0x800;
1003       
1004       if (gloTrg->PairLikeLpt())    gloTrigPat|= 0x1000;
1005       if (gloTrg->PairLikeHpt())    gloTrigPat|= 0x2000;
1006       if (gloTrg->PairLikeApt())    gloTrigPat|= 0x4000;
1007     }
1008
1009     // local trigger for tracking 
1010     localTrigger = fMUONData->LocalTrigger();    
1011     Int_t nlocals = (Int_t) (localTrigger->GetEntries());
1012
1013     Float_t z11 = AliMUONConstants::DefaultChamberZ(10);
1014     Float_t z21 = AliMUONConstants::DefaultChamberZ(12);
1015
1016     for (Int_t i=0; i<nlocals; i++) { // loop on Local Trigger
1017       locTrg = (AliMUONLocalTrigger*)localTrigger->UncheckedAt(i);      
1018       circuit = &(pMUON->TriggerCircuit(locTrg->LoCircuit()));
1019       Float_t y11 = circuit->GetY11Pos(locTrg->LoStripX()); 
1020       Int_t stripX21 = locTrg->LoStripX()+locTrg->LoDev()+1;
1021       Float_t y21 = circuit->GetY21Pos(stripX21);       
1022       Float_t x11 = circuit->GetX11Pos(locTrg->LoStripY());
1023       Float_t thetax = TMath::ATan2( x11 , z11 );
1024       Float_t thetay = TMath::ATan2( (y21-y11) , (z21-z11) );
1025
1026       recTriggerTrack = new AliMUONTriggerTrack(x11,y11,thetax,thetay,gloTrigPat);
1027       
1028       // since static statement does not work, set gloTrigPat for each track
1029
1030       fMUONData->AddRecTriggerTrack(*recTriggerTrack);
1031       delete recTriggerTrack;
1032     } // end of loop on Local Trigger
1033     return kTRUE;    
1034 }
1035
1036   //__________________________________________________________________________
1037 Int_t AliMUONEventReconstructor::MakeTrackCandidatesWithTwoSegments(AliMUONSegment *BegSegment)
1038 {
1039   // To make track candidates with two segments in stations(1..) 4 and 5,
1040   // the first segment being pointed to by "BegSegment".
1041   // Returns the number of such track candidates.
1042   Int_t endStation, iEndSegment, nbCan2Seg;
1043   AliMUONSegment *endSegment;
1044   AliMUONSegment *extrapSegment = NULL;
1045   AliMUONTrack *recTrack;
1046   Double_t mcsFactor;
1047   AliDebug(1,"Enter MakeTrackCandidatesWithTwoSegments");
1048   // Station for the end segment
1049   endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
1050   // multiple scattering factor corresponding to one chamber
1051   mcsFactor = 0.0136 /
1052     GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
1053   mcsFactor     = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1054   // linear extrapolation to end station
1055   // number of candidates with 2 segments to 0
1056   nbCan2Seg = 0;
1057   // Loop over segments in the end station
1058   for (iEndSegment = 0; iEndSegment < fNSegments[endStation]; iEndSegment++) {
1059     // pointer to segment
1060     endSegment = (AliMUONSegment*) ((*fSegmentsPtr[endStation])[iEndSegment]);
1061     // test compatibility between current segment and "extrapSegment"
1062     // 4 because 4 quantities in chi2
1063     extrapSegment =
1064       BegSegment->CreateSegmentFromLinearExtrapToStation(endSegment->GetZ(), mcsFactor);
1065     if ((endSegment->
1066          NormalizedChi2WithSegment(extrapSegment,
1067                                    fMaxSigma2Distance)) <= 4.0) {
1068       // both segments compatible:
1069       // make track candidate from "begSegment" and "endSegment"
1070       AliDebug(2,Form("TrackCandidate with Segment %d in Station(0..) %d", iEndSegment, endStation));
1071       // flag for both segments in one track:
1072       // to be done in track constructor ????
1073       BegSegment->SetInTrack(kTRUE);
1074       endSegment->SetInTrack(kTRUE);
1075       recTrack = new ((*fRecTracksPtr)[fNRecTracks])
1076         AliMUONTrack(BegSegment, endSegment, this);
1077       fNRecTracks++;
1078       if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump();
1079       // increment number of track candidates with 2 segments
1080       nbCan2Seg++;
1081     }
1082     delete extrapSegment; // should not delete HitForRec's it points to !!!!
1083   } // for (iEndSegment = 0;...
1084   return nbCan2Seg;
1085 }
1086
1087   //__________________________________________________________________________
1088 Int_t AliMUONEventReconstructor::MakeTrackCandidatesWithOneSegmentAndOnePoint(AliMUONSegment *BegSegment)
1089 {
1090   // To make track candidates with one segment and one point
1091   // in stations(1..) 4 and 5,
1092   // the segment being pointed to by "BegSegment".
1093   Int_t ch, ch1, ch2, endStation, iHit, iHitMax, iHitMin, nbCan1Seg1Hit;
1094   AliMUONHitForRec *extrapHitForRec= NULL;
1095   AliMUONHitForRec *hit;
1096   AliMUONTrack *recTrack;
1097   Double_t mcsFactor;
1098   AliDebug(1,"Enter MakeTrackCandidatesWithOneSegmentAndOnePoint");
1099   // station for the end point
1100   endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
1101   // multiple scattering factor corresponding to one chamber
1102   mcsFactor = 0.0136 /
1103     GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
1104   mcsFactor     = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1105   // first and second chambers(0..) in the end station
1106   ch1 = 2 * endStation;
1107   ch2 = ch1 + 1;
1108   // number of candidates to 0
1109   nbCan1Seg1Hit = 0;
1110   // Loop over chambers of the end station
1111   for (ch = ch2; ch >= ch1; ch--) {
1112     // limits for the hit index in the loop
1113     iHitMin = fIndexOfFirstHitForRecPerChamber[ch];
1114     iHitMax = iHitMin + fNHitsForRecPerChamber[ch];
1115     // Loop over HitForRec's in the chamber
1116     for (iHit = iHitMin; iHit < iHitMax; iHit++) {
1117       // pointer to HitForRec
1118       hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1119       // test compatibility between current HitForRec and "extrapHitForRec"
1120       // 2 because 2 quantities in chi2
1121       // linear extrapolation to chamber
1122       extrapHitForRec =
1123         BegSegment->CreateHitForRecFromLinearExtrapToChamber( hit->GetZ(), mcsFactor);
1124       if ((hit->
1125            NormalizedChi2WithHitForRec(extrapHitForRec,
1126                                        fMaxSigma2Distance)) <= 2.0) {
1127         // both HitForRec's compatible:
1128         // make track candidate from begSegment and current HitForRec
1129         AliDebug(2, Form("TrackCandidate with HitForRec  %d in Chamber(0..) %d", iHit, ch));
1130         // flag for beginning segments in one track:
1131         // to be done in track constructor ????
1132         BegSegment->SetInTrack(kTRUE);
1133         recTrack = new ((*fRecTracksPtr)[fNRecTracks])
1134           AliMUONTrack(BegSegment, hit, this);
1135         // the right place to eliminate "double counting" ???? how ????
1136         fNRecTracks++;
1137         if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump();
1138         // increment number of track candidates
1139         nbCan1Seg1Hit++;
1140       }
1141       delete extrapHitForRec;
1142     } // for (iHit = iHitMin;...
1143   } // for (ch = ch2;...
1144   return nbCan1Seg1Hit;
1145 }
1146
1147   //__________________________________________________________________________
1148 void AliMUONEventReconstructor::MakeTrackCandidates(void)
1149 {
1150   // To make track candidates
1151   // with at least 3 aligned points in stations(1..) 4 and 5
1152   // (two Segment's or one Segment and one HitForRec)
1153   Int_t begStation, iBegSegment, nbCan1Seg1Hit, nbCan2Seg;
1154   AliMUONSegment *begSegment;
1155   AliDebug(1,"Enter MakeTrackCandidates");
1156   // Loop over stations(1..) 5 and 4 for the beginning segment
1157   for (begStation = 4; begStation > 2; begStation--) {
1158     // Loop over segments in the beginning station
1159     for (iBegSegment = 0; iBegSegment < fNSegments[begStation]; iBegSegment++) {
1160       // pointer to segment
1161       begSegment = (AliMUONSegment*) ((*fSegmentsPtr[begStation])[iBegSegment]);
1162       AliDebug(2,Form("Look for TrackCandidate's with Segment %d  in Station(0..) %d", iBegSegment, begStation));
1163       // Look for track candidates with two segments,
1164       // "begSegment" and all compatible segments in other station.
1165       // Only for beginning station(1..) 5
1166       // because candidates with 2 segments have to looked for only once.
1167       if (begStation == 4)
1168         nbCan2Seg = MakeTrackCandidatesWithTwoSegments(begSegment);
1169       // Look for track candidates with one segment and one point,
1170       // "begSegment" and all compatible HitForRec's in other station.
1171       // Only if "begSegment" does not belong already to a track candidate.
1172       // Is that a too strong condition ????
1173       if (!(begSegment->GetInTrack()))
1174         nbCan1Seg1Hit = MakeTrackCandidatesWithOneSegmentAndOnePoint(begSegment);
1175     } // for (iBegSegment = 0;...
1176   } // for (begStation = 4;...
1177   return;
1178 }
1179
1180   //__________________________________________________________________________
1181 void AliMUONEventReconstructor::FollowTracks(void)
1182 {
1183   // Follow tracks in stations(1..) 3, 2 and 1
1184   // too long: should be made more modular !!!!
1185   AliMUONHitForRec *bestHit, *extrapHit, *hit;
1186   AliMUONSegment *bestSegment, *extrapSegment, *segment;
1187   AliMUONTrack *track, *nextTrack;
1188   AliMUONTrackParam *trackParam1, trackParam[2], trackParamVertex;
1189   // -1 to avoid compilation warnings
1190   Int_t ch = -1, chInStation, chBestHit = -1, iHit, iSegment, station, trackIndex; 
1191   Double_t bestChi2, chi2, dZ1, dZ2, dZ3, maxSigma2Distance, mcsFactor;
1192   Double_t bendingMomentum, chi2Norm = 0.;
1193   AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
1194   // local maxSigma2Distance, for easy increase in testing
1195   maxSigma2Distance = fMaxSigma2Distance;
1196   AliDebug(2,"Enter FollowTracks");
1197   // Loop over track candidates
1198   track = (AliMUONTrack*) fRecTracksPtr->First();
1199   trackIndex = -1;
1200   while (track) {
1201     // Follow function for each track candidate ????
1202     trackIndex++;
1203     nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track); // prepare next track
1204     AliDebug(2,Form("FollowTracks: track candidate(0..): %d", trackIndex));
1205     // Fit track candidate
1206     track->SetFitMCS(0); // without multiple Coulomb scattering
1207     track->SetFitNParam(3); // with 3 parameters (X = Y = 0)
1208     track->SetFitStart(0); // from parameters at vertex
1209     track->Fit();
1210     if (AliLog::GetGlobalDebugLevel()> 2) {
1211       cout << "FollowTracks: track candidate(0..): " << trackIndex
1212            << " after fit in stations(0..) 3 and 4" << endl;
1213       track->RecursiveDump();
1214     }
1215     // Loop over stations(1..) 3, 2 and 1
1216     // something SPECIAL for stations 2 and 1 for majority 3 coincidence ????
1217     // otherwise: majority coincidence 2 !!!!
1218     for (station = 2; station >= 0; station--) {
1219       // Track parameters at first track hit (smallest Z)
1220       trackParam1 = ((AliMUONTrackHit*)
1221                      (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1222       // extrapolation to station
1223       trackParam1->ExtrapToStation(station, trackParam);
1224       extrapSegment = new AliMUONSegment(); //  empty segment
1225       // multiple scattering factor corresponding to one chamber
1226       // and momentum in bending plane (not total)
1227       mcsFactor = 0.0136 * trackParam1->GetInverseBendingMomentum();
1228       mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1229       // Z difference from previous station
1230       dZ1 = (&(pMUON->Chamber(2 * station)))->Z() -
1231         (&(pMUON->Chamber(2 * station + 2)))->Z();
1232       // Z difference between the two previous stations
1233       dZ2 = (&(pMUON->Chamber(2 * station + 2)))->Z() -
1234         (&(pMUON->Chamber(2 * station + 4)))->Z();
1235       // Z difference between the two chambers in the previous station
1236       dZ3 = (&(pMUON->Chamber(2 * station)))->Z() -
1237         (&(pMUON->Chamber(2 * station + 1)))->Z();
1238       extrapSegment->SetBendingCoorReso2(fBendingResolution * fBendingResolution);
1239       extrapSegment->
1240         SetNonBendingCoorReso2(fNonBendingResolution * fNonBendingResolution);
1241       extrapSegment->UpdateFromStationTrackParam
1242         (trackParam, mcsFactor, dZ1, dZ2, dZ3, station,
1243          trackParam1->GetInverseBendingMomentum());
1244       bestChi2 = 5.0;
1245       bestSegment = NULL;
1246       if (AliLog::GetGlobalDebugLevel() > 2) {
1247         cout << "FollowTracks: track candidate(0..): " << trackIndex
1248              << " Look for segment in station(0..): " << station << endl;
1249       }
1250       // Loop over segments in station
1251       for (iSegment = 0; iSegment < fNSegments[station]; iSegment++) {
1252         // Look for best compatible Segment in station
1253         // should consider all possibilities ????
1254         // multiple scattering ????
1255         // separation in 2 functions: Segment and HitForRec ????
1256         segment = (AliMUONSegment*) ((*fSegmentsPtr[station])[iSegment]);
1257         // correction of corrected segment (fBendingCoor and fNonBendingCoor)
1258         // according to real Z value of "segment" and slopes of "extrapSegment"
1259         (&(trackParam[0]))->ExtrapToZ(segment->GetZ());
1260         (&(trackParam[1]))->ExtrapToZ(segment->GetZ());
1261         extrapSegment->SetBendingCoor((&(trackParam[0]))->GetBendingCoor());
1262         extrapSegment->SetNonBendingCoor((&(trackParam[0]))->GetNonBendingCoor());
1263         extrapSegment->SetBendingSlope((&(trackParam[0]))->GetBendingSlope());
1264         extrapSegment->SetNonBendingSlope((&(trackParam[0]))->GetNonBendingSlope());
1265         chi2 = segment->
1266           NormalizedChi2WithSegment(extrapSegment, maxSigma2Distance);
1267         if (chi2 < bestChi2) {
1268           // update best Chi2 and Segment if better found
1269           bestSegment = segment;
1270           bestChi2 = chi2;
1271         }
1272       }
1273       if (bestSegment) {
1274         // best segment found: add it to track candidate
1275         (&(trackParam[0]))->ExtrapToZ(bestSegment->GetZ());
1276         (&(trackParam[1]))->ExtrapToZ(bestSegment->GetZ());
1277         track->AddSegment(bestSegment);
1278         // set track parameters at these two TrakHit's
1279         track->SetTrackParamAtHit(track->GetNTrackHits() - 2, &(trackParam[0]));
1280         track->SetTrackParamAtHit(track->GetNTrackHits() - 1, &(trackParam[1]));
1281         AliDebug(3, Form("FollowTracks: track candidate(0..): %d  Added segment in station(0..): %d", trackIndex, station));
1282         if (AliLog::GetGlobalDebugLevel()>2) track->RecursiveDump();
1283       }
1284       else {
1285         // No best segment found:
1286         // Look for best compatible HitForRec in station:
1287         // should consider all possibilities ????
1288         // multiple scattering ???? do about like for extrapSegment !!!!
1289         extrapHit = new AliMUONHitForRec(); //  empty hit
1290         bestChi2 = 3.0;
1291         bestHit = NULL;
1292         AliDebug(3, Form("FollowTracks: track candidate(0..): %d Segment not found, look for hit in station(0..): %d ", 
1293                          trackIndex, station));
1294         
1295         // Loop over chambers of the station
1296         for (chInStation = 0; chInStation < 2; chInStation++) {
1297           ch = 2 * station + chInStation;
1298           for (iHit = fIndexOfFirstHitForRecPerChamber[ch];
1299                iHit < fIndexOfFirstHitForRecPerChamber[ch] +
1300                  fNHitsForRecPerChamber[ch];
1301                iHit++) {
1302             hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1303             // coordinates of extrapolated hit
1304             (&(trackParam[chInStation]))->ExtrapToZ(hit->GetZ());
1305             extrapHit->
1306               SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor());
1307             extrapHit->
1308               SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor());
1309             // resolutions from "extrapSegment"
1310             extrapHit->SetBendingReso2(extrapSegment->GetBendingCoorReso2());
1311             extrapHit->SetNonBendingReso2(extrapSegment->GetNonBendingCoorReso2());
1312             // Loop over hits in the chamber
1313             // condition for hit not already in segment ????
1314             chi2 = hit->NormalizedChi2WithHitForRec(extrapHit, maxSigma2Distance);
1315             if (chi2 < bestChi2) {
1316               // update best Chi2 and HitForRec if better found
1317               bestHit = hit;
1318               bestChi2 = chi2;
1319               chBestHit = chInStation;
1320             }
1321           }
1322         }
1323         if (bestHit) {
1324           // best hit found: add it to track candidate
1325           (&(trackParam[chBestHit]))->ExtrapToZ(bestHit->GetZ());
1326           track->AddHitForRec(bestHit);
1327           // set track parameters at this TrackHit
1328           track->SetTrackParamAtHit(track->GetNTrackHits() - 1,
1329                                     &(trackParam[chBestHit]));
1330           if (AliLog::GetGlobalDebugLevel() > 2) {
1331             cout << "FollowTracks: track candidate(0..): " << trackIndex
1332                  << " Added hit in station(0..): " << station << endl;
1333             track->RecursiveDump();
1334           }
1335         }
1336         else {
1337           // Remove current track candidate
1338           // and corresponding TrackHit's, ...
1339           track->Remove();
1340           delete extrapSegment;
1341           delete extrapHit;
1342           break; // stop the search for this candidate:
1343           // exit from the loop over station
1344         }
1345         delete extrapHit;
1346       }
1347       delete extrapSegment;
1348       // Sort track hits according to increasing Z
1349       track->GetTrackHitsPtr()->Sort();
1350       // Update track parameters at first track hit (smallest Z)
1351       trackParam1 = ((AliMUONTrackHit*)
1352                      (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1353       bendingMomentum = 0.;
1354       if (TMath::Abs(trackParam1->GetInverseBendingMomentum()) > 0.)
1355         bendingMomentum = TMath::Abs(1/(trackParam1->GetInverseBendingMomentum()));
1356       // Track removed if bendingMomentum not in window [min, max]
1357       if ((bendingMomentum < fMinBendingMomentum) || (bendingMomentum > fMaxBendingMomentum)) {
1358         track->Remove();
1359         break; // stop the search for this candidate:
1360         // exit from the loop over station 
1361       }
1362       // Track fit
1363       // with multiple Coulomb scattering if all stations
1364       if (station == 0) track->SetFitMCS(1);
1365       // without multiple Coulomb scattering if not all stations
1366       else track->SetFitMCS(0);
1367       track->SetFitNParam(5);  // with 5 parameters (momentum and position)
1368       track->SetFitStart(1);  // from parameters at first hit
1369       track->Fit();
1370       Double_t numberOfDegFree = (2.0 * track->GetNTrackHits() - 5);
1371       if (numberOfDegFree > 0) {
1372         chi2Norm =  track->GetFitFMin() / numberOfDegFree;
1373       } else {
1374         chi2Norm = 1.e10;
1375       }
1376       // Track removed if normalized chi2 too high
1377       if (chi2Norm > fMaxChi2) {
1378         track->Remove();
1379         break; // stop the search for this candidate:
1380         // exit from the loop over station 
1381       }
1382       if (AliLog::GetGlobalDebugLevel() > 2) {
1383         cout << "FollowTracks: track candidate(0..): " << trackIndex
1384              << " after fit from station(0..): " << station << " to 4" << endl;
1385         track->RecursiveDump();
1386       }
1387       // Track extrapolation to the vertex through the absorber (Branson)
1388       // after going through the first station
1389       if (station == 0) {
1390         trackParamVertex = *trackParam1;
1391         (&trackParamVertex)->ExtrapToVertex(0.,0.,0.);
1392         track->SetTrackParamAtVertex(&trackParamVertex);
1393         if (AliLog::GetGlobalDebugLevel() > 0) {
1394           cout << "FollowTracks: track candidate(0..): " << trackIndex
1395                << " after extrapolation to vertex" << endl;
1396           track->RecursiveDump();
1397         }
1398       }
1399     } // for (station = 2;...
1400     // go really to next track
1401     track = nextTrack;
1402   } // while (track)
1403   // Compression of track array (necessary after Remove ????)
1404   fRecTracksPtr->Compress();
1405   return;
1406 }
1407
1408   //__________________________________________________________________________
1409 void AliMUONEventReconstructor::RemoveDoubleTracks(void)
1410 {
1411   // To remove double tracks.
1412   // Tracks are considered identical
1413   // if they have at least half of their hits in common.
1414   // Among two identical tracks, one keeps the track with the larger number of hits
1415   // or, if these numbers are equal, the track with the minimum Chi2.
1416   AliMUONTrack *track1, *track2, *trackToRemove;
1417   Bool_t identicalTracks;
1418   Int_t hitsInCommon, nHits1, nHits2;
1419   identicalTracks = kTRUE;
1420   while (identicalTracks) {
1421     identicalTracks = kFALSE;
1422     // Loop over first track of the pair
1423     track1 = (AliMUONTrack*) fRecTracksPtr->First();
1424     while (track1 && (!identicalTracks)) {
1425       nHits1 = track1->GetNTrackHits();
1426       // Loop over second track of the pair
1427       track2 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1428       while (track2 && (!identicalTracks)) {
1429         nHits2 = track2->GetNTrackHits();
1430         // number of hits in common between two tracks
1431         hitsInCommon = track1->HitsInCommon(track2);
1432         // check for identical tracks
1433         if ((4 * hitsInCommon) >= (nHits1 + nHits2)) {
1434           identicalTracks = kTRUE;
1435           // decide which track to remove
1436           if (nHits1 > nHits2) trackToRemove = track2;
1437           else if (nHits1 < nHits2) trackToRemove = track1;
1438           else if ((track1->GetFitFMin()) < (track2->GetFitFMin()))
1439             trackToRemove = track2;
1440           else trackToRemove = track1;
1441           // remove it
1442           trackToRemove->Remove();
1443         }
1444         track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
1445       } // track2
1446       track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1447     } // track1
1448   }
1449   return;
1450 }
1451
1452   //__________________________________________________________________________
1453 void AliMUONEventReconstructor::UpdateTrackParamAtHit()
1454 {
1455   // Set track parameters after track fitting. Fill fTrackParamAtHit of AliMUONTrack's
1456   AliMUONTrack *track;
1457   AliMUONTrackHit *trackHit;
1458   AliMUONTrackParam *trackParam;
1459   track = (AliMUONTrack*) fRecTracksPtr->First();
1460   while (track) {
1461     trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->First();
1462     while (trackHit) {
1463       trackParam = trackHit->GetTrackParam();
1464       track->AddTrackParamAtHit(trackParam);
1465       trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->After(trackHit); 
1466     } // trackHit    
1467     track = (AliMUONTrack*) fRecTracksPtr->After(track);
1468   } // track
1469   return;
1470 }
1471
1472   //__________________________________________________________________________
1473 void AliMUONEventReconstructor::UpdateHitForRecAtHit()
1474 {
1475   // Set cluster parameterss after track fitting. Fill fHitForRecAtHit of AliMUONTrack's
1476   AliMUONTrack *track;
1477   AliMUONTrackHit *trackHit;
1478   AliMUONHitForRec *hitForRec;
1479   track = (AliMUONTrack*) fRecTracksPtr->First();
1480   while (track) {
1481     trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->First();
1482     while (trackHit) {
1483       hitForRec = trackHit->GetHitForRecPtr();
1484       track->AddHitForRecAtHit(hitForRec);
1485       trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->After(trackHit); 
1486     } // trackHit    
1487     track = (AliMUONTrack*) fRecTracksPtr->After(track);
1488   } // track
1489   return;
1490 }
1491
1492   //__________________________________________________________________________
1493 void AliMUONEventReconstructor::FillMUONTrack()
1494 {
1495   // Set track parameters at hits for Kalman track. Fill fTrackParamAtHit of AliMUONTrack's
1496   AliMUONTrackK *track;
1497   track = (AliMUONTrackK*) fRecTracksPtr->First();
1498   while (track) {
1499     track->FillMUONTrack();
1500     track = (AliMUONTrackK*) fRecTracksPtr->After(track);
1501   } 
1502   return;
1503 }
1504
1505   //__________________________________________________________________________
1506 void AliMUONEventReconstructor::EventDump(void)
1507 {
1508   // Dump reconstructed event (track parameters at vertex and at first hit),
1509   // and the particle parameters
1510
1511   AliMUONTrack *track;
1512   AliMUONTrackParam *trackParam, *trackParam1;
1513   Double_t bendingSlope, nonBendingSlope, pYZ;
1514   Double_t pX, pY, pZ, x, y, z, c;
1515   Int_t np, trackIndex, nTrackHits;
1516  
1517   AliDebug(1,"****** enter EventDump ******");
1518   AliDebug(1, Form("Number of Reconstructed tracks : %d", fNRecTracks)); 
1519   
1520   fRecTracksPtr->Compress(); // for simple loop without "Next" since no hole
1521   // Loop over reconstructed tracks
1522   for (trackIndex = 0; trackIndex < fNRecTracks; trackIndex++) {
1523     if (fTrackMethod != 1) continue; //AZ - skip the rest for now
1524     AliDebug(1, Form("track number: %d", trackIndex));
1525     // function for each track for modularity ????
1526     track = (AliMUONTrack*) ((*fRecTracksPtr)[trackIndex]);
1527     nTrackHits = track->GetNTrackHits();
1528     AliDebug(1, Form("Number of track hits: %d ", nTrackHits));
1529     // track parameters at Vertex
1530     trackParam = track->GetTrackParamAtVertex();
1531     x = trackParam->GetNonBendingCoor();
1532     y = trackParam->GetBendingCoor();
1533     z = trackParam->GetZ();
1534     bendingSlope = trackParam->GetBendingSlope();
1535     nonBendingSlope = trackParam->GetNonBendingSlope();
1536     pYZ = 1/TMath::Abs(trackParam->GetInverseBendingMomentum());
1537     pZ = pYZ/TMath::Sqrt(1+bendingSlope*bendingSlope);
1538     pX = pZ * nonBendingSlope;
1539     pY = pZ * bendingSlope;
1540     c = TMath::Sign(1.0, trackParam->GetInverseBendingMomentum());
1541     AliDebug(1, Form("Track parameters at Vertex z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1542                      z, x, y, pX, pY, pZ, c));
1543
1544     // track parameters at first hit
1545     trackParam1 = ((AliMUONTrackHit*)
1546                    (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1547     x = trackParam1->GetNonBendingCoor();
1548     y = trackParam1->GetBendingCoor();
1549     z = trackParam1->GetZ();
1550     bendingSlope = trackParam1->GetBendingSlope();
1551     nonBendingSlope = trackParam1->GetNonBendingSlope();
1552     pYZ = 1/TMath::Abs(trackParam1->GetInverseBendingMomentum());
1553     pZ = pYZ/TMath::Sqrt(1.0 + bendingSlope * bendingSlope);
1554     pX = pZ * nonBendingSlope;
1555     pY = pZ * bendingSlope;
1556     c = TMath::Sign(1.0, trackParam1->GetInverseBendingMomentum());
1557     AliDebug(1, Form("track parameters at z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1558                      z, x, y, pX, pY, pZ, c));
1559   }
1560   // informations about generated particles
1561   np = gAlice->GetMCApp()->GetNtrack();
1562   printf(" **** number of generated particles: %d  \n", np);
1563   
1564 //    for (Int_t iPart = 0; iPart < np; iPart++) {
1565 //      p = gAlice->Particle(iPart);
1566 //      printf(" particle %d: type= %d px= %f py= %f pz= %f pdg= %d\n",
1567 //         iPart, p->GetPdgCode(), p->Px(), p->Py(), p->Pz(), p->GetPdgCode());    
1568 //    }
1569   return;
1570 }
1571
1572
1573 //__________________________________________________________________________
1574 void AliMUONEventReconstructor::EventDumpTrigger(void)
1575 {
1576   // Dump reconstructed trigger event 
1577   // and the particle parameters
1578     
1579   AliMUONTriggerTrack *triggertrack ;
1580   Int_t nTriggerTracks = fMUONData->RecTriggerTracks()->GetEntriesFast();
1581  
1582   AliDebug(1, "****** enter EventDumpTrigger ******");
1583   AliDebug(1, Form("Number of Reconstructed tracks : %d ",  nTriggerTracks));
1584   
1585   // Loop over reconstructed tracks
1586   for (Int_t trackIndex = 0; trackIndex < nTriggerTracks; trackIndex++) {
1587     triggertrack = (AliMUONTriggerTrack*)fMUONData->RecTriggerTracks()->At(trackIndex);
1588       printf(" trigger track number %i x11=%f y11=%f thetax=%f thetay=%f \n",
1589              trackIndex,
1590              triggertrack->GetX11(),triggertrack->GetY11(),
1591              triggertrack->GetThetax(),triggertrack->GetThetay());      
1592   } 
1593 }
1594
1595 //__________________________________________________________________________
1596 void AliMUONEventReconstructor::MakeTrackCandidatesK(void)
1597 {
1598   // To make initial tracks for Kalman filter from the list of segments
1599   Int_t istat, iseg;
1600   AliMUONSegment *segment;
1601   AliMUONTrackK *trackK;
1602
1603   AliDebug(1,"Enter MakeTrackCandidatesK");
1604   // Reset the TClonesArray of reconstructed tracks
1605   if (fRecTracksPtr) fRecTracksPtr->Delete();
1606   // Delete in order that the Track destructors are called,
1607   // hence the space for the TClonesArray of pointers to TrackHit's is freed
1608   fNRecTracks = 0;
1609
1610   AliMUONTrackK a(this, fHitsForRecPtr); // bad idea ???
1611   // Loop over stations(1...) 5 and 4
1612   for (istat=4; istat>=3; istat--) {
1613     // Loop over segments in the station
1614     for (iseg=0; iseg<fNSegments[istat]; iseg++) {
1615       // Transform segments to tracks and evaluate covariance matrix
1616       segment = (AliMUONSegment*) ((*fSegmentsPtr[istat])[iseg]);
1617       trackK = new ((*fRecTracksPtr)[fNRecTracks]) AliMUONTrackK(segment);
1618       fNRecTracks++;
1619     } // for (iseg=0;...)
1620   } // for (istat=4;...)
1621   return;
1622 }
1623
1624 //__________________________________________________________________________
1625 void AliMUONEventReconstructor::FollowTracksK(void)
1626 {
1627   // Follow tracks using Kalman filter
1628   Bool_t ok;
1629   Int_t icand, ichamBeg = 0, ichamEnd, chamBits;
1630   Double_t zDipole1, zDipole2;
1631   AliMUONTrackK *trackK;
1632   AliMUONHitForRec *hit;
1633   AliMUONRawCluster *clus;
1634   TClonesArray *rawclusters;
1635   clus = 0; rawclusters = 0;
1636
1637   //AZ(z->-z) zDipole1 = GetSimpleBPosition() - GetSimpleBLength()/2;
1638   //AZ(z->-z) zDipole2 = zDipole1 + GetSimpleBLength();
1639   zDipole1 = GetSimpleBPosition() + GetSimpleBLength()/2;
1640   zDipole2 = zDipole1 - GetSimpleBLength();
1641
1642   // Print hits
1643   trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[0]);
1644   if (trackK->DebugLevel() > 0) {
1645     for (Int_t i1=0; i1<fNHitsForRec; i1++) {
1646       hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[i1]);
1647       //if (hit->GetTTRTrack() > 1 || hit->GetTrackRefSignal() == 0) continue;
1648       printf(" Hit # %d %10.4f %10.4f %10.4f",
1649              hit->GetChamberNumber(), hit->GetBendingCoor(),
1650              hit->GetNonBendingCoor(), hit->GetZ());
1651       if (fRecTrackRefHits) {
1652         // from track ref hits
1653         printf(" %3d %3d \n", hit->GetTrackRefSignal(), hit->GetTTRTrack());
1654       } else {
1655         // from raw clusters
1656         rawclusters = fMUONData->RawClusters(hit->GetChamberNumber());
1657         clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(hit->
1658                                                  GetHitNumber());
1659         printf("%3d", clus->GetTrack(1)-1);
1660         if (clus->GetTrack(2) != 0) printf("%3d \n", clus->GetTrack(2)-1);
1661         else printf("\n");
1662       } // if (fRecTrackRefHits)
1663     }
1664   } // if (trackK->DebugLevel() > 0)
1665
1666   icand = -1;
1667   Int_t nSeeds;
1668   nSeeds = fNRecTracks; // starting number of seeds
1669   // Loop over track candidates
1670   while (icand < fNRecTracks-1) {
1671     icand ++;
1672     if (trackK->DebugLevel()>0) cout << " *** Kalman track candidate No. " << icand << endl;
1673     trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]);
1674     if (trackK->GetRecover() < 0) continue; // failed track
1675
1676     // Discard candidate which will produce the double track
1677     /*
1678     if (icand > 0) {
1679       ok = CheckCandidateK(icand,nSeeds);
1680       if (!ok) {
1681         trackK->SetRecover(-1); // mark candidate to be removed
1682         continue;
1683       }
1684     }
1685     */
1686
1687     ok = kTRUE;
1688     if (trackK->GetRecover() == 0) hit = (AliMUONHitForRec*) 
1689                                    trackK->GetHitOnTrack()->Last(); // last hit
1690     //else hit = (AliMUONHitForRec*) (*trackK->GetHitOnTrack())[1]; // 2'nd hit
1691     else hit = trackK->GetHitLastOk(); // hit where track stopped
1692     if (hit) ichamBeg = hit->GetChamberNumber();
1693     ichamEnd = 0;
1694     // Check propagation direction
1695     if (hit == NULL) ichamBeg = ichamEnd;
1696     //AZ(z->-z) else if (trackK->GetTrackDir() > 0) {
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         //AZ trackK->SetTrackDir(-1);
1706         trackK->SetTrackDir(1);
1707         ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kTRUE,zDipole1,zDipole2);
1708         ichamBeg = ichamEnd;
1709         ichamEnd = 0;
1710       }
1711     } else {
1712       if (trackK->GetBPFlag()) {
1713         // backpropagation
1714         ichamEnd = 6; // backward propagation
1715         // Change weight matrix and zero fChi2 for backpropagation
1716         trackK->StartBack();
1717         ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kTRUE,zDipole1,zDipole2);
1718         ichamBeg = ichamEnd;
1719         ichamEnd = 0;
1720       }
1721     }
1722
1723     if (ok) {
1724       //AZ trackK->SetTrackDir(-1);
1725       trackK->SetTrackDir(1);
1726       trackK->SetBPFlag(kFALSE);
1727       ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kFALSE,zDipole1,zDipole2);
1728     }
1729     if (!ok) { trackK->SetRecover(-1); continue; } // mark candidate to be removed
1730
1731     // Apply smoother
1732     if (trackK->GetRecover() >= 0) {
1733       ok = trackK->Smooth();
1734       if (!ok) trackK->SetRecover(-1); // mark candidate to be removed
1735     }
1736
1737     // Majority 3 of 4 in first 2 stations
1738     if (!ok) continue;
1739     chamBits = 0;
1740     Double_t chi2_max = 0;
1741     for (Int_t i=0; i<trackK->GetNTrackHits(); i++) {
1742       hit = (AliMUONHitForRec*) (*trackK->GetHitOnTrack())[i];
1743       chamBits |= BIT(hit->GetChamberNumber());
1744       if (trackK->GetChi2PerPoint(i) > chi2_max) chi2_max = trackK->GetChi2PerPoint(i);
1745     }
1746     if (!((chamBits&3)==3 || (chamBits>>2&3)==3) && chi2_max > 25) {
1747       //trackK->Recover();
1748       trackK->SetRecover(-1); //mark candidate to be removed
1749       continue;
1750     }
1751     if (ok) trackK->SetTrackQuality(0); // compute "track quality"
1752   } // while
1753
1754   for (Int_t i=0; i<fNRecTracks; i++) {
1755     trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[i]);
1756     if (trackK->GetRecover() < 0) fRecTracksPtr->RemoveAt(i);
1757   }
1758
1759   // Compress TClonesArray
1760   fRecTracksPtr->Compress();
1761   fNRecTracks = fRecTracksPtr->GetEntriesFast();
1762   return;
1763 }
1764
1765 //__________________________________________________________________________
1766 Bool_t AliMUONEventReconstructor::CheckCandidateK(Int_t icand, Int_t nSeeds) const
1767 {
1768   // Discards track candidate if it will produce the double track (having
1769   // the same seed segment hits as hits of a good track found before)
1770   AliMUONTrackK *track1, *track2;
1771   AliMUONHitForRec *hit1, *hit2, *hit;
1772
1773   track1 = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]);
1774   hit1 = (AliMUONHitForRec*) (*track1->GetHitOnTrack())[0]; // 1'st hit
1775   hit2 = (AliMUONHitForRec*) (*track1->GetHitOnTrack())[1]; // 2'nd hit
1776
1777   for (Int_t i=0; i<icand; i++) {
1778     track2 = (AliMUONTrackK*) ((*fRecTracksPtr)[i]);
1779     //if (track2->GetRecover() < 0) continue;
1780     if (track2->GetRecover() < 0 && icand >= nSeeds) continue;
1781
1782     if (track1->GetStartSegment() == track2->GetStartSegment()) {
1783       return kFALSE;
1784     } else {
1785       Int_t nSame = 0;
1786       for (Int_t j=0; j<track2->GetNTrackHits(); j++) {
1787         hit = (AliMUONHitForRec*) (*track2->GetHitOnTrack())[j];
1788         if (hit == hit1 || hit == hit2) {
1789           nSame++;
1790           if (nSame == 2) return kFALSE;
1791         }
1792       } // for (Int_t j=0;
1793     }
1794   } // for (Int_t i=0;
1795   return kTRUE;
1796 }
1797
1798 //__________________________________________________________________________
1799 void AliMUONEventReconstructor::RemoveDoubleTracksK(void)
1800 {
1801   // Removes double tracks (sharing more than half of their hits). Keeps
1802   // the track with higher quality
1803   AliMUONTrackK *track1, *track2, *trackToKill;
1804
1805   // Sort tracks according to their quality
1806   fRecTracksPtr->Sort();
1807
1808   // Loop over first track of the pair
1809   track1 = (AliMUONTrackK*) fRecTracksPtr->First();
1810   Int_t debug = track1->DebugLevel();
1811   while (track1) {
1812     // Loop over second track of the pair
1813     track2 = (AliMUONTrackK*) fRecTracksPtr->After(track1);
1814     while (track2) {
1815       // Check whether or not to keep track2
1816       if (!track2->KeepTrack(track1)) {
1817         if (debug >= 0) cout << " Killed track: " << 1/(*track2->GetTrackParameters())(4,0) <<
1818           " " << track2->GetTrackQuality() << endl;
1819         trackToKill = track2;
1820         track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2);
1821         trackToKill->Kill();
1822         fRecTracksPtr->Compress();
1823       } else track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2);
1824     } // track2
1825     track1 = (AliMUONTrackK*) fRecTracksPtr->After(track1);
1826   } // track1
1827
1828   fNRecTracks = fRecTracksPtr->GetEntriesFast();
1829   if (debug >= 0) cout << " Number of Kalman tracks: " << fNRecTracks << endl;
1830 }
1831
1832 //__________________________________________________________________________
1833 void AliMUONEventReconstructor::GoToVertex(void)
1834 {
1835   // Propagates track to the vertex thru absorber
1836   // (using Branson correction for now)
1837
1838   Double_t zVertex;
1839   zVertex = 0;
1840   for (Int_t i=0; i<fNRecTracks; i++) {
1841     //((AliMUONTrackK*)(*fRecTracksPtr)[i])->Branson();
1842     ((AliMUONTrackK*)(*fRecTracksPtr)[i])->SetTrackQuality(1); // compute Chi2
1843     //((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToZ(zVertex); // w/out absorber
1844     ((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToVertex(1); // with absorber
1845   }
1846 }
1847
1848 //__________________________________________________________________________
1849 void AliMUONEventReconstructor::SetTrackMethod(Int_t iTrackMethod)
1850 {
1851   // Set track method and recreate track container if necessary
1852   
1853   fTrackMethod = iTrackMethod == 2 ? 2 : 1;
1854   if (fTrackMethod == 2) {
1855     if (fRecTracksPtr) delete fRecTracksPtr;
1856     fRecTracksPtr = new TClonesArray("AliMUONTrackK", 10);
1857     cout << " *** Tracking with the Kalman filter *** " << endl;
1858   } else cout << " *** Traditional tracking *** " << endl;
1859
1860 }