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