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