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