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