1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 ////////////////////////////////////
18 // MUON track reconstructor in ALICE (class renamed from AliMUONEventReconstructor)
20 // This class contains as data:
21 // * the parameters for the track reconstruction
22 // * a pointer to the array of hits to be reconstructed (the event)
23 // * a pointer to the array of segments made with these hits inside each station
24 // * a pointer to the array of reconstructed tracks
26 // It contains as methods, among others:
27 // * MakeEventToBeReconstructed to build the array of hits to be reconstructed
28 // * MakeSegments to build the segments
29 // * MakeTracks to build the tracks
31 ////////////////////////////////////
34 #include <Riostream.h>
35 #include <TDirectory.h>
37 #include <TMatrixD.h> //AZ
39 #include "AliMUONTrackReconstructor.h"
41 #include "AliMUONConstants.h"
42 #include "AliMUONHitForRec.h"
43 #include "AliMUONTriggerTrack.h"
44 #include "AliMUONTriggerCircuit.h"
45 #include "AliMUONRawCluster.h"
46 #include "AliMUONLocalTrigger.h"
47 #include "AliMUONGlobalTrigger.h"
48 #include "AliMUONSegment.h"
49 #include "AliMUONTrack.h"
50 #include "AliMUONTrackHit.h"
52 #include "AliRun.h" // for gAlice
53 #include "AliRunLoader.h"
54 #include "AliLoader.h"
55 #include "AliMUONTrackK.h" //AZ
58 #include "AliTrackReference.h"
60 //************* Defaults parameters for reconstruction
61 const Double_t AliMUONTrackReconstructor::fgkDefaultMinBendingMomentum = 3.0;
62 const Double_t AliMUONTrackReconstructor::fgkDefaultMaxBendingMomentum = 3000.0;
63 const Double_t AliMUONTrackReconstructor::fgkDefaultMaxChi2 = 100.0;
64 const Double_t AliMUONTrackReconstructor::fgkDefaultMaxSigma2Distance = 16.0;
65 const Double_t AliMUONTrackReconstructor::fgkDefaultBendingResolution = 0.01;
66 const Double_t AliMUONTrackReconstructor::fgkDefaultNonBendingResolution = 0.144;
67 const Double_t AliMUONTrackReconstructor::fgkDefaultChamberThicknessInX0 = 0.03;
68 // Simple magnetic field:
69 // Value taken from macro MUONtracking.C: 0.7 T, hence 7 kG
70 // Length and Position from reco_muon.F, with opposite sign:
71 // Length = ZMAGEND-ZCOIL
72 // Position = (ZMAGEND+ZCOIL)/2
73 // to be ajusted differently from real magnetic field ????
74 const Double_t AliMUONTrackReconstructor::fgkDefaultSimpleBValue = 7.0;
75 const Double_t AliMUONTrackReconstructor::fgkDefaultSimpleBLength = 428.0;
76 const Double_t AliMUONTrackReconstructor::fgkDefaultSimpleBPosition = 1019.0;
77 const Int_t AliMUONTrackReconstructor::fgkDefaultRecTrackRefHits = 0;
78 const Double_t AliMUONTrackReconstructor::fgkDefaultEfficiency = 0.95;
80 ClassImp(AliMUONTrackReconstructor) // Class implementation in ROOT context
82 //__________________________________________________________________________
83 AliMUONTrackReconstructor::AliMUONTrackReconstructor(AliLoader* loader)
86 // Constructor for class AliMUONTrackReconstructor
87 SetReconstructionParametersToDefaults();
88 fTrackMethod = 1; //AZ - tracking method (1-default, 2-Kalman)
89 // Memory allocation for the TClonesArray of hits for reconstruction
90 // Is 10000 the right size ????
91 fHitsForRecPtr = new TClonesArray("AliMUONHitForRec", 10000);
92 fNHitsForRec = 0; // really needed or GetEntriesFast sufficient ????
93 // Memory allocation for the TClonesArray's of segments in stations
94 // Is 2000 the right size ????
95 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) {
96 fSegmentsPtr[st] = new TClonesArray("AliMUONSegment", 2000);
97 fNSegments[st] = 0; // really needed or GetEntriesFast sufficient ????
99 // Memory allocation for the TClonesArray of reconstructed tracks
100 // Is 10 the right size ????
101 fRecTracksPtr = new TClonesArray("AliMUONTrack", 10);
102 fNRecTracks = 0; // really needed or GetEntriesFast sufficient ????
103 // Memory allocation for the TClonesArray of hits on reconstructed tracks
104 // Is 100 the right size ????
105 //AZ fRecTrackHitsPtr = new TClonesArray("AliMUONTrack", 100);
106 fRecTrackHitsPtr = new TClonesArray("AliMUONTrackHit", 100);
107 fNRecTrackHits = 0; // really needed or GetEntriesFast sufficient ????
109 // Sign of fSimpleBValue according to sign of Bx value at (50,50,-950).
111 x[0] = 50.; x[1] = 50.; x[2] = -950.;
112 gAlice->Field()->Field(x, b);
113 fSimpleBValue = TMath::Sign(fSimpleBValue,(Double_t) b[0]);
114 fSimpleBPosition = TMath::Sign(fSimpleBPosition,(Double_t) x[2]);
115 // See how to get fSimple(BValue, BLength, BPosition)
116 // automatically calculated from the actual magnetic field ????
118 AliDebug(1,"AliMUONTrackReconstructor constructed with defaults");
119 if ( AliLog::GetGlobalDebugLevel()>0) Dump();
120 AliDebug(1,"Magnetic field from root file:");
121 if ( AliLog::GetGlobalDebugLevel()>0) gAlice->Field()->Dump();
124 // Initializions for track ref. background events
125 fBkgTrackRefFile = 0;
127 fBkgTrackRefParticles = 0;
129 fBkgTrackRefEventNumber = -1;
131 // initialize loader's
134 // initialize container
135 fMUONData = new AliMUONData(fLoader,"MUON","MUON");
139 //__________________________________________________________________________
140 AliMUONTrackReconstructor::AliMUONTrackReconstructor (const AliMUONTrackReconstructor& rhs)
143 // Protected copy constructor
145 AliFatal("Not implemented.");
148 AliMUONTrackReconstructor &
149 AliMUONTrackReconstructor::operator=(const AliMUONTrackReconstructor& rhs)
151 // Protected assignement operator
153 if (this == &rhs) return *this;
155 AliFatal("Not implemented.");
160 //__________________________________________________________________________
161 AliMUONTrackReconstructor::~AliMUONTrackReconstructor(void)
163 // Destructor for class AliMUONTrackReconstructor
164 delete fHitsForRecPtr; // Correct destruction of everything ???? or delete [] ????
165 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++)
166 delete fSegmentsPtr[st]; // Correct destruction of everything ????
170 //__________________________________________________________________________
171 void AliMUONTrackReconstructor::SetReconstructionParametersToDefaults(void)
173 // Set reconstruction parameters to default values
174 // Would be much more convenient with a structure (or class) ????
175 fMinBendingMomentum = fgkDefaultMinBendingMomentum;
176 fMaxBendingMomentum = fgkDefaultMaxBendingMomentum;
177 fMaxChi2 = fgkDefaultMaxChi2;
178 fMaxSigma2Distance = fgkDefaultMaxSigma2Distance;
180 // ******** Parameters for making HitsForRec
182 // like in TRACKF_STAT:
183 // 2 degrees for stations 1 and 2, or ch(0...) from 0 to 3;
184 // 30 cm for stations 3 to 5, or ch(0...) from 4 to 9
185 for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
186 if (ch < 4) fRMin[ch] = TMath::Abs(AliMUONConstants::DefaultChamberZ(ch)) *
187 2.0 * TMath::Pi() / 180.0;
188 else fRMin[ch] = 30.0;
189 // maximum radius at 10 degrees and Z of chamber
190 fRMax[ch] = TMath::Abs(AliMUONConstants::DefaultChamberZ(ch)) *
191 10.0 * TMath::Pi() / 180.0;
194 // ******** Parameters for making segments
195 // should be parametrized ????
196 // according to interval between chambers in a station ????
197 // Maximum distance in non bending plane
198 // 5 * 0.22 just to remember the way it was made in TRACKF_STAT
200 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++)
201 fSegmentMaxDistNonBending[st] = 5. * 0.22;
202 // Maximum distance in bending plane:
203 // values from TRACKF_STAT, corresponding to (J psi 20cm),
204 // scaled to the real distance between chambers in a station
205 fSegmentMaxDistBending[0] = TMath::Abs( 1.5 *
206 (AliMUONConstants::DefaultChamberZ(1) - AliMUONConstants::DefaultChamberZ(0)) / 20.0);
207 fSegmentMaxDistBending[1] = TMath::Abs( 1.5 *
208 (AliMUONConstants::DefaultChamberZ(3) - AliMUONConstants::DefaultChamberZ(2)) / 20.0);
209 fSegmentMaxDistBending[2] = TMath::Abs( 3.0 *
210 (AliMUONConstants::DefaultChamberZ(5) - AliMUONConstants::DefaultChamberZ(4)) / 20.0);
211 fSegmentMaxDistBending[3] = TMath::Abs( 6.0 *
212 (AliMUONConstants::DefaultChamberZ(7) - AliMUONConstants::DefaultChamberZ(6)) / 20.0);
213 fSegmentMaxDistBending[4] = TMath::Abs( 6.0 *
214 (AliMUONConstants::DefaultChamberZ(9) - AliMUONConstants::DefaultChamberZ(8)) / 20.0);
217 fBendingResolution = fgkDefaultBendingResolution;
218 fNonBendingResolution = fgkDefaultNonBendingResolution;
219 fChamberThicknessInX0 = fgkDefaultChamberThicknessInX0;
220 fSimpleBValue = fgkDefaultSimpleBValue;
221 fSimpleBLength = fgkDefaultSimpleBLength;
222 fSimpleBPosition = fgkDefaultSimpleBPosition;
223 fRecTrackRefHits = fgkDefaultRecTrackRefHits;
224 fEfficiency = fgkDefaultEfficiency;
228 //__________________________________________________________________________
229 Double_t AliMUONTrackReconstructor::GetImpactParamFromBendingMomentum(Double_t BendingMomentum) const
231 // Returns impact parameter at vertex in bending plane (cm),
232 // from the signed bending momentum "BendingMomentum" in bending plane (GeV/c),
233 // using simple values for dipole magnetic field.
234 // The sign of "BendingMomentum" is the sign of the charge.
235 return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
239 //__________________________________________________________________________
240 Double_t AliMUONTrackReconstructor::GetBendingMomentumFromImpactParam(Double_t ImpactParam) const
242 // Returns signed bending momentum in bending plane (GeV/c),
243 // the sign being the sign of the charge for particles moving forward in Z,
244 // from the impact parameter "ImpactParam" at vertex in bending plane (cm),
245 // using simple values for dipole magnetic field.
246 return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
250 //__________________________________________________________________________
251 void AliMUONTrackReconstructor::SetBkgTrackRefFile(Text_t *BkgTrackRefFileName)
253 // Set background file ... for track ref. hits
254 // Must be called after having loaded the firts signal event
255 AliDebug(1,Form("Enter SetBkgTrackRefFile with BkgTrackRefFileName %s",BkgTrackRefFileName));
257 if (strlen(BkgTrackRefFileName)) {
258 // BkgTrackRefFileName not empty: try to open the file
260 if(AliLog::GetGlobalDebugLevel()>1) {
261 cout << "Before File(Bkg)" << endl; gDirectory->Dump();
263 fBkgTrackRefFile = new TFile(BkgTrackRefFileName);
264 if(AliLog::GetGlobalDebugLevel()>1) {
265 cout << "After File(Bkg)" << endl; gDirectory->Dump();
267 if (fBkgTrackRefFile-> IsOpen()) {
268 if(AliLog::GetGlobalDebugLevel()>0) {
269 cout << "Background for Track ref. hits in file: ``" << BkgTrackRefFileName
270 << "'' successfully opened" << endl;}
273 cout << "Background for Track Ref. hits in file: " << BkgTrackRefFileName << endl;
274 cout << "NOT FOUND: EXIT" << endl;
275 exit(0); // right instruction for exit ????
277 // Arrays for "particles" and "hits"
278 fBkgTrackRefParticles = new TClonesArray("TParticle", 200);
279 // Event number to -1 for initialization
280 fBkgTrackRefEventNumber = -1;
281 // Back to the signal file:
282 // first signal event must have been loaded previously,
283 // otherwise, Segmentation violation at the next instruction
284 // How is it possible to do smething better ????
285 ((gAlice->TreeK())->GetCurrentFile())->cd();
286 if(AliLog::GetGlobalDebugLevel()>1) cout << "After cd(gAlice)" << endl; gDirectory->Dump();
291 //__________________________________________________________________________
292 void AliMUONTrackReconstructor::NextBkgTrackRefEvent(void)
294 // Get next event in background file for track ref. hits
295 // Goes back to event number 0 when end of file is reached
298 AliDebug(1,"Enter NextBkgTrackRefEvent");
299 // Clean previous event
300 if(fBkgTrackRefTK) delete fBkgTrackRefTK;
301 fBkgTrackRefTK = NULL;
302 if(fBkgTrackRefParticles) fBkgTrackRefParticles->Clear();
303 if(fBkgTrackRefTTR) delete fBkgTrackRefTTR;
304 fBkgTrackRefTTR = NULL;
305 // Increment event number
306 fBkgTrackRefEventNumber++;
307 // Get access to Particles and Hits for event from background file
308 if (AliLog::GetGlobalDebugLevel()>1) cout << "Before cd(Bkg)" << endl; gDirectory->Dump();
309 fBkgTrackRefFile->cd();
310 if (AliLog::GetGlobalDebugLevel()>1) cout << "After cd(Bkg)" << endl; gDirectory->Dump();
311 // Particles: TreeK for event and branch "Particles"
312 sprintf(treeName, "TreeK%d", fBkgTrackRefEventNumber);
313 fBkgTrackRefTK = (TTree*)gDirectory->Get(treeName);
314 if (!fBkgTrackRefTK) {
316 AliDebug(1,Form("Cannot find Kine Tree for background event: %d",fBkgTrackRefEventNumber));
317 AliDebug(1,"Goes back to event 0");
319 fBkgTrackRefEventNumber = 0;
320 sprintf(treeName, "TreeK%d", fBkgTrackRefEventNumber);
321 fBkgTrackRefTK = (TTree*)gDirectory->Get(treeName);
322 if (!fBkgTrackRefTK) {
323 AliError(Form("cannot find Kine Tree for background event: %d",fBkgTrackRefEventNumber));
328 fBkgTrackRefTK->SetBranchAddress("Particles", &fBkgTrackRefParticles);
329 fBkgTrackRefTK->GetEvent(0); // why event 0 ???? necessary ????
330 // Hits: TreeH for event and branch "MUON"
331 sprintf(treeName, "TreeTR%d", fBkgTrackRefEventNumber);
332 fBkgTrackRefTTR = (TTree*)gDirectory->Get(treeName);
333 if (!fBkgTrackRefTTR) {
334 AliError(Form("cannot find Hits Tree for background event: %d",fBkgTrackRefEventNumber));
337 fBkgTrackRefTTR->GetEntries(); // necessary ????
338 // Back to the signal file
339 ((gAlice->TreeK())->GetCurrentFile())->cd();
340 if (AliLog::GetGlobalDebugLevel()>1)
341 cout << "After cd(gAlice)" << endl; gDirectory->Dump();
345 //__________________________________________________________________________
346 void AliMUONTrackReconstructor::EventReconstruct(void)
348 // To reconstruct one event
349 AliDebug(1,"Enter EventReconstruct");
350 MakeEventToBeReconstructed();
353 if (fMUONData->IsTriggerTrackBranchesInTree())
354 ValidateTracksWithTrigger();
356 // Add tracks to MUON data container
357 for(Int_t i=0; i<GetNRecTracks(); i++) {
358 AliMUONTrack * track = (AliMUONTrack*) GetRecTracksPtr()->At(i);
359 fMUONData->AddRecTrack(*track);
365 //__________________________________________________________________________
366 void AliMUONTrackReconstructor::EventReconstructTrigger(void)
368 // To reconstruct one event
369 AliDebug(1,"Enter EventReconstructTrigger");
374 //__________________________________________________________________________
375 void AliMUONTrackReconstructor::ResetHitsForRec(void)
377 // To reset the array and the number of HitsForRec,
378 // and also the number of HitsForRec
379 // and the index of the first HitForRec per chamber
380 if (fHitsForRecPtr) fHitsForRecPtr->Clear();
382 for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++)
383 fNHitsForRecPerChamber[ch] = fIndexOfFirstHitForRecPerChamber[ch] = 0;
387 //__________________________________________________________________________
388 void AliMUONTrackReconstructor::ResetSegments(void)
390 // To reset the TClonesArray of segments and the number of Segments
392 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) {
393 if (fSegmentsPtr[st]) fSegmentsPtr[st]->Clear();
399 //__________________________________________________________________________
400 void AliMUONTrackReconstructor::ResetTracks(void)
402 // To reset the TClonesArray of reconstructed tracks
403 if (fRecTracksPtr) fRecTracksPtr->Delete();
404 // Delete in order that the Track destructors are called,
405 // hence the space for the TClonesArray of pointers to TrackHit's is freed
411 //__________________________________________________________________________
412 void AliMUONTrackReconstructor::ResetTrackHits(void)
414 // To reset the TClonesArray of hits on reconstructed tracks
415 if (fRecTrackHitsPtr) fRecTrackHitsPtr->Clear();
420 //__________________________________________________________________________
421 void AliMUONTrackReconstructor::MakeEventToBeReconstructed(void)
423 // To make the list of hits to be reconstructed,
424 // either from the track ref. hits or from the raw clusters
425 // according to the parameter set for the reconstructor
426 // TString evfoldname = AliConfig::GetDefaultEventFolderName();//to be interfaced properly
428 // AliRunLoader* rl = AliRunLoader::GetRunLoader(evfoldname);
431 // Error("MakeEventToBeReconstructed",
432 // "Can not find Run Loader in Event Folder named %s.",
433 // evfoldname.Data());
436 // AliLoader* gime = rl->GetLoader("MUONLoader");
439 // Error("MakeEventToBeReconstructed","Can not get MUON Loader from Run Loader.");
442 AliRunLoader *runLoader = fLoader->GetRunLoader();
444 AliDebug(1,"Enter MakeEventToBeReconstructed");
446 if (fRecTrackRefHits == 1) {
447 // Reconstruction from track ref. hits
448 // Back to the signal file
449 TTree* treeTR = runLoader->TreeTR();
452 Int_t retval = runLoader->LoadTrackRefs("READ");
455 AliError("Error occured while loading hits.");
458 treeTR = runLoader->TreeTR();
461 AliError("Can not get TreeTR");
466 AddHitsForRecFromTrackRef(treeTR,1);
469 AddHitsForRecFromTrackRef(fBkgTrackRefTTR,0);
470 // Sort HitsForRec in increasing order with respect to chamber number
471 SortHitsForRecWithIncreasingChamber();
474 // Reconstruction from raw clusters
475 // AliMUON *MUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
476 // Security on MUON ????
477 // TreeR assumed to be be "prepared" in calling function
478 // by "MUON->GetTreeR(nev)" ????
479 TTree *treeR = fLoader->TreeR();
480 fMUONData->SetTreeAddress("RC");
481 AddHitsForRecFromRawClusters(treeR);
482 // No sorting: it is done automatically in the previous function
485 AliDebug(1,"End of MakeEventToBeReconstructed");
486 if (AliLog::GetGlobalDebugLevel() > 0) {
487 AliDebug(1, Form("NHitsForRec: %d",fNHitsForRec));
488 for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
489 AliDebug(1, Form("Chamber(0...): %d",ch));
490 AliDebug(1, Form("NHitsForRec: %d", fNHitsForRecPerChamber[ch]));
491 AliDebug(1, Form("Index(first HitForRec): %d", fIndexOfFirstHitForRecPerChamber[ch]));
492 for (Int_t hit = fIndexOfFirstHitForRecPerChamber[ch];
493 hit < fIndexOfFirstHitForRecPerChamber[ch] + fNHitsForRecPerChamber[ch];
495 AliDebug(1, Form("HitForRec index(0...): %d",hit));
496 ((*fHitsForRecPtr)[hit])->Dump();
503 //__________________________________________________________________________
504 void AliMUONTrackReconstructor::AddHitsForRecFromTrackRef(TTree *TTR, Int_t signal)
506 // To add to the list of hits for reconstruction
507 // the signal hits from a track reference tree TreeTR.
508 TClonesArray *listOfTrackRefs = NULL;
509 AliTrackReference *trackRef;
512 AliDebug(2,Form("Enter AddHitsForRecFromTrackRef with TTR: %d", TTR));
513 if (TTR == NULL) return;
515 listOfTrackRefs = CleanTrackRefs(TTR);
517 Int_t ntracks = listOfTrackRefs->GetEntriesFast();
519 AliDebug(2,Form("ntracks: %d", ntracks));
521 for (Int_t index = 0; index < ntracks; index++) {
522 trackRef = (AliTrackReference*) listOfTrackRefs->At(index);
523 track = trackRef->GetTrack();
525 NewHitForRecFromTrackRef(trackRef,track,signal);
526 } // end of track ref.
528 listOfTrackRefs->Delete();
529 delete listOfTrackRefs;
534 //__________________________________________________________________________
535 AliMUONHitForRec* AliMUONTrackReconstructor::NewHitForRecFromTrackRef(AliTrackReference* Hit, Int_t TrackNumber, Int_t Signal)
537 // To make a new hit for reconstruction from a track ref. hit pointed to by "Hit",
538 // with the track numbered "TrackNumber",
539 // either from signal ("Signal" = 1) or background ("Signal" = 0) event.
540 // Selects hits in tracking (not trigger) chambers.
541 // Takes into account the efficiency (fEfficiency)
542 // and the smearing from resolution (fBendingResolution and fNonBendingResolution).
543 // Adds a condition on the radius between RMin and RMax
544 // to better simulate the real chambers.
545 // Returns the pointer to the new hit for reconstruction,
546 // or NULL in case of inefficiency or non tracking chamber or bad radius.
547 // No condition on at most 20 cm from a muon from a resonance
548 // like in Fortran TRACKF_STAT.
549 AliMUONHitForRec* hitForRec;
550 Double_t bendCoor, nonBendCoor, radius;
551 Int_t chamber = AliMUONConstants::ChamberNumber(Hit->Z()); // chamber(0...)
552 if (chamber < 0) return NULL;
553 // only in tracking chambers (fChamber starts at 1)
554 if (chamber >= AliMUONConstants::NTrackingCh()) return NULL;
555 // only if hit is efficient (keep track for checking ????)
556 if (gRandom->Rndm() > fEfficiency) return NULL;
557 // only if radius between RMin and RMax
559 nonBendCoor = Hit->X();
560 radius = TMath::Sqrt((bendCoor * bendCoor) + (nonBendCoor * nonBendCoor));
561 // This cut is not needed with a realistic chamber geometry !!!!
562 // if ((radius < fRMin[chamber]) || (radius > fRMax[chamber])) return NULL;
563 // new AliMUONHitForRec from track ref. hit and increment number of AliMUONHitForRec's
564 hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(Hit);
566 // add smearing from resolution
567 hitForRec->SetBendingCoor(bendCoor + gRandom->Gaus(0., fBendingResolution));
568 hitForRec->SetNonBendingCoor(nonBendCoor
569 + gRandom->Gaus(0., fNonBendingResolution));
570 // // !!!! without smearing
571 // hitForRec->SetBendingCoor(bendCoor);
572 // hitForRec->SetNonBendingCoor(nonBendCoor);
573 // more information into HitForRec
574 // resolution: angular effect to be added here ????
575 hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
576 hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
578 hitForRec->SetTTRTrack(TrackNumber);
579 hitForRec->SetTrackRefSignal(Signal);
580 if (AliLog::GetGlobalDebugLevel()> 1) {
581 AliDebug(2,Form("Track: %d", TrackNumber));
583 cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
588 //__________________________________________________________________________
589 TClonesArray* AliMUONTrackReconstructor::CleanTrackRefs(TTree *treeTR)
591 // Make hits from track ref..
592 // Re-calculate hits parameters because two AliTrackReferences are recorded for
593 // each chamber (one when particle is entering + one when particle is leaving
594 // the sensitive volume)
596 AliTrackReference *trackReference;
597 Float_t x1, y1, z1, pX1, pY1, pZ1;
598 Float_t x2, y2, z2, pX2, pY2, pZ2;
599 Int_t track1, track2;
601 Float_t maxGasGap = 1.; // cm
605 AliTrackReference *trackReferenceNew = new AliTrackReference();
607 TClonesArray* trackRefs = new TClonesArray("AliTrackReference", 10);
608 TClonesArray* cleanTrackRefs = new TClonesArray("AliTrackReference", 10);
610 if (treeTR == NULL) return NULL;
611 TBranch* branch = treeTR->GetBranch("MUON");
612 if (branch == NULL) return NULL;
613 branch->SetAddress(&trackRefs);
615 Int_t nTrack = (Int_t)treeTR->GetEntries();
616 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
617 treeTR->GetEntry(iTrack);
620 while (iHit1 < trackRefs->GetEntries()) {
621 trackReference = (AliTrackReference*)trackRefs->At(iHit1);
622 x1 = trackReference->X();
623 y1 = trackReference->Y();
624 z1 = trackReference->Z();
625 pX1 = trackReference->Px();
626 pY1 = trackReference->Py();
627 pZ1 = trackReference->Pz();
628 track1 = trackReference->GetTrack();
631 for (Int_t iHit2 = iHit1+1; iHit2 < trackRefs->GetEntries(); iHit2++) {
632 trackReference = (AliTrackReference*)trackRefs->At(iHit2);
633 x2 = trackReference->X();
634 y2 = trackReference->Y();
635 z2 = trackReference->Z();
636 pX2 = trackReference->Px();
637 pY2 = trackReference->Py();
638 pZ2 = trackReference->Pz();
639 track2 = trackReference->GetTrack();
640 if (track2 == track1 && TMath::Abs(z2-z1) < maxGasGap ) {
655 pX1 /= (Float_t)nRec;
656 pY1 /= (Float_t)nRec;
657 pZ1 /= (Float_t)nRec;
658 trackReferenceNew->SetPosition(x1,y1,z1);
659 trackReferenceNew->SetMomentum(pX1,pY1,pZ1);
660 trackReferenceNew->SetTrack(track1);
661 {new ((*cleanTrackRefs)[cleanTrackRefs->GetEntriesFast()]) AliTrackReference(*trackReferenceNew);}
668 delete trackReferenceNew;
669 return cleanTrackRefs;
672 //__________________________________________________________________________
673 void AliMUONTrackReconstructor::SortHitsForRecWithIncreasingChamber()
675 // Sort HitsForRec's in increasing order with respect to chamber number.
676 // Uses the function "Compare".
677 // Update the information for HitsForRec per chamber too.
678 Int_t ch, nhits, prevch;
679 fHitsForRecPtr->Sort();
680 for (ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
681 fNHitsForRecPerChamber[ch] = 0;
682 fIndexOfFirstHitForRecPerChamber[ch] = 0;
684 prevch = 0; // previous chamber
685 nhits = 0; // number of hits in current chamber
686 // Loop over HitsForRec
687 for (Int_t hit = 0; hit < fNHitsForRec; hit++) {
688 // chamber number (0...)
689 ch = ((AliMUONHitForRec*) ((*fHitsForRecPtr)[hit]))->GetChamberNumber();
690 // increment number of hits in current chamber
691 (fNHitsForRecPerChamber[ch])++;
692 // update index of first HitForRec in current chamber
693 // if chamber number different from previous one
695 fIndexOfFirstHitForRecPerChamber[ch] = hit;
702 //__________________________________________________________________________
703 void AliMUONTrackReconstructor::AddHitsForRecFromRawClusters(TTree* TR)
705 // To add to the list of hits for reconstruction all the raw clusters
706 // No condition added, like in Fortran TRACKF_STAT,
707 // on the radius between RMin and RMax.
708 AliMUONHitForRec *hitForRec;
709 AliMUONRawCluster *clus;
710 Int_t iclus, nclus, nTRentries;
711 TClonesArray *rawclusters;
712 AliDebug(1,"Enter AddHitsForRecFromRawClusters");
714 nTRentries = Int_t(TR->GetEntries());
715 if (nTRentries != 1) {
716 AliError(Form("nTRentries = %d not equal to 1 ",nTRentries));
719 fLoader->TreeR()->GetEvent(0); // only one entry
721 // Loop over tracking chambers
722 for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
723 // number of HitsForRec to 0 for the chamber
724 fNHitsForRecPerChamber[ch] = 0;
725 // index of first HitForRec for the chamber
726 if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0;
727 else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec;
728 rawclusters =fMUONData->RawClusters(ch);
729 // pMUON->ResetRawClusters();
730 // TR->GetEvent((Int_t) (TR->GetEntries()) - 1); // to be checked ????
731 nclus = (Int_t) (rawclusters->GetEntries());
732 // Loop over (cathode correlated) raw clusters
733 for (iclus = 0; iclus < nclus; iclus++) {
734 clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(iclus);
735 // new AliMUONHitForRec from raw cluster
736 // and increment number of AliMUONHitForRec's (total and in chamber)
737 hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(clus);
739 (fNHitsForRecPerChamber[ch])++;
740 // more information into HitForRec
741 // resolution: info should be already in raw cluster and taken from it ????
742 hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
743 hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
744 // original raw cluster
745 hitForRec->SetChamberNumber(ch);
746 hitForRec->SetHitNumber(iclus);
747 // Z coordinate of the raw cluster (cm)
748 hitForRec->SetZ(clus->GetZ(0));
749 if (AliLog::GetGlobalDebugLevel() > 1) {
750 cout << "chamber (0...): " << ch <<
751 " raw cluster (0...): " << iclus << endl;
753 cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
755 } // end of cluster loop
756 } // end of chamber loop
757 SortHitsForRecWithIncreasingChamber(); //AZ
761 //__________________________________________________________________________
762 void AliMUONTrackReconstructor::MakeSegments(void)
764 // To make the list of segments in all stations,
765 // from the list of hits to be reconstructed
766 AliDebug(1,"Enter MakeSegments");
768 // Loop over stations
769 Int_t nb = (fTrackMethod == 2) ? 3 : 0; //AZ
770 //AZ for (Int_t st = 0; st < fgkMaxMuonTrackingStations; st++)
771 for (Int_t st = nb; st < AliMUONConstants::NTrackingCh()/2; st++) //AZ
772 MakeSegmentsPerStation(st);
773 if (AliLog::GetGlobalDebugLevel() > 1) {
774 cout << "end of MakeSegments" << endl;
775 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) {
776 cout << "station(0...): " << st
777 << " Segments: " << fNSegments[st]
780 seg < fNSegments[st];
782 cout << "Segment index(0...): " << seg << endl;
783 ((*fSegmentsPtr[st])[seg])->Dump();
790 //__________________________________________________________________________
791 void AliMUONTrackReconstructor::MakeSegmentsPerStation(Int_t Station)
793 // To make the list of segments in station number "Station" (0...)
794 // from the list of hits to be reconstructed.
795 // Updates "fNSegments"[Station].
796 // Segments in stations 4 and 5 are sorted
797 // according to increasing absolute value of "impact parameter"
798 AliMUONHitForRec *hit1Ptr, *hit2Ptr;
799 AliMUONSegment *segment;
801 Double_t bendingSlope, distBend, distNonBend, extBendCoor, extNonBendCoor,
802 impactParam = 0., maxImpactParam = 0., minImpactParam = 0.; // =0 to avoid compilation warnings.
803 AliDebug(1,Form("Enter MakeSegmentsPerStation (0...) %d",Station));
804 // first and second chambers (0...) in the station
805 Int_t ch1 = 2 * Station;
807 // variable true for stations downstream of the dipole:
808 // Station(0..4) equal to 3 or 4
809 if ((Station == 3) || (Station == 4)) {
811 // maximum impact parameter (cm) according to fMinBendingMomentum...
813 TMath::Abs(GetImpactParamFromBendingMomentum(fMinBendingMomentum));
814 // minimum impact parameter (cm) according to fMaxBendingMomentum...
816 TMath::Abs(GetImpactParamFromBendingMomentum(fMaxBendingMomentum));
818 else last2st = kFALSE;
819 // extrapolation factor from Z of first chamber to Z of second chamber
820 // dZ to be changed to take into account fine structure of chambers ????
822 // index for current segment
823 Int_t segmentIndex = 0;
824 // Loop over HitsForRec in the first chamber of the station
825 for (Int_t hit1 = fIndexOfFirstHitForRecPerChamber[ch1];
826 hit1 < fIndexOfFirstHitForRecPerChamber[ch1] + fNHitsForRecPerChamber[ch1];
828 // pointer to the HitForRec
829 hit1Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit1]);
831 // on the straight line joining the HitForRec to the vertex (0,0,0),
832 // to the Z of the second chamber of the station
833 // Loop over HitsForRec in the second chamber of the station
834 for (Int_t hit2 = fIndexOfFirstHitForRecPerChamber[ch2];
835 hit2 < fIndexOfFirstHitForRecPerChamber[ch2] + fNHitsForRecPerChamber[ch2];
837 // pointer to the HitForRec
838 hit2Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit2]);
839 // absolute values of distances, in bending and non bending planes,
840 // between the HitForRec in the second chamber
841 // and the previous extrapolation
842 extrapFact = hit2Ptr->GetZ()/ hit1Ptr->GetZ();
843 extBendCoor = extrapFact * hit1Ptr->GetBendingCoor();
844 extNonBendCoor = extrapFact * hit1Ptr->GetNonBendingCoor();
845 distBend = TMath::Abs(hit2Ptr->GetBendingCoor() - extBendCoor);
846 distNonBend = TMath::Abs(hit2Ptr->GetNonBendingCoor() - extNonBendCoor);
849 bendingSlope = (hit1Ptr->GetBendingCoor() - hit2Ptr->GetBendingCoor()) /
850 (hit1Ptr->GetZ() - hit2Ptr->GetZ());
851 // absolute value of impact parameter
853 TMath::Abs(hit1Ptr->GetBendingCoor() - hit1Ptr->GetZ() * bendingSlope);
855 // check for distances not too large,
856 // and impact parameter not too big if stations downstream of the dipole.
857 // Conditions "distBend" and "impactParam" correlated for these stations ????
858 if ((distBend < fSegmentMaxDistBending[Station]) &&
859 (distNonBend < fSegmentMaxDistNonBending[Station]) &&
860 (!last2st || (impactParam < maxImpactParam)) &&
861 (!last2st || (impactParam > minImpactParam))) {
863 segment = new ((*fSegmentsPtr[Station])[segmentIndex])
864 AliMUONSegment(hit1Ptr, hit2Ptr);
865 // update "link" to this segment from the hit in the first chamber
866 if (hit1Ptr->GetNSegments() == 0)
867 hit1Ptr->SetIndexOfFirstSegment(segmentIndex);
868 hit1Ptr->SetNSegments(hit1Ptr->GetNSegments() + 1);
869 if (AliLog::GetGlobalDebugLevel() > 1) {
870 cout << "segmentIndex(0...): " << segmentIndex
871 << " distBend: " << distBend
872 << " distNonBend: " << distNonBend
875 cout << "HitForRec in first chamber" << endl;
877 cout << "HitForRec in second chamber" << endl;
880 // increment index for current segment
884 } // for (Int_t hit1...
885 fNSegments[Station] = segmentIndex;
886 // Sorting according to "impact parameter" if station(1..5) 4 or 5,
887 // i.e. Station(0..4) 3 or 4, using the function "Compare".
888 // After this sorting, it is impossible to use
889 // the "fNSegments" and "fIndexOfFirstSegment"
890 // of the HitForRec in the first chamber to explore all segments formed with it.
891 // Is this sorting really needed ????
892 if ((Station == 3) || (Station == 4)) (fSegmentsPtr[Station])->Sort();
893 AliDebug(1,Form("Station: %d NSegments: %d ", Station, fNSegments[Station]));
897 //__________________________________________________________________________
898 void AliMUONTrackReconstructor::MakeTracks(void)
900 // To make the tracks,
901 // from the list of segments and points in all stations
902 AliDebug(1,"Enter MakeTracks");
903 // The order may be important for the following Reset's
906 if (fTrackMethod == 2) { //AZ - Kalman filter
907 MakeTrackCandidatesK();
908 if (fRecTracksPtr->GetEntriesFast() == 0) return;
909 // Follow tracks in stations(1..) 3, 2 and 1
911 // Remove double tracks
912 RemoveDoubleTracksK();
913 // Propagate tracks to the vertex thru absorber
915 // Fill AliMUONTrack data members
918 // Look for candidates from at least 3 aligned points in stations(1..) 4 and 5
919 MakeTrackCandidates();
920 // Follow tracks in stations(1..) 3, 2 and 1
922 // Remove double tracks
923 RemoveDoubleTracks();
924 UpdateTrackParamAtHit();
925 UpdateHitForRecAtHit();
930 //__________________________________________________________________________
931 void AliMUONTrackReconstructor::ValidateTracksWithTrigger(void)
933 // Try to match track from tracking system with trigger track
935 TClonesArray *recTriggerTracks;
937 fMUONData->SetTreeAddress("RL");
938 fMUONData->GetRecTriggerTracks();
939 recTriggerTracks = fMUONData->RecTriggerTracks();
941 track = (AliMUONTrack*) fRecTracksPtr->First();
943 track->MatchTriggerTrack(recTriggerTracks);
944 track = (AliMUONTrack*) fRecTracksPtr->After(track);
949 //__________________________________________________________________________
950 Bool_t AliMUONTrackReconstructor::MakeTriggerTracks(void)
952 // To make the trigger tracks from Local Trigger
953 AliDebug(1, "Enter MakeTriggerTracks");
957 TClonesArray *localTrigger;
958 TClonesArray *globalTrigger;
959 AliMUONLocalTrigger *locTrg;
960 AliMUONGlobalTrigger *gloTrg;
961 AliMUONTriggerCircuit *circuit;
962 AliMUONTriggerTrack *recTriggerTrack = 0;
964 TTree* treeR = fLoader->TreeR();
966 // Loading MUON subsystem
967 AliMUON * pMUON = (AliMUON *) gAlice->GetDetector("MUON");
969 nTRentries = Int_t(treeR->GetEntries());
971 treeR->GetEvent(0); // only one entry
973 if (!(fMUONData->IsTriggerBranchesInTree())) {
974 AliWarning(Form("Trigger information is not avalaible, nTRentries = %d not equal to 1",nTRentries));
978 fMUONData->SetTreeAddress("TC");
979 fMUONData->GetTrigger();
981 // global trigger for trigger pattern
983 globalTrigger = fMUONData->GlobalTrigger();
984 gloTrg = (AliMUONGlobalTrigger*)globalTrigger->UncheckedAt(0);
986 if (gloTrg->SinglePlusLpt()) gloTrigPat|= 0x1;
987 if (gloTrg->SinglePlusHpt()) gloTrigPat|= 0x2;
988 if (gloTrg->SinglePlusApt()) gloTrigPat|= 0x4;
990 if (gloTrg->SingleMinusLpt()) gloTrigPat|= 0x8;
991 if (gloTrg->SingleMinusHpt()) gloTrigPat|= 0x10;
992 if (gloTrg->SingleMinusApt()) gloTrigPat|= 0x20;
994 if (gloTrg->SingleUndefLpt()) gloTrigPat|= 0x40;
995 if (gloTrg->SingleUndefHpt()) gloTrigPat|= 0x80;
996 if (gloTrg->SingleUndefApt()) gloTrigPat|= 0x100;
998 if (gloTrg->PairUnlikeLpt()) gloTrigPat|= 0x200;
999 if (gloTrg->PairUnlikeHpt()) gloTrigPat|= 0x400;
1000 if (gloTrg->PairUnlikeApt()) gloTrigPat|= 0x800;
1002 if (gloTrg->PairLikeLpt()) gloTrigPat|= 0x1000;
1003 if (gloTrg->PairLikeHpt()) gloTrigPat|= 0x2000;
1004 if (gloTrg->PairLikeApt()) gloTrigPat|= 0x4000;
1007 // local trigger for tracking
1008 localTrigger = fMUONData->LocalTrigger();
1009 Int_t nlocals = (Int_t) (localTrigger->GetEntries());
1011 Float_t z11 = AliMUONConstants::DefaultChamberZ(10);
1012 Float_t z21 = AliMUONConstants::DefaultChamberZ(12);
1014 for (Int_t i=0; i<nlocals; i++) { // loop on Local Trigger
1015 locTrg = (AliMUONLocalTrigger*)localTrigger->UncheckedAt(i);
1016 circuit = &(pMUON->TriggerCircuit(locTrg->LoCircuit()));
1017 Float_t y11 = circuit->GetY11Pos(locTrg->LoStripX());
1018 Int_t stripX21 = locTrg->LoStripX()+locTrg->LoDev()+1;
1019 Float_t y21 = circuit->GetY21Pos(stripX21);
1020 Float_t x11 = circuit->GetX11Pos(locTrg->LoStripY());
1021 Float_t thetax = TMath::ATan2( x11 , z11 );
1022 Float_t thetay = TMath::ATan2( (y21-y11) , (z21-z11) );
1024 recTriggerTrack = new AliMUONTriggerTrack(x11,y11,thetax,thetay,gloTrigPat);
1026 // since static statement does not work, set gloTrigPat for each track
1028 fMUONData->AddRecTriggerTrack(*recTriggerTrack);
1029 delete recTriggerTrack;
1030 } // end of loop on Local Trigger
1034 //__________________________________________________________________________
1035 Int_t AliMUONTrackReconstructor::MakeTrackCandidatesWithTwoSegments(AliMUONSegment *BegSegment)
1037 // To make track candidates with two segments in stations(1..) 4 and 5,
1038 // the first segment being pointed to by "BegSegment".
1039 // Returns the number of such track candidates.
1040 Int_t endStation, iEndSegment, nbCan2Seg;
1041 AliMUONSegment *endSegment;
1042 AliMUONSegment *extrapSegment = NULL;
1043 AliMUONTrack *recTrack;
1045 AliDebug(1,"Enter MakeTrackCandidatesWithTwoSegments");
1046 // Station for the end segment
1047 endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
1048 // multiple scattering factor corresponding to one chamber
1049 mcsFactor = 0.0136 /
1050 GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
1051 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1052 // linear extrapolation to end station
1053 // number of candidates with 2 segments to 0
1055 // Loop over segments in the end station
1056 for (iEndSegment = 0; iEndSegment < fNSegments[endStation]; iEndSegment++) {
1057 // pointer to segment
1058 endSegment = (AliMUONSegment*) ((*fSegmentsPtr[endStation])[iEndSegment]);
1059 // test compatibility between current segment and "extrapSegment"
1060 // 4 because 4 quantities in chi2
1062 BegSegment->CreateSegmentFromLinearExtrapToStation(endSegment->GetZ(), mcsFactor);
1064 NormalizedChi2WithSegment(extrapSegment,
1065 fMaxSigma2Distance)) <= 4.0) {
1066 // both segments compatible:
1067 // make track candidate from "begSegment" and "endSegment"
1068 AliDebug(2,Form("TrackCandidate with Segment %d in Station(0..) %d", iEndSegment, endStation));
1069 // flag for both segments in one track:
1070 // to be done in track constructor ????
1071 BegSegment->SetInTrack(kTRUE);
1072 endSegment->SetInTrack(kTRUE);
1073 recTrack = new ((*fRecTracksPtr)[fNRecTracks])
1074 AliMUONTrack(BegSegment, endSegment, this);
1076 if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump();
1077 // increment number of track candidates with 2 segments
1080 delete extrapSegment; // should not delete HitForRec's it points to !!!!
1081 } // for (iEndSegment = 0;...
1085 //__________________________________________________________________________
1086 Int_t AliMUONTrackReconstructor::MakeTrackCandidatesWithOneSegmentAndOnePoint(AliMUONSegment *BegSegment)
1088 // To make track candidates with one segment and one point
1089 // in stations(1..) 4 and 5,
1090 // the segment being pointed to by "BegSegment".
1091 Int_t ch, ch1, ch2, endStation, iHit, iHitMax, iHitMin, nbCan1Seg1Hit;
1092 AliMUONHitForRec *extrapHitForRec= NULL;
1093 AliMUONHitForRec *hit;
1094 AliMUONTrack *recTrack;
1096 AliDebug(1,"Enter MakeTrackCandidatesWithOneSegmentAndOnePoint");
1097 // station for the end point
1098 endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
1099 // multiple scattering factor corresponding to one chamber
1100 mcsFactor = 0.0136 /
1101 GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
1102 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1103 // first and second chambers(0..) in the end station
1104 ch1 = 2 * endStation;
1106 // number of candidates to 0
1108 // Loop over chambers of the end station
1109 for (ch = ch2; ch >= ch1; ch--) {
1110 // limits for the hit index in the loop
1111 iHitMin = fIndexOfFirstHitForRecPerChamber[ch];
1112 iHitMax = iHitMin + fNHitsForRecPerChamber[ch];
1113 // Loop over HitForRec's in the chamber
1114 for (iHit = iHitMin; iHit < iHitMax; iHit++) {
1115 // pointer to HitForRec
1116 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1117 // test compatibility between current HitForRec and "extrapHitForRec"
1118 // 2 because 2 quantities in chi2
1119 // linear extrapolation to chamber
1121 BegSegment->CreateHitForRecFromLinearExtrapToChamber( hit->GetZ(), mcsFactor);
1123 NormalizedChi2WithHitForRec(extrapHitForRec,
1124 fMaxSigma2Distance)) <= 2.0) {
1125 // both HitForRec's compatible:
1126 // make track candidate from begSegment and current HitForRec
1127 AliDebug(2, Form("TrackCandidate with HitForRec %d in Chamber(0..) %d", iHit, ch));
1128 // flag for beginning segments in one track:
1129 // to be done in track constructor ????
1130 BegSegment->SetInTrack(kTRUE);
1131 recTrack = new ((*fRecTracksPtr)[fNRecTracks])
1132 AliMUONTrack(BegSegment, hit, this);
1133 // the right place to eliminate "double counting" ???? how ????
1135 if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump();
1136 // increment number of track candidates
1139 delete extrapHitForRec;
1140 } // for (iHit = iHitMin;...
1141 } // for (ch = ch2;...
1142 return nbCan1Seg1Hit;
1145 //__________________________________________________________________________
1146 void AliMUONTrackReconstructor::MakeTrackCandidates(void)
1148 // To make track candidates
1149 // with at least 3 aligned points in stations(1..) 4 and 5
1150 // (two Segment's or one Segment and one HitForRec)
1151 Int_t begStation, iBegSegment, nbCan1Seg1Hit, nbCan2Seg;
1152 AliMUONSegment *begSegment;
1153 AliDebug(1,"Enter MakeTrackCandidates");
1154 // Loop over stations(1..) 5 and 4 for the beginning segment
1155 for (begStation = 4; begStation > 2; begStation--) {
1156 // Loop over segments in the beginning station
1157 for (iBegSegment = 0; iBegSegment < fNSegments[begStation]; iBegSegment++) {
1158 // pointer to segment
1159 begSegment = (AliMUONSegment*) ((*fSegmentsPtr[begStation])[iBegSegment]);
1160 AliDebug(2,Form("Look for TrackCandidate's with Segment %d in Station(0..) %d", iBegSegment, begStation));
1161 // Look for track candidates with two segments,
1162 // "begSegment" and all compatible segments in other station.
1163 // Only for beginning station(1..) 5
1164 // because candidates with 2 segments have to looked for only once.
1165 if (begStation == 4)
1166 nbCan2Seg = MakeTrackCandidatesWithTwoSegments(begSegment);
1167 // Look for track candidates with one segment and one point,
1168 // "begSegment" and all compatible HitForRec's in other station.
1169 // Only if "begSegment" does not belong already to a track candidate.
1170 // Is that a too strong condition ????
1171 if (!(begSegment->GetInTrack()))
1172 nbCan1Seg1Hit = MakeTrackCandidatesWithOneSegmentAndOnePoint(begSegment);
1173 } // for (iBegSegment = 0;...
1174 } // for (begStation = 4;...
1178 //__________________________________________________________________________
1179 void AliMUONTrackReconstructor::FollowTracks(void)
1181 // Follow tracks in stations(1..) 3, 2 and 1
1182 // too long: should be made more modular !!!!
1183 AliMUONHitForRec *bestHit, *extrapHit, *hit;
1184 AliMUONSegment *bestSegment, *extrapSegment, *segment;
1185 AliMUONTrack *track, *nextTrack;
1186 AliMUONTrackParam *trackParam1, trackParam[2], trackParamVertex;
1187 // -1 to avoid compilation warnings
1188 Int_t ch = -1, chInStation, chBestHit = -1, iHit, iSegment, station, trackIndex;
1189 Double_t bestChi2, chi2, dZ1, dZ2, dZ3, maxSigma2Distance, mcsFactor;
1190 Double_t bendingMomentum, chi2Norm = 0.;
1191 // local maxSigma2Distance, for easy increase in testing
1192 maxSigma2Distance = fMaxSigma2Distance;
1193 AliDebug(2,"Enter FollowTracks");
1194 // Loop over track candidates
1195 track = (AliMUONTrack*) fRecTracksPtr->First();
1198 // Follow function for each track candidate ????
1200 nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track); // prepare next track
1201 AliDebug(2,Form("FollowTracks: track candidate(0..): %d", trackIndex));
1202 // Fit track candidate
1203 track->SetFitMCS(0); // without multiple Coulomb scattering
1204 track->SetFitNParam(3); // with 3 parameters (X = Y = 0)
1205 track->SetFitStart(0); // from parameters at vertex
1207 if (AliLog::GetGlobalDebugLevel()> 2) {
1208 cout << "FollowTracks: track candidate(0..): " << trackIndex
1209 << " after fit in stations(0..) 3 and 4" << endl;
1210 track->RecursiveDump();
1212 // Loop over stations(1..) 3, 2 and 1
1213 // something SPECIAL for stations 2 and 1 for majority 3 coincidence ????
1214 // otherwise: majority coincidence 2 !!!!
1215 for (station = 2; station >= 0; station--) {
1216 // Track parameters at first track hit (smallest Z)
1217 trackParam1 = ((AliMUONTrackHit*)
1218 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1219 // extrapolation to station
1220 trackParam1->ExtrapToStation(station, trackParam);
1221 extrapSegment = new AliMUONSegment(); // empty segment
1222 // multiple scattering factor corresponding to one chamber
1223 // and momentum in bending plane (not total)
1224 mcsFactor = 0.0136 * trackParam1->GetInverseBendingMomentum();
1225 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1226 // Z difference from previous station
1227 dZ1 = AliMUONConstants::DefaultChamberZ(2 * station) -
1228 AliMUONConstants::DefaultChamberZ(2 * station + 2);
1229 // Z difference between the two previous stations
1230 dZ2 = AliMUONConstants::DefaultChamberZ(2 * station + 2) -
1231 AliMUONConstants::DefaultChamberZ(2 * station + 4);
1232 // Z difference between the two chambers in the previous station
1233 dZ3 = AliMUONConstants::DefaultChamberZ(2 * station) -
1234 AliMUONConstants::DefaultChamberZ(2 * station + 1);
1235 extrapSegment->SetBendingCoorReso2(fBendingResolution * fBendingResolution);
1237 SetNonBendingCoorReso2(fNonBendingResolution * fNonBendingResolution);
1238 extrapSegment->UpdateFromStationTrackParam
1239 (trackParam, mcsFactor, dZ1, dZ2, dZ3, station,
1240 trackParam1->GetInverseBendingMomentum());
1243 if (AliLog::GetGlobalDebugLevel() > 2) {
1244 cout << "FollowTracks: track candidate(0..): " << trackIndex
1245 << " Look for segment in station(0..): " << station << endl;
1247 // Loop over segments in station
1248 for (iSegment = 0; iSegment < fNSegments[station]; iSegment++) {
1249 // Look for best compatible Segment in station
1250 // should consider all possibilities ????
1251 // multiple scattering ????
1252 // separation in 2 functions: Segment and HitForRec ????
1253 segment = (AliMUONSegment*) ((*fSegmentsPtr[station])[iSegment]);
1254 // correction of corrected segment (fBendingCoor and fNonBendingCoor)
1255 // according to real Z value of "segment" and slopes of "extrapSegment"
1256 (&(trackParam[0]))->ExtrapToZ(segment->GetZ());
1257 (&(trackParam[1]))->ExtrapToZ(segment->GetZ());
1258 extrapSegment->SetBendingCoor((&(trackParam[0]))->GetBendingCoor());
1259 extrapSegment->SetNonBendingCoor((&(trackParam[0]))->GetNonBendingCoor());
1260 extrapSegment->SetBendingSlope((&(trackParam[0]))->GetBendingSlope());
1261 extrapSegment->SetNonBendingSlope((&(trackParam[0]))->GetNonBendingSlope());
1263 NormalizedChi2WithSegment(extrapSegment, maxSigma2Distance);
1264 if (chi2 < bestChi2) {
1265 // update best Chi2 and Segment if better found
1266 bestSegment = segment;
1271 // best segment found: add it to track candidate
1272 (&(trackParam[0]))->ExtrapToZ(bestSegment->GetZ());
1273 (&(trackParam[1]))->ExtrapToZ(bestSegment->GetZ());
1274 track->AddSegment(bestSegment);
1275 // set track parameters at these two TrakHit's
1276 track->SetTrackParamAtHit(track->GetNTrackHits() - 2, &(trackParam[0]));
1277 track->SetTrackParamAtHit(track->GetNTrackHits() - 1, &(trackParam[1]));
1278 AliDebug(3, Form("FollowTracks: track candidate(0..): %d Added segment in station(0..): %d", trackIndex, station));
1279 if (AliLog::GetGlobalDebugLevel()>2) track->RecursiveDump();
1282 // No best segment found:
1283 // Look for best compatible HitForRec in station:
1284 // should consider all possibilities ????
1285 // multiple scattering ???? do about like for extrapSegment !!!!
1286 extrapHit = new AliMUONHitForRec(); // empty hit
1289 AliDebug(3, Form("FollowTracks: track candidate(0..): %d Segment not found, look for hit in station(0..): %d ",
1290 trackIndex, station));
1292 // Loop over chambers of the station
1293 for (chInStation = 0; chInStation < 2; chInStation++) {
1294 ch = 2 * station + chInStation;
1295 for (iHit = fIndexOfFirstHitForRecPerChamber[ch];
1296 iHit < fIndexOfFirstHitForRecPerChamber[ch] +
1297 fNHitsForRecPerChamber[ch];
1299 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1300 // coordinates of extrapolated hit
1301 (&(trackParam[chInStation]))->ExtrapToZ(hit->GetZ());
1303 SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor());
1305 SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor());
1306 // resolutions from "extrapSegment"
1307 extrapHit->SetBendingReso2(extrapSegment->GetBendingCoorReso2());
1308 extrapHit->SetNonBendingReso2(extrapSegment->GetNonBendingCoorReso2());
1309 // Loop over hits in the chamber
1310 // condition for hit not already in segment ????
1311 chi2 = hit->NormalizedChi2WithHitForRec(extrapHit, maxSigma2Distance);
1312 if (chi2 < bestChi2) {
1313 // update best Chi2 and HitForRec if better found
1316 chBestHit = chInStation;
1321 // best hit found: add it to track candidate
1322 (&(trackParam[chBestHit]))->ExtrapToZ(bestHit->GetZ());
1323 track->AddHitForRec(bestHit);
1324 // set track parameters at this TrackHit
1325 track->SetTrackParamAtHit(track->GetNTrackHits() - 1,
1326 &(trackParam[chBestHit]));
1327 if (AliLog::GetGlobalDebugLevel() > 2) {
1328 cout << "FollowTracks: track candidate(0..): " << trackIndex
1329 << " Added hit in station(0..): " << station << endl;
1330 track->RecursiveDump();
1334 // Remove current track candidate
1335 // and corresponding TrackHit's, ...
1337 delete extrapSegment;
1339 break; // stop the search for this candidate:
1340 // exit from the loop over station
1344 delete extrapSegment;
1345 // Sort track hits according to increasing Z
1346 track->GetTrackHitsPtr()->Sort();
1347 // Update track parameters at first track hit (smallest Z)
1348 trackParam1 = ((AliMUONTrackHit*)
1349 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1350 bendingMomentum = 0.;
1351 if (TMath::Abs(trackParam1->GetInverseBendingMomentum()) > 0.)
1352 bendingMomentum = TMath::Abs(1/(trackParam1->GetInverseBendingMomentum()));
1353 // Track removed if bendingMomentum not in window [min, max]
1354 if ((bendingMomentum < fMinBendingMomentum) || (bendingMomentum > fMaxBendingMomentum)) {
1356 break; // stop the search for this candidate:
1357 // exit from the loop over station
1360 // with multiple Coulomb scattering if all stations
1361 if (station == 0) track->SetFitMCS(1);
1362 // without multiple Coulomb scattering if not all stations
1363 else track->SetFitMCS(0);
1364 track->SetFitNParam(5); // with 5 parameters (momentum and position)
1365 track->SetFitStart(1); // from parameters at first hit
1367 Double_t numberOfDegFree = (2.0 * track->GetNTrackHits() - 5);
1368 if (numberOfDegFree > 0) {
1369 chi2Norm = track->GetFitFMin() / numberOfDegFree;
1373 // Track removed if normalized chi2 too high
1374 if (chi2Norm > fMaxChi2) {
1376 break; // stop the search for this candidate:
1377 // exit from the loop over station
1379 if (AliLog::GetGlobalDebugLevel() > 2) {
1380 cout << "FollowTracks: track candidate(0..): " << trackIndex
1381 << " after fit from station(0..): " << station << " to 4" << endl;
1382 track->RecursiveDump();
1384 // Track extrapolation to the vertex through the absorber (Branson)
1385 // after going through the first station
1387 trackParamVertex = *trackParam1;
1388 (&trackParamVertex)->ExtrapToVertex(0.,0.,0.);
1389 track->SetTrackParamAtVertex(&trackParamVertex);
1390 if (AliLog::GetGlobalDebugLevel() > 0) {
1391 cout << "FollowTracks: track candidate(0..): " << trackIndex
1392 << " after extrapolation to vertex" << endl;
1393 track->RecursiveDump();
1396 } // for (station = 2;...
1397 // go really to next track
1400 // Compression of track array (necessary after Remove ????)
1401 fRecTracksPtr->Compress();
1405 //__________________________________________________________________________
1406 void AliMUONTrackReconstructor::RemoveDoubleTracks(void)
1408 // To remove double tracks.
1409 // Tracks are considered identical
1410 // if they have at least half of their hits in common.
1411 // Among two identical tracks, one keeps the track with the larger number of hits
1412 // or, if these numbers are equal, the track with the minimum Chi2.
1413 AliMUONTrack *track1, *track2, *trackToRemove;
1414 Bool_t identicalTracks;
1415 Int_t hitsInCommon, nHits1, nHits2;
1416 identicalTracks = kTRUE;
1417 while (identicalTracks) {
1418 identicalTracks = kFALSE;
1419 // Loop over first track of the pair
1420 track1 = (AliMUONTrack*) fRecTracksPtr->First();
1421 while (track1 && (!identicalTracks)) {
1422 nHits1 = track1->GetNTrackHits();
1423 // Loop over second track of the pair
1424 track2 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1425 while (track2 && (!identicalTracks)) {
1426 nHits2 = track2->GetNTrackHits();
1427 // number of hits in common between two tracks
1428 hitsInCommon = track1->HitsInCommon(track2);
1429 // check for identical tracks
1430 if ((4 * hitsInCommon) >= (nHits1 + nHits2)) {
1431 identicalTracks = kTRUE;
1432 // decide which track to remove
1433 if (nHits1 > nHits2) trackToRemove = track2;
1434 else if (nHits1 < nHits2) trackToRemove = track1;
1435 else if ((track1->GetFitFMin()) < (track2->GetFitFMin()))
1436 trackToRemove = track2;
1437 else trackToRemove = track1;
1439 trackToRemove->Remove();
1441 track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
1443 track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1449 //__________________________________________________________________________
1450 void AliMUONTrackReconstructor::UpdateTrackParamAtHit()
1452 // Set track parameters after track fitting. Fill fTrackParamAtHit of AliMUONTrack's
1453 AliMUONTrack *track;
1454 AliMUONTrackHit *trackHit;
1455 AliMUONTrackParam *trackParam;
1456 track = (AliMUONTrack*) fRecTracksPtr->First();
1458 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->First();
1460 trackParam = trackHit->GetTrackParam();
1461 track->AddTrackParamAtHit(trackParam);
1462 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->After(trackHit);
1464 track = (AliMUONTrack*) fRecTracksPtr->After(track);
1469 //__________________________________________________________________________
1470 void AliMUONTrackReconstructor::UpdateHitForRecAtHit()
1472 // Set cluster parameterss after track fitting. Fill fHitForRecAtHit of AliMUONTrack's
1473 AliMUONTrack *track;
1474 AliMUONTrackHit *trackHit;
1475 AliMUONHitForRec *hitForRec;
1476 track = (AliMUONTrack*) fRecTracksPtr->First();
1478 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->First();
1480 hitForRec = trackHit->GetHitForRecPtr();
1481 track->AddHitForRecAtHit(hitForRec);
1482 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->After(trackHit);
1484 track = (AliMUONTrack*) fRecTracksPtr->After(track);
1489 //__________________________________________________________________________
1490 void AliMUONTrackReconstructor::FillMUONTrack()
1492 // Set track parameters at hits for Kalman track. Fill fTrackParamAtHit of AliMUONTrack's
1493 AliMUONTrackK *track;
1494 track = (AliMUONTrackK*) fRecTracksPtr->First();
1496 track->FillMUONTrack();
1497 track = (AliMUONTrackK*) fRecTracksPtr->After(track);
1502 //__________________________________________________________________________
1503 void AliMUONTrackReconstructor::EventDump(void)
1505 // Dump reconstructed event (track parameters at vertex and at first hit),
1506 // and the particle parameters
1508 AliMUONTrack *track;
1509 AliMUONTrackParam *trackParam, *trackParam1;
1510 Double_t bendingSlope, nonBendingSlope, pYZ;
1511 Double_t pX, pY, pZ, x, y, z, c;
1512 Int_t np, trackIndex, nTrackHits;
1514 AliDebug(1,"****** enter EventDump ******");
1515 AliDebug(1, Form("Number of Reconstructed tracks : %d", fNRecTracks));
1517 fRecTracksPtr->Compress(); // for simple loop without "Next" since no hole
1518 // Loop over reconstructed tracks
1519 for (trackIndex = 0; trackIndex < fNRecTracks; trackIndex++) {
1520 if (fTrackMethod != 1) continue; //AZ - skip the rest for now
1521 AliDebug(1, Form("track number: %d", trackIndex));
1522 // function for each track for modularity ????
1523 track = (AliMUONTrack*) ((*fRecTracksPtr)[trackIndex]);
1524 nTrackHits = track->GetNTrackHits();
1525 AliDebug(1, Form("Number of track hits: %d ", nTrackHits));
1526 // track parameters at Vertex
1527 trackParam = track->GetTrackParamAtVertex();
1528 x = trackParam->GetNonBendingCoor();
1529 y = trackParam->GetBendingCoor();
1530 z = trackParam->GetZ();
1531 bendingSlope = trackParam->GetBendingSlope();
1532 nonBendingSlope = trackParam->GetNonBendingSlope();
1533 pYZ = 1/TMath::Abs(trackParam->GetInverseBendingMomentum());
1534 pZ = pYZ/TMath::Sqrt(1+bendingSlope*bendingSlope);
1535 pX = pZ * nonBendingSlope;
1536 pY = pZ * bendingSlope;
1537 c = TMath::Sign(1.0, trackParam->GetInverseBendingMomentum());
1538 AliDebug(1, Form("Track parameters at Vertex z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1539 z, x, y, pX, pY, pZ, c));
1541 // track parameters at first hit
1542 trackParam1 = ((AliMUONTrackHit*)
1543 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1544 x = trackParam1->GetNonBendingCoor();
1545 y = trackParam1->GetBendingCoor();
1546 z = trackParam1->GetZ();
1547 bendingSlope = trackParam1->GetBendingSlope();
1548 nonBendingSlope = trackParam1->GetNonBendingSlope();
1549 pYZ = 1/TMath::Abs(trackParam1->GetInverseBendingMomentum());
1550 pZ = pYZ/TMath::Sqrt(1.0 + bendingSlope * bendingSlope);
1551 pX = pZ * nonBendingSlope;
1552 pY = pZ * bendingSlope;
1553 c = TMath::Sign(1.0, trackParam1->GetInverseBendingMomentum());
1554 AliDebug(1, Form("track parameters at z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1555 z, x, y, pX, pY, pZ, c));
1557 // informations about generated particles
1558 np = gAlice->GetMCApp()->GetNtrack();
1559 printf(" **** number of generated particles: %d \n", np);
1561 // for (Int_t iPart = 0; iPart < np; iPart++) {
1562 // p = gAlice->Particle(iPart);
1563 // printf(" particle %d: type= %d px= %f py= %f pz= %f pdg= %d\n",
1564 // iPart, p->GetPdgCode(), p->Px(), p->Py(), p->Pz(), p->GetPdgCode());
1570 //__________________________________________________________________________
1571 void AliMUONTrackReconstructor::EventDumpTrigger(void)
1573 // Dump reconstructed trigger event
1574 // and the particle parameters
1576 AliMUONTriggerTrack *triggertrack ;
1577 Int_t nTriggerTracks = fMUONData->RecTriggerTracks()->GetEntriesFast();
1579 AliDebug(1, "****** enter EventDumpTrigger ******");
1580 AliDebug(1, Form("Number of Reconstructed tracks : %d ", nTriggerTracks));
1582 // Loop over reconstructed tracks
1583 for (Int_t trackIndex = 0; trackIndex < nTriggerTracks; trackIndex++) {
1584 triggertrack = (AliMUONTriggerTrack*)fMUONData->RecTriggerTracks()->At(trackIndex);
1585 printf(" trigger track number %i x11=%f y11=%f thetax=%f thetay=%f \n",
1587 triggertrack->GetX11(),triggertrack->GetY11(),
1588 triggertrack->GetThetax(),triggertrack->GetThetay());
1592 //__________________________________________________________________________
1593 void AliMUONTrackReconstructor::MakeTrackCandidatesK(void)
1595 // To make initial tracks for Kalman filter from the list of segments
1597 AliMUONSegment *segment;
1598 AliMUONTrackK *trackK;
1600 AliDebug(1,"Enter MakeTrackCandidatesK");
1601 // Reset the TClonesArray of reconstructed tracks
1602 if (fRecTracksPtr) fRecTracksPtr->Delete();
1603 // Delete in order that the Track destructors are called,
1604 // hence the space for the TClonesArray of pointers to TrackHit's is freed
1607 AliMUONTrackK a(this, fHitsForRecPtr); // bad idea ???
1608 // Loop over stations(1...) 5 and 4
1609 for (istat=4; istat>=3; istat--) {
1610 // Loop over segments in the station
1611 for (iseg=0; iseg<fNSegments[istat]; iseg++) {
1612 // Transform segments to tracks and evaluate covariance matrix
1613 segment = (AliMUONSegment*) ((*fSegmentsPtr[istat])[iseg]);
1614 trackK = new ((*fRecTracksPtr)[fNRecTracks]) AliMUONTrackK(segment);
1616 } // for (iseg=0;...)
1617 } // for (istat=4;...)
1621 //__________________________________________________________________________
1622 void AliMUONTrackReconstructor::FollowTracksK(void)
1624 // Follow tracks using Kalman filter
1626 Int_t icand, ichamBeg = 0, ichamEnd, chamBits;
1627 Double_t zDipole1, zDipole2;
1628 AliMUONTrackK *trackK;
1629 AliMUONHitForRec *hit;
1630 AliMUONRawCluster *clus;
1631 TClonesArray *rawclusters;
1632 clus = 0; rawclusters = 0;
1634 //AZ(z->-z) zDipole1 = GetSimpleBPosition() - GetSimpleBLength()/2;
1635 //AZ(z->-z) zDipole2 = zDipole1 + GetSimpleBLength();
1636 zDipole1 = GetSimpleBPosition() + GetSimpleBLength()/2;
1637 zDipole2 = zDipole1 - GetSimpleBLength();
1640 trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[0]);
1641 if (trackK->DebugLevel() > 0) {
1642 for (Int_t i1=0; i1<fNHitsForRec; i1++) {
1643 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[i1]);
1644 //if (hit->GetTTRTrack() > 1 || hit->GetTrackRefSignal() == 0) continue;
1645 printf(" Hit # %d %10.4f %10.4f %10.4f",
1646 hit->GetChamberNumber(), hit->GetBendingCoor(),
1647 hit->GetNonBendingCoor(), hit->GetZ());
1648 if (fRecTrackRefHits) {
1649 // from track ref hits
1650 printf(" %3d %3d \n", hit->GetTrackRefSignal(), hit->GetTTRTrack());
1652 // from raw clusters
1653 rawclusters = fMUONData->RawClusters(hit->GetChamberNumber());
1654 clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(hit->
1656 printf("%3d", clus->GetTrack(1)-1);
1657 if (clus->GetTrack(2) != 0) printf("%3d \n", clus->GetTrack(2)-1);
1659 } // if (fRecTrackRefHits)
1661 } // if (trackK->DebugLevel() > 0)
1665 nSeeds = fNRecTracks; // starting number of seeds
1666 // Loop over track candidates
1667 while (icand < fNRecTracks-1) {
1669 if (trackK->DebugLevel()>0) cout << " *** Kalman track candidate No. " << icand << endl;
1670 trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]);
1671 if (trackK->GetRecover() < 0) continue; // failed track
1673 // Discard candidate which will produce the double track
1676 ok = CheckCandidateK(icand,nSeeds);
1678 trackK->SetRecover(-1); // mark candidate to be removed
1685 if (trackK->GetRecover() == 0) hit = (AliMUONHitForRec*)
1686 trackK->GetHitOnTrack()->Last(); // last hit
1687 //else hit = (AliMUONHitForRec*) (*trackK->GetHitOnTrack())[1]; // 2'nd hit
1688 else hit = trackK->GetHitLastOk(); // hit where track stopped
1689 if (hit) ichamBeg = hit->GetChamberNumber();
1691 // Check propagation direction
1692 if (hit == NULL) ichamBeg = ichamEnd;
1693 //AZ(z->-z) else if (trackK->GetTrackDir() > 0) {
1694 else if (trackK->GetTrackDir() < 0) {
1695 ichamEnd = 9; // forward propagation
1696 ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kFALSE,zDipole1,zDipole2);
1698 ichamBeg = ichamEnd;
1699 ichamEnd = 6; // backward propagation
1700 // Change weight matrix and zero fChi2 for backpropagation
1701 trackK->StartBack();
1702 //AZ trackK->SetTrackDir(-1);
1703 trackK->SetTrackDir(1);
1704 ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kTRUE,zDipole1,zDipole2);
1705 ichamBeg = ichamEnd;
1709 if (trackK->GetBPFlag()) {
1711 ichamEnd = 6; // backward propagation
1712 // Change weight matrix and zero fChi2 for backpropagation
1713 trackK->StartBack();
1714 ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kTRUE,zDipole1,zDipole2);
1715 ichamBeg = ichamEnd;
1721 //AZ trackK->SetTrackDir(-1);
1722 trackK->SetTrackDir(1);
1723 trackK->SetBPFlag(kFALSE);
1724 ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kFALSE,zDipole1,zDipole2);
1726 if (!ok) { trackK->SetRecover(-1); continue; } // mark candidate to be removed
1729 if (trackK->GetRecover() >= 0) {
1730 ok = trackK->Smooth();
1731 if (!ok) trackK->SetRecover(-1); // mark candidate to be removed
1734 // Majority 3 of 4 in first 2 stations
1737 Double_t chi2max = 0;
1738 for (Int_t i=0; i<trackK->GetNTrackHits(); i++) {
1739 hit = (AliMUONHitForRec*) (*trackK->GetHitOnTrack())[i];
1740 chamBits |= BIT(hit->GetChamberNumber());
1741 if (trackK->GetChi2PerPoint(i) > chi2max) chi2max = trackK->GetChi2PerPoint(i);
1743 if (!((chamBits&3)==3 || (chamBits>>2&3)==3) && chi2max > 25) {
1744 //trackK->Recover();
1745 trackK->SetRecover(-1); //mark candidate to be removed
1748 if (ok) trackK->SetTrackQuality(0); // compute "track quality"
1751 for (Int_t i=0; i<fNRecTracks; i++) {
1752 trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[i]);
1753 if (trackK->GetRecover() < 0) fRecTracksPtr->RemoveAt(i);
1756 // Compress TClonesArray
1757 fRecTracksPtr->Compress();
1758 fNRecTracks = fRecTracksPtr->GetEntriesFast();
1762 //__________________________________________________________________________
1763 Bool_t AliMUONTrackReconstructor::CheckCandidateK(Int_t icand, Int_t nSeeds) const
1765 // Discards track candidate if it will produce the double track (having
1766 // the same seed segment hits as hits of a good track found before)
1767 AliMUONTrackK *track1, *track2;
1768 AliMUONHitForRec *hit1, *hit2, *hit;
1770 track1 = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]);
1771 hit1 = (AliMUONHitForRec*) (*track1->GetHitOnTrack())[0]; // 1'st hit
1772 hit2 = (AliMUONHitForRec*) (*track1->GetHitOnTrack())[1]; // 2'nd hit
1774 for (Int_t i=0; i<icand; i++) {
1775 track2 = (AliMUONTrackK*) ((*fRecTracksPtr)[i]);
1776 //if (track2->GetRecover() < 0) continue;
1777 if (track2->GetRecover() < 0 && icand >= nSeeds) continue;
1779 if (track1->GetStartSegment() == track2->GetStartSegment()) {
1783 for (Int_t j=0; j<track2->GetNTrackHits(); j++) {
1784 hit = (AliMUONHitForRec*) (*track2->GetHitOnTrack())[j];
1785 if (hit == hit1 || hit == hit2) {
1787 if (nSame == 2) return kFALSE;
1789 } // for (Int_t j=0;
1791 } // for (Int_t i=0;
1795 //__________________________________________________________________________
1796 void AliMUONTrackReconstructor::RemoveDoubleTracksK(void)
1798 // Removes double tracks (sharing more than half of their hits). Keeps
1799 // the track with higher quality
1800 AliMUONTrackK *track1, *track2, *trackToKill;
1802 // Sort tracks according to their quality
1803 fRecTracksPtr->Sort();
1805 // Loop over first track of the pair
1806 track1 = (AliMUONTrackK*) fRecTracksPtr->First();
1807 Int_t debug = track1->DebugLevel();
1809 // Loop over second track of the pair
1810 track2 = (AliMUONTrackK*) fRecTracksPtr->After(track1);
1812 // Check whether or not to keep track2
1813 if (!track2->KeepTrack(track1)) {
1814 if (debug >= 0) cout << " Killed track: " << 1/(*track2->GetTrackParameters())(4,0) <<
1815 " " << track2->GetTrackQuality() << endl;
1816 trackToKill = track2;
1817 track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2);
1818 trackToKill->Kill();
1819 fRecTracksPtr->Compress();
1820 } else track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2);
1822 track1 = (AliMUONTrackK*) fRecTracksPtr->After(track1);
1825 fNRecTracks = fRecTracksPtr->GetEntriesFast();
1826 if (debug >= 0) cout << " Number of Kalman tracks: " << fNRecTracks << endl;
1829 //__________________________________________________________________________
1830 void AliMUONTrackReconstructor::GoToVertex(void)
1832 // Propagates track to the vertex thru absorber
1833 // (using Branson correction for now)
1837 for (Int_t i=0; i<fNRecTracks; i++) {
1838 //((AliMUONTrackK*)(*fRecTracksPtr)[i])->Branson();
1839 ((AliMUONTrackK*)(*fRecTracksPtr)[i])->SetTrackQuality(1); // compute Chi2
1840 //((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToZ(zVertex); // w/out absorber
1841 ((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToVertex(1); // with absorber
1845 //__________________________________________________________________________
1846 void AliMUONTrackReconstructor::SetTrackMethod(Int_t iTrackMethod)
1848 // Set track method and recreate track container if necessary
1850 fTrackMethod = iTrackMethod == 2 ? 2 : 1;
1851 if (fTrackMethod == 2) {
1852 if (fRecTracksPtr) delete fRecTracksPtr;
1853 fRecTracksPtr = new TClonesArray("AliMUONTrackK", 10);
1854 cout << " *** Tracking with the Kalman filter *** " << endl;
1855 } else cout << " *** Traditional tracking *** " << endl;