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