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 **************************************************************************/
18 ////////////////////////////////////
20 // MUON track reconstructor in ALICE (class renamed from AliMUONEventReconstructor)
22 // This class contains as data:
23 // * the parameters for the track 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
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
33 ////////////////////////////////////
36 #include <Riostream.h>
37 #include <TDirectory.h>
41 #include "AliMUONTrackReconstructor.h"
43 #include "AliMUONConstants.h"
44 #include "AliMUONHitForRec.h"
45 #include "AliMUONTriggerTrack.h"
46 #include "AliMUONTriggerCircuit.h"
47 #include "AliMUONTriggerCircuitNew.h"
48 #include "AliMUONRawCluster.h"
49 #include "AliMUONLocalTrigger.h"
50 #include "AliMUONGlobalTrigger.h"
51 #include "AliMUONSegment.h"
52 #include "AliMUONTrack.h"
53 #include "AliMUONTrackHit.h"
55 #include "AliRun.h" // for gAlice
56 #include "AliRunLoader.h"
57 #include "AliLoader.h"
58 #include "AliMUONTrackK.h"
61 #include "AliTrackReference.h"
63 //************* Defaults parameters for reconstruction
64 const Double_t AliMUONTrackReconstructor::fgkDefaultMinBendingMomentum = 3.0;
65 const Double_t AliMUONTrackReconstructor::fgkDefaultMaxBendingMomentum = 3000.0;
66 const Double_t AliMUONTrackReconstructor::fgkDefaultMaxChi2 = 100.0;
67 const Double_t AliMUONTrackReconstructor::fgkDefaultMaxSigma2Distance = 16.0;
68 const Double_t AliMUONTrackReconstructor::fgkDefaultBendingResolution = 0.01;
69 const Double_t AliMUONTrackReconstructor::fgkDefaultNonBendingResolution = 0.144;
70 const Double_t AliMUONTrackReconstructor::fgkDefaultChamberThicknessInX0 = 0.03;
71 // Simple magnetic field:
72 // Value taken from macro MUONtracking.C: 0.7 T, hence 7 kG
73 // Length and Position from reco_muon.F, with opposite sign:
74 // Length = ZMAGEND-ZCOIL
75 // Position = (ZMAGEND+ZCOIL)/2
76 // to be ajusted differently from real magnetic field ????
77 const Double_t AliMUONTrackReconstructor::fgkDefaultSimpleBValue = 7.0;
78 const Double_t AliMUONTrackReconstructor::fgkDefaultSimpleBLength = 428.0;
79 const Double_t AliMUONTrackReconstructor::fgkDefaultSimpleBPosition = 1019.0;
80 const Int_t AliMUONTrackReconstructor::fgkDefaultRecTrackRefHits = 0;
81 const Double_t AliMUONTrackReconstructor::fgkDefaultEfficiency = 0.95;
83 ClassImp(AliMUONTrackReconstructor) // Class implementation in ROOT context
85 //__________________________________________________________________________
86 AliMUONTrackReconstructor::AliMUONTrackReconstructor(AliLoader* loader, AliMUONData* data)
89 // Constructor for class AliMUONTrackReconstructor
90 SetReconstructionParametersToDefaults();
91 fTrackMethod = 1; //AZ - tracking method (1-default, 2-Kalman)
92 // Memory allocation for the TClonesArray of hits for reconstruction
93 // Is 10000 the right size ????
94 fHitsForRecPtr = new TClonesArray("AliMUONHitForRec", 10000);
95 fNHitsForRec = 0; // really needed or GetEntriesFast sufficient ????
96 // Memory allocation for the TClonesArray's of segments in stations
97 // Is 2000 the right size ????
98 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) {
99 fSegmentsPtr[st] = new TClonesArray("AliMUONSegment", 2000);
100 fNSegments[st] = 0; // really needed or GetEntriesFast sufficient ????
102 // Memory allocation for the TClonesArray of reconstructed tracks
103 // Is 10 the right size ????
104 fRecTracksPtr = new TClonesArray("AliMUONTrack", 10);
105 fNRecTracks = 0; // really needed or GetEntriesFast sufficient ????
106 // Memory allocation for the TClonesArray of hits on reconstructed tracks
107 // Is 100 the right size ????
108 fRecTrackHitsPtr = new TClonesArray("AliMUONTrackHit", 100);
109 fNRecTrackHits = 0; // really needed or GetEntriesFast sufficient ????
111 // Sign of fSimpleBValue according to sign of Bx value at (50,50,-950).
113 x[0] = 50.; x[1] = 50.; x[2] = -950.;
114 gAlice->Field()->Field(x, b);
115 fSimpleBValue = TMath::Sign(fSimpleBValue,(Double_t) b[0]);
116 fSimpleBPosition = TMath::Sign(fSimpleBPosition,(Double_t) x[2]);
117 // See how to get fSimple(BValue, BLength, BPosition)
118 // automatically calculated from the actual magnetic field ????
120 AliDebug(1,"AliMUONTrackReconstructor constructed with defaults");
121 if ( AliLog::GetGlobalDebugLevel()>0) Dump();
122 AliDebug(1,"Magnetic field from root file:");
123 if ( AliLog::GetGlobalDebugLevel()>0) gAlice->Field()->Dump();
126 // Initializions for track ref. background events
127 fBkgTrackRefFile = 0;
129 fBkgTrackRefParticles = 0;
131 fBkgTrackRefEventNumber = -1;
133 // initialize loader's
136 // initialize container
137 // fMUONData = new AliMUONData(fLoader,"MUON","MUON");
142 //__________________________________________________________________________
143 AliMUONTrackReconstructor::AliMUONTrackReconstructor (const AliMUONTrackReconstructor& rhs)
146 // Protected copy constructor
148 AliFatal("Not implemented.");
151 AliMUONTrackReconstructor &
152 AliMUONTrackReconstructor::operator=(const AliMUONTrackReconstructor& rhs)
154 // Protected assignement operator
156 if (this == &rhs) return *this;
158 AliFatal("Not implemented.");
163 //__________________________________________________________________________
164 AliMUONTrackReconstructor::~AliMUONTrackReconstructor(void)
166 // Destructor for class AliMUONTrackReconstructor
167 delete fHitsForRecPtr; // Correct destruction of everything ???? or delete [] ????
168 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++)
169 delete fSegmentsPtr[st]; // Correct destruction of everything ????
173 //__________________________________________________________________________
174 void AliMUONTrackReconstructor::SetReconstructionParametersToDefaults(void)
176 // Set reconstruction parameters to default values
177 // Would be much more convenient with a structure (or class) ????
178 fMinBendingMomentum = fgkDefaultMinBendingMomentum;
179 fMaxBendingMomentum = fgkDefaultMaxBendingMomentum;
180 fMaxChi2 = fgkDefaultMaxChi2;
181 fMaxSigma2Distance = fgkDefaultMaxSigma2Distance;
183 // ******** Parameters for making HitsForRec
185 // like in TRACKF_STAT:
186 // 2 degrees for stations 1 and 2, or ch(0...) from 0 to 3;
187 // 30 cm for stations 3 to 5, or ch(0...) from 4 to 9
188 for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
189 if (ch < 4) fRMin[ch] = TMath::Abs(AliMUONConstants::DefaultChamberZ(ch)) *
190 2.0 * TMath::Pi() / 180.0;
191 else fRMin[ch] = 30.0;
192 // maximum radius at 10 degrees and Z of chamber
193 fRMax[ch] = TMath::Abs(AliMUONConstants::DefaultChamberZ(ch)) *
194 10.0 * TMath::Pi() / 180.0;
197 // ******** Parameters for making segments
198 // should be parametrized ????
199 // according to interval between chambers in a station ????
200 // Maximum distance in non bending plane
201 // 5 * 0.22 just to remember the way it was made in TRACKF_STAT
203 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++)
204 fSegmentMaxDistNonBending[st] = 5. * 0.22;
205 // Maximum distance in bending plane:
206 // values from TRACKF_STAT, corresponding to (J psi 20cm),
207 // scaled to the real distance between chambers in a station
208 fSegmentMaxDistBending[0] = TMath::Abs( 1.5 *
209 (AliMUONConstants::DefaultChamberZ(1) - AliMUONConstants::DefaultChamberZ(0)) / 20.0);
210 fSegmentMaxDistBending[1] = TMath::Abs( 1.5 *
211 (AliMUONConstants::DefaultChamberZ(3) - AliMUONConstants::DefaultChamberZ(2)) / 20.0);
212 fSegmentMaxDistBending[2] = TMath::Abs( 3.0 *
213 (AliMUONConstants::DefaultChamberZ(5) - AliMUONConstants::DefaultChamberZ(4)) / 20.0);
214 fSegmentMaxDistBending[3] = TMath::Abs( 6.0 *
215 (AliMUONConstants::DefaultChamberZ(7) - AliMUONConstants::DefaultChamberZ(6)) / 20.0);
216 fSegmentMaxDistBending[4] = TMath::Abs( 6.0 *
217 (AliMUONConstants::DefaultChamberZ(9) - AliMUONConstants::DefaultChamberZ(8)) / 20.0);
220 fBendingResolution = fgkDefaultBendingResolution;
221 fNonBendingResolution = fgkDefaultNonBendingResolution;
222 fChamberThicknessInX0 = fgkDefaultChamberThicknessInX0;
223 fSimpleBValue = fgkDefaultSimpleBValue;
224 fSimpleBLength = fgkDefaultSimpleBLength;
225 fSimpleBPosition = fgkDefaultSimpleBPosition;
226 fRecTrackRefHits = fgkDefaultRecTrackRefHits;
227 fEfficiency = fgkDefaultEfficiency;
231 //__________________________________________________________________________
232 Double_t AliMUONTrackReconstructor::GetImpactParamFromBendingMomentum(Double_t BendingMomentum) const
234 // Returns impact parameter at vertex in bending plane (cm),
235 // from the signed bending momentum "BendingMomentum" in bending plane (GeV/c),
236 // using simple values for dipole magnetic field.
237 // The sign of "BendingMomentum" is the sign of the charge.
238 return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
242 //__________________________________________________________________________
243 Double_t AliMUONTrackReconstructor::GetBendingMomentumFromImpactParam(Double_t ImpactParam) const
245 // Returns signed bending momentum in bending plane (GeV/c),
246 // the sign being the sign of the charge for particles moving forward in Z,
247 // from the impact parameter "ImpactParam" at vertex in bending plane (cm),
248 // using simple values for dipole magnetic field.
249 return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
253 //__________________________________________________________________________
254 void AliMUONTrackReconstructor::SetBkgTrackRefFile(Text_t *BkgTrackRefFileName)
256 // Set background file ... for track ref. hits
257 // Must be called after having loaded the firts signal event
258 AliDebug(1,Form("Enter SetBkgTrackRefFile with BkgTrackRefFileName %s",BkgTrackRefFileName));
260 if (strlen(BkgTrackRefFileName)) {
261 // BkgTrackRefFileName not empty: try to open the file
263 if(AliLog::GetGlobalDebugLevel()>1) {
264 cout << "Before File(Bkg)" << endl; gDirectory->Dump();
266 fBkgTrackRefFile = new TFile(BkgTrackRefFileName);
267 if(AliLog::GetGlobalDebugLevel()>1) {
268 cout << "After File(Bkg)" << endl; gDirectory->Dump();
270 if (fBkgTrackRefFile-> IsOpen()) {
271 if(AliLog::GetGlobalDebugLevel()>0) {
272 cout << "Background for Track ref. hits in file: ``" << BkgTrackRefFileName
273 << "'' successfully opened" << endl;}
276 cout << "Background for Track Ref. hits in file: " << BkgTrackRefFileName << endl;
277 cout << "NOT FOUND: EXIT" << endl;
278 exit(0); // right instruction for exit ????
280 // Arrays for "particles" and "hits"
281 fBkgTrackRefParticles = new TClonesArray("TParticle", 200);
282 // Event number to -1 for initialization
283 fBkgTrackRefEventNumber = -1;
284 // Back to the signal file:
285 // first signal event must have been loaded previously,
286 // otherwise, Segmentation violation at the next instruction
287 // How is it possible to do smething better ????
288 ((gAlice->TreeK())->GetCurrentFile())->cd();
289 if(AliLog::GetGlobalDebugLevel()>1) cout << "After cd(gAlice)" << endl; gDirectory->Dump();
294 //__________________________________________________________________________
295 void AliMUONTrackReconstructor::NextBkgTrackRefEvent(void)
297 // Get next event in background file for track ref. hits
298 // Goes back to event number 0 when end of file is reached
301 AliDebug(1,"Enter NextBkgTrackRefEvent");
302 // Clean previous event
303 if(fBkgTrackRefTK) delete fBkgTrackRefTK;
304 fBkgTrackRefTK = NULL;
305 if(fBkgTrackRefParticles) fBkgTrackRefParticles->Clear();
306 if(fBkgTrackRefTTR) delete fBkgTrackRefTTR;
307 fBkgTrackRefTTR = NULL;
308 // Increment event number
309 fBkgTrackRefEventNumber++;
310 // Get access to Particles and Hits for event from background file
311 if (AliLog::GetGlobalDebugLevel()>1) cout << "Before cd(Bkg)" << endl; gDirectory->Dump();
312 fBkgTrackRefFile->cd();
313 if (AliLog::GetGlobalDebugLevel()>1) cout << "After cd(Bkg)" << endl; gDirectory->Dump();
314 // Particles: TreeK for event and branch "Particles"
315 sprintf(treeName, "TreeK%d", fBkgTrackRefEventNumber);
316 fBkgTrackRefTK = (TTree*)gDirectory->Get(treeName);
317 if (!fBkgTrackRefTK) {
319 AliDebug(1,Form("Cannot find Kine Tree for background event: %d",fBkgTrackRefEventNumber));
320 AliDebug(1,"Goes back to event 0");
322 fBkgTrackRefEventNumber = 0;
323 sprintf(treeName, "TreeK%d", fBkgTrackRefEventNumber);
324 fBkgTrackRefTK = (TTree*)gDirectory->Get(treeName);
325 if (!fBkgTrackRefTK) {
326 AliError(Form("cannot find Kine Tree for background event: %d",fBkgTrackRefEventNumber));
331 fBkgTrackRefTK->SetBranchAddress("Particles", &fBkgTrackRefParticles);
332 fBkgTrackRefTK->GetEvent(0); // why event 0 ???? necessary ????
333 // Hits: TreeH for event and branch "MUON"
334 sprintf(treeName, "TreeTR%d", fBkgTrackRefEventNumber);
335 fBkgTrackRefTTR = (TTree*)gDirectory->Get(treeName);
336 if (!fBkgTrackRefTTR) {
337 AliError(Form("cannot find Hits Tree for background event: %d",fBkgTrackRefEventNumber));
340 fBkgTrackRefTTR->GetEntries(); // necessary ????
341 // Back to the signal file
342 ((gAlice->TreeK())->GetCurrentFile())->cd();
343 if (AliLog::GetGlobalDebugLevel()>1)
344 cout << "After cd(gAlice)" << endl; gDirectory->Dump();
348 //__________________________________________________________________________
349 void AliMUONTrackReconstructor::EventReconstruct(void)
351 // To reconstruct one event
352 AliDebug(1,"Enter EventReconstruct");
354 ResetTrackHits(); //AZ
355 ResetSegments(); //AZ
356 ResetHitsForRec(); //AZ
357 MakeEventToBeReconstructed();
360 if (fMUONData->IsTriggerTrackBranchesInTree())
361 ValidateTracksWithTrigger();
363 // Add tracks to MUON data container
364 for(Int_t i=0; i<GetNRecTracks(); i++) {
365 AliMUONTrack * track = (AliMUONTrack*) GetRecTracksPtr()->At(i);
366 fMUONData->AddRecTrack(*track);
372 //__________________________________________________________________________
373 void AliMUONTrackReconstructor::EventReconstructTrigger(void)
375 // To reconstruct one event
376 AliDebug(1,"Enter EventReconstructTrigger");
381 //__________________________________________________________________________
382 void AliMUONTrackReconstructor::ResetHitsForRec(void)
384 // To reset the array and the number of HitsForRec,
385 // and also the number of HitsForRec
386 // and the index of the first HitForRec per chamber
387 if (fHitsForRecPtr) fHitsForRecPtr->Delete();
389 for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++)
390 fNHitsForRecPerChamber[ch] = fIndexOfFirstHitForRecPerChamber[ch] = 0;
394 //__________________________________________________________________________
395 void AliMUONTrackReconstructor::ResetSegments(void)
397 // To reset the TClonesArray of segments and the number of Segments
399 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) {
400 if (fSegmentsPtr[st]) fSegmentsPtr[st]->Clear();
406 //__________________________________________________________________________
407 void AliMUONTrackReconstructor::ResetTracks(void)
409 // To reset the TClonesArray of reconstructed tracks
410 if (fRecTracksPtr) fRecTracksPtr->Delete();
411 // Delete in order that the Track destructors are called,
412 // hence the space for the TClonesArray of pointers to TrackHit's is freed
418 //__________________________________________________________________________
419 void AliMUONTrackReconstructor::ResetTrackHits(void)
421 // To reset the TClonesArray of hits on reconstructed tracks
422 if (fRecTrackHitsPtr) fRecTrackHitsPtr->Clear();
427 //__________________________________________________________________________
428 void AliMUONTrackReconstructor::MakeEventToBeReconstructed(void)
430 // To make the list of hits to be reconstructed,
431 // either from the track ref. hits or from the raw clusters
432 // according to the parameter set for the reconstructor
433 // TString evfoldname = AliConfig::GetDefaultEventFolderName();//to be interfaced properly
435 // AliRunLoader* rl = AliRunLoader::GetRunLoader(evfoldname);
438 // Error("MakeEventToBeReconstructed",
439 // "Can not find Run Loader in Event Folder named %s.",
440 // evfoldname.Data());
443 // AliLoader* gime = rl->GetLoader("MUONLoader");
446 // Error("MakeEventToBeReconstructed","Can not get MUON Loader from Run Loader.");
449 AliRunLoader *runLoader = fLoader->GetRunLoader();
451 AliDebug(1,"Enter MakeEventToBeReconstructed");
452 //AZ ResetHitsForRec();
453 if (fRecTrackRefHits == 1) {
454 // Reconstruction from track ref. hits
455 // Back to the signal file
456 TTree* treeTR = runLoader->TreeTR();
459 Int_t retval = runLoader->LoadTrackRefs("READ");
462 AliError("Error occured while loading hits.");
465 treeTR = runLoader->TreeTR();
468 AliError("Can not get TreeTR");
473 AddHitsForRecFromTrackRef(treeTR,1);
476 AddHitsForRecFromTrackRef(fBkgTrackRefTTR,0);
477 // Sort HitsForRec in increasing order with respect to chamber number
478 SortHitsForRecWithIncreasingChamber();
481 // Reconstruction from raw clusters
482 // AliMUON *MUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
483 // Security on MUON ????
484 // TreeR assumed to be be "prepared" in calling function
485 // by "MUON->GetTreeR(nev)" ????
486 TTree *treeR = fLoader->TreeR();
487 //AZ? fMUONData->SetTreeAddress("RC");
488 AddHitsForRecFromRawClusters(treeR);
489 // No sorting: it is done automatically in the previous function
492 AliDebug(1,"End of MakeEventToBeReconstructed");
493 if (AliLog::GetGlobalDebugLevel() > 0) {
494 AliDebug(1, Form("NHitsForRec: %d",fNHitsForRec));
495 for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
496 AliDebug(1, Form("Chamber(0...): %d",ch));
497 AliDebug(1, Form("NHitsForRec: %d", fNHitsForRecPerChamber[ch]));
498 AliDebug(1, Form("Index(first HitForRec): %d", fIndexOfFirstHitForRecPerChamber[ch]));
499 for (Int_t hit = fIndexOfFirstHitForRecPerChamber[ch];
500 hit < fIndexOfFirstHitForRecPerChamber[ch] + fNHitsForRecPerChamber[ch];
502 AliDebug(1, Form("HitForRec index(0...): %d",hit));
503 ((*fHitsForRecPtr)[hit])->Dump();
510 //__________________________________________________________________________
511 void AliMUONTrackReconstructor::AddHitsForRecFromTrackRef(TTree *TTR, Int_t signal)
513 // To add to the list of hits for reconstruction
514 // the signal hits from a track reference tree TreeTR.
515 TClonesArray *listOfTrackRefs = NULL;
516 AliTrackReference *trackRef;
519 AliDebug(2,Form("Enter AddHitsForRecFromTrackRef with TTR: %d", TTR));
520 if (TTR == NULL) return;
522 listOfTrackRefs = CleanTrackRefs(TTR);
524 Int_t ntracks = listOfTrackRefs->GetEntriesFast();
526 AliDebug(2,Form("ntracks: %d", ntracks));
528 for (Int_t index = 0; index < ntracks; index++) {
529 trackRef = (AliTrackReference*) listOfTrackRefs->At(index);
530 track = trackRef->GetTrack();
532 NewHitForRecFromTrackRef(trackRef,track,signal);
533 } // end of track ref.
535 listOfTrackRefs->Delete();
536 delete listOfTrackRefs;
541 //__________________________________________________________________________
542 AliMUONHitForRec* AliMUONTrackReconstructor::NewHitForRecFromTrackRef(AliTrackReference* Hit, Int_t TrackNumber, Int_t Signal)
544 // To make a new hit for reconstruction from a track ref. hit pointed to by "Hit",
545 // with the track numbered "TrackNumber",
546 // either from signal ("Signal" = 1) or background ("Signal" = 0) event.
547 // Selects hits in tracking (not trigger) chambers.
548 // Takes into account the efficiency (fEfficiency)
549 // and the smearing from resolution (fBendingResolution and fNonBendingResolution).
550 // Adds a condition on the radius between RMin and RMax
551 // to better simulate the real chambers.
552 // Returns the pointer to the new hit for reconstruction,
553 // or NULL in case of inefficiency or non tracking chamber or bad radius.
554 // No condition on at most 20 cm from a muon from a resonance
555 // like in Fortran TRACKF_STAT.
556 AliMUONHitForRec* hitForRec;
557 Double_t bendCoor, nonBendCoor, radius;
558 Int_t chamber = AliMUONConstants::ChamberNumber(Hit->Z()); // chamber(0...)
559 if (chamber < 0) return NULL;
560 // only in tracking chambers (fChamber starts at 1)
561 if (chamber >= AliMUONConstants::NTrackingCh()) return NULL;
562 // only if hit is efficient (keep track for checking ????)
563 if (gRandom->Rndm() > fEfficiency) return NULL;
564 // only if radius between RMin and RMax
566 nonBendCoor = Hit->X();
567 radius = TMath::Sqrt((bendCoor * bendCoor) + (nonBendCoor * nonBendCoor));
568 // This cut is not needed with a realistic chamber geometry !!!!
569 // if ((radius < fRMin[chamber]) || (radius > fRMax[chamber])) return NULL;
570 // new AliMUONHitForRec from track ref. hit and increment number of AliMUONHitForRec's
571 hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(Hit);
573 // add smearing from resolution
574 hitForRec->SetBendingCoor(bendCoor + gRandom->Gaus(0., fBendingResolution));
575 hitForRec->SetNonBendingCoor(nonBendCoor
576 + gRandom->Gaus(0., fNonBendingResolution));
577 // // !!!! without smearing
578 // hitForRec->SetBendingCoor(bendCoor);
579 // hitForRec->SetNonBendingCoor(nonBendCoor);
580 // more information into HitForRec
581 // resolution: angular effect to be added here ????
582 hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
583 hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
585 hitForRec->SetTTRTrack(TrackNumber);
586 hitForRec->SetTrackRefSignal(Signal);
587 if (AliLog::GetGlobalDebugLevel()> 1) {
588 AliDebug(2,Form("Track: %d", TrackNumber));
590 cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
595 //__________________________________________________________________________
596 TClonesArray* AliMUONTrackReconstructor::CleanTrackRefs(TTree *treeTR)
598 // Make hits from track ref..
599 // Re-calculate hits parameters because two AliTrackReferences are recorded for
600 // each chamber (one when particle is entering + one when particle is leaving
601 // the sensitive volume)
603 AliTrackReference *trackReference;
604 Float_t x1, y1, z1, pX1, pY1, pZ1;
605 Float_t x2, y2, z2, pX2, pY2, pZ2;
606 Int_t track1, track2;
608 Float_t maxGasGap = 1.; // cm
612 AliTrackReference *trackReferenceNew = new AliTrackReference();
614 TClonesArray* trackRefs = new TClonesArray("AliTrackReference", 10);
615 TClonesArray* cleanTrackRefs = new TClonesArray("AliTrackReference", 10);
617 if (treeTR == NULL) return NULL;
618 TBranch* branch = treeTR->GetBranch("MUON");
619 if (branch == NULL) return NULL;
620 branch->SetAddress(&trackRefs);
622 Int_t nTrack = (Int_t)treeTR->GetEntries();
623 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
624 treeTR->GetEntry(iTrack);
627 while (iHit1 < trackRefs->GetEntries()) {
628 trackReference = (AliTrackReference*)trackRefs->At(iHit1);
629 x1 = trackReference->X();
630 y1 = trackReference->Y();
631 z1 = trackReference->Z();
632 pX1 = trackReference->Px();
633 pY1 = trackReference->Py();
634 pZ1 = trackReference->Pz();
635 track1 = trackReference->GetTrack();
638 for (Int_t iHit2 = iHit1+1; iHit2 < trackRefs->GetEntries(); iHit2++) {
639 trackReference = (AliTrackReference*)trackRefs->At(iHit2);
640 x2 = trackReference->X();
641 y2 = trackReference->Y();
642 z2 = trackReference->Z();
643 pX2 = trackReference->Px();
644 pY2 = trackReference->Py();
645 pZ2 = trackReference->Pz();
646 track2 = trackReference->GetTrack();
647 if (track2 == track1 && TMath::Abs(z2-z1) < maxGasGap ) {
662 pX1 /= (Float_t)nRec;
663 pY1 /= (Float_t)nRec;
664 pZ1 /= (Float_t)nRec;
665 trackReferenceNew->SetPosition(x1,y1,z1);
666 trackReferenceNew->SetMomentum(pX1,pY1,pZ1);
667 trackReferenceNew->SetTrack(track1);
668 {new ((*cleanTrackRefs)[cleanTrackRefs->GetEntriesFast()]) AliTrackReference(*trackReferenceNew);}
675 delete trackReferenceNew;
676 return cleanTrackRefs;
679 //__________________________________________________________________________
680 void AliMUONTrackReconstructor::SortHitsForRecWithIncreasingChamber()
682 // Sort HitsForRec's in increasing order with respect to chamber number.
683 // Uses the function "Compare".
684 // Update the information for HitsForRec per chamber too.
685 Int_t ch, nhits, prevch;
686 fHitsForRecPtr->Sort();
687 for (ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
688 fNHitsForRecPerChamber[ch] = 0;
689 fIndexOfFirstHitForRecPerChamber[ch] = 0;
691 prevch = 0; // previous chamber
692 nhits = 0; // number of hits in current chamber
693 // Loop over HitsForRec
694 for (Int_t hit = 0; hit < fNHitsForRec; hit++) {
695 // chamber number (0...)
696 ch = ((AliMUONHitForRec*) ((*fHitsForRecPtr)[hit]))->GetChamberNumber();
697 // increment number of hits in current chamber
698 (fNHitsForRecPerChamber[ch])++;
699 // update index of first HitForRec in current chamber
700 // if chamber number different from previous one
702 fIndexOfFirstHitForRecPerChamber[ch] = hit;
709 //__________________________________________________________________________
710 void AliMUONTrackReconstructor::AddHitsForRecFromRawClusters(TTree* TR)
712 // To add to the list of hits for reconstruction all the raw clusters
713 // No condition added, like in Fortran TRACKF_STAT,
714 // on the radius between RMin and RMax.
715 AliMUONHitForRec *hitForRec;
716 AliMUONRawCluster *clus;
717 Int_t iclus, nclus, nTRentries;
718 TClonesArray *rawclusters;
719 AliDebug(1,"Enter AddHitsForRecFromRawClusters");
721 if (fTrackMethod != 3) { //AZ
722 fMUONData->SetTreeAddress("RC"); //AZ
723 nTRentries = Int_t(TR->GetEntries());
724 if (nTRentries != 1) {
725 AliError(Form("nTRentries = %d not equal to 1 ",nTRentries));
728 fLoader->TreeR()->GetEvent(0); // only one entry
731 // Loop over tracking chambers
732 for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
733 // number of HitsForRec to 0 for the chamber
734 fNHitsForRecPerChamber[ch] = 0;
735 // index of first HitForRec for the chamber
736 if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0;
737 else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec;
738 rawclusters =fMUONData->RawClusters(ch);
739 // pMUON->ResetRawClusters();
740 // TR->GetEvent((Int_t) (TR->GetEntries()) - 1); // to be checked ????
741 nclus = (Int_t) (rawclusters->GetEntries());
742 // Loop over (cathode correlated) raw clusters
743 for (iclus = 0; iclus < nclus; iclus++) {
744 clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(iclus);
745 // new AliMUONHitForRec from raw cluster
746 // and increment number of AliMUONHitForRec's (total and in chamber)
747 hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(clus);
749 (fNHitsForRecPerChamber[ch])++;
750 // more information into HitForRec
751 // resolution: info should be already in raw cluster and taken from it ????
752 //hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
753 //hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
754 hitForRec->SetBendingReso2(clus->GetErrY() * clus->GetErrY());
755 hitForRec->SetNonBendingReso2(clus->GetErrX() * clus->GetErrX());
756 // original raw cluster
757 hitForRec->SetChamberNumber(ch);
758 hitForRec->SetHitNumber(iclus);
759 // Z coordinate of the raw cluster (cm)
760 hitForRec->SetZ(clus->GetZ(0));
761 if (AliLog::GetGlobalDebugLevel() > 1) {
762 cout << "chamber (0...): " << ch <<
763 " raw cluster (0...): " << iclus << endl;
765 cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
767 } // end of cluster loop
768 } // end of chamber loop
769 SortHitsForRecWithIncreasingChamber();
773 //__________________________________________________________________________
774 void AliMUONTrackReconstructor::MakeSegments(void)
776 // To make the list of segments in all stations,
777 // from the list of hits to be reconstructed
778 AliDebug(1,"Enter MakeSegments");
779 //AZ ResetSegments();
780 // Loop over stations
781 Int_t nb = (fTrackMethod != 1) ? 3 : 0; //AZ
782 for (Int_t st = nb; st < AliMUONConstants::NTrackingCh()/2; st++)
783 MakeSegmentsPerStation(st);
784 if (AliLog::GetGlobalDebugLevel() > 1) {
785 cout << "end of MakeSegments" << endl;
786 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) {
787 cout << "station(0...): " << st
788 << " Segments: " << fNSegments[st]
791 seg < fNSegments[st];
793 cout << "Segment index(0...): " << seg << endl;
794 ((*fSegmentsPtr[st])[seg])->Dump();
801 //__________________________________________________________________________
802 void AliMUONTrackReconstructor::MakeSegmentsPerStation(Int_t Station)
804 // To make the list of segments in station number "Station" (0...)
805 // from the list of hits to be reconstructed.
806 // Updates "fNSegments"[Station].
807 // Segments in stations 4 and 5 are sorted
808 // according to increasing absolute value of "impact parameter"
809 AliMUONHitForRec *hit1Ptr, *hit2Ptr;
810 AliMUONSegment *segment;
812 Double_t bendingSlope, distBend, distNonBend, extBendCoor, extNonBendCoor,
813 impactParam = 0., maxImpactParam = 0., minImpactParam = 0.; // =0 to avoid compilation warnings.
814 AliDebug(1,Form("Enter MakeSegmentsPerStation (0...) %d",Station));
815 // first and second chambers (0...) in the station
816 Int_t ch1 = 2 * Station;
818 // variable true for stations downstream of the dipole:
819 // Station(0..4) equal to 3 or 4
820 if ((Station == 3) || (Station == 4)) {
822 // maximum impact parameter (cm) according to fMinBendingMomentum...
824 TMath::Abs(GetImpactParamFromBendingMomentum(fMinBendingMomentum));
825 // minimum impact parameter (cm) according to fMaxBendingMomentum...
827 TMath::Abs(GetImpactParamFromBendingMomentum(fMaxBendingMomentum));
829 else last2st = kFALSE;
830 // extrapolation factor from Z of first chamber to Z of second chamber
831 // dZ to be changed to take into account fine structure of chambers ????
833 // index for current segment
834 Int_t segmentIndex = 0;
835 // Loop over HitsForRec in the first chamber of the station
836 for (Int_t hit1 = fIndexOfFirstHitForRecPerChamber[ch1];
837 hit1 < fIndexOfFirstHitForRecPerChamber[ch1] + fNHitsForRecPerChamber[ch1];
839 // pointer to the HitForRec
840 hit1Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit1]);
842 // on the straight line joining the HitForRec to the vertex (0,0,0),
843 // to the Z of the second chamber of the station
844 // Loop over HitsForRec in the second chamber of the station
845 for (Int_t hit2 = fIndexOfFirstHitForRecPerChamber[ch2];
846 hit2 < fIndexOfFirstHitForRecPerChamber[ch2] + fNHitsForRecPerChamber[ch2];
848 // pointer to the HitForRec
849 hit2Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit2]);
850 // absolute values of distances, in bending and non bending planes,
851 // between the HitForRec in the second chamber
852 // and the previous extrapolation
853 extrapFact = hit2Ptr->GetZ()/ hit1Ptr->GetZ();
854 extBendCoor = extrapFact * hit1Ptr->GetBendingCoor();
855 extNonBendCoor = extrapFact * hit1Ptr->GetNonBendingCoor();
856 distBend = TMath::Abs(hit2Ptr->GetBendingCoor() - extBendCoor);
857 distNonBend = TMath::Abs(hit2Ptr->GetNonBendingCoor() - extNonBendCoor);
860 if ( hit1Ptr->GetZ() - hit2Ptr->GetZ() != 0.0 ) {
861 bendingSlope = (hit1Ptr->GetBendingCoor() - hit2Ptr->GetBendingCoor()) /
862 (hit1Ptr->GetZ() - hit2Ptr->GetZ());
863 // absolute value of impact parameter
865 TMath::Abs(hit1Ptr->GetBendingCoor() - hit1Ptr->GetZ() * bendingSlope);
868 AliWarning("hit1Ptr->GetZ() = hit2Ptr->GetZ(): impactParam set to maxImpactParam");
869 impactParam = maxImpactParam;
872 // check for distances not too large,
873 // and impact parameter not too big if stations downstream of the dipole.
874 // Conditions "distBend" and "impactParam" correlated for these stations ????
875 if ((distBend < fSegmentMaxDistBending[Station]) &&
876 (distNonBend < fSegmentMaxDistNonBending[Station]) &&
877 (!last2st || (impactParam < maxImpactParam)) &&
878 (!last2st || (impactParam > minImpactParam))) {
880 segment = new ((*fSegmentsPtr[Station])[segmentIndex])
881 AliMUONSegment(hit1Ptr, hit2Ptr);
882 // update "link" to this segment from the hit in the first chamber
883 if (hit1Ptr->GetNSegments() == 0)
884 hit1Ptr->SetIndexOfFirstSegment(segmentIndex);
885 hit1Ptr->SetNSegments(hit1Ptr->GetNSegments() + 1);
886 if (AliLog::GetGlobalDebugLevel() > 1) {
887 cout << "segmentIndex(0...): " << segmentIndex
888 << " distBend: " << distBend
889 << " distNonBend: " << distNonBend
892 cout << "HitForRec in first chamber" << endl;
894 cout << "HitForRec in second chamber" << endl;
897 // increment index for current segment
901 } // for (Int_t hit1...
902 fNSegments[Station] = segmentIndex;
903 // Sorting according to "impact parameter" if station(1..5) 4 or 5,
904 // i.e. Station(0..4) 3 or 4, using the function "Compare".
905 // After this sorting, it is impossible to use
906 // the "fNSegments" and "fIndexOfFirstSegment"
907 // of the HitForRec in the first chamber to explore all segments formed with it.
908 // Is this sorting really needed ????
909 if ((Station == 3) || (Station == 4)) (fSegmentsPtr[Station])->Sort();
910 AliDebug(1,Form("Station: %d NSegments: %d ", Station, fNSegments[Station]));
914 //__________________________________________________________________________
915 void AliMUONTrackReconstructor::MakeTracks(void)
917 // To make the tracks,
918 // from the list of segments and points in all stations
919 AliDebug(1,"Enter MakeTracks");
920 // The order may be important for the following Reset's
922 //AZ ResetTrackHits();
923 if (fTrackMethod != 1) { //AZ - Kalman filter
924 MakeTrackCandidatesK();
925 if (fRecTracksPtr->GetEntriesFast() == 0) return;
926 // Follow tracks in stations(1..) 3, 2 and 1
928 // Remove double tracks
929 RemoveDoubleTracksK();
930 // Propagate tracks to the vertex thru absorber
932 // Fill AliMUONTrack data members
935 // Look for candidates from at least 3 aligned points in stations(1..) 4 and 5
936 MakeTrackCandidates();
937 // Follow tracks in stations(1..) 3, 2 and 1
939 // Remove double tracks
940 RemoveDoubleTracks();
941 UpdateTrackParamAtHit();
942 UpdateHitForRecAtHit();
947 //__________________________________________________________________________
948 void AliMUONTrackReconstructor::ValidateTracksWithTrigger(void)
950 // Try to match track from tracking system with trigger track
952 TClonesArray *recTriggerTracks;
954 fMUONData->SetTreeAddress("RL");
955 fMUONData->GetRecTriggerTracks();
956 recTriggerTracks = fMUONData->RecTriggerTracks();
958 track = (AliMUONTrack*) fRecTracksPtr->First();
960 track->MatchTriggerTrack(recTriggerTracks);
961 track = (AliMUONTrack*) fRecTracksPtr->After(track);
966 //__________________________________________________________________________
967 Bool_t AliMUONTrackReconstructor::MakeTriggerTracks(void)
969 // To make the trigger tracks from Local Trigger
970 AliDebug(1, "Enter MakeTriggerTracks");
974 TClonesArray *localTrigger;
975 TClonesArray *globalTrigger;
976 AliMUONLocalTrigger *locTrg;
977 AliMUONGlobalTrigger *gloTrg;
979 AliMUONTriggerTrack *recTriggerTrack = 0;
981 TTree* treeR = fLoader->TreeR();
983 // Loading MUON subsystem
984 AliMUON * pMUON = (AliMUON *) gAlice->GetDetector("MUON");
985 // Not really clean, but for the moment we must check whether the
986 // trigger uses new or old TriggerCircuit
987 Bool_t newTrigger=kFALSE;
988 if ( pMUON->DigitizerType().Contains("NewTrigger") ) newTrigger = kTRUE;
990 nTRentries = Int_t(treeR->GetEntries());
992 treeR->GetEvent(0); // only one entry
994 if (!(fMUONData->IsTriggerBranchesInTree())) {
995 AliWarning(Form("Trigger information is not avalaible, nTRentries = %d not equal to 1",nTRentries));
999 fMUONData->SetTreeAddress("TC");
1000 fMUONData->GetTrigger();
1002 // global trigger for trigger pattern
1004 globalTrigger = fMUONData->GlobalTrigger();
1005 gloTrg = (AliMUONGlobalTrigger*)globalTrigger->UncheckedAt(0);
1007 if (gloTrg->SinglePlusLpt()) gloTrigPat|= 0x1;
1008 if (gloTrg->SinglePlusHpt()) gloTrigPat|= 0x2;
1009 if (gloTrg->SinglePlusApt()) gloTrigPat|= 0x4;
1011 if (gloTrg->SingleMinusLpt()) gloTrigPat|= 0x8;
1012 if (gloTrg->SingleMinusHpt()) gloTrigPat|= 0x10;
1013 if (gloTrg->SingleMinusApt()) gloTrigPat|= 0x20;
1015 if (gloTrg->SingleUndefLpt()) gloTrigPat|= 0x40;
1016 if (gloTrg->SingleUndefHpt()) gloTrigPat|= 0x80;
1017 if (gloTrg->SingleUndefApt()) gloTrigPat|= 0x100;
1019 if (gloTrg->PairUnlikeLpt()) gloTrigPat|= 0x200;
1020 if (gloTrg->PairUnlikeHpt()) gloTrigPat|= 0x400;
1021 if (gloTrg->PairUnlikeApt()) gloTrigPat|= 0x800;
1023 if (gloTrg->PairLikeLpt()) gloTrigPat|= 0x1000;
1024 if (gloTrg->PairLikeHpt()) gloTrigPat|= 0x2000;
1025 if (gloTrg->PairLikeApt()) gloTrigPat|= 0x4000;
1028 // local trigger for tracking
1029 localTrigger = fMUONData->LocalTrigger();
1030 Int_t nlocals = (Int_t) (localTrigger->GetEntries());
1032 Float_t z11 = AliMUONConstants::DefaultChamberZ(10);
1033 Float_t z21 = AliMUONConstants::DefaultChamberZ(12);
1040 for (Int_t i=0; i<nlocals; i++) { // loop on Local Trigger
1041 locTrg = (AliMUONLocalTrigger*)localTrigger->UncheckedAt(i);
1043 if (!newTrigger) { // old trigger
1044 // printf("AliMUONTrackReconstructor::MakeTriggerTrack using OLD trigger \n");
1045 AliMUONTriggerCircuit * circuit = &(pMUON->TriggerCircuit(locTrg->LoCircuit()));
1046 y11 = circuit->GetY11Pos(locTrg->LoStripX());
1047 stripX21 = locTrg->LoStripX()+locTrg->LoDev()+1;
1048 y21 = circuit->GetY21Pos(stripX21);
1049 x11 = circuit->GetX11Pos(locTrg->LoStripY());
1050 } else { // new trigger
1051 // printf("AliMUONTrackReconstructor::MakeTriggerTrack using NEW trigger \n");
1052 AliMUONTriggerCircuitNew * circuit =
1053 &(pMUON->TriggerCircuitNew(locTrg->LoCircuit()-1)); // -1 !!!
1054 y11 = circuit->GetY11Pos(locTrg->LoStripX());
1055 stripX21 = locTrg->LoStripX()+locTrg->LoDev()+1;
1056 y21 = circuit->GetY21Pos(stripX21);
1057 x11 = circuit->GetX11Pos(locTrg->LoStripY());
1059 // printf(" MakeTriggerTrack %d %d %d %d %d %f %f %f \n",i,locTrg->LoCircuit(),locTrg->LoStripX(),locTrg->LoStripX()+locTrg->LoDev()+1,locTrg->LoStripY(),y11, y21, x11);
1061 Float_t thetax = TMath::ATan2( x11 , z11 );
1062 Float_t thetay = TMath::ATan2( (y21-y11) , (z21-z11) );
1064 recTriggerTrack = new AliMUONTriggerTrack(x11,y11,thetax,thetay,gloTrigPat);
1066 // since static statement does not work, set gloTrigPat for each track
1068 fMUONData->AddRecTriggerTrack(*recTriggerTrack);
1069 delete recTriggerTrack;
1070 } // end of loop on Local Trigger
1074 //__________________________________________________________________________
1075 Int_t AliMUONTrackReconstructor::MakeTrackCandidatesWithTwoSegments(AliMUONSegment *BegSegment)
1077 // To make track candidates with two segments in stations(1..) 4 and 5,
1078 // the first segment being pointed to by "BegSegment".
1079 // Returns the number of such track candidates.
1080 Int_t endStation, iEndSegment, nbCan2Seg;
1081 AliMUONSegment *endSegment;
1082 AliMUONSegment *extrapSegment = NULL;
1083 AliMUONTrack *recTrack;
1085 AliDebug(1,"Enter MakeTrackCandidatesWithTwoSegments");
1086 // Station for the end segment
1087 endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
1088 // multiple scattering factor corresponding to one chamber
1089 mcsFactor = 0.0136 /
1090 GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
1091 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1092 // linear extrapolation to end station
1093 // number of candidates with 2 segments to 0
1095 // Loop over segments in the end station
1096 for (iEndSegment = 0; iEndSegment < fNSegments[endStation]; iEndSegment++) {
1097 // pointer to segment
1098 endSegment = (AliMUONSegment*) ((*fSegmentsPtr[endStation])[iEndSegment]);
1099 // test compatibility between current segment and "extrapSegment"
1100 // 4 because 4 quantities in chi2
1102 BegSegment->CreateSegmentFromLinearExtrapToStation(endSegment->GetZ(), mcsFactor);
1104 NormalizedChi2WithSegment(extrapSegment,
1105 fMaxSigma2Distance)) <= 4.0) {
1106 // both segments compatible:
1107 // make track candidate from "begSegment" and "endSegment"
1108 AliDebug(2,Form("TrackCandidate with Segment %d in Station(0..) %d", iEndSegment, endStation));
1109 // flag for both segments in one track:
1110 // to be done in track constructor ????
1111 BegSegment->SetInTrack(kTRUE);
1112 endSegment->SetInTrack(kTRUE);
1113 recTrack = new ((*fRecTracksPtr)[fNRecTracks])
1114 AliMUONTrack(BegSegment, endSegment, this);
1116 if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump();
1117 // increment number of track candidates with 2 segments
1120 delete extrapSegment; // should not delete HitForRec's it points to !!!!
1121 } // for (iEndSegment = 0;...
1125 //__________________________________________________________________________
1126 Int_t AliMUONTrackReconstructor::MakeTrackCandidatesWithOneSegmentAndOnePoint(AliMUONSegment *BegSegment)
1128 // To make track candidates with one segment and one point
1129 // in stations(1..) 4 and 5,
1130 // the segment being pointed to by "BegSegment".
1131 Int_t ch, ch1, ch2, endStation, iHit, iHitMax, iHitMin, nbCan1Seg1Hit;
1132 AliMUONHitForRec *extrapHitForRec= NULL;
1133 AliMUONHitForRec *hit;
1134 AliMUONTrack *recTrack;
1136 AliDebug(1,"Enter MakeTrackCandidatesWithOneSegmentAndOnePoint");
1137 // station for the end point
1138 endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
1139 // multiple scattering factor corresponding to one chamber
1140 mcsFactor = 0.0136 /
1141 GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
1142 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1143 // first and second chambers(0..) in the end station
1144 ch1 = 2 * endStation;
1146 // number of candidates to 0
1148 // Loop over chambers of the end station
1149 for (ch = ch2; ch >= ch1; ch--) {
1150 // limits for the hit index in the loop
1151 iHitMin = fIndexOfFirstHitForRecPerChamber[ch];
1152 iHitMax = iHitMin + fNHitsForRecPerChamber[ch];
1153 // Loop over HitForRec's in the chamber
1154 for (iHit = iHitMin; iHit < iHitMax; iHit++) {
1155 // pointer to HitForRec
1156 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1157 // test compatibility between current HitForRec and "extrapHitForRec"
1158 // 2 because 2 quantities in chi2
1159 // linear extrapolation to chamber
1161 BegSegment->CreateHitForRecFromLinearExtrapToChamber( hit->GetZ(), mcsFactor);
1163 NormalizedChi2WithHitForRec(extrapHitForRec,
1164 fMaxSigma2Distance)) <= 2.0) {
1165 // both HitForRec's compatible:
1166 // make track candidate from begSegment and current HitForRec
1167 AliDebug(2, Form("TrackCandidate with HitForRec %d in Chamber(0..) %d", iHit, ch));
1168 // flag for beginning segments in one track:
1169 // to be done in track constructor ????
1170 BegSegment->SetInTrack(kTRUE);
1171 recTrack = new ((*fRecTracksPtr)[fNRecTracks])
1172 AliMUONTrack(BegSegment, hit, this);
1173 // the right place to eliminate "double counting" ???? how ????
1175 if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump();
1176 // increment number of track candidates
1179 delete extrapHitForRec;
1180 } // for (iHit = iHitMin;...
1181 } // for (ch = ch2;...
1182 return nbCan1Seg1Hit;
1185 //__________________________________________________________________________
1186 void AliMUONTrackReconstructor::MakeTrackCandidates(void)
1188 // To make track candidates
1189 // with at least 3 aligned points in stations(1..) 4 and 5
1190 // (two Segment's or one Segment and one HitForRec)
1191 Int_t begStation, iBegSegment, nbCan1Seg1Hit, nbCan2Seg;
1192 AliMUONSegment *begSegment;
1193 AliDebug(1,"Enter MakeTrackCandidates");
1194 // Loop over stations(1..) 5 and 4 for the beginning segment
1195 for (begStation = 4; begStation > 2; begStation--) {
1196 // Loop over segments in the beginning station
1197 for (iBegSegment = 0; iBegSegment < fNSegments[begStation]; iBegSegment++) {
1198 // pointer to segment
1199 begSegment = (AliMUONSegment*) ((*fSegmentsPtr[begStation])[iBegSegment]);
1200 AliDebug(2,Form("Look for TrackCandidate's with Segment %d in Station(0..) %d", iBegSegment, begStation));
1201 // Look for track candidates with two segments,
1202 // "begSegment" and all compatible segments in other station.
1203 // Only for beginning station(1..) 5
1204 // because candidates with 2 segments have to looked for only once.
1205 if (begStation == 4)
1206 nbCan2Seg = MakeTrackCandidatesWithTwoSegments(begSegment);
1207 // Look for track candidates with one segment and one point,
1208 // "begSegment" and all compatible HitForRec's in other station.
1209 // Only if "begSegment" does not belong already to a track candidate.
1210 // Is that a too strong condition ????
1211 if (!(begSegment->GetInTrack()))
1212 nbCan1Seg1Hit = MakeTrackCandidatesWithOneSegmentAndOnePoint(begSegment);
1213 } // for (iBegSegment = 0;...
1214 } // for (begStation = 4;...
1218 //__________________________________________________________________________
1219 void AliMUONTrackReconstructor::FollowTracks(void)
1221 // Follow tracks in stations(1..) 3, 2 and 1
1222 // too long: should be made more modular !!!!
1223 AliMUONHitForRec *bestHit, *extrapHit, *hit;
1224 AliMUONSegment *bestSegment, *extrapSegment, *segment;
1225 AliMUONTrack *track, *nextTrack;
1226 AliMUONTrackParam *trackParam1, trackParam[2], trackParamVertex;
1227 // -1 to avoid compilation warnings
1228 Int_t ch = -1, chInStation, chBestHit = -1, iHit, iSegment, station, trackIndex;
1229 Double_t bestChi2, chi2, dZ1, dZ2, dZ3, maxSigma2Distance, mcsFactor;
1230 Double_t bendingMomentum, chi2Norm = 0.;
1232 // local maxSigma2Distance, for easy increase in testing
1233 maxSigma2Distance = fMaxSigma2Distance;
1234 AliDebug(2,"Enter FollowTracks");
1235 // Loop over track candidates
1236 track = (AliMUONTrack*) fRecTracksPtr->First();
1239 // Follow function for each track candidate ????
1241 nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track); // prepare next track
1242 AliDebug(2,Form("FollowTracks: track candidate(0..): %d", trackIndex));
1243 // Fit track candidate
1244 track->SetFitMCS(0); // without multiple Coulomb scattering
1245 track->SetFitNParam(3); // with 3 parameters (X = Y = 0)
1246 track->SetFitStart(0); // from parameters at vertex
1248 if (AliLog::GetGlobalDebugLevel()> 2) {
1249 cout << "FollowTracks: track candidate(0..): " << trackIndex
1250 << " after fit in stations(0..) 3 and 4" << endl;
1251 track->RecursiveDump();
1253 // Loop over stations(1..) 3, 2 and 1
1254 // something SPECIAL for stations 2 and 1 for majority 3 coincidence ????
1255 // otherwise: majority coincidence 2 !!!!
1256 for (station = 2; station >= 0; station--) {
1257 // Track parameters at first track hit (smallest Z)
1258 trackParam1 = ((AliMUONTrackHit*)
1259 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1260 // extrapolation to station
1261 trackParam1->ExtrapToStation(station, trackParam);
1262 extrapSegment = new AliMUONSegment(); // empty segment
1263 // multiple scattering factor corresponding to one chamber
1264 // and momentum in bending plane (not total)
1265 mcsFactor = 0.0136 * trackParam1->GetInverseBendingMomentum();
1266 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1267 // Z difference from previous station
1268 dZ1 = AliMUONConstants::DefaultChamberZ(2 * station) -
1269 AliMUONConstants::DefaultChamberZ(2 * station + 2);
1270 // Z difference between the two previous stations
1271 dZ2 = AliMUONConstants::DefaultChamberZ(2 * station + 2) -
1272 AliMUONConstants::DefaultChamberZ(2 * station + 4);
1273 // Z difference between the two chambers in the previous station
1274 dZ3 = AliMUONConstants::DefaultChamberZ(2 * station) -
1275 AliMUONConstants::DefaultChamberZ(2 * station + 1);
1276 extrapSegment->SetBendingCoorReso2(fBendingResolution * fBendingResolution);
1278 SetNonBendingCoorReso2(fNonBendingResolution * fNonBendingResolution);
1279 extrapSegment->UpdateFromStationTrackParam
1280 (trackParam, mcsFactor, dZ1, dZ2, dZ3, station,
1281 trackParam1->GetInverseBendingMomentum());
1284 if (AliLog::GetGlobalDebugLevel() > 2) {
1285 cout << "FollowTracks: track candidate(0..): " << trackIndex
1286 << " Look for segment in station(0..): " << station << endl;
1288 // Loop over segments in station
1289 for (iSegment = 0; iSegment < fNSegments[station]; iSegment++) {
1290 // Look for best compatible Segment in station
1291 // should consider all possibilities ????
1292 // multiple scattering ????
1293 // separation in 2 functions: Segment and HitForRec ????
1294 segment = (AliMUONSegment*) ((*fSegmentsPtr[station])[iSegment]);
1295 // correction of corrected segment (fBendingCoor and fNonBendingCoor)
1296 // according to real Z value of "segment" and slopes of "extrapSegment"
1297 (&(trackParam[0]))->ExtrapToZ(segment->GetZ());
1298 (&(trackParam[1]))->ExtrapToZ(segment->GetZ());
1299 extrapSegment->SetBendingCoor((&(trackParam[0]))->GetBendingCoor());
1300 extrapSegment->SetNonBendingCoor((&(trackParam[0]))->GetNonBendingCoor());
1301 extrapSegment->SetBendingSlope((&(trackParam[0]))->GetBendingSlope());
1302 extrapSegment->SetNonBendingSlope((&(trackParam[0]))->GetNonBendingSlope());
1304 NormalizedChi2WithSegment(extrapSegment, maxSigma2Distance);
1305 if (chi2 < bestChi2) {
1306 // update best Chi2 and Segment if better found
1307 bestSegment = segment;
1312 // best segment found: add it to track candidate
1313 (&(trackParam[0]))->ExtrapToZ(bestSegment->GetZ());
1314 (&(trackParam[1]))->ExtrapToZ(bestSegment->GetZ());
1315 track->AddSegment(bestSegment);
1316 // set track parameters at these two TrakHit's
1317 track->SetTrackParamAtHit(track->GetNTrackHits() - 2, &(trackParam[0]));
1318 track->SetTrackParamAtHit(track->GetNTrackHits() - 1, &(trackParam[1]));
1319 AliDebug(3, Form("FollowTracks: track candidate(0..): %d Added segment in station(0..): %d", trackIndex, station));
1320 if (AliLog::GetGlobalDebugLevel()>2) track->RecursiveDump();
1323 // No best segment found:
1324 // Look for best compatible HitForRec in station:
1325 // should consider all possibilities ????
1326 // multiple scattering ???? do about like for extrapSegment !!!!
1327 extrapHit = new AliMUONHitForRec(); // empty hit
1330 AliDebug(3, Form("FollowTracks: track candidate(0..): %d Segment not found, look for hit in station(0..): %d ",
1331 trackIndex, station));
1333 // Loop over chambers of the station
1334 for (chInStation = 0; chInStation < 2; chInStation++) {
1335 ch = 2 * station + chInStation;
1336 for (iHit = fIndexOfFirstHitForRecPerChamber[ch];
1337 iHit < fIndexOfFirstHitForRecPerChamber[ch] +
1338 fNHitsForRecPerChamber[ch];
1340 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1341 // coordinates of extrapolated hit
1342 (&(trackParam[chInStation]))->ExtrapToZ(hit->GetZ());
1344 SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor());
1346 SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor());
1347 // resolutions from "extrapSegment"
1348 extrapHit->SetBendingReso2(extrapSegment->GetBendingCoorReso2());
1349 extrapHit->SetNonBendingReso2(extrapSegment->GetNonBendingCoorReso2());
1350 // Loop over hits in the chamber
1351 // condition for hit not already in segment ????
1352 chi2 = hit->NormalizedChi2WithHitForRec(extrapHit, maxSigma2Distance);
1353 if (chi2 < bestChi2) {
1354 // update best Chi2 and HitForRec if better found
1357 chBestHit = chInStation;
1362 // best hit found: add it to track candidate
1363 (&(trackParam[chBestHit]))->ExtrapToZ(bestHit->GetZ());
1364 track->AddHitForRec(bestHit);
1365 // set track parameters at this TrackHit
1366 track->SetTrackParamAtHit(track->GetNTrackHits() - 1,
1367 &(trackParam[chBestHit]));
1368 if (AliLog::GetGlobalDebugLevel() > 2) {
1369 cout << "FollowTracks: track candidate(0..): " << trackIndex
1370 << " Added hit in station(0..): " << station << endl;
1371 track->RecursiveDump();
1375 // Remove current track candidate
1376 // and corresponding TrackHit's, ...
1378 delete extrapSegment;
1380 break; // stop the search for this candidate:
1381 // exit from the loop over station
1385 delete extrapSegment;
1386 // Sort track hits according to increasing Z
1387 track->GetTrackHitsPtr()->Sort();
1388 // Update track parameters at first track hit (smallest Z)
1389 trackParam1 = ((AliMUONTrackHit*)
1390 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1391 bendingMomentum = 0.;
1392 if (TMath::Abs(trackParam1->GetInverseBendingMomentum()) > 0.)
1393 bendingMomentum = TMath::Abs(1/(trackParam1->GetInverseBendingMomentum()));
1394 // Track removed if bendingMomentum not in window [min, max]
1395 if ((bendingMomentum < fMinBendingMomentum) || (bendingMomentum > fMaxBendingMomentum)) {
1397 break; // stop the search for this candidate:
1398 // exit from the loop over station
1401 // with multiple Coulomb scattering if all stations
1402 if (station == 0) track->SetFitMCS(1);
1403 // without multiple Coulomb scattering if not all stations
1404 else track->SetFitMCS(0);
1405 track->SetFitNParam(5); // with 5 parameters (momentum and position)
1406 track->SetFitStart(1); // from parameters at first hit
1408 Double_t numberOfDegFree = (2.0 * track->GetNTrackHits() - 5);
1409 if (numberOfDegFree > 0) {
1410 chi2Norm = track->GetFitFMin() / numberOfDegFree;
1414 // Track removed if normalized chi2 too high
1415 if (chi2Norm > fMaxChi2) {
1417 break; // stop the search for this candidate:
1418 // exit from the loop over station
1420 if (AliLog::GetGlobalDebugLevel() > 2) {
1421 cout << "FollowTracks: track candidate(0..): " << trackIndex
1422 << " after fit from station(0..): " << station << " to 4" << endl;
1423 track->RecursiveDump();
1425 // Track extrapolation to the vertex through the absorber (Branson)
1426 // after going through the first station
1428 trackParamVertex = *trackParam1;
1429 (&trackParamVertex)->ExtrapToVertex(0.,0.,0.);
1430 track->SetTrackParamAtVertex(&trackParamVertex);
1431 if (AliLog::GetGlobalDebugLevel() > 0) {
1432 cout << "FollowTracks: track candidate(0..): " << trackIndex
1433 << " after extrapolation to vertex" << endl;
1434 track->RecursiveDump();
1437 } // for (station = 2;...
1438 // go really to next track
1441 // Compression of track array (necessary after Remove ????)
1442 fRecTracksPtr->Compress();
1446 //__________________________________________________________________________
1447 void AliMUONTrackReconstructor::RemoveDoubleTracks(void)
1449 // To remove double tracks.
1450 // Tracks are considered identical
1451 // if they have at least half of their hits in common.
1452 // Among two identical tracks, one keeps the track with the larger number of hits
1453 // or, if these numbers are equal, the track with the minimum Chi2.
1454 AliMUONTrack *track1, *track2, *trackToRemove;
1455 Bool_t identicalTracks;
1456 Int_t hitsInCommon, nHits1, nHits2;
1457 identicalTracks = kTRUE;
1458 while (identicalTracks) {
1459 identicalTracks = kFALSE;
1460 // Loop over first track of the pair
1461 track1 = (AliMUONTrack*) fRecTracksPtr->First();
1462 while (track1 && (!identicalTracks)) {
1463 nHits1 = track1->GetNTrackHits();
1464 // Loop over second track of the pair
1465 track2 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1466 while (track2 && (!identicalTracks)) {
1467 nHits2 = track2->GetNTrackHits();
1468 // number of hits in common between two tracks
1469 hitsInCommon = track1->HitsInCommon(track2);
1470 // check for identical tracks
1471 if ((4 * hitsInCommon) >= (nHits1 + nHits2)) {
1472 identicalTracks = kTRUE;
1473 // decide which track to remove
1474 if (nHits1 > nHits2) trackToRemove = track2;
1475 else if (nHits1 < nHits2) trackToRemove = track1;
1476 else if ((track1->GetFitFMin()) < (track2->GetFitFMin()))
1477 trackToRemove = track2;
1478 else trackToRemove = track1;
1480 trackToRemove->Remove();
1482 track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
1484 track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1490 //__________________________________________________________________________
1491 void AliMUONTrackReconstructor::UpdateTrackParamAtHit()
1493 // Set track parameters after track fitting. Fill fTrackParamAtHit of AliMUONTrack's
1494 AliMUONTrack *track;
1495 AliMUONTrackHit *trackHit;
1496 AliMUONTrackParam *trackParam;
1497 track = (AliMUONTrack*) fRecTracksPtr->First();
1499 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->First();
1501 trackParam = trackHit->GetTrackParam();
1502 track->AddTrackParamAtHit(trackParam);
1503 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->After(trackHit);
1505 track = (AliMUONTrack*) fRecTracksPtr->After(track);
1510 //__________________________________________________________________________
1511 void AliMUONTrackReconstructor::UpdateHitForRecAtHit()
1513 // Set cluster parameterss after track fitting. Fill fHitForRecAtHit of AliMUONTrack's
1514 AliMUONTrack *track;
1515 AliMUONTrackHit *trackHit;
1516 AliMUONHitForRec *hitForRec;
1517 track = (AliMUONTrack*) fRecTracksPtr->First();
1519 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->First();
1521 hitForRec = trackHit->GetHitForRecPtr();
1522 track->AddHitForRecAtHit(hitForRec);
1523 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->After(trackHit);
1525 track = (AliMUONTrack*) fRecTracksPtr->After(track);
1530 //__________________________________________________________________________
1531 void AliMUONTrackReconstructor::FillMUONTrack()
1533 // Set track parameters at hits for Kalman track. Fill fTrackParamAtHit of AliMUONTrack's
1534 AliMUONTrackK *track;
1535 track = (AliMUONTrackK*) fRecTracksPtr->First();
1537 track->FillMUONTrack();
1538 track = (AliMUONTrackK*) fRecTracksPtr->After(track);
1543 //__________________________________________________________________________
1544 void AliMUONTrackReconstructor::EventDump(void)
1546 // Dump reconstructed event (track parameters at vertex and at first hit),
1547 // and the particle parameters
1549 AliMUONTrack *track;
1550 AliMUONTrackParam *trackParam, *trackParam1;
1551 Double_t bendingSlope, nonBendingSlope, pYZ;
1552 Double_t pX, pY, pZ, x, y, z, c;
1553 Int_t np, trackIndex, nTrackHits;
1555 AliDebug(1,"****** enter EventDump ******");
1556 AliDebug(1, Form("Number of Reconstructed tracks : %d", fNRecTracks));
1558 fRecTracksPtr->Compress(); // for simple loop without "Next" since no hole
1559 // Loop over reconstructed tracks
1560 for (trackIndex = 0; trackIndex < fNRecTracks; trackIndex++) {
1561 if (fTrackMethod != 1) continue; //AZ - skip the rest for now
1562 AliDebug(1, Form("track number: %d", trackIndex));
1563 // function for each track for modularity ????
1564 track = (AliMUONTrack*) ((*fRecTracksPtr)[trackIndex]);
1565 nTrackHits = track->GetNTrackHits();
1566 AliDebug(1, Form("Number of track hits: %d ", nTrackHits));
1567 // track parameters at Vertex
1568 trackParam = track->GetTrackParamAtVertex();
1569 x = trackParam->GetNonBendingCoor();
1570 y = trackParam->GetBendingCoor();
1571 z = trackParam->GetZ();
1572 bendingSlope = trackParam->GetBendingSlope();
1573 nonBendingSlope = trackParam->GetNonBendingSlope();
1574 pYZ = 1/TMath::Abs(trackParam->GetInverseBendingMomentum());
1575 pZ = pYZ/TMath::Sqrt(1+bendingSlope*bendingSlope);
1576 pX = pZ * nonBendingSlope;
1577 pY = pZ * bendingSlope;
1578 c = TMath::Sign(1.0, trackParam->GetInverseBendingMomentum());
1579 AliDebug(1, Form("Track parameters at Vertex z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1580 z, x, y, pX, pY, pZ, c));
1582 // track parameters at first hit
1583 trackParam1 = ((AliMUONTrackHit*)
1584 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1585 x = trackParam1->GetNonBendingCoor();
1586 y = trackParam1->GetBendingCoor();
1587 z = trackParam1->GetZ();
1588 bendingSlope = trackParam1->GetBendingSlope();
1589 nonBendingSlope = trackParam1->GetNonBendingSlope();
1590 pYZ = 1/TMath::Abs(trackParam1->GetInverseBendingMomentum());
1591 pZ = pYZ/TMath::Sqrt(1.0 + bendingSlope * bendingSlope);
1592 pX = pZ * nonBendingSlope;
1593 pY = pZ * bendingSlope;
1594 c = TMath::Sign(1.0, trackParam1->GetInverseBendingMomentum());
1595 AliDebug(1, Form("track parameters at z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1596 z, x, y, pX, pY, pZ, c));
1598 // informations about generated particles
1599 np = gAlice->GetMCApp()->GetNtrack();
1600 printf(" **** number of generated particles: %d \n", np);
1602 // for (Int_t iPart = 0; iPart < np; iPart++) {
1603 // p = gAlice->Particle(iPart);
1604 // printf(" particle %d: type= %d px= %f py= %f pz= %f pdg= %d\n",
1605 // iPart, p->GetPdgCode(), p->Px(), p->Py(), p->Pz(), p->GetPdgCode());
1611 //__________________________________________________________________________
1612 void AliMUONTrackReconstructor::EventDumpTrigger(void)
1614 // Dump reconstructed trigger event
1615 // and the particle parameters
1617 AliMUONTriggerTrack *triggertrack ;
1618 Int_t nTriggerTracks = fMUONData->RecTriggerTracks()->GetEntriesFast();
1620 AliDebug(1, "****** enter EventDumpTrigger ******");
1621 AliDebug(1, Form("Number of Reconstructed tracks : %d ", nTriggerTracks));
1623 // Loop over reconstructed tracks
1624 for (Int_t trackIndex = 0; trackIndex < nTriggerTracks; trackIndex++) {
1625 triggertrack = (AliMUONTriggerTrack*)fMUONData->RecTriggerTracks()->At(trackIndex);
1626 printf(" trigger track number %i x11=%f y11=%f thetax=%f thetay=%f \n",
1628 triggertrack->GetX11(),triggertrack->GetY11(),
1629 triggertrack->GetThetax(),triggertrack->GetThetay());
1633 //__________________________________________________________________________
1634 void AliMUONTrackReconstructor::MakeTrackCandidatesK(void)
1636 // To make initial tracks for Kalman filter from the list of segments
1638 AliMUONSegment *segment;
1639 AliMUONTrackK *trackK;
1641 AliDebug(1,"Enter MakeTrackCandidatesK");
1643 AliMUONTrackK a(this, fHitsForRecPtr);
1644 // Loop over stations(1...) 5 and 4
1645 for (istat=4; istat>=3; istat--) {
1646 // Loop over segments in the station
1647 for (iseg=0; iseg<fNSegments[istat]; iseg++) {
1648 // Transform segments to tracks and evaluate covariance matrix
1649 segment = (AliMUONSegment*) ((*fSegmentsPtr[istat])[iseg]);
1650 trackK = new ((*fRecTracksPtr)[fNRecTracks++]) AliMUONTrackK(segment);
1651 } // for (iseg=0;...)
1652 } // for (istat=4;...)
1656 //__________________________________________________________________________
1657 void AliMUONTrackReconstructor::FollowTracksK(void)
1659 // Follow tracks using Kalman filter
1661 Int_t icand, ichamBeg = 0, ichamEnd, chamBits;
1662 Double_t zDipole1, zDipole2;
1663 AliMUONTrackK *trackK;
1664 AliMUONHitForRec *hit;
1665 AliMUONRawCluster *clus;
1666 TClonesArray *rawclusters;
1667 clus = 0; rawclusters = 0;
1669 zDipole1 = GetSimpleBPosition() + GetSimpleBLength()/2;
1670 zDipole2 = zDipole1 - GetSimpleBLength();
1673 trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[0]);
1674 if (trackK->DebugLevel() > 0) {
1675 for (Int_t i1=0; i1<fNHitsForRec; i1++) {
1676 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[i1]);
1677 //if (hit->GetTTRTrack() > 1 || hit->GetTrackRefSignal() == 0) continue;
1678 printf(" Hit # %d %10.4f %10.4f %10.4f",
1679 hit->GetChamberNumber(), hit->GetBendingCoor(),
1680 hit->GetNonBendingCoor(), hit->GetZ());
1681 if (fRecTrackRefHits) {
1682 // from track ref hits
1683 printf(" %3d %3d \n", hit->GetTrackRefSignal(), hit->GetTTRTrack());
1685 // from raw clusters
1686 rawclusters = fMUONData->RawClusters(hit->GetChamberNumber());
1687 clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(hit->
1689 printf(" %d", clus->GetTrack(1));
1690 if (clus->GetTrack(2) != -1) printf(" %d \n", clus->GetTrack(2));
1692 } // if (fRecTrackRefHits)
1694 } // if (trackK->DebugLevel() > 0)
1698 nSeeds = fNRecTracks; // starting number of seeds
1699 // Loop over track candidates
1700 while (icand < fNRecTracks-1) {
1702 if (trackK->DebugLevel()>0) cout << " *** Kalman track candidate No. " << icand << endl;
1703 trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]);
1704 if (trackK->GetRecover() < 0) continue; // failed track
1706 // Discard candidate which will produce the double track
1709 ok = CheckCandidateK(icand,nSeeds);
1711 trackK->SetRecover(-1); // mark candidate to be removed
1718 if (trackK->GetRecover() == 0) hit = (AliMUONHitForRec*)
1719 trackK->GetTrackHits()->Last(); // last hit
1720 else hit = trackK->GetHitLastOk(); // hit where track stopped
1721 if (hit) ichamBeg = hit->GetChamberNumber();
1723 // Check propagation direction
1724 if (!hit) { ichamBeg = ichamEnd; AliFatal(" ??? "); }
1725 else if (trackK->GetTrackDir() < 0) {
1726 ichamEnd = 9; // forward propagation
1727 ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kFALSE,zDipole1,zDipole2);
1729 ichamBeg = ichamEnd;
1730 ichamEnd = 6; // backward propagation
1731 // Change weight matrix and zero fChi2 for backpropagation
1732 trackK->StartBack();
1733 trackK->SetTrackDir(1);
1734 ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kTRUE,zDipole1,zDipole2);
1735 ichamBeg = ichamEnd;
1739 if (trackK->GetBPFlag()) {
1741 ichamEnd = 6; // backward propagation
1742 // Change weight matrix and zero fChi2 for backpropagation
1743 trackK->StartBack();
1744 ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kTRUE,zDipole1,zDipole2);
1745 ichamBeg = ichamEnd;
1751 trackK->SetTrackDir(1);
1752 trackK->SetBPFlag(kFALSE);
1753 ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kFALSE,zDipole1,zDipole2);
1755 if (!ok) { trackK->SetRecover(-1); continue; } // mark candidate to be removed
1758 if (trackK->GetRecover() >= 0) {
1759 ok = trackK->Smooth();
1760 if (!ok) trackK->SetRecover(-1); // mark candidate to be removed
1763 // Majority 3 of 4 in first 2 stations
1766 Double_t chi2max = 0;
1767 for (Int_t i=0; i<trackK->GetNTrackHits(); i++) {
1768 hit = (AliMUONHitForRec*) (*trackK->GetTrackHits())[i];
1769 chamBits |= BIT(hit->GetChamberNumber());
1770 if (trackK->GetChi2PerPoint(i) > chi2max) chi2max = trackK->GetChi2PerPoint(i);
1772 if (!((chamBits&3)==3 || (chamBits>>2&3)==3) && chi2max > 25) {
1773 //trackK->Recover();
1774 trackK->SetRecover(-1); //mark candidate to be removed
1777 if (ok) trackK->SetTrackQuality(0); // compute "track quality"
1780 for (Int_t i=0; i<fNRecTracks; i++) {
1781 trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[i]);
1782 if (trackK->GetRecover() < 0) fRecTracksPtr->RemoveAt(i);
1785 // Compress TClonesArray
1786 fRecTracksPtr->Compress();
1787 fNRecTracks = fRecTracksPtr->GetEntriesFast();
1791 //__________________________________________________________________________
1792 Bool_t AliMUONTrackReconstructor::CheckCandidateK(Int_t icand, Int_t nSeeds) const
1794 // Discards track candidate if it will produce the double track (having
1795 // the same seed segment hits as hits of a good track found before)
1796 AliMUONTrackK *track1, *track2;
1797 AliMUONHitForRec *hit1, *hit2, *hit;
1799 track1 = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]);
1800 hit1 = (AliMUONHitForRec*) (*track1->GetTrackHits())[0]; // 1'st hit
1801 hit2 = (AliMUONHitForRec*) (*track1->GetTrackHits())[1]; // 2'nd hit
1803 for (Int_t i=0; i<icand; i++) {
1804 track2 = (AliMUONTrackK*) ((*fRecTracksPtr)[i]);
1805 //if (track2->GetRecover() < 0) continue;
1806 if (track2->GetRecover() < 0 && icand >= nSeeds) continue;
1808 if (track1->GetStartSegment() == track2->GetStartSegment()) {
1812 for (Int_t j=0; j<track2->GetNTrackHits(); j++) {
1813 hit = (AliMUONHitForRec*) (*track2->GetTrackHits())[j];
1814 if (hit == hit1 || hit == hit2) {
1816 if (nSame == 2) return kFALSE;
1818 } // for (Int_t j=0;
1820 } // for (Int_t i=0;
1824 //__________________________________________________________________________
1825 void AliMUONTrackReconstructor::RemoveDoubleTracksK(void)
1827 // Removes double tracks (sharing more than half of their hits). Keeps
1828 // the track with higher quality
1829 AliMUONTrackK *track1, *track2, *trackToKill;
1831 // Sort tracks according to their quality
1832 fRecTracksPtr->Sort();
1834 // Loop over first track of the pair
1835 track1 = (AliMUONTrackK*) fRecTracksPtr->First();
1836 Int_t debug = track1->DebugLevel();
1838 // Loop over second track of the pair
1839 track2 = (AliMUONTrackK*) fRecTracksPtr->After(track1);
1841 // Check whether or not to keep track2
1842 if (!track2->KeepTrack(track1)) {
1843 if (debug >= 0) cout << " Killed track: " << 1/(*track2->GetTrackParameters())(4,0) <<
1844 " " << track2->GetTrackQuality() << endl;
1845 trackToKill = track2;
1846 track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2);
1847 trackToKill->Kill();
1848 fRecTracksPtr->Compress();
1849 } else track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2);
1851 track1 = (AliMUONTrackK*) fRecTracksPtr->After(track1);
1854 fNRecTracks = fRecTracksPtr->GetEntriesFast();
1855 if (debug >= 0) cout << " Number of Kalman tracks: " << fNRecTracks << endl;
1858 //__________________________________________________________________________
1859 void AliMUONTrackReconstructor::GoToVertex(void)
1861 // Propagates track to the vertex thru absorber
1862 // (using Branson correction for now)
1866 for (Int_t i=0; i<fNRecTracks; i++) {
1867 //((AliMUONTrackK*)(*fRecTracksPtr)[i])->Branson();
1868 ((AliMUONTrackK*)(*fRecTracksPtr)[i])->SetTrackQuality(1); // compute Chi2
1869 //((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToZ(zVertex); // w/out absorber
1870 ((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToVertex(1); // with absorber
1874 //__________________________________________________________________________
1875 void AliMUONTrackReconstructor::SetTrackMethod(Int_t iTrackMethod)
1877 // Set track method and recreate track container if necessary
1879 fTrackMethod = TMath::Min (iTrackMethod, 3);
1880 fTrackMethod = TMath::Max (fTrackMethod, 1);
1881 if (fTrackMethod != 1) {
1882 if (fRecTracksPtr) delete fRecTracksPtr;
1883 fRecTracksPtr = new TClonesArray("AliMUONTrackK", 10);
1884 if (fTrackMethod == 2) cout << " *** Tracking with the Kalman filter *** " << endl;
1885 else cout << " *** Combined cluster / track finder ***" << endl;
1886 } else cout << " *** Original tracking *** " << endl;