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