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