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