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