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