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