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