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