]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONEventReconstructor.cxx
AliMUONEventReconstructor package:
[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 /*
17 $Log$
18 Revision 1.7  2000/07/03 12:28:06  gosset
19 Printout at the right place after extrapolation to vertex
20
21 Revision 1.6  2000/06/30 12:01:06  gosset
22 Correction for hit search in the right chamber (JPC)
23
24 Revision 1.5  2000/06/30 10:15:48  gosset
25 Changes to EventReconstructor...:
26 precision fit with multiple Coulomb scattering;
27 extrapolation to vertex with Branson correction in absorber (JPC)
28
29 Revision 1.4  2000/06/27 14:11:36  gosset
30 Corrections against violations of coding conventions
31
32 Revision 1.3  2000/06/16 07:27:08  gosset
33 To remove problem in running RuleChecker, like in MUON-dev
34
35 Revision 1.1.2.5  2000/06/16 07:00:26  gosset
36 To remove problem in running RuleChecker
37
38 Revision 1.1.2.4  2000/06/12 08:00:07  morsch
39 Dummy streamer to solve CINT compilation problem (to be investigated !)
40
41 Revision 1.1.2.3  2000/06/09 20:59:57  morsch
42 Make includes consistent with new file structure.
43
44 Revision 1.1.2.2  2000/06/09 12:58:05  gosset
45 Removed comment beginnings in Log sections of .cxx files
46 Suppressed most violations of coding rules
47
48 Revision 1.1.2.1  2000/06/07 14:44:53  gosset
49 Addition of files for track reconstruction in C++
50 */
51
52 //__________________________________________________________________________
53 //
54 // MUON event reconstructor in ALICE
55 //
56 // This class contains as data:
57 // * the parameters for the event reconstruction
58 // * a pointer to the array of hits to be reconstructed (the event)
59 // * a pointer to the array of segments made with these hits inside each station
60 // * a pointer to the array of reconstructed tracks
61 //
62 // It contains as methods, among others:
63 // * MakeEventToBeReconstructed to build the array of hits to be reconstructed
64 // * MakeSegments to build the segments
65 // * MakeTracks to build the tracks
66 //__________________________________________________________________________
67
68 #include <iostream.h>
69
70 #include <TRandom.h>
71 #include <TFile.h>
72
73 #include "AliCallf77.h" 
74 #include "AliMUONEventReconstructor.h"
75 #include "AliMUON.h"
76 #include "AliMUONHitForRec.h"
77 #include "AliMUONSegment.h"
78 #include "AliMUONHit.h"
79 #include "AliMUONRawCluster.h"
80 #include "AliMUONTrack.h"
81 #include "AliMUONChamber.h"
82 #include "AliMUONTrackHit.h"
83 #include "AliRun.h"
84 #include "TParticle.h"
85
86 #ifndef WIN32 
87 # define initfield initfield_
88 # define reco_gufld reco_gufld_
89 #else 
90 # define initfield INITFIELD
91 # define reco_gufld RECO_GUFLD
92 #endif 
93
94 extern "C"
95 {
96 void type_of_call initfield();
97 void type_of_call reco_gufld(Double_t *Coor, Double_t *Field);
98 }
99
100 //************* Defaults parameters for reconstruction
101 static const Double_t kDefaultMinBendingMomentum = 3.0;
102 static const Double_t kDefaultMaxSigma2Distance = 16.0;
103 static const Double_t kDefaultBendingResolution = 0.01;
104 static const Double_t kDefaultNonBendingResolution = 0.144;
105 static const Double_t kDefaultChamberThicknessInX0 = 0.03;
106 // Simple magnetic field:
107 // Value taken from macro MUONtracking.C: 0.7 T, hence 7 kG
108 // Length and Position from reco_muon.F, with opposite sign:
109 // Length = ZMAGEND-ZCOIL
110 // Position = (ZMAGEND+ZCOIL)/2
111 // to be ajusted differently from real magnetic field ????
112 static const Double_t kDefaultSimpleBValue = 7.0;
113 static const Double_t kDefaultSimpleBLength = 428.0;
114 static const Double_t kDefaultSimpleBPosition = 1019.0;
115 static const Int_t kDefaultRecGeantHits = 0;
116 static const Double_t kDefaultEfficiency = 0.95;
117
118 static const Int_t kDefaultPrintLevel = 0;
119
120 ClassImp(AliMUONEventReconstructor) // Class implementation in ROOT context
121
122   //__________________________________________________________________________
123 AliMUONEventReconstructor::AliMUONEventReconstructor(void)
124 {
125   // Constructor for class AliMUONEventReconstructor
126   SetReconstructionParametersToDefaults();
127   // Memory allocation for the TClonesArray of hits for reconstruction
128   // Is 10000 the right size ????
129   fHitsForRecPtr = new TClonesArray("AliMUONHitForRec", 10000);
130   fNHitsForRec = 0; // really needed or GetEntriesFast sufficient ????
131   // Memory allocation for the TClonesArray's of segments in stations
132   // Is 2000 the right size ????
133   for (Int_t st = 0; st < kMaxMuonTrackingStations; st++) {
134     fSegmentsPtr[st] = new TClonesArray("AliMUONSegment", 2000);
135     fNSegments[st] = 0; // really needed or GetEntriesFast sufficient ????
136   }
137   // Memory allocation for the TClonesArray of reconstructed tracks
138   // Is 10 the right size ????
139   fRecTracksPtr = new TClonesArray("AliMUONTrack", 10);
140   fNRecTracks = 0; // really needed or GetEntriesFast sufficient ????
141
142   // Initialize magnetic field
143   // using Fortran subroutine INITFIELD in "reco_muon.F".
144   // Should rather use AliMagF ???? and remove prototyping ...
145   initfield();
146   // Impression de quelques valeurs
147   Double_t coor[3], field[3];
148   coor[0] = 50.0;
149   coor[1] = 50.0;
150   coor[2] = 950.0;
151   reco_gufld(coor, field);
152   cout << "coor: " << coor[0] << ", " << coor[1] << ", " << coor[2] << endl;
153   cout << "field: " << field[0] << ", " << field[1] << ", " << field[2] << endl;
154   coor[2] = -950.0;
155   reco_gufld(coor, field);
156   cout << "coor: " << coor[0] << ", " << coor[1] << ", " << coor[2] << endl;
157   cout << "field: " << field[0] << ", " << field[1] << ", " << field[2] << endl;
158   coor[2] = -950.0;
159
160   if (fPrintLevel >= 0) {
161     cout << "AliMUONEventReconstructor constructed with defaults" << endl; Dump();}
162   return;
163 }
164
165 AliMUONEventReconstructor::AliMUONEventReconstructor (const AliMUONEventReconstructor& Reconstructor)
166 {
167   // Dummy copy constructor
168 }
169
170 AliMUONEventReconstructor & AliMUONEventReconstructor::operator=(const AliMUONEventReconstructor& Reconstructor)
171 {
172   // Dummy assignment operator
173     return *this;
174 }
175
176   //__________________________________________________________________________
177 AliMUONEventReconstructor::~AliMUONEventReconstructor(void)
178 {
179   // Destructor for class AliMUONEventReconstructor
180   delete fHitsForRecPtr; // Correct destruction of everything ???? or delete [] ????
181   for (Int_t st = 0; st < kMaxMuonTrackingStations; st++)
182     delete fSegmentsPtr[st]; // Correct destruction of everything ????
183   return;
184 }
185
186   //__________________________________________________________________________
187 void AliMUONEventReconstructor::SetReconstructionParametersToDefaults(void)
188 {
189   // Set reconstruction parameters to default values
190   // Would be much more convenient with a structure (or class) ????
191   fMinBendingMomentum = kDefaultMinBendingMomentum;
192   fMaxSigma2Distance = kDefaultMaxSigma2Distance;
193
194   AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON");
195   // ******** Parameters for making HitsForRec
196   // minimum radius,
197   // like in TRACKF_STAT:
198   // 2 degrees for stations 1 and 2, or ch(0...) from 0 to 3;
199   // 30 cm for stations 3 to 5, or ch(0...) from 4 to 9
200   for (Int_t ch = 0; ch < kMaxMuonTrackingChambers; ch++) {
201     if (ch < 4) fRMin[ch] = TMath::Abs((&(pMUON->Chamber(ch)))->Z()) *
202                   2.0 * TMath::Pi() / 180.0;
203     else fRMin[ch] = 30.0;
204   }
205   // maximum radius
206   // like in TRACKF_STAT (10 degrees ????)
207   fRMax[0] = fRMax[1] = 91.5;
208   fRMax[2] = fRMax[3] = 122.5;
209   fRMax[4] = fRMax[5] = 158.3;
210   fRMax[6] = fRMax[7] = 260.0;
211   fRMax[8] = fRMax[9] = 260.0;
212
213   // ******** Parameters for making segments
214   // should be parametrized ????
215   // according to interval between chambers in a station ????
216   // Maximum distance in non bending plane
217   // 5 * 0.22 just to remember the way it was made in TRACKF_STAT
218   // SIGCUT*DYMAX(IZ)
219   for (Int_t st = 0; st < kMaxMuonTrackingStations; st++)
220     fSegmentMaxDistNonBending[st] = 5. * 0.22;
221   // Maximum distance in bending plane
222   // values from TRACKF_STAT corresponding to (J psi 20cm)
223   fSegmentMaxDistBending[0] = 1.5;
224   fSegmentMaxDistBending[1] = 1.5;
225   fSegmentMaxDistBending[2] = 3.0;
226   fSegmentMaxDistBending[3] = 6.0;
227   fSegmentMaxDistBending[4] = 6.0;
228   
229   fBendingResolution = kDefaultBendingResolution;
230   fNonBendingResolution = kDefaultNonBendingResolution;
231   fChamberThicknessInX0 = kDefaultChamberThicknessInX0;
232   fSimpleBValue = kDefaultSimpleBValue;
233   fSimpleBLength = kDefaultSimpleBLength;
234   fSimpleBPosition = kDefaultSimpleBPosition;
235   fRecGeantHits = kDefaultRecGeantHits;
236   fEfficiency = kDefaultEfficiency;
237   fPrintLevel = kDefaultPrintLevel;
238   return;
239 }
240
241 //__________________________________________________________________________
242 Double_t AliMUONEventReconstructor::GetImpactParamFromBendingMomentum(Double_t BendingMomentum)
243 {
244   // Returns impact parameter at vertex in bending plane (cm),
245   // from the signed bending momentum "BendingMomentum" in bending plane (GeV/c),
246   // using simple values for dipole magnetic field.
247   // The sign is the sign of the charge.
248   return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
249           BendingMomentum);
250 }
251
252 //__________________________________________________________________________
253 Double_t AliMUONEventReconstructor::GetBendingMomentumFromImpactParam(Double_t ImpactParam)
254 {
255   // Returns signed bending momentum in bending plane (GeV/c),
256   // from the impact parameter "ImpactParam" at vertex in bending plane (cm),
257   // using simple values for dipole magnetic field.
258   // The sign is the sign of the charge.
259   return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
260           ImpactParam);
261 }
262
263 //__________________________________________________________________________
264 void AliMUONEventReconstructor::SetBkgGeantFile(Text_t *BkgGeantFileName)
265 {
266   // Set background file ... for GEANT hits
267   // Must be called after having loaded the firts signal event
268   if (fPrintLevel >= 0) {
269     cout << "Enter SetBkgGeantFile with BkgGeantFileName ``"
270          << BkgGeantFileName << "''" << endl;}
271   if (strlen(BkgGeantFileName)) {
272     // BkgGeantFileName not empty: try to open the file
273     if (fPrintLevel >= 2) {cout << "Before File(Bkg)" << endl; gDirectory->Dump();}
274     fBkgGeantFile = new TFile(BkgGeantFileName);
275     if (fPrintLevel >= 2) {cout << "After File(Bkg)" << endl; gDirectory->Dump();}
276     if (fBkgGeantFile-> IsOpen()) {
277       if (fPrintLevel >= 0) {
278         cout << "Background for GEANT hits in file: ``" << BkgGeantFileName
279              << "'' successfully opened" << endl;}
280     }
281     else {
282       cout << "Background for GEANT hits in file: " << BkgGeantFileName << endl;
283       cout << "NOT FOUND: EXIT" << endl;
284       exit(0); // right instruction for exit ????
285     }
286     // Arrays for "particles" and "hits"
287     fBkgGeantParticles = new TClonesArray("TParticle", 200);
288     fBkgGeantHits = new TClonesArray("AliMUONHit", 2000);
289     // Event number to -1 for initialization
290     fBkgGeantEventNumber = -1;
291     // Back to the signal file:
292     // first signal event must have been loaded previously,
293     // otherwise, Segmentation violation at the next instruction
294     // How is it possible to do smething better ????
295     ((gAlice->TreeK())->GetCurrentFile())->cd();
296     if (fPrintLevel >= 2) {cout << "After cd(gAlice)" << endl; gDirectory->Dump();}
297   }
298   return;
299 }
300
301 //__________________________________________________________________________
302 void AliMUONEventReconstructor::NextBkgGeantEvent(void)
303 {
304   // Get next event in background file for GEANT hits
305   // Goes back to event number 0 when end of file is reached
306   char treeName[20];
307   TBranch *branch;
308   if (fPrintLevel >= 0) {
309     cout << "Enter NextBkgGeantEvent" << endl;}
310   // Clean previous event
311   if(fBkgGeantTK) delete fBkgGeantTK;
312   fBkgGeantTK = NULL;
313   if(fBkgGeantParticles) fBkgGeantParticles->Clear();
314   if(fBkgGeantTH) delete fBkgGeantTH;
315   fBkgGeantTH = NULL;
316   if(fBkgGeantHits) fBkgGeantHits->Clear();
317   // Increment event number
318   fBkgGeantEventNumber++;
319   // Get access to Particles and Hits for event from background file
320   if (fPrintLevel >= 2) {cout << "Before cd(Bkg)" << endl; gDirectory->Dump();}
321   fBkgGeantFile->cd();
322   if (fPrintLevel >= 2) {cout << "After cd(Bkg)" << endl; gDirectory->Dump();}
323   // Particles: TreeK for event and branch "Particles"
324   sprintf(treeName, "TreeK%d", fBkgGeantEventNumber);
325   fBkgGeantTK = (TTree*)gDirectory->Get(treeName);
326   if (!fBkgGeantTK) {
327     if (fPrintLevel >= 0) {
328       cout << "Cannot find Kine Tree for background event: " <<
329         fBkgGeantEventNumber << endl;
330       cout << "Goes back to event 0" << endl;
331     }
332     fBkgGeantEventNumber = 0;
333     sprintf(treeName, "TreeK%d", fBkgGeantEventNumber);
334     fBkgGeantTK = (TTree*)gDirectory->Get(treeName);
335     if (!fBkgGeantTK) {
336       cout << "ERROR: cannot find Kine Tree for background event: " <<
337         fBkgGeantEventNumber << endl;
338       exit(0);
339     }
340   }
341   if (fBkgGeantTK) 
342     fBkgGeantTK->SetBranchAddress("Particles", &fBkgGeantParticles);
343   fBkgGeantTK->GetEvent(0); // why event 0 ???? necessary ????
344   // Hits: TreeH for event and branch "MUON"
345   sprintf(treeName, "TreeH%d", fBkgGeantEventNumber);
346   fBkgGeantTH = (TTree*)gDirectory->Get(treeName);
347   if (!fBkgGeantTH) {
348     cout << "ERROR: cannot find Hits Tree for background event: " <<
349       fBkgGeantEventNumber << endl;
350       exit(0);
351   }
352   if (fBkgGeantTH && fBkgGeantHits) {
353     branch = fBkgGeantTH->GetBranch("MUON");
354     if (branch) branch->SetAddress(&fBkgGeantHits);
355   }
356   fBkgGeantTH->GetEntries(); // necessary ????
357   // Back to the signal file
358   ((gAlice->TreeK())->GetCurrentFile())->cd();
359   if (fPrintLevel >= 2) {cout << "After cd(gAlice)" << endl; gDirectory->Dump();}
360   return;
361 }
362
363 //__________________________________________________________________________
364 void AliMUONEventReconstructor::EventReconstruct(void)
365 {
366   // To reconstruct one event
367   if (fPrintLevel >= 1) cout << "enter EventReconstruct" << endl;
368   MakeEventToBeReconstructed();
369   MakeSegments();
370   MakeTracks();
371   return;
372 }
373
374   //__________________________________________________________________________
375 void AliMUONEventReconstructor::ResetHitsForRec(void)
376 {
377   // To reset the array and the number of HitsForRec,
378   // and also the number of HitsForRec
379   // and the index of the first HitForRec per chamber
380   if (fHitsForRecPtr) fHitsForRecPtr->Clear();
381   fNHitsForRec = 0;
382   for (Int_t ch = 0; ch < kMaxMuonTrackingChambers; ch++)
383     fNHitsForRecPerChamber[ch] = fIndexOfFirstHitForRecPerChamber[ch] = 0;
384   return;
385 }
386
387   //__________________________________________________________________________
388 void AliMUONEventReconstructor::ResetSegments(void)
389 {
390   // To reset the TClonesArray of segments and the number of Segments
391   // for all stations
392   for (Int_t st = 0; st < kMaxMuonTrackingStations; st++) {
393     if (fSegmentsPtr[st]) fSegmentsPtr[st]->Clear();
394     fNSegments[st] = 0;
395   }
396   return;
397 }
398
399   //__________________________________________________________________________
400 void AliMUONEventReconstructor::ResetTracks(void)
401 {
402   // To reset the TClonesArray of reconstructed tracks
403   if (fRecTracksPtr) fRecTracksPtr->Clear();
404   fNRecTracks = 0;
405   return;
406 }
407
408   //__________________________________________________________________________
409 void AliMUONEventReconstructor::MakeEventToBeReconstructed(void)
410 {
411   // To make the list of hits to be reconstructed,
412   // either from the GEANT hits or from the raw clusters
413   // according to the parameter set for the reconstructor
414   if (fPrintLevel >= 1) cout << "enter MakeEventToBeReconstructed" << endl;
415   ResetHitsForRec();
416   if (fRecGeantHits == 1) {
417     // Reconstruction from GEANT hits
418     // Back to the signal file
419     ((gAlice->TreeK())->GetCurrentFile())->cd();
420     // Signal hits
421     // AliMUON *MUON  = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
422     // Security on MUON ????
423     AddHitsForRecFromGEANT(gAlice->TreeH());
424     // Background hits
425     AddHitsForRecFromBkgGEANT(fBkgGeantTH, fBkgGeantHits);
426     // Sort HitsForRec in increasing order with respect to chamber number
427     SortHitsForRecWithIncreasingChamber();
428   }
429   else {
430     // Reconstruction from raw clusters
431     // AliMUON *MUON  = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
432     // Security on MUON ????
433     // TreeR assumed to be be "prepared" in calling function
434     // by "MUON->GetTreeR(nev)" ????
435     TTree *treeR = gAlice->TreeR();
436     AddHitsForRecFromRawClusters(treeR);
437     // No sorting: it is done automatically in the previous function
438   }
439   if (fPrintLevel >= 10) {
440     cout << "end of MakeEventToBeReconstructed" << endl;
441     cout << "NHitsForRec: " << fNHitsForRec << endl;
442     for (Int_t ch = 0; ch < kMaxMuonTrackingChambers; ch++) {
443       cout << "chamber(0...): " << ch
444            << "  NHitsForRec: " << fNHitsForRecPerChamber[ch]
445            << "  index(first HitForRec): " << fIndexOfFirstHitForRecPerChamber[ch]
446            << endl;
447       for (Int_t hit = fIndexOfFirstHitForRecPerChamber[ch];
448            hit < fIndexOfFirstHitForRecPerChamber[ch] + fNHitsForRecPerChamber[ch];
449            hit++) {
450         cout << "HitForRec index(0...): " << hit << endl;
451         ((*fHitsForRecPtr)[hit])->Dump();
452       }
453     }
454   }
455   return;
456 }
457
458   //__________________________________________________________________________
459 void AliMUONEventReconstructor::AddHitsForRecFromGEANT(TTree *TH)
460 {
461   // To add to the list of hits for reconstruction
462   // the GEANT signal hits from a hit tree TH.
463   if (fPrintLevel >= 2)
464     cout << "enter AddHitsForRecFromGEANT with TH: " << TH << endl;
465   if (TH == NULL) return;
466   AliMUON *pMUON  = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
467   // Security on MUON ????
468   // See whether it could be the same for signal and background ????
469   // Loop over tracks in tree
470   Int_t ntracks = (Int_t) TH->GetEntries();
471   if (fPrintLevel >= 2)
472     cout << "ntracks: " << ntracks << endl;
473   for (Int_t track = 0; track < ntracks; track++) {
474     gAlice->ResetHits();
475     TH->GetEvent(track);
476     // Loop over hits
477     Int_t hit = 0;
478     for (AliMUONHit* mHit = (AliMUONHit*) pMUON->FirstHit(-1); 
479          mHit;
480          mHit = (AliMUONHit*) pMUON->NextHit(), hit++) {
481       NewHitForRecFromGEANT(mHit,track, hit, 1);
482     } // end of hit loop
483   } // end of track loop
484   return;
485 }
486
487   //__________________________________________________________________________
488 void AliMUONEventReconstructor::AddHitsForRecFromBkgGEANT(TTree *TH, TClonesArray *Hits)
489 {
490   // To add to the list of hits for reconstruction
491   // the GEANT background hits from a hit tree TH and a pointer Hits to a hit list.
492   // How to have only one function "AddHitsForRecFromGEANT" ????
493   if (fPrintLevel >= 2)
494     cout << "enter AddHitsForRecFromBkgGEANT with TH: " << TH << endl;
495   if (TH == NULL) return;
496   // Loop over tracks in tree
497   Int_t ntracks = (Int_t) TH->GetEntries();
498   if (fPrintLevel >= 2)
499     cout << "ntracks: " << ntracks << endl;
500   for (Int_t track = 0; track < ntracks; track++) {
501     if (Hits) Hits->Clear();
502     TH->GetEvent(track);
503     // Loop over hits
504     for (Int_t hit = 0; hit < Hits->GetEntriesFast(); hit++) {
505       NewHitForRecFromGEANT((AliMUONHit*) (*Hits)[hit], track, hit, 0);
506     } // end of hit loop
507   } // end of track loop
508   return;
509 }
510
511   //__________________________________________________________________________
512 AliMUONHitForRec* AliMUONEventReconstructor::NewHitForRecFromGEANT(AliMUONHit* Hit, Int_t TrackNumber, Int_t HitNumber, Int_t Signal)
513 {
514   // To make a new hit for reconstruction from a GEANT hit pointed to by "Hit",
515   // with hit number "HitNumber" in the track numbered "TrackNumber",
516   // either from signal ("Signal" = 1) or background ("Signal" = 0) event.
517   // Selects hits in tracking (not trigger) chambers.
518   // Takes into account the efficiency (fEfficiency)
519   // and the smearing from resolution (fBendingResolution and fNonBendingResolution).
520   // Adds a condition on the radius between RMin and RMax
521   // to better simulate the real chambers.
522   // Returns the pointer to the new hit for reconstruction,
523   // or NULL in case of inefficiency or non tracking chamber or bad radius.
524   // No condition on at most 20 cm from a muon from a resonance
525   // like in Fortran TRACKF_STAT.
526   AliMUONHitForRec* hitForRec;
527   Double_t bendCoor, nonBendCoor, radius;
528   Int_t chamber = Hit->fChamber - 1; // chamber(0...)
529   // only in tracking chambers (fChamber starts at 1)
530   if (chamber >= kMaxMuonTrackingChambers) return NULL;
531   // only if hit is efficient (keep track for checking ????)
532   if (gRandom->Rndm() > fEfficiency) return NULL;
533   // only if radius between RMin and RMax
534   bendCoor = Hit->fY;
535   nonBendCoor = Hit->fX;
536   radius = TMath::Sqrt((bendCoor * bendCoor) + (nonBendCoor * nonBendCoor));
537   if ((radius < fRMin[chamber]) || (radius > fRMax[chamber])) return NULL;
538   // new AliMUONHitForRec from GEANT hit and increment number of AliMUONHitForRec's
539   hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(Hit);
540   fNHitsForRec++;
541   // add smearing from resolution
542   hitForRec->SetBendingCoor(bendCoor + gRandom->Gaus(0., fBendingResolution));
543   hitForRec->SetNonBendingCoor(nonBendCoor
544                                + gRandom->Gaus(0., fNonBendingResolution));
545   // more information into HitForRec
546   //  resolution: angular effect to be added here ????
547   hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
548   hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
549   //  GEANT track info
550   hitForRec->SetHitNumber(HitNumber);
551   hitForRec->SetTHTrack(TrackNumber);
552   hitForRec->SetGeantSignal(Signal);
553   if (fPrintLevel >= 10) {
554     cout << "track: " << TrackNumber << " hit: " << HitNumber << endl;
555     Hit->Dump();
556     cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
557     hitForRec->Dump();}
558   return hitForRec;
559 }
560
561   //__________________________________________________________________________
562 void AliMUONEventReconstructor::SortHitsForRecWithIncreasingChamber()
563 {
564   // Sort HitsForRec's in increasing order with respect to chamber number.
565   // Uses the function "Compare".
566   // Update the information for HitsForRec per chamber too.
567   Int_t ch, nhits, prevch;
568   fHitsForRecPtr->Sort();
569   for (ch = 0; ch < kMaxMuonTrackingChambers; ch++) {
570     fNHitsForRecPerChamber[ch] = 0;
571     fIndexOfFirstHitForRecPerChamber[ch] = 0;
572   }
573   prevch = 0; // previous chamber
574   nhits = 0; // number of hits in current chamber
575   // Loop over HitsForRec
576   for (Int_t hit = 0; hit < fNHitsForRec; hit++) {
577     // chamber number (0...)
578     ch = ((AliMUONHitForRec*)  ((*fHitsForRecPtr)[hit]))->GetChamberNumber();
579     // increment number of hits in current chamber
580     (fNHitsForRecPerChamber[ch])++;
581     // update index of first HitForRec in current chamber
582     // if chamber number different from previous one
583     if (ch != prevch) {
584       fIndexOfFirstHitForRecPerChamber[ch] = hit;
585       prevch = ch;
586     }
587   }
588   return;
589 }
590
591 //   //__________________________________________________________________________
592 // void AliMUONEventReconstructor::AddHitsForRecFromCathodeCorrelations(TTree* TC)
593 // {
594 //   // OLD VERSION WHEN ONE ONE WAS USING SO CALLED CATHODE CORRELATIONS
595 //   // To add to the list of hits for reconstruction
596 //   // the (cathode correlated) raw clusters
597 //   // No condition added, like in Fortran TRACKF_STAT,
598 //   // on the radius between RMin and RMax.
599 //   AliMUONHitForRec *hitForRec;
600 //   if (fPrintLevel >= 1) cout << "enter AddHitsForRecFromCathodeCorrelations" << endl;
601 //   AliMUON *MUON  = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
602 //   // Security on MUON ????
603 //   // Loop over tracking chambers
604 //   for (Int_t ch = 0; ch < kMaxMuonTrackingChambers; ch++) {
605 //     // number of HitsForRec to 0 for the chamber
606 //     fNHitsForRecPerChamber[ch] = 0;
607 //     // index of first HitForRec for the chamber
608 //     if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0;
609 //     else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec;
610 //     TClonesArray *reconst_hits  = MUON->ReconstHitsAddress(ch);
611 //     MUON->ResetReconstHits();
612 //     TC->GetEvent();
613 //     Int_t ncor = (Int_t)reconst_hits->GetEntries();
614 //     // Loop over (cathode correlated) raw clusters
615 //     for (Int_t cor = 0; cor < ncor; cor++) {
616 //       AliMUONReconstHit * mCor = 
617 //      (AliMUONReconstHit*) reconst_hits->UncheckedAt(cor);
618 //       // new AliMUONHitForRec from (cathode correlated) raw cluster
619 //       // and increment number of AliMUONHitForRec's (total and in chamber)
620 //       hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(mCor);
621 //       fNHitsForRec++;
622 //       (fNHitsForRecPerChamber[ch])++;
623 //       // more information into HitForRec
624 //       hitForRec->SetChamberNumber(ch);
625 //       hitForRec->SetHitNumber(cor);
626 //       // Z coordinate of the chamber (cm) with sign opposite to GEANT sign
627 //       // could (should) be more exact from chamber geometry ???? 
628 //       hitForRec->SetZ(-(&(MUON->Chamber(ch)))->Z());
629 //       if (fPrintLevel >= 10) {
630 //      cout << "chamber (0...): " << ch <<
631 //        " cathcorrel (0...): " << cor << endl;
632 //      mCor->Dump();
633 //      cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
634 //      hitForRec->Dump();}
635 //     } // end of cluster loop
636 //   } // end of chamber loop
637 //   return;
638 // }
639
640   //__________________________________________________________________________
641 void AliMUONEventReconstructor::AddHitsForRecFromRawClusters(TTree* TR)
642 {
643   // To add to the list of hits for reconstruction all the raw clusters
644   // No condition added, like in Fortran TRACKF_STAT,
645   // on the radius between RMin and RMax.
646   AliMUONHitForRec *hitForRec;
647   AliMUONRawCluster *clus;
648   Int_t iclus, nclus;
649   TClonesArray *rawclusters;
650   if (fPrintLevel >= 1) cout << "enter AddHitsForRecFromRawClusters" << endl;
651   AliMUON *pMUON  = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
652   // Security on MUON ????
653   // Loop over tracking chambers
654   for (Int_t ch = 0; ch < kMaxMuonTrackingChambers; ch++) {
655     // number of HitsForRec to 0 for the chamber
656     fNHitsForRecPerChamber[ch] = 0;
657     // index of first HitForRec for the chamber
658     if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0;
659     else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec;
660     rawclusters = pMUON->RawClustAddress(ch);
661     pMUON->ResetRawClusters();
662     TR->GetEvent((Int_t) (TR->GetEntries()) - 1); // to be checked ????
663     nclus = (Int_t) (rawclusters->GetEntries());
664     // Loop over (cathode correlated) raw clusters
665     for (iclus = 0; iclus < nclus; iclus++) {
666       clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(iclus);
667       // new AliMUONHitForRec from raw cluster
668       // and increment number of AliMUONHitForRec's (total and in chamber)
669       hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(clus);
670       fNHitsForRec++;
671       (fNHitsForRecPerChamber[ch])++;
672       // more information into HitForRec
673       //  resolution: info should be already in raw cluster and taken from it ????
674       hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
675       hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
676       //  original raw cluster
677       hitForRec->SetChamberNumber(ch);
678       hitForRec->SetHitNumber(iclus);
679       // Z coordinate of the chamber (cm)
680       // could (should) be more exact from chamber geometry ???? 
681       hitForRec->SetZ((&(pMUON->Chamber(ch)))->Z());
682       if (fPrintLevel >= 10) {
683         cout << "chamber (0...): " << ch <<
684           " raw cluster (0...): " << iclus << endl;
685         clus->Dump();
686         cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
687         hitForRec->Dump();}
688     } // end of cluster loop
689   } // end of chamber loop
690   return;
691 }
692
693   //__________________________________________________________________________
694 void AliMUONEventReconstructor::MakeSegments(void)
695 {
696   // To make the list of segments in all stations,
697   // from the list of hits to be reconstructed
698   if (fPrintLevel >= 1) cout << "enter MakeSegments" << endl;
699   ResetSegments();
700   // Loop over stations
701   for (Int_t st = 0; st < kMaxMuonTrackingStations; st++)
702     MakeSegmentsPerStation(st); 
703   if (fPrintLevel >= 10) {
704     cout << "end of MakeSegments" << endl;
705     for (Int_t st = 0; st < kMaxMuonTrackingStations; st++) {
706       cout << "station(0...): " << st
707            << "  Segments: " << fNSegments[st]
708            << endl;
709       for (Int_t seg = 0;
710            seg < fNSegments[st];
711            seg++) {
712         cout << "Segment index(0...): " << seg << endl;
713         ((*fSegmentsPtr[st])[seg])->Dump();
714       }
715     }
716   }
717   return;
718 }
719
720   //__________________________________________________________________________
721 void AliMUONEventReconstructor::MakeSegmentsPerStation(Int_t Station)
722 {
723   // To make the list of segments in station number "Station" (0...)
724   // from the list of hits to be reconstructed.
725   // Updates "fNSegments"[Station].
726   // Segments in stations 4 and 5 are sorted
727   // according to increasing absolute value of "impact parameter"
728   AliMUONHitForRec *hit1Ptr, *hit2Ptr;
729   AliMUONSegment *segment;
730   Bool_t last2st;
731   Double_t bendingSlope, distBend, distNonBend, extBendCoor, extNonBendCoor,
732     impactParam, maxImpactParam;
733   AliMUON *pMUON  = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
734   if (fPrintLevel >= 1)
735     cout << "enter MakeSegmentsPerStation (0...) " << Station << endl;
736   // first and second chambers (0...) in the station
737   Int_t ch1 = 2 * Station;
738   Int_t ch2 = ch1 + 1;
739   // variable true for stations downstream of the dipole:
740   // Station(0..4) equal to 3 or 4
741   if ((Station == 3) || (Station == 4)) {
742     last2st = kTRUE;
743     // maximum impact parameter (cm) according to fMinBendingMomentum...
744     maxImpactParam =
745       TMath::Abs(GetImpactParamFromBendingMomentum(fMinBendingMomentum));
746   }
747   else last2st = kFALSE;
748   // extrapolation factor from Z of first chamber to Z of second chamber
749   // dZ to be changed to take into account fine structure of chambers ????
750   Double_t extrapFact =
751     (&(pMUON->Chamber(ch2)))->Z() / (&(pMUON->Chamber(ch1)))->Z();
752   // index for current segment
753   Int_t segmentIndex = 0;
754   // Loop over HitsForRec in the first chamber of the station
755   for (Int_t hit1 = fIndexOfFirstHitForRecPerChamber[ch1];
756        hit1 < fIndexOfFirstHitForRecPerChamber[ch1] + fNHitsForRecPerChamber[ch1];
757        hit1++) {
758     // pointer to the HitForRec
759     hit1Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit1]);
760     // extrapolation,
761     // on the straight line joining the HitForRec to the vertex (0,0,0),
762     // to the Z of the second chamber of the station
763     extBendCoor = extrapFact * hit1Ptr->GetBendingCoor();
764     extNonBendCoor = extrapFact * hit1Ptr->GetNonBendingCoor();
765     // Loop over HitsForRec in the second chamber of the station
766     for (Int_t hit2 = fIndexOfFirstHitForRecPerChamber[ch2];
767          hit2 < fIndexOfFirstHitForRecPerChamber[ch2] + fNHitsForRecPerChamber[ch2];
768          hit2++) {
769       // pointer to the HitForRec
770       hit2Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit2]);
771       // absolute values of distances, in bending and non bending planes,
772       // between the HitForRec in the second chamber
773       // and the previous extrapolation
774       distBend = TMath::Abs(hit2Ptr->GetBendingCoor() - extBendCoor);
775       distNonBend = TMath::Abs(hit2Ptr->GetNonBendingCoor() - extNonBendCoor);
776       if (last2st) {
777         // bending slope
778         bendingSlope = (hit1Ptr->GetBendingCoor() - hit2Ptr->GetBendingCoor()) /
779           (hit1Ptr->GetZ() - hit2Ptr->GetZ());
780         // absolute value of impact parameter
781         impactParam =
782           TMath::Abs(hit1Ptr->GetBendingCoor() - hit2Ptr->GetZ() * bendingSlope);
783       }
784       // check for distances not too large,
785       // and impact parameter not too big if stations downstream of the dipole.
786       // Conditions "distBend" and "impactParam" correlated for these stations ????
787       if ((distBend < fSegmentMaxDistBending[Station]) &&
788           (distNonBend < fSegmentMaxDistNonBending[Station]) &&
789           (!last2st || (impactParam < maxImpactParam))) {
790         // make new segment
791         segment = new ((*fSegmentsPtr[Station])[segmentIndex])
792           AliMUONSegment(hit1Ptr, hit2Ptr);
793         // update "link" to this segment from the hit in the first chamber
794         if (hit1Ptr->GetNSegments() == 0)
795           hit1Ptr->SetIndexOfFirstSegment(segmentIndex);
796         hit1Ptr->SetNSegments(hit1Ptr->GetNSegments() + 1);
797         if (fPrintLevel >= 10) {
798           cout << "segmentIndex(0...): " << segmentIndex
799                << "  distBend: " << distBend
800                << "  distNonBend: " << distNonBend
801                << endl;
802           segment->Dump();
803           cout << "HitForRec in first chamber" << endl;
804           hit1Ptr->Dump();
805           cout << "HitForRec in second chamber" << endl;
806           hit2Ptr->Dump();
807         };
808         // increment index for current segment
809         segmentIndex++;
810       }
811     } //for (Int_t hit2
812   } // for (Int_t hit1...
813   fNSegments[Station] = segmentIndex;
814   // Sorting according to "impact parameter" if station(1..5) 4 or 5,
815   // i.e. Station(0..4) 3 or 4, using the function "Compare".
816   // After this sorting, it is impossible to use
817   // the "fNSegments" and "fIndexOfFirstSegment"
818   // of the HitForRec in the first chamber to explore all segments formed with it.
819   // Is this sorting really needed ????
820   if ((Station == 3) || (Station == 4)) (fSegmentsPtr[Station])->Sort();
821   if (fPrintLevel >= 1) cout << "Station: " << Station << "  NSegments: "
822                              << fNSegments[Station] << endl;
823   return;
824 }
825
826   //__________________________________________________________________________
827 void AliMUONEventReconstructor::MakeTracks(void)
828 {
829   // To make the tracks,
830   // from the list of segments and points in all stations
831   if (fPrintLevel >= 1) cout << "enter MakeTracks" << endl;
832   ResetTracks();
833   // Look for candidates from at least 3 aligned points in stations(1..) 4 and 5
834   MakeTrackCandidates();
835   // Follow tracks in stations(1..) 3, 2 and 1
836   FollowTracks();
837   return;
838 }
839
840   //__________________________________________________________________________
841 Int_t AliMUONEventReconstructor::MakeTrackCandidatesWithTwoSegments(AliMUONSegment *BegSegment)
842 {
843   // To make track candidates with two segments in stations(1..) 4 and 5,
844   // the first segment being pointed to by "BegSegment".
845   // Returns the number of such track candidates.
846   Int_t endStation, iEndSegment, nbCan2Seg;
847   AliMUONSegment *endSegment, *extrapSegment;
848   AliMUONTrack *recTrack;
849   Double_t mcsFactor;
850   if (fPrintLevel >= 1) cout << "enter MakeTrackCandidatesWithTwoSegments" << endl;
851   // Station for the end segment
852   endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
853   // multiple scattering factor corresponding to one chamber
854   mcsFactor = 0.0136 /
855     GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
856   mcsFactor     = fChamberThicknessInX0 * mcsFactor * mcsFactor;
857   // linear extrapolation to end station
858   extrapSegment =
859     BegSegment->CreateSegmentFromLinearExtrapToStation(endStation, mcsFactor);
860   // number of candidates with 2 segments to 0
861   nbCan2Seg = 0;
862   // Loop over segments in the end station
863   for (iEndSegment = 0; iEndSegment < fNSegments[endStation]; iEndSegment++) {
864     // pointer to segment
865     endSegment = (AliMUONSegment*) ((*fSegmentsPtr[endStation])[iEndSegment]);
866     // test compatibility between current segment and "extrapSegment"
867     // 4 because 4 quantities in chi2
868     if ((endSegment->
869          NormalizedChi2WithSegment(extrapSegment,
870                                    fMaxSigma2Distance)) <= 4.0) {
871       // both segments compatible:
872       // make track candidate from "begSegment" and "endSegment"
873       if (fPrintLevel >= 2)
874         cout << "TrackCandidate with Segment " << iEndSegment <<
875           " in Station(0..) " << endStation << endl;
876       // flag for both segments in one track:
877       // to be done in track constructor ????
878       BegSegment->SetInTrack(kTRUE);
879       endSegment->SetInTrack(kTRUE);
880       recTrack = new ((*fRecTracksPtr)[fNRecTracks])
881         AliMUONTrack(BegSegment, endSegment, this);
882       fNRecTracks++;
883       if (fPrintLevel >= 10) recTrack->RecursiveDump();
884       // increment number of track candidates with 2 segments
885       nbCan2Seg++;
886     }
887   } // for (iEndSegment = 0;...
888   delete extrapSegment; // should not delete HitForRec's it points to !!!!
889   return nbCan2Seg;
890 }
891
892   //__________________________________________________________________________
893 Int_t AliMUONEventReconstructor::MakeTrackCandidatesWithOneSegmentAndOnePoint(AliMUONSegment *BegSegment)
894 {
895   // To make track candidates with one segment and one point
896   // in stations(1..) 4 and 5,
897   // the segment being pointed to by "BegSegment".
898   Int_t ch, ch1, ch2, endStation, iHit, iHitMax, iHitMin, nbCan1Seg1Hit;
899   AliMUONHitForRec *extrapHitForRec, *hit;
900   AliMUONTrack *recTrack;
901   Double_t mcsFactor;
902   if (fPrintLevel >= 1)
903     cout << "enter MakeTrackCandidatesWithOneSegmentAndOnePoint" << endl;
904   // station for the end point
905   endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
906   // multiple scattering factor corresponding to one chamber
907   mcsFactor = 0.0136 /
908     GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
909   mcsFactor     = fChamberThicknessInX0 * mcsFactor * mcsFactor;
910   // first and second chambers(0..) in the end station
911   ch1 = 2 * endStation;
912   ch2 = ch1 + 1;
913   // number of candidates to 0
914   nbCan1Seg1Hit = 0;
915   // Loop over chambers of the end station
916   for (ch = ch2; ch >= ch1; ch--) {
917     // linear extrapolation to chamber
918     extrapHitForRec =
919       BegSegment->CreateHitForRecFromLinearExtrapToChamber(ch, mcsFactor);
920     // limits for the hit index in the loop
921     iHitMin = fIndexOfFirstHitForRecPerChamber[ch];
922     iHitMax = iHitMin + fNHitsForRecPerChamber[ch];
923     // Loop over HitForRec's in the chamber
924     for (iHit = iHitMin; iHit < iHitMax; iHit++) {
925       // pointer to HitForRec
926       hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
927       // test compatibility between current HitForRec and "extrapHitForRec"
928       // 2 because 2 quantities in chi2
929       if ((hit->
930            NormalizedChi2WithHitForRec(extrapHitForRec,
931                                        fMaxSigma2Distance)) <= 2.0) {
932         // both HitForRec's compatible:
933         // make track candidate from begSegment and current HitForRec
934         if (fPrintLevel >= 2)
935           cout << "TrackCandidate with HitForRec " << iHit <<
936             " in Chamber(0..) " << ch << endl;
937         // flag for beginning segments in one track:
938         // to be done in track constructor ????
939         BegSegment->SetInTrack(kTRUE);
940         recTrack = new ((*fRecTracksPtr)[fNRecTracks])
941           AliMUONTrack(BegSegment, hit, this);
942         // the right place to eliminate "double counting" ???? how ????
943         fNRecTracks++;
944         if (fPrintLevel >= 10) recTrack->RecursiveDump();
945         // increment number of track candidates
946         nbCan1Seg1Hit++;
947       }
948     } // for (iHit = iHitMin;...
949     delete extrapHitForRec;
950   } // for (ch = ch2;...
951   return nbCan1Seg1Hit;
952 }
953
954   //__________________________________________________________________________
955 void AliMUONEventReconstructor::MakeTrackCandidates(void)
956 {
957   // To make track candidates
958   // with at least 3 aligned points in stations(1..) 4 and 5
959   // (two Segment's or one Segment and one HitForRec)
960   Int_t begStation, iBegSegment, nbCan1Seg1Hit, nbCan2Seg;
961   AliMUONSegment *begSegment;
962   if (fPrintLevel >= 1) cout << "enter MakeTrackCandidates" << endl;
963   // Loop over stations(1..) 5 and 4 for the beginning segment
964   for (begStation = 4; begStation > 2; begStation--) {
965     // Loop over segments in the beginning station
966     for (iBegSegment = 0; iBegSegment < fNSegments[begStation]; iBegSegment++) {
967       // pointer to segment
968       begSegment = (AliMUONSegment*) ((*fSegmentsPtr[begStation])[iBegSegment]);
969       if (fPrintLevel >= 2)
970         cout << "look for TrackCandidate's with Segment " << iBegSegment <<
971           " in Station(0..) " << begStation << endl;
972       // Look for track candidates with two segments,
973       // "begSegment" and all compatible segments in other station.
974       // Only for beginning station(1..) 5
975       // because candidates with 2 segments have to looked for only once.
976       if (begStation == 4)
977         nbCan2Seg = MakeTrackCandidatesWithTwoSegments(begSegment);
978       // Look for track candidates with one segment and one point,
979       // "begSegment" and all compatible HitForRec's in other station.
980       // Only if "begSegment" does not belong already to a track candidate.
981       // Is that a too strong condition ????
982       if (!(begSegment->GetInTrack()))
983         nbCan1Seg1Hit = MakeTrackCandidatesWithOneSegmentAndOnePoint(begSegment);
984     } // for (iBegSegment = 0;...
985   } // for (begStation = 4;...
986   return;
987 }
988
989   //__________________________________________________________________________
990 void AliMUONEventReconstructor::FollowTracks(void)
991 {
992   // Follow tracks in stations(1..) 3, 2 and 1
993   // too long: should be made more modular !!!!
994   AliMUONHitForRec *bestHit, *extrapHit, *hit;
995   AliMUONSegment *bestSegment, *extrapSegment, *segment;
996   AliMUONTrack *track, *nextTrack;
997   AliMUONTrackParam *trackParam1, trackParam[2], trackParamVertex;
998   Int_t ch, chInStation, chBestHit, iHit, iSegment, station, trackIndex;
999   Double_t bestChi2, chi2, dZ1, dZ2, dZ3, maxSigma2Distance, mcsFactor;
1000   AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
1001   // local maxSigma2Distance, for easy increase in testing
1002   maxSigma2Distance = fMaxSigma2Distance;
1003   if (fPrintLevel >= 2)
1004     cout << "enter FollowTracks" << endl;
1005   // Loop over track candidates
1006   track = (AliMUONTrack*) fRecTracksPtr->First();
1007   trackIndex = -1;
1008   while (track) {
1009     // Follow function for each track candidate ????
1010     trackIndex++;
1011     nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track); // prepare next track
1012     if (fPrintLevel >= 2)
1013       cout << "FollowTracks: track candidate(0..): " << trackIndex << endl;
1014     // Fit track candidate
1015     track->SetFitMCS(0); // without multiple Coulomb scattering
1016     track->SetFitNParam(3); // with 3 parameters (X = Y = 0)
1017     track->SetFitStart(0); // from parameters at vertex
1018     track->Fit();
1019     if (fPrintLevel >= 10) {
1020       cout << "FollowTracks: track candidate(0..): " << trackIndex
1021            << " after fit in stations(0..) 3 and 4" << endl;
1022       track->RecursiveDump();
1023     }
1024     // Loop over stations(1..) 3, 2 and 1
1025     // something SPECIAL for stations 2 and 1 for majority 3 coincidence ????
1026     // otherwise: majority coincidence 2 !!!!
1027     for (station = 2; station >= 0; station--) {
1028       // Track parameters at first track hit (smallest Z)
1029       trackParam1 = ((AliMUONTrackHit*)
1030                      (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1031       // extrapolation to station
1032       trackParam1->ExtrapToStation(station, trackParam);
1033       extrapSegment = new AliMUONSegment(); //  empty segment
1034       // multiple scattering factor corresponding to one chamber
1035       // and momentum in bending plane (not total)
1036       mcsFactor = 0.0136 * trackParam1->GetInverseBendingMomentum();
1037       mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1038       // Z difference from previous station
1039       dZ1 = (&(pMUON->Chamber(2 * station)))->Z() -
1040         (&(pMUON->Chamber(2 * station + 2)))->Z();
1041       // Z difference between the two previous stations
1042       dZ2 = (&(pMUON->Chamber(2 * station + 2)))->Z() -
1043         (&(pMUON->Chamber(2 * station + 4)))->Z();
1044       // Z difference between the two chambers in the previous station
1045       dZ3 = (&(pMUON->Chamber(2 * station)))->Z() -
1046         (&(pMUON->Chamber(2 * station + 1)))->Z();
1047       extrapSegment->SetBendingCoorReso2(fBendingResolution * fBendingResolution);
1048       extrapSegment->
1049         SetNonBendingCoorReso2(fNonBendingResolution * fNonBendingResolution);
1050       extrapSegment->UpdateFromStationTrackParam
1051         (trackParam, mcsFactor, dZ1, dZ2, dZ3, station,
1052          trackParam1->GetInverseBendingMomentum());
1053       bestChi2 = 5.0;
1054       bestSegment = NULL;
1055       if (fPrintLevel >= 10) {
1056         cout << "FollowTracks: track candidate(0..): " << trackIndex
1057              << " Look for segment in station(0..): " << station << endl;
1058       }
1059       // Loop over segments in station
1060       for (iSegment = 0; iSegment < fNSegments[station]; iSegment++) {
1061         // Look for best compatible Segment in station
1062         // should consider all possibilities ????
1063         // multiple scattering ????
1064         // separation in 2 functions: Segment and HitForRec ????
1065         segment = (AliMUONSegment*) ((*fSegmentsPtr[station])[iSegment]);
1066         chi2 = segment->NormalizedChi2WithSegment(extrapSegment, maxSigma2Distance);
1067         if (chi2 < bestChi2) {
1068           // update best Chi2 and Segment if better found
1069           bestSegment = segment;
1070           bestChi2 = chi2;
1071         }
1072       }
1073       if (bestSegment) {
1074         // best segment found: add it to track candidate
1075         track->AddSegment(bestSegment);
1076         // set track parameters at these two TrakHit's
1077         track->SetTrackParamAtHit(track->GetNTrackHits() - 2, &(trackParam[0]));
1078         track->SetTrackParamAtHit(track->GetNTrackHits() - 1, &(trackParam[1]));
1079         if (fPrintLevel >= 10) {
1080           cout << "FollowTracks: track candidate(0..): " << trackIndex
1081                << " Added segment in station(0..): " << station << endl;
1082           track->RecursiveDump();
1083         }
1084       }
1085       else {
1086         // No best segment found:
1087         // Look for best compatible HitForRec in station:
1088         // should consider all possibilities ????
1089         // multiple scattering ???? do about like for extrapSegment !!!!
1090         extrapHit = new AliMUONHitForRec(); //  empty hit
1091         bestChi2 = 3.0;
1092         bestHit = NULL;
1093         if (fPrintLevel >= 10) {
1094           cout << "FollowTracks: track candidate(0..): " << trackIndex
1095                << " Segment not found, look for hit in station(0..): " << station
1096                << endl;
1097         }
1098         // Loop over chambers of the station
1099         for (chInStation = 0; chInStation < 2; chInStation++) {
1100           // coordinates of extrapolated hit
1101           extrapHit->
1102             SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor());
1103           extrapHit->
1104             SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor());
1105           // resolutions from "extrapSegment"
1106           extrapHit->SetBendingReso2(extrapSegment->GetBendingCoorReso2());
1107           extrapHit->SetNonBendingReso2(extrapSegment->GetNonBendingCoorReso2());
1108           // Loop over hits in the chamber
1109           ch = 2 * station + chInStation;
1110           for (iHit = fIndexOfFirstHitForRecPerChamber[ch];
1111                iHit < fIndexOfFirstHitForRecPerChamber[ch] +
1112                  fNHitsForRecPerChamber[ch];
1113                iHit++) {
1114             hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1115             // condition for hit not already in segment ????
1116             chi2 = hit->NormalizedChi2WithHitForRec(extrapHit, maxSigma2Distance);
1117             if (chi2 < bestChi2) {
1118               // update best Chi2 and HitForRec if better found
1119               bestHit = hit;
1120               bestChi2 = chi2;
1121               chBestHit = chInStation;
1122             }
1123           }
1124         }
1125         if (bestHit) {
1126           // best hit found: add it to track candidate
1127           track->AddHitForRec(bestHit);
1128           // set track parameters at this TrackHit
1129           track->SetTrackParamAtHit(track->GetNTrackHits() - 1,
1130                                     &(trackParam[chBestHit]));
1131           if (fPrintLevel >= 10) {
1132             cout << "FollowTracks: track candidate(0..): " << trackIndex
1133                  << " Added hit in station(0..): " << station << endl;
1134             track->RecursiveDump();
1135           }
1136         }
1137         else {
1138           // Remove current track candidate and update fNRecTracks
1139           // To be checked: recursive delete of TrackHit's !!!!
1140           // For cleaner implementation: call track->Remove()
1141           // to be coded, with all cleanings,
1142           // including links between HitForRec's and TrackHit's !!!!
1143           fRecTracksPtr->Remove(track);
1144           fNRecTracks--;
1145           delete extrapSegment;
1146           break; // stop the search for this candidate:
1147           // exit from the loop over station
1148         }
1149       }
1150       delete extrapSegment;
1151       // Sort track hits according to increasing Z
1152       track->GetTrackHitsPtr()->Sort();
1153       // Update track parameters at first track hit (smallest Z)
1154       trackParam1 = ((AliMUONTrackHit*)
1155                      (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1156       // Track fit
1157       // with multiple Coulomb scattering if all stations
1158       if (station == 0) track->SetFitMCS(1);
1159       // without multiple Coulomb scattering if not all stations
1160       else track->SetFitMCS(0);
1161       track->SetFitNParam(5);  // with 5 parameters (momentum and position)
1162       track->SetFitStart(1);  // from parameters at first hit
1163       track->Fit();
1164       if (fPrintLevel >= 10) {
1165         cout << "FollowTracks: track candidate(0..): " << trackIndex
1166              << " after fit from station(0..): " << station << " to 4" << endl;
1167         track->RecursiveDump();
1168       }
1169       // Track extrapolation to the vertex through the absorber (Branson)
1170       // after going through the first station
1171       if (station == 0) {
1172         trackParamVertex = *trackParam1;
1173         (&trackParamVertex)->ExtrapToVertex();
1174         track->SetTrackParamAtVertex(&trackParamVertex);
1175         if (fPrintLevel >= 1) {
1176           cout << "FollowTracks: track candidate(0..): " << trackIndex
1177                << " after extrapolation to vertex" << endl;
1178           track->RecursiveDump();
1179         }
1180       }
1181     } // for (station = 2;...
1182     // go really to next track
1183     track = nextTrack;
1184   } // while (track)
1185   // Compression of track array (necessary after Remove ????)
1186   fRecTracksPtr->Compress();
1187   return;
1188 }
1189
1190   //__________________________________________________________________________
1191 void AliMUONEventReconstructor::EventDump(void)
1192 {
1193   // Dump reconstructed event (track parameters at vertex and at first hit),
1194   // and the particle parameters
1195
1196   AliMUONTrack *track;
1197   AliMUONTrackParam *trackParam, *trackParam1;
1198   TClonesArray *particles; // pointer to the particle list
1199   TParticle *p;
1200   Double_t bendingSlope, nonBendingSlope, pYZ;
1201   Double_t pX, pY, pZ, x, y, z, c;
1202   Int_t np, trackIndex, nTrackHits;
1203  
1204   if (fPrintLevel >= 1) cout << "****** enter EventDump ******" << endl;
1205   if (fPrintLevel >= 1) {
1206     cout << " Number of Reconstructed tracks :" <<  fNRecTracks << endl; 
1207   }
1208   fRecTracksPtr->Compress(); // for simple loop without "Next" since no hole
1209   // Loop over reconstructed tracks
1210   for (trackIndex = 0; trackIndex < fNRecTracks; trackIndex++) {
1211     if (fPrintLevel >= 1)
1212       cout << " track number: " << trackIndex << endl;
1213     // function for each track for modularity ????
1214     track = (AliMUONTrack*) ((*fRecTracksPtr)[trackIndex]);
1215     nTrackHits = track->GetNTrackHits();
1216     if (fPrintLevel >= 1)
1217       cout << " number of track hits: " << nTrackHits << endl;
1218     // track parameters at Vertex
1219     trackParam = track->GetTrackParamAtVertex();
1220     x = trackParam->GetNonBendingCoor();
1221     y = trackParam->GetBendingCoor();
1222     z = trackParam->GetZ();
1223     bendingSlope = trackParam->GetBendingSlope();
1224     nonBendingSlope = trackParam->GetNonBendingSlope();
1225     pYZ = 1/TMath::Abs(trackParam->GetInverseBendingMomentum());
1226     pZ = pYZ/TMath::Sqrt(1+bendingSlope*bendingSlope);
1227     pX = pZ * nonBendingSlope;
1228     pY = pZ * bendingSlope;
1229     c = TMath::Sign(1.0, trackParam->GetInverseBendingMomentum());
1230     if (fPrintLevel >= 1)
1231       printf(" track parameters at Vertex z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1232              z, x, y, pX, pY, pZ, c);
1233
1234     // track parameters at first hit
1235     trackParam1 = ((AliMUONTrackHit*)
1236                    (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1237     x = trackParam1->GetNonBendingCoor();
1238     y = trackParam1->GetBendingCoor();
1239     z = trackParam1->GetZ();
1240     bendingSlope = trackParam1->GetBendingSlope();
1241     nonBendingSlope = trackParam1->GetNonBendingSlope();
1242     pYZ = 1/TMath::Abs(trackParam1->GetInverseBendingMomentum());
1243     pZ = pYZ/TMath::Sqrt(1.0 + bendingSlope * bendingSlope);
1244     pX = pZ * nonBendingSlope;
1245     pY = pZ * bendingSlope;
1246     c = TMath::Sign(1.0, trackParam1->GetInverseBendingMomentum());
1247     if (fPrintLevel >= 1)
1248       printf(" track parameters at z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1249              z, x, y, pX, pY, pZ, c);
1250   }
1251   // informations about generated particles
1252   particles = gAlice->Particles();
1253   np = particles->GetEntriesFast();
1254   printf(" **** number of generated particles: %d  \n", np);
1255   
1256   for (Int_t iPart = 0; iPart < np; iPart++) {
1257     p = (TParticle*) particles->UncheckedAt(iPart);
1258     printf(" particle %d: type= %d px= %f py= %f pz= %f pdg= %d\n",
1259            iPart, p->GetPdgCode(), p->Px(), p->Py(), p->Pz(), p->GetPdgCode());    
1260   }
1261   return;
1262 }
1263
1264 void AliMUONEventReconstructor::Streamer(TBuffer &R__b)
1265 {
1266     ;
1267 }