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