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