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