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