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