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