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