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