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 event reconstructor in ALICE
22 // This class contains as data:
23 // * the parameters for the event reconstruction
24 // * a pointer to the array of hits to be reconstructed (the event)
25 // * a pointer to the array of segments made with these hits inside each station
26 // * a pointer to the array of reconstructed tracks
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>
39 #include <TMatrixD.h> //AZ
41 #include "AliMUONEventReconstructor.h"
43 #include "AliMUONConstants.h"
44 #include "AliMUONHitForRec.h"
45 #include "AliMUONTriggerTrack.h"
46 #include "AliMUONTriggerCircuit.h"
47 #include "AliMUONRawCluster.h"
48 #include "AliMUONLocalTrigger.h"
49 #include "AliMUONGlobalTrigger.h"
50 #include "AliMUONRecoEvent.h"
51 #include "AliMUONSegment.h"
52 #include "AliMUONTrack.h"
53 #include "AliMUONTrackHit.h"
55 #include "AliRun.h" // for gAlice
56 #include "AliRunLoader.h"
57 #include "AliLoader.h"
58 #include "AliMUONTrackK.h" //AZ
61 #include "AliTrackReference.h"
63 //************* Defaults parameters for reconstruction
64 const Double_t AliMUONEventReconstructor::fgkDefaultMinBendingMomentum = 3.0;
65 const Double_t AliMUONEventReconstructor::fgkDefaultMaxBendingMomentum = 2000.0;
66 const Double_t AliMUONEventReconstructor::fgkDefaultMaxChi2 = 100.0;
67 const Double_t AliMUONEventReconstructor::fgkDefaultMaxSigma2Distance = 16.0;
68 const Double_t AliMUONEventReconstructor::fgkDefaultBendingResolution = 0.01;
69 const Double_t AliMUONEventReconstructor::fgkDefaultNonBendingResolution = 0.144;
70 const Double_t AliMUONEventReconstructor::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 AliMUONEventReconstructor::fgkDefaultSimpleBValue = 7.0;
78 const Double_t AliMUONEventReconstructor::fgkDefaultSimpleBLength = 428.0;
79 const Double_t AliMUONEventReconstructor::fgkDefaultSimpleBPosition = 1019.0;
80 const Int_t AliMUONEventReconstructor::fgkDefaultRecTrackRefHits = 0;
81 const Double_t AliMUONEventReconstructor::fgkDefaultEfficiency = 0.95;
83 ClassImp(AliMUONEventReconstructor) // Class implementation in ROOT context
85 //__________________________________________________________________________
86 AliMUONEventReconstructor::AliMUONEventReconstructor(AliLoader* loader)
89 // Constructor for class AliMUONEventReconstructor
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 ????
107 fRecTriggerTracksPtr = new TClonesArray("AliMUONTriggerTrack", 10);
108 fNRecTriggerTracks = 0; // really needed or GetEntriesFast sufficient ????
109 // Memory allocation for the TClonesArray of hits on reconstructed tracks
110 // Is 100 the right size ????
111 //AZ fRecTrackHitsPtr = new TClonesArray("AliMUONTrack", 100);
112 fRecTrackHitsPtr = new TClonesArray("AliMUONTrackHit", 100);
113 fNRecTrackHits = 0; // really needed or GetEntriesFast sufficient ????
115 // Sign of fSimpleBValue according to sign of Bx value at (50,50,-950).
117 x[0] = 50.; x[1] = 50.; x[2] = -950.;
118 gAlice->Field()->Field(x, b);
119 fSimpleBValue = TMath::Sign(fSimpleBValue,(Double_t) b[0]);
120 fSimpleBPosition = TMath::Sign(fSimpleBPosition,(Double_t) x[2]);
121 // See how to get fSimple(BValue, BLength, BPosition)
122 // automatically calculated from the actual magnetic field ????
124 AliDebug(1,"AliMUONEventReconstructor constructed with defaults");
125 if ( AliLog::GetGlobalDebugLevel()>0) Dump();
126 AliDebug(1,"Magnetic field from root file:");
127 if ( AliLog::GetGlobalDebugLevel()>0) gAlice->Field()->Dump();
130 // Initializions for track ref. background events
131 fBkgTrackRefFile = 0;
133 fBkgTrackRefParticles = 0;
135 fBkgTrackRefEventNumber = -1;
137 // Initialize to 0 pointers to RecoEvent, tree and tree file
142 // initialize loader's
145 // initialize container
146 fMUONData = new AliMUONData(fLoader,"MUON","MUON");
150 //__________________________________________________________________________
151 AliMUONEventReconstructor::AliMUONEventReconstructor (const AliMUONEventReconstructor& rhs)
154 // Protected copy constructor
156 AliFatal("Not implemented.");
159 AliMUONEventReconstructor &
160 AliMUONEventReconstructor::operator=(const AliMUONEventReconstructor& rhs)
162 // Protected assignement operator
164 if (this == &rhs) return *this;
166 AliFatal("Not implemented.");
171 //__________________________________________________________________________
172 AliMUONEventReconstructor::~AliMUONEventReconstructor(void)
174 // Destructor for class AliMUONEventReconstructor
179 // if (fEventTree) delete fEventTree;
180 if (fRecoEvent) delete fRecoEvent;
181 delete fHitsForRecPtr; // Correct destruction of everything ???? or delete [] ????
182 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++)
183 delete fSegmentsPtr[st]; // Correct destruction of everything ????
187 //__________________________________________________________________________
188 void AliMUONEventReconstructor::SetReconstructionParametersToDefaults(void)
190 // Set reconstruction parameters to default values
191 // Would be much more convenient with a structure (or class) ????
192 fMinBendingMomentum = fgkDefaultMinBendingMomentum;
193 fMaxBendingMomentum = fgkDefaultMaxBendingMomentum;
194 fMaxChi2 = fgkDefaultMaxChi2;
195 fMaxSigma2Distance = fgkDefaultMaxSigma2Distance;
197 // ******** Parameters for making HitsForRec
199 // like in TRACKF_STAT:
200 // 2 degrees for stations 1 and 2, or ch(0...) from 0 to 3;
201 // 30 cm for stations 3 to 5, or ch(0...) from 4 to 9
202 for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
203 if (ch < 4) fRMin[ch] = TMath::Abs(AliMUONConstants::DefaultChamberZ(ch)) *
204 2.0 * TMath::Pi() / 180.0;
205 else fRMin[ch] = 30.0;
206 // maximum radius at 10 degrees and Z of chamber
207 fRMax[ch] = TMath::Abs(AliMUONConstants::DefaultChamberZ(ch)) *
208 10.0 * TMath::Pi() / 180.0;
211 // ******** Parameters for making segments
212 // should be parametrized ????
213 // according to interval between chambers in a station ????
214 // Maximum distance in non bending plane
215 // 5 * 0.22 just to remember the way it was made in TRACKF_STAT
217 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++)
218 fSegmentMaxDistNonBending[st] = 5. * 0.22;
219 // Maximum distance in bending plane:
220 // values from TRACKF_STAT, corresponding to (J psi 20cm),
221 // scaled to the real distance between chambers in a station
222 fSegmentMaxDistBending[0] = TMath::Abs( 1.5 *
223 (AliMUONConstants::DefaultChamberZ(1) - AliMUONConstants::DefaultChamberZ(0)) / 20.0);
224 fSegmentMaxDistBending[1] = TMath::Abs( 1.5 *
225 (AliMUONConstants::DefaultChamberZ(3) - AliMUONConstants::DefaultChamberZ(2)) / 20.0);
226 fSegmentMaxDistBending[2] = TMath::Abs( 3.0 *
227 (AliMUONConstants::DefaultChamberZ(5) - AliMUONConstants::DefaultChamberZ(4)) / 20.0);
228 fSegmentMaxDistBending[3] = TMath::Abs( 6.0 *
229 (AliMUONConstants::DefaultChamberZ(7) - AliMUONConstants::DefaultChamberZ(6)) / 20.0);
230 fSegmentMaxDistBending[4] = TMath::Abs( 6.0 *
231 (AliMUONConstants::DefaultChamberZ(9) - AliMUONConstants::DefaultChamberZ(8)) / 20.0);
234 fBendingResolution = fgkDefaultBendingResolution;
235 fNonBendingResolution = fgkDefaultNonBendingResolution;
236 fChamberThicknessInX0 = fgkDefaultChamberThicknessInX0;
237 fSimpleBValue = fgkDefaultSimpleBValue;
238 fSimpleBLength = fgkDefaultSimpleBLength;
239 fSimpleBPosition = fgkDefaultSimpleBPosition;
240 fRecTrackRefHits = fgkDefaultRecTrackRefHits;
241 fEfficiency = fgkDefaultEfficiency;
245 //__________________________________________________________________________
246 Double_t AliMUONEventReconstructor::GetImpactParamFromBendingMomentum(Double_t BendingMomentum) const
248 // Returns impact parameter at vertex in bending plane (cm),
249 // from the signed bending momentum "BendingMomentum" in bending plane (GeV/c),
250 // using simple values for dipole magnetic field.
251 // The sign of "BendingMomentum" is the sign of the charge.
252 return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
256 //__________________________________________________________________________
257 Double_t AliMUONEventReconstructor::GetBendingMomentumFromImpactParam(Double_t ImpactParam) const
259 // Returns signed bending momentum in bending plane (GeV/c),
260 // the sign being the sign of the charge for particles moving forward in Z,
261 // from the impact parameter "ImpactParam" at vertex in bending plane (cm),
262 // using simple values for dipole magnetic field.
263 return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
267 //__________________________________________________________________________
268 void AliMUONEventReconstructor::SetBkgTrackRefFile(Text_t *BkgTrackRefFileName)
270 // Set background file ... for track ref. hits
271 // Must be called after having loaded the firts signal event
272 AliDebug(1,Form("Enter SetBkgTrackRefFile with BkgTrackRefFileName %s",BkgTrackRefFileName));
274 if (strlen(BkgTrackRefFileName)) {
275 // BkgTrackRefFileName not empty: try to open the file
277 if(AliLog::GetGlobalDebugLevel()>1) {
278 cout << "Before File(Bkg)" << endl; gDirectory->Dump();
280 fBkgTrackRefFile = new TFile(BkgTrackRefFileName);
281 if(AliLog::GetGlobalDebugLevel()>1) {
282 cout << "After File(Bkg)" << endl; gDirectory->Dump();
284 if (fBkgTrackRefFile-> IsOpen()) {
285 if(AliLog::GetGlobalDebugLevel()>0) {
286 cout << "Background for Track ref. hits in file: ``" << BkgTrackRefFileName
287 << "'' successfully opened" << endl;}
290 cout << "Background for Track Ref. hits in file: " << BkgTrackRefFileName << endl;
291 cout << "NOT FOUND: EXIT" << endl;
292 exit(0); // right instruction for exit ????
294 // Arrays for "particles" and "hits"
295 fBkgTrackRefParticles = new TClonesArray("TParticle", 200);
296 // Event number to -1 for initialization
297 fBkgTrackRefEventNumber = -1;
298 // Back to the signal file:
299 // first signal event must have been loaded previously,
300 // otherwise, Segmentation violation at the next instruction
301 // How is it possible to do smething better ????
302 ((gAlice->TreeK())->GetCurrentFile())->cd();
303 if(AliLog::GetGlobalDebugLevel()>1) cout << "After cd(gAlice)" << endl; gDirectory->Dump();
308 //__________________________________________________________________________
309 void AliMUONEventReconstructor::NextBkgTrackRefEvent(void)
311 // Get next event in background file for track ref. hits
312 // Goes back to event number 0 when end of file is reached
315 AliDebug(1,"Enter NextBkgTrackRefEvent");
316 // Clean previous event
317 if(fBkgTrackRefTK) delete fBkgTrackRefTK;
318 fBkgTrackRefTK = NULL;
319 if(fBkgTrackRefParticles) fBkgTrackRefParticles->Clear();
320 if(fBkgTrackRefTTR) delete fBkgTrackRefTTR;
321 fBkgTrackRefTTR = NULL;
322 // Increment event number
323 fBkgTrackRefEventNumber++;
324 // Get access to Particles and Hits for event from background file
325 if (AliLog::GetGlobalDebugLevel()>1) cout << "Before cd(Bkg)" << endl; gDirectory->Dump();
326 fBkgTrackRefFile->cd();
327 if (AliLog::GetGlobalDebugLevel()>1) cout << "After cd(Bkg)" << endl; gDirectory->Dump();
328 // Particles: TreeK for event and branch "Particles"
329 sprintf(treeName, "TreeK%d", fBkgTrackRefEventNumber);
330 fBkgTrackRefTK = (TTree*)gDirectory->Get(treeName);
331 if (!fBkgTrackRefTK) {
333 AliDebug(1,Form("Cannot find Kine Tree for background event: %d",fBkgTrackRefEventNumber));
334 AliDebug(1,"Goes back to event 0");
336 fBkgTrackRefEventNumber = 0;
337 sprintf(treeName, "TreeK%d", fBkgTrackRefEventNumber);
338 fBkgTrackRefTK = (TTree*)gDirectory->Get(treeName);
339 if (!fBkgTrackRefTK) {
340 AliError(Form("cannot find Kine Tree for background event: %d",fBkgTrackRefEventNumber));
345 fBkgTrackRefTK->SetBranchAddress("Particles", &fBkgTrackRefParticles);
346 fBkgTrackRefTK->GetEvent(0); // why event 0 ???? necessary ????
347 // Hits: TreeH for event and branch "MUON"
348 sprintf(treeName, "TreeTR%d", fBkgTrackRefEventNumber);
349 fBkgTrackRefTTR = (TTree*)gDirectory->Get(treeName);
350 if (!fBkgTrackRefTTR) {
351 AliError(Form("cannot find Hits Tree for background event: %d",fBkgTrackRefEventNumber));
354 fBkgTrackRefTTR->GetEntries(); // necessary ????
355 // Back to the signal file
356 ((gAlice->TreeK())->GetCurrentFile())->cd();
357 if (AliLog::GetGlobalDebugLevel()>1)
358 cout << "After cd(gAlice)" << endl; gDirectory->Dump();
362 //__________________________________________________________________________
363 void AliMUONEventReconstructor::EventReconstruct(void)
365 // To reconstruct one event
366 AliDebug(1,"Enter EventReconstruct");
367 MakeEventToBeReconstructed();
370 if (fMUONData->IsTriggerTrackBranchesInTree())
371 ValidateTracksWithTrigger();
373 // Add tracks to MUON data container
374 for(Int_t i=0; i<GetNRecTracks(); i++) {
375 AliMUONTrack * track = (AliMUONTrack*) GetRecTracksPtr()->At(i);
376 fMUONData->AddRecTrack(*track);
382 //__________________________________________________________________________
383 void AliMUONEventReconstructor::EventReconstructTrigger(void)
385 // To reconstruct one event
386 AliDebug(1,"Enter EventReconstructTrigger");
391 //__________________________________________________________________________
392 void AliMUONEventReconstructor::ResetHitsForRec(void)
394 // To reset the array and the number of HitsForRec,
395 // and also the number of HitsForRec
396 // and the index of the first HitForRec per chamber
397 if (fHitsForRecPtr) fHitsForRecPtr->Clear();
399 for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++)
400 fNHitsForRecPerChamber[ch] = fIndexOfFirstHitForRecPerChamber[ch] = 0;
404 //__________________________________________________________________________
405 void AliMUONEventReconstructor::ResetSegments(void)
407 // To reset the TClonesArray of segments and the number of Segments
409 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) {
410 if (fSegmentsPtr[st]) fSegmentsPtr[st]->Clear();
416 //__________________________________________________________________________
417 void AliMUONEventReconstructor::ResetTracks(void)
419 // To reset the TClonesArray of reconstructed tracks
420 if (fRecTracksPtr) fRecTracksPtr->Delete();
421 // Delete in order that the Track destructors are called,
422 // hence the space for the TClonesArray of pointers to TrackHit's is freed
427 //__________________________________________________________________________
428 void AliMUONEventReconstructor::ResetTriggerTracks(void)
430 // To reset the TClonesArray of reconstructed trigger tracks
431 if (fRecTriggerTracksPtr) fRecTriggerTracksPtr->Delete();
432 // Delete in order that the Track destructors are called,
433 // hence the space for the TClonesArray of pointers to TrackHit's is freed
434 fNRecTriggerTracks = 0;
438 //__________________________________________________________________________
439 void AliMUONEventReconstructor::ResetTrackHits(void)
441 // To reset the TClonesArray of hits on reconstructed tracks
442 if (fRecTrackHitsPtr) fRecTrackHitsPtr->Clear();
447 //__________________________________________________________________________
448 void AliMUONEventReconstructor::MakeEventToBeReconstructed(void)
450 // To make the list of hits to be reconstructed,
451 // either from the track ref. hits or from the raw clusters
452 // according to the parameter set for the reconstructor
453 // TString evfoldname = AliConfig::GetDefaultEventFolderName();//to be interfaced properly
455 // AliRunLoader* rl = AliRunLoader::GetRunLoader(evfoldname);
458 // Error("MakeEventToBeReconstructed",
459 // "Can not find Run Loader in Event Folder named %s.",
460 // evfoldname.Data());
463 // AliLoader* gime = rl->GetLoader("MUONLoader");
466 // Error("MakeEventToBeReconstructed","Can not get MUON Loader from Run Loader.");
469 AliRunLoader *runLoader = fLoader->GetRunLoader();
471 AliDebug(1,"Enter MakeEventToBeReconstructed");
473 if (fRecTrackRefHits == 1) {
474 // Reconstruction from track ref. hits
475 // Back to the signal file
476 TTree* treeTR = runLoader->TreeTR();
479 Int_t retval = runLoader->LoadTrackRefs("READ");
482 AliError("Error occured while loading hits.");
485 treeTR = runLoader->TreeTR();
488 AliError("Can not get TreeTR");
493 AddHitsForRecFromTrackRef(treeTR,1);
496 AddHitsForRecFromTrackRef(fBkgTrackRefTTR,0);
497 // Sort HitsForRec in increasing order with respect to chamber number
498 SortHitsForRecWithIncreasingChamber();
501 // Reconstruction from raw clusters
502 // AliMUON *MUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
503 // Security on MUON ????
504 // TreeR assumed to be be "prepared" in calling function
505 // by "MUON->GetTreeR(nev)" ????
506 TTree *treeR = fLoader->TreeR();
507 fMUONData->SetTreeAddress("RC");
508 AddHitsForRecFromRawClusters(treeR);
509 // No sorting: it is done automatically in the previous function
512 AliDebug(1,"End of MakeEventToBeReconstructed");
513 if (AliLog::GetGlobalDebugLevel() > 0) {
514 AliDebug(1, Form("NHitsForRec: %d",fNHitsForRec));
515 for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
516 AliDebug(1, Form("Chamber(0...): %d",ch));
517 AliDebug(1, Form("NHitsForRec: %d", fNHitsForRecPerChamber[ch]));
518 AliDebug(1, Form("Index(first HitForRec): %d", fIndexOfFirstHitForRecPerChamber[ch]));
519 for (Int_t hit = fIndexOfFirstHitForRecPerChamber[ch];
520 hit < fIndexOfFirstHitForRecPerChamber[ch] + fNHitsForRecPerChamber[ch];
522 AliDebug(1, Form("HitForRec index(0...): %d",hit));
523 ((*fHitsForRecPtr)[hit])->Dump();
530 //__________________________________________________________________________
531 void AliMUONEventReconstructor::AddHitsForRecFromTrackRef(TTree *TTR, Int_t signal)
533 // To add to the list of hits for reconstruction
534 // the signal hits from a track reference tree TreeTR.
535 TClonesArray *listOfTrackRefs = NULL;
536 AliTrackReference *trackRef;
539 AliDebug(2,Form("Enter AddHitsForRecFromTrackRef with TTR: %d", TTR));
540 if (TTR == NULL) return;
542 listOfTrackRefs = CleanTrackRefs(TTR);
544 Int_t ntracks = listOfTrackRefs->GetEntriesFast();
546 AliDebug(2,Form("ntracks: %d", ntracks));
548 for (Int_t index = 0; index < ntracks; index++) {
549 trackRef = (AliTrackReference*) listOfTrackRefs->At(index);
550 track = trackRef->GetTrack();
552 NewHitForRecFromTrackRef(trackRef,track,signal);
553 } // end of track ref.
555 listOfTrackRefs->Delete();
556 delete listOfTrackRefs;
561 //__________________________________________________________________________
562 AliMUONHitForRec* AliMUONEventReconstructor::NewHitForRecFromTrackRef(AliTrackReference* Hit, Int_t TrackNumber, Int_t Signal)
564 // To make a new hit for reconstruction from a track ref. hit pointed to by "Hit",
565 // with the track numbered "TrackNumber",
566 // either from signal ("Signal" = 1) or background ("Signal" = 0) event.
567 // Selects hits in tracking (not trigger) chambers.
568 // Takes into account the efficiency (fEfficiency)
569 // and the smearing from resolution (fBendingResolution and fNonBendingResolution).
570 // Adds a condition on the radius between RMin and RMax
571 // to better simulate the real chambers.
572 // Returns the pointer to the new hit for reconstruction,
573 // or NULL in case of inefficiency or non tracking chamber or bad radius.
574 // No condition on at most 20 cm from a muon from a resonance
575 // like in Fortran TRACKF_STAT.
576 AliMUONHitForRec* hitForRec;
577 Double_t bendCoor, nonBendCoor, radius;
578 Int_t chamber = AliMUONConstants::ChamberNumber(Hit->Z()); // chamber(0...)
579 if (chamber < 0) return NULL;
580 // only in tracking chambers (fChamber starts at 1)
581 if (chamber >= AliMUONConstants::NTrackingCh()) return NULL;
582 // only if hit is efficient (keep track for checking ????)
583 if (gRandom->Rndm() > fEfficiency) return NULL;
584 // only if radius between RMin and RMax
586 nonBendCoor = Hit->X();
587 radius = TMath::Sqrt((bendCoor * bendCoor) + (nonBendCoor * nonBendCoor));
588 // This cut is not needed with a realistic chamber geometry !!!!
589 // if ((radius < fRMin[chamber]) || (radius > fRMax[chamber])) return NULL;
590 // new AliMUONHitForRec from track ref. hit and increment number of AliMUONHitForRec's
591 hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(Hit);
593 // add smearing from resolution
594 hitForRec->SetBendingCoor(bendCoor + gRandom->Gaus(0., fBendingResolution));
595 hitForRec->SetNonBendingCoor(nonBendCoor
596 + gRandom->Gaus(0., fNonBendingResolution));
597 // // !!!! without smearing
598 // hitForRec->SetBendingCoor(bendCoor);
599 // hitForRec->SetNonBendingCoor(nonBendCoor);
600 // more information into HitForRec
601 // resolution: angular effect to be added here ????
602 hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
603 hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
605 hitForRec->SetTTRTrack(TrackNumber);
606 hitForRec->SetTrackRefSignal(Signal);
607 if (AliLog::GetGlobalDebugLevel()> 1) {
608 AliDebug(2,Form("Track: %d", TrackNumber));
610 cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
615 //__________________________________________________________________________
616 TClonesArray* AliMUONEventReconstructor::CleanTrackRefs(TTree *treeTR)
618 // Make hits from track ref..
619 // Re-calculate hits parameters because two AliTrackReferences are recorded for
620 // each chamber (one when particle is entering + one when particle is leaving
621 // the sensitive volume)
623 AliTrackReference *trackReference;
624 Float_t x1, y1, z1, pX1, pY1, pZ1;
625 Float_t x2, y2, z2, pX2, pY2, pZ2;
626 Int_t track1, track2;
628 Float_t maxGasGap = 1.; // cm
632 AliTrackReference *trackReferenceNew = new AliTrackReference();
634 TClonesArray* trackRefs = new TClonesArray("AliTrackReference", 10);
635 TClonesArray* cleanTrackRefs = new TClonesArray("AliTrackReference", 10);
637 if (treeTR == NULL) return NULL;
638 TBranch* branch = treeTR->GetBranch("MUON");
639 if (branch == NULL) return NULL;
640 branch->SetAddress(&trackRefs);
642 Int_t nTrack = (Int_t)treeTR->GetEntries();
643 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
644 treeTR->GetEntry(iTrack);
647 while (iHit1 < trackRefs->GetEntries()) {
648 trackReference = (AliTrackReference*)trackRefs->At(iHit1);
649 x1 = trackReference->X();
650 y1 = trackReference->Y();
651 z1 = trackReference->Z();
652 pX1 = trackReference->Px();
653 pY1 = trackReference->Py();
654 pZ1 = trackReference->Pz();
655 track1 = trackReference->GetTrack();
658 for (Int_t iHit2 = iHit1+1; iHit2 < trackRefs->GetEntries(); iHit2++) {
659 trackReference = (AliTrackReference*)trackRefs->At(iHit2);
660 x2 = trackReference->X();
661 y2 = trackReference->Y();
662 z2 = trackReference->Z();
663 pX2 = trackReference->Px();
664 pY2 = trackReference->Py();
665 pZ2 = trackReference->Pz();
666 track2 = trackReference->GetTrack();
667 if (track2 == track1 && TMath::Abs(z2-z1) < maxGasGap ) {
682 pX1 /= (Float_t)nRec;
683 pY1 /= (Float_t)nRec;
684 pZ1 /= (Float_t)nRec;
685 trackReferenceNew->SetPosition(x1,y1,z1);
686 trackReferenceNew->SetMomentum(pX1,pY1,pZ1);
687 trackReferenceNew->SetTrack(track1);
688 {new ((*cleanTrackRefs)[cleanTrackRefs->GetEntriesFast()]) AliTrackReference(*trackReferenceNew);}
695 delete trackReferenceNew;
696 return cleanTrackRefs;
699 //__________________________________________________________________________
700 void AliMUONEventReconstructor::SortHitsForRecWithIncreasingChamber()
702 // Sort HitsForRec's in increasing order with respect to chamber number.
703 // Uses the function "Compare".
704 // Update the information for HitsForRec per chamber too.
705 Int_t ch, nhits, prevch;
706 fHitsForRecPtr->Sort();
707 for (ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
708 fNHitsForRecPerChamber[ch] = 0;
709 fIndexOfFirstHitForRecPerChamber[ch] = 0;
711 prevch = 0; // previous chamber
712 nhits = 0; // number of hits in current chamber
713 // Loop over HitsForRec
714 for (Int_t hit = 0; hit < fNHitsForRec; hit++) {
715 // chamber number (0...)
716 ch = ((AliMUONHitForRec*) ((*fHitsForRecPtr)[hit]))->GetChamberNumber();
717 // increment number of hits in current chamber
718 (fNHitsForRecPerChamber[ch])++;
719 // update index of first HitForRec in current chamber
720 // if chamber number different from previous one
722 fIndexOfFirstHitForRecPerChamber[ch] = hit;
729 //__________________________________________________________________________
730 void AliMUONEventReconstructor::AddHitsForRecFromRawClusters(TTree* TR)
732 // To add to the list of hits for reconstruction all the raw clusters
733 // No condition added, like in Fortran TRACKF_STAT,
734 // on the radius between RMin and RMax.
735 AliMUONHitForRec *hitForRec;
736 AliMUONRawCluster *clus;
737 Int_t iclus, nclus, nTRentries;
738 TClonesArray *rawclusters;
739 AliDebug(1,"Enter AddHitsForRecFromRawClusters");
741 nTRentries = Int_t(TR->GetEntries());
742 if (nTRentries != 1) {
743 AliError(Form("nTRentries = %d not equal to 1 ",nTRentries));
746 fLoader->TreeR()->GetEvent(0); // only one entry
748 // Loop over tracking chambers
749 for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
750 // number of HitsForRec to 0 for the chamber
751 fNHitsForRecPerChamber[ch] = 0;
752 // index of first HitForRec for the chamber
753 if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0;
754 else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec;
755 rawclusters =fMUONData->RawClusters(ch);
756 // pMUON->ResetRawClusters();
757 // TR->GetEvent((Int_t) (TR->GetEntries()) - 1); // to be checked ????
758 nclus = (Int_t) (rawclusters->GetEntries());
759 // Loop over (cathode correlated) raw clusters
760 for (iclus = 0; iclus < nclus; iclus++) {
761 clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(iclus);
762 // new AliMUONHitForRec from raw cluster
763 // and increment number of AliMUONHitForRec's (total and in chamber)
764 hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(clus);
766 (fNHitsForRecPerChamber[ch])++;
767 // more information into HitForRec
768 // resolution: info should be already in raw cluster and taken from it ????
769 hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
770 hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
771 // original raw cluster
772 hitForRec->SetChamberNumber(ch);
773 hitForRec->SetHitNumber(iclus);
774 // Z coordinate of the raw cluster (cm)
775 hitForRec->SetZ(clus->GetZ(0));
776 if (AliLog::GetGlobalDebugLevel() > 1) {
777 cout << "chamber (0...): " << ch <<
778 " raw cluster (0...): " << iclus << endl;
780 cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
782 } // end of cluster loop
783 } // end of chamber loop
784 SortHitsForRecWithIncreasingChamber(); //AZ
788 //__________________________________________________________________________
789 void AliMUONEventReconstructor::MakeSegments(void)
791 // To make the list of segments in all stations,
792 // from the list of hits to be reconstructed
793 AliDebug(1,"Enter MakeSegments");
795 // Loop over stations
796 Int_t nb = (fTrackMethod == 2) ? 3 : 0; //AZ
797 //AZ for (Int_t st = 0; st < fgkMaxMuonTrackingStations; st++)
798 for (Int_t st = nb; st < AliMUONConstants::NTrackingCh()/2; st++) //AZ
799 MakeSegmentsPerStation(st);
800 if (AliLog::GetGlobalDebugLevel() > 1) {
801 cout << "end of MakeSegments" << endl;
802 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) {
803 cout << "station(0...): " << st
804 << " Segments: " << fNSegments[st]
807 seg < fNSegments[st];
809 cout << "Segment index(0...): " << seg << endl;
810 ((*fSegmentsPtr[st])[seg])->Dump();
817 //__________________________________________________________________________
818 void AliMUONEventReconstructor::MakeSegmentsPerStation(Int_t Station)
820 // To make the list of segments in station number "Station" (0...)
821 // from the list of hits to be reconstructed.
822 // Updates "fNSegments"[Station].
823 // Segments in stations 4 and 5 are sorted
824 // according to increasing absolute value of "impact parameter"
825 AliMUONHitForRec *hit1Ptr, *hit2Ptr;
826 AliMUONSegment *segment;
828 Double_t bendingSlope, distBend, distNonBend, extBendCoor, extNonBendCoor,
829 impactParam = 0., maxImpactParam = 0., minImpactParam = 0.; // =0 to avoid compilation warnings.
830 AliDebug(1,Form("Enter MakeSegmentsPerStation (0...) %d",Station));
831 // first and second chambers (0...) in the station
832 Int_t ch1 = 2 * Station;
834 // variable true for stations downstream of the dipole:
835 // Station(0..4) equal to 3 or 4
836 if ((Station == 3) || (Station == 4)) {
838 // maximum impact parameter (cm) according to fMinBendingMomentum...
840 TMath::Abs(GetImpactParamFromBendingMomentum(fMinBendingMomentum));
841 // minimum impact parameter (cm) according to fMaxBendingMomentum...
843 TMath::Abs(GetImpactParamFromBendingMomentum(fMaxBendingMomentum));
845 else last2st = kFALSE;
846 // extrapolation factor from Z of first chamber to Z of second chamber
847 // dZ to be changed to take into account fine structure of chambers ????
849 // index for current segment
850 Int_t segmentIndex = 0;
851 // Loop over HitsForRec in the first chamber of the station
852 for (Int_t hit1 = fIndexOfFirstHitForRecPerChamber[ch1];
853 hit1 < fIndexOfFirstHitForRecPerChamber[ch1] + fNHitsForRecPerChamber[ch1];
855 // pointer to the HitForRec
856 hit1Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit1]);
858 // on the straight line joining the HitForRec to the vertex (0,0,0),
859 // to the Z of the second chamber of the station
860 // Loop over HitsForRec in the second chamber of the station
861 for (Int_t hit2 = fIndexOfFirstHitForRecPerChamber[ch2];
862 hit2 < fIndexOfFirstHitForRecPerChamber[ch2] + fNHitsForRecPerChamber[ch2];
864 // pointer to the HitForRec
865 hit2Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit2]);
866 // absolute values of distances, in bending and non bending planes,
867 // between the HitForRec in the second chamber
868 // and the previous extrapolation
869 extrapFact = hit2Ptr->GetZ()/ hit1Ptr->GetZ();
870 extBendCoor = extrapFact * hit1Ptr->GetBendingCoor();
871 extNonBendCoor = extrapFact * hit1Ptr->GetNonBendingCoor();
872 distBend = TMath::Abs(hit2Ptr->GetBendingCoor() - extBendCoor);
873 distNonBend = TMath::Abs(hit2Ptr->GetNonBendingCoor() - extNonBendCoor);
876 bendingSlope = (hit1Ptr->GetBendingCoor() - hit2Ptr->GetBendingCoor()) /
877 (hit1Ptr->GetZ() - hit2Ptr->GetZ());
878 // absolute value of impact parameter
880 TMath::Abs(hit1Ptr->GetBendingCoor() - hit1Ptr->GetZ() * bendingSlope);
882 // check for distances not too large,
883 // and impact parameter not too big if stations downstream of the dipole.
884 // Conditions "distBend" and "impactParam" correlated for these stations ????
885 if ((distBend < fSegmentMaxDistBending[Station]) &&
886 (distNonBend < fSegmentMaxDistNonBending[Station]) &&
887 (!last2st || (impactParam < maxImpactParam)) &&
888 (!last2st || (impactParam > minImpactParam))) {
890 segment = new ((*fSegmentsPtr[Station])[segmentIndex])
891 AliMUONSegment(hit1Ptr, hit2Ptr);
892 // update "link" to this segment from the hit in the first chamber
893 if (hit1Ptr->GetNSegments() == 0)
894 hit1Ptr->SetIndexOfFirstSegment(segmentIndex);
895 hit1Ptr->SetNSegments(hit1Ptr->GetNSegments() + 1);
896 if (AliLog::GetGlobalDebugLevel() > 1) {
897 cout << "segmentIndex(0...): " << segmentIndex
898 << " distBend: " << distBend
899 << " distNonBend: " << distNonBend
902 cout << "HitForRec in first chamber" << endl;
904 cout << "HitForRec in second chamber" << endl;
907 // increment index for current segment
911 } // for (Int_t hit1...
912 fNSegments[Station] = segmentIndex;
913 // Sorting according to "impact parameter" if station(1..5) 4 or 5,
914 // i.e. Station(0..4) 3 or 4, using the function "Compare".
915 // After this sorting, it is impossible to use
916 // the "fNSegments" and "fIndexOfFirstSegment"
917 // of the HitForRec in the first chamber to explore all segments formed with it.
918 // Is this sorting really needed ????
919 if ((Station == 3) || (Station == 4)) (fSegmentsPtr[Station])->Sort();
920 AliDebug(1,Form("Station: %d NSegments: %d ", Station, fNSegments[Station]));
924 //__________________________________________________________________________
925 void AliMUONEventReconstructor::MakeTracks(void)
927 // To make the tracks,
928 // from the list of segments and points in all stations
929 AliDebug(1,"Enter MakeTracks");
930 // The order may be important for the following Reset's
933 if (fTrackMethod == 2) { //AZ - Kalman filter
934 MakeTrackCandidatesK();
935 if (fRecTracksPtr->GetEntriesFast() == 0) return;
936 // Follow tracks in stations(1..) 3, 2 and 1
938 // Remove double tracks
939 RemoveDoubleTracksK();
940 // Propagate tracks to the vertex thru absorber
942 // Fill AliMUONTrack data members
945 // Look for candidates from at least 3 aligned points in stations(1..) 4 and 5
946 MakeTrackCandidates();
947 // Follow tracks in stations(1..) 3, 2 and 1
949 // Remove double tracks
950 RemoveDoubleTracks();
951 UpdateTrackParamAtHit();
952 UpdateHitForRecAtHit();
957 //__________________________________________________________________________
958 void AliMUONEventReconstructor::ValidateTracksWithTrigger(void)
960 // Try to match track from tracking system with trigger track
962 TClonesArray *recTriggerTracks;
964 fMUONData->ResetRecTriggerTracks();
965 fMUONData->SetTreeAddress("RL");
966 fMUONData->GetRecTriggerTracks();
967 recTriggerTracks = fMUONData->RecTriggerTracks();
969 track = (AliMUONTrack*) fRecTracksPtr->First();
971 track->MatchTriggerTrack(recTriggerTracks);
972 track = (AliMUONTrack*) fRecTracksPtr->After(track);
977 //__________________________________________________________________________
978 Bool_t AliMUONEventReconstructor::MakeTriggerTracks(void)
980 // To make the trigger tracks from Local Trigger
981 AliDebug(1, "Enter MakeTriggerTracks");
982 // ResetTriggerTracks();
986 TClonesArray *localTrigger;
987 TClonesArray *globalTrigger;
988 AliMUONLocalTrigger *locTrg;
989 AliMUONGlobalTrigger *gloTrg;
990 AliMUONTriggerCircuit *circuit;
991 AliMUONTriggerTrack *recTriggerTrack = 0;
993 TTree* treeR = fLoader->TreeR();
995 // Loading MUON subsystem
996 AliMUON * pMUON = (AliMUON *) gAlice->GetDetector("MUON");
998 nTRentries = Int_t(treeR->GetEntries());
1000 treeR->GetEvent(0); // only one entry
1002 if (!(fMUONData->IsTriggerBranchesInTree())) {
1003 AliWarning(Form("Trigger information is not avalaible, nTRentries = %d not equal to 1",nTRentries));
1007 fMUONData->SetTreeAddress("TC");
1008 fMUONData->GetTrigger();
1010 // global trigger for trigger pattern
1012 globalTrigger = fMUONData->GlobalTrigger();
1013 gloTrg = (AliMUONGlobalTrigger*)globalTrigger->UncheckedAt(0);
1015 if (gloTrg->SinglePlusLpt()) gloTrigPat|= 0x1;
1016 if (gloTrg->SinglePlusHpt()) gloTrigPat|= 0x2;
1017 if (gloTrg->SinglePlusApt()) gloTrigPat|= 0x4;
1019 if (gloTrg->SingleMinusLpt()) gloTrigPat|= 0x8;
1020 if (gloTrg->SingleMinusHpt()) gloTrigPat|= 0x10;
1021 if (gloTrg->SingleMinusApt()) gloTrigPat|= 0x20;
1023 if (gloTrg->SingleUndefLpt()) gloTrigPat|= 0x40;
1024 if (gloTrg->SingleUndefHpt()) gloTrigPat|= 0x80;
1025 if (gloTrg->SingleUndefApt()) gloTrigPat|= 0x100;
1027 if (gloTrg->PairUnlikeLpt()) gloTrigPat|= 0x200;
1028 if (gloTrg->PairUnlikeHpt()) gloTrigPat|= 0x400;
1029 if (gloTrg->PairUnlikeApt()) gloTrigPat|= 0x800;
1031 if (gloTrg->PairLikeLpt()) gloTrigPat|= 0x1000;
1032 if (gloTrg->PairLikeHpt()) gloTrigPat|= 0x2000;
1033 if (gloTrg->PairLikeApt()) gloTrigPat|= 0x4000;
1036 // local trigger for tracking
1037 localTrigger = fMUONData->LocalTrigger();
1038 Int_t nlocals = (Int_t) (localTrigger->GetEntries());
1040 Float_t z11 = AliMUONConstants::DefaultChamberZ(10);
1041 Float_t z21 = AliMUONConstants::DefaultChamberZ(12);
1043 for (Int_t i=0; i<nlocals; i++) { // loop on Local Trigger
1044 locTrg = (AliMUONLocalTrigger*)localTrigger->UncheckedAt(i);
1045 circuit = &(pMUON->TriggerCircuit(locTrg->LoCircuit()));
1046 Float_t y11 = circuit->GetY11Pos(locTrg->LoStripX());
1047 Int_t stripX21 = locTrg->LoStripX()+locTrg->LoDev()+1;
1048 Float_t y21 = circuit->GetY21Pos(stripX21);
1049 Float_t x11 = circuit->GetX11Pos(locTrg->LoStripY());
1050 Float_t thetax = TMath::ATan2( x11 , z11 );
1051 Float_t thetay = TMath::ATan2( (y21-y11) , (z21-z11) );
1053 recTriggerTrack = new AliMUONTriggerTrack(x11,y11,thetax,thetay,gloTrigPat);
1054 // since static statement does not work, set gloTrigPat for each track
1056 fMUONData->AddRecTriggerTrack(*recTriggerTrack);
1057 } // end of loop on Local Trigger
1061 //__________________________________________________________________________
1062 Int_t AliMUONEventReconstructor::MakeTrackCandidatesWithTwoSegments(AliMUONSegment *BegSegment)
1064 // To make track candidates with two segments in stations(1..) 4 and 5,
1065 // the first segment being pointed to by "BegSegment".
1066 // Returns the number of such track candidates.
1067 Int_t endStation, iEndSegment, nbCan2Seg;
1068 AliMUONSegment *endSegment;
1069 AliMUONSegment *extrapSegment = NULL;
1070 AliMUONTrack *recTrack;
1072 AliDebug(1,"Enter MakeTrackCandidatesWithTwoSegments");
1073 // Station for the end segment
1074 endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
1075 // multiple scattering factor corresponding to one chamber
1076 mcsFactor = 0.0136 /
1077 GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
1078 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1079 // linear extrapolation to end station
1080 // number of candidates with 2 segments to 0
1082 // Loop over segments in the end station
1083 for (iEndSegment = 0; iEndSegment < fNSegments[endStation]; iEndSegment++) {
1084 // pointer to segment
1085 endSegment = (AliMUONSegment*) ((*fSegmentsPtr[endStation])[iEndSegment]);
1086 // test compatibility between current segment and "extrapSegment"
1087 // 4 because 4 quantities in chi2
1089 BegSegment->CreateSegmentFromLinearExtrapToStation(endSegment->GetZ(), mcsFactor);
1091 NormalizedChi2WithSegment(extrapSegment,
1092 fMaxSigma2Distance)) <= 4.0) {
1093 // both segments compatible:
1094 // make track candidate from "begSegment" and "endSegment"
1095 AliDebug(2,Form("TrackCandidate with Segment %d in Station(0..) %d", iEndSegment, endStation));
1096 // flag for both segments in one track:
1097 // to be done in track constructor ????
1098 BegSegment->SetInTrack(kTRUE);
1099 endSegment->SetInTrack(kTRUE);
1100 recTrack = new ((*fRecTracksPtr)[fNRecTracks])
1101 AliMUONTrack(BegSegment, endSegment, this);
1103 if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump();
1104 // increment number of track candidates with 2 segments
1107 delete extrapSegment; // should not delete HitForRec's it points to !!!!
1108 } // for (iEndSegment = 0;...
1112 //__________________________________________________________________________
1113 Int_t AliMUONEventReconstructor::MakeTrackCandidatesWithOneSegmentAndOnePoint(AliMUONSegment *BegSegment)
1115 // To make track candidates with one segment and one point
1116 // in stations(1..) 4 and 5,
1117 // the segment being pointed to by "BegSegment".
1118 Int_t ch, ch1, ch2, endStation, iHit, iHitMax, iHitMin, nbCan1Seg1Hit;
1119 AliMUONHitForRec *extrapHitForRec= NULL;
1120 AliMUONHitForRec *hit;
1121 AliMUONTrack *recTrack;
1123 AliDebug(1,"Enter MakeTrackCandidatesWithOneSegmentAndOnePoint");
1124 // station for the end point
1125 endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
1126 // multiple scattering factor corresponding to one chamber
1127 mcsFactor = 0.0136 /
1128 GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
1129 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1130 // first and second chambers(0..) in the end station
1131 ch1 = 2 * endStation;
1133 // number of candidates to 0
1135 // Loop over chambers of the end station
1136 for (ch = ch2; ch >= ch1; ch--) {
1137 // limits for the hit index in the loop
1138 iHitMin = fIndexOfFirstHitForRecPerChamber[ch];
1139 iHitMax = iHitMin + fNHitsForRecPerChamber[ch];
1140 // Loop over HitForRec's in the chamber
1141 for (iHit = iHitMin; iHit < iHitMax; iHit++) {
1142 // pointer to HitForRec
1143 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1144 // test compatibility between current HitForRec and "extrapHitForRec"
1145 // 2 because 2 quantities in chi2
1146 // linear extrapolation to chamber
1148 BegSegment->CreateHitForRecFromLinearExtrapToChamber( hit->GetZ(), mcsFactor);
1150 NormalizedChi2WithHitForRec(extrapHitForRec,
1151 fMaxSigma2Distance)) <= 2.0) {
1152 // both HitForRec's compatible:
1153 // make track candidate from begSegment and current HitForRec
1154 AliDebug(2, Form("TrackCandidate with HitForRec %d in Chamber(0..) %d", iHit, ch));
1155 // flag for beginning segments in one track:
1156 // to be done in track constructor ????
1157 BegSegment->SetInTrack(kTRUE);
1158 recTrack = new ((*fRecTracksPtr)[fNRecTracks])
1159 AliMUONTrack(BegSegment, hit, this);
1160 // the right place to eliminate "double counting" ???? how ????
1162 if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump();
1163 // increment number of track candidates
1166 delete extrapHitForRec;
1167 } // for (iHit = iHitMin;...
1168 } // for (ch = ch2;...
1169 return nbCan1Seg1Hit;
1172 //__________________________________________________________________________
1173 void AliMUONEventReconstructor::MakeTrackCandidates(void)
1175 // To make track candidates
1176 // with at least 3 aligned points in stations(1..) 4 and 5
1177 // (two Segment's or one Segment and one HitForRec)
1178 Int_t begStation, iBegSegment, nbCan1Seg1Hit, nbCan2Seg;
1179 AliMUONSegment *begSegment;
1180 AliDebug(1,"Enter MakeTrackCandidates");
1181 // Loop over stations(1..) 5 and 4 for the beginning segment
1182 for (begStation = 4; begStation > 2; begStation--) {
1183 // Loop over segments in the beginning station
1184 for (iBegSegment = 0; iBegSegment < fNSegments[begStation]; iBegSegment++) {
1185 // pointer to segment
1186 begSegment = (AliMUONSegment*) ((*fSegmentsPtr[begStation])[iBegSegment]);
1187 AliDebug(2,Form("Look for TrackCandidate's with Segment %d in Station(0..) %d", iBegSegment, begStation));
1188 // Look for track candidates with two segments,
1189 // "begSegment" and all compatible segments in other station.
1190 // Only for beginning station(1..) 5
1191 // because candidates with 2 segments have to looked for only once.
1192 if (begStation == 4)
1193 nbCan2Seg = MakeTrackCandidatesWithTwoSegments(begSegment);
1194 // Look for track candidates with one segment and one point,
1195 // "begSegment" and all compatible HitForRec's in other station.
1196 // Only if "begSegment" does not belong already to a track candidate.
1197 // Is that a too strong condition ????
1198 if (!(begSegment->GetInTrack()))
1199 nbCan1Seg1Hit = MakeTrackCandidatesWithOneSegmentAndOnePoint(begSegment);
1200 } // for (iBegSegment = 0;...
1201 } // for (begStation = 4;...
1205 //__________________________________________________________________________
1206 void AliMUONEventReconstructor::FollowTracks(void)
1208 // Follow tracks in stations(1..) 3, 2 and 1
1209 // too long: should be made more modular !!!!
1210 AliMUONHitForRec *bestHit, *extrapHit, *hit;
1211 AliMUONSegment *bestSegment, *extrapSegment, *segment;
1212 AliMUONTrack *track, *nextTrack;
1213 AliMUONTrackParam *trackParam1, trackParam[2], trackParamVertex;
1214 // -1 to avoid compilation warnings
1215 Int_t ch = -1, chInStation, chBestHit = -1, iHit, iSegment, station, trackIndex;
1216 Double_t bestChi2, chi2, dZ1, dZ2, dZ3, maxSigma2Distance, mcsFactor;
1217 Double_t bendingMomentum, chi2Norm = 0.;
1218 AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
1219 // local maxSigma2Distance, for easy increase in testing
1220 maxSigma2Distance = fMaxSigma2Distance;
1221 AliDebug(2,"Enter FollowTracks");
1222 // Loop over track candidates
1223 track = (AliMUONTrack*) fRecTracksPtr->First();
1226 // Follow function for each track candidate ????
1228 nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track); // prepare next track
1229 AliDebug(2,Form("FollowTracks: track candidate(0..): %d", trackIndex));
1230 // Fit track candidate
1231 track->SetFitMCS(0); // without multiple Coulomb scattering
1232 track->SetFitNParam(3); // with 3 parameters (X = Y = 0)
1233 track->SetFitStart(0); // from parameters at vertex
1235 if (AliLog::GetGlobalDebugLevel()> 2) {
1236 cout << "FollowTracks: track candidate(0..): " << trackIndex
1237 << " after fit in stations(0..) 3 and 4" << endl;
1238 track->RecursiveDump();
1240 // Loop over stations(1..) 3, 2 and 1
1241 // something SPECIAL for stations 2 and 1 for majority 3 coincidence ????
1242 // otherwise: majority coincidence 2 !!!!
1243 for (station = 2; station >= 0; station--) {
1244 // Track parameters at first track hit (smallest Z)
1245 trackParam1 = ((AliMUONTrackHit*)
1246 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1247 // extrapolation to station
1248 trackParam1->ExtrapToStation(station, trackParam);
1249 extrapSegment = new AliMUONSegment(); // empty segment
1250 // multiple scattering factor corresponding to one chamber
1251 // and momentum in bending plane (not total)
1252 mcsFactor = 0.0136 * trackParam1->GetInverseBendingMomentum();
1253 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1254 // Z difference from previous station
1255 dZ1 = (&(pMUON->Chamber(2 * station)))->Z() -
1256 (&(pMUON->Chamber(2 * station + 2)))->Z();
1257 // Z difference between the two previous stations
1258 dZ2 = (&(pMUON->Chamber(2 * station + 2)))->Z() -
1259 (&(pMUON->Chamber(2 * station + 4)))->Z();
1260 // Z difference between the two chambers in the previous station
1261 dZ3 = (&(pMUON->Chamber(2 * station)))->Z() -
1262 (&(pMUON->Chamber(2 * station + 1)))->Z();
1263 extrapSegment->SetBendingCoorReso2(fBendingResolution * fBendingResolution);
1265 SetNonBendingCoorReso2(fNonBendingResolution * fNonBendingResolution);
1266 extrapSegment->UpdateFromStationTrackParam
1267 (trackParam, mcsFactor, dZ1, dZ2, dZ3, station,
1268 trackParam1->GetInverseBendingMomentum());
1271 if (AliLog::GetGlobalDebugLevel() > 2) {
1272 cout << "FollowTracks: track candidate(0..): " << trackIndex
1273 << " Look for segment in station(0..): " << station << endl;
1275 // Loop over segments in station
1276 for (iSegment = 0; iSegment < fNSegments[station]; iSegment++) {
1277 // Look for best compatible Segment in station
1278 // should consider all possibilities ????
1279 // multiple scattering ????
1280 // separation in 2 functions: Segment and HitForRec ????
1281 segment = (AliMUONSegment*) ((*fSegmentsPtr[station])[iSegment]);
1282 // correction of corrected segment (fBendingCoor and fNonBendingCoor)
1283 // according to real Z value of "segment" and slopes of "extrapSegment"
1284 (&(trackParam[0]))->ExtrapToZ(segment->GetZ());
1285 (&(trackParam[1]))->ExtrapToZ(segment->GetZ());
1286 extrapSegment->SetBendingCoor((&(trackParam[0]))->GetBendingCoor());
1287 extrapSegment->SetNonBendingCoor((&(trackParam[0]))->GetNonBendingCoor());
1288 extrapSegment->SetBendingSlope((&(trackParam[0]))->GetBendingSlope());
1289 extrapSegment->SetNonBendingSlope((&(trackParam[0]))->GetNonBendingSlope());
1291 NormalizedChi2WithSegment(extrapSegment, maxSigma2Distance);
1292 if (chi2 < bestChi2) {
1293 // update best Chi2 and Segment if better found
1294 bestSegment = segment;
1299 // best segment found: add it to track candidate
1300 (&(trackParam[0]))->ExtrapToZ(bestSegment->GetZ());
1301 (&(trackParam[1]))->ExtrapToZ(bestSegment->GetZ());
1302 track->AddSegment(bestSegment);
1303 // set track parameters at these two TrakHit's
1304 track->SetTrackParamAtHit(track->GetNTrackHits() - 2, &(trackParam[0]));
1305 track->SetTrackParamAtHit(track->GetNTrackHits() - 1, &(trackParam[1]));
1306 AliDebug(3, Form("FollowTracks: track candidate(0..): %d Added segment in station(0..): %d", trackIndex, station));
1307 if (AliLog::GetGlobalDebugLevel()>2) track->RecursiveDump();
1310 // No best segment found:
1311 // Look for best compatible HitForRec in station:
1312 // should consider all possibilities ????
1313 // multiple scattering ???? do about like for extrapSegment !!!!
1314 extrapHit = new AliMUONHitForRec(); // empty hit
1317 AliDebug(3, Form("FollowTracks: track candidate(0..): %d Segment not found, look for hit in station(0..): %d ",
1318 trackIndex, station));
1320 // Loop over chambers of the station
1321 for (chInStation = 0; chInStation < 2; chInStation++) {
1322 ch = 2 * station + chInStation;
1323 for (iHit = fIndexOfFirstHitForRecPerChamber[ch];
1324 iHit < fIndexOfFirstHitForRecPerChamber[ch] +
1325 fNHitsForRecPerChamber[ch];
1327 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1328 // coordinates of extrapolated hit
1329 (&(trackParam[chInStation]))->ExtrapToZ(hit->GetZ());
1331 SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor());
1333 SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor());
1334 // resolutions from "extrapSegment"
1335 extrapHit->SetBendingReso2(extrapSegment->GetBendingCoorReso2());
1336 extrapHit->SetNonBendingReso2(extrapSegment->GetNonBendingCoorReso2());
1337 // Loop over hits in the chamber
1338 // condition for hit not already in segment ????
1339 chi2 = hit->NormalizedChi2WithHitForRec(extrapHit, maxSigma2Distance);
1340 if (chi2 < bestChi2) {
1341 // update best Chi2 and HitForRec if better found
1344 chBestHit = chInStation;
1349 // best hit found: add it to track candidate
1350 (&(trackParam[chBestHit]))->ExtrapToZ(bestHit->GetZ());
1351 track->AddHitForRec(bestHit);
1352 // set track parameters at this TrackHit
1353 track->SetTrackParamAtHit(track->GetNTrackHits() - 1,
1354 &(trackParam[chBestHit]));
1355 if (AliLog::GetGlobalDebugLevel() > 2) {
1356 cout << "FollowTracks: track candidate(0..): " << trackIndex
1357 << " Added hit in station(0..): " << station << endl;
1358 track->RecursiveDump();
1362 // Remove current track candidate
1363 // and corresponding TrackHit's, ...
1365 delete extrapSegment;
1367 break; // stop the search for this candidate:
1368 // exit from the loop over station
1372 delete extrapSegment;
1373 // Sort track hits according to increasing Z
1374 track->GetTrackHitsPtr()->Sort();
1375 // Update track parameters at first track hit (smallest Z)
1376 trackParam1 = ((AliMUONTrackHit*)
1377 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1378 bendingMomentum = 0.;
1379 if (TMath::Abs(trackParam1->GetInverseBendingMomentum()) > 0.)
1380 bendingMomentum = TMath::Abs(1/(trackParam1->GetInverseBendingMomentum()));
1381 // Track removed if bendingMomentum not in window [min, max]
1382 if ((bendingMomentum < fMinBendingMomentum) || (bendingMomentum > fMaxBendingMomentum)) {
1384 break; // stop the search for this candidate:
1385 // exit from the loop over station
1388 // with multiple Coulomb scattering if all stations
1389 if (station == 0) track->SetFitMCS(1);
1390 // without multiple Coulomb scattering if not all stations
1391 else track->SetFitMCS(0);
1392 track->SetFitNParam(5); // with 5 parameters (momentum and position)
1393 track->SetFitStart(1); // from parameters at first hit
1395 Double_t numberOfDegFree = (2.0 * track->GetNTrackHits() - 5);
1396 if (numberOfDegFree > 0) {
1397 chi2Norm = track->GetFitFMin() / numberOfDegFree;
1401 // Track removed if normalized chi2 too high
1402 if (chi2Norm > fMaxChi2) {
1404 break; // stop the search for this candidate:
1405 // exit from the loop over station
1407 if (AliLog::GetGlobalDebugLevel() > 2) {
1408 cout << "FollowTracks: track candidate(0..): " << trackIndex
1409 << " after fit from station(0..): " << station << " to 4" << endl;
1410 track->RecursiveDump();
1412 // Track extrapolation to the vertex through the absorber (Branson)
1413 // after going through the first station
1415 trackParamVertex = *trackParam1;
1416 (&trackParamVertex)->ExtrapToVertex(0.,0.,0.);
1417 track->SetTrackParamAtVertex(&trackParamVertex);
1418 if (AliLog::GetGlobalDebugLevel() > 0) {
1419 cout << "FollowTracks: track candidate(0..): " << trackIndex
1420 << " after extrapolation to vertex" << endl;
1421 track->RecursiveDump();
1424 } // for (station = 2;...
1425 // go really to next track
1428 // Compression of track array (necessary after Remove ????)
1429 fRecTracksPtr->Compress();
1433 //__________________________________________________________________________
1434 void AliMUONEventReconstructor::RemoveDoubleTracks(void)
1436 // To remove double tracks.
1437 // Tracks are considered identical
1438 // if they have at least half of their hits in common.
1439 // Among two identical tracks, one keeps the track with the larger number of hits
1440 // or, if these numbers are equal, the track with the minimum Chi2.
1441 AliMUONTrack *track1, *track2, *trackToRemove;
1442 Bool_t identicalTracks;
1443 Int_t hitsInCommon, nHits1, nHits2;
1444 identicalTracks = kTRUE;
1445 while (identicalTracks) {
1446 identicalTracks = kFALSE;
1447 // Loop over first track of the pair
1448 track1 = (AliMUONTrack*) fRecTracksPtr->First();
1449 while (track1 && (!identicalTracks)) {
1450 nHits1 = track1->GetNTrackHits();
1451 // Loop over second track of the pair
1452 track2 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1453 while (track2 && (!identicalTracks)) {
1454 nHits2 = track2->GetNTrackHits();
1455 // number of hits in common between two tracks
1456 hitsInCommon = track1->HitsInCommon(track2);
1457 // check for identical tracks
1458 if ((4 * hitsInCommon) >= (nHits1 + nHits2)) {
1459 identicalTracks = kTRUE;
1460 // decide which track to remove
1461 if (nHits1 > nHits2) trackToRemove = track2;
1462 else if (nHits1 < nHits2) trackToRemove = track1;
1463 else if ((track1->GetFitFMin()) < (track2->GetFitFMin()))
1464 trackToRemove = track2;
1465 else trackToRemove = track1;
1467 trackToRemove->Remove();
1469 track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
1471 track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1477 //__________________________________________________________________________
1478 void AliMUONEventReconstructor::UpdateTrackParamAtHit()
1480 // Set track parameters after track fitting. Fill fTrackParamAtHit of AliMUONTrack's
1481 AliMUONTrack *track;
1482 AliMUONTrackHit *trackHit;
1483 AliMUONTrackParam *trackParam;
1484 track = (AliMUONTrack*) fRecTracksPtr->First();
1486 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->First();
1488 trackParam = trackHit->GetTrackParam();
1489 track->AddTrackParamAtHit(trackParam);
1490 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->After(trackHit);
1492 track = (AliMUONTrack*) fRecTracksPtr->After(track);
1497 //__________________________________________________________________________
1498 void AliMUONEventReconstructor::UpdateHitForRecAtHit()
1500 // Set cluster parameterss after track fitting. Fill fHitForRecAtHit of AliMUONTrack's
1501 AliMUONTrack *track;
1502 AliMUONTrackHit *trackHit;
1503 AliMUONHitForRec *hitForRec;
1504 track = (AliMUONTrack*) fRecTracksPtr->First();
1506 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->First();
1508 hitForRec = trackHit->GetHitForRecPtr();
1509 track->AddHitForRecAtHit(hitForRec);
1510 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->After(trackHit);
1512 track = (AliMUONTrack*) fRecTracksPtr->After(track);
1517 //__________________________________________________________________________
1518 void AliMUONEventReconstructor::FillMUONTrack()
1520 // Set track parameters at hits for Kalman track. Fill fTrackParamAtHit of AliMUONTrack's
1521 AliMUONTrackK *track;
1522 track = (AliMUONTrackK*) fRecTracksPtr->First();
1524 track->FillMUONTrack();
1525 track = (AliMUONTrackK*) fRecTracksPtr->After(track);
1530 //__________________________________________________________________________
1531 void AliMUONEventReconstructor::EventDump(void)
1533 // Dump reconstructed event (track parameters at vertex and at first hit),
1534 // and the particle parameters
1536 AliMUONTrack *track;
1537 AliMUONTrackParam *trackParam, *trackParam1;
1538 Double_t bendingSlope, nonBendingSlope, pYZ;
1539 Double_t pX, pY, pZ, x, y, z, c;
1540 Int_t np, trackIndex, nTrackHits;
1542 AliDebug(1,"****** enter EventDump ******");
1543 AliDebug(1, Form("Number of Reconstructed tracks : %d", fNRecTracks));
1545 fRecTracksPtr->Compress(); // for simple loop without "Next" since no hole
1546 // Loop over reconstructed tracks
1547 for (trackIndex = 0; trackIndex < fNRecTracks; trackIndex++) {
1548 if (fTrackMethod != 1) continue; //AZ - skip the rest for now
1549 AliDebug(1, Form("track number: %d", trackIndex));
1550 // function for each track for modularity ????
1551 track = (AliMUONTrack*) ((*fRecTracksPtr)[trackIndex]);
1552 nTrackHits = track->GetNTrackHits();
1553 AliDebug(1, Form("Number of track hits: %d ", nTrackHits));
1554 // track parameters at Vertex
1555 trackParam = track->GetTrackParamAtVertex();
1556 x = trackParam->GetNonBendingCoor();
1557 y = trackParam->GetBendingCoor();
1558 z = trackParam->GetZ();
1559 bendingSlope = trackParam->GetBendingSlope();
1560 nonBendingSlope = trackParam->GetNonBendingSlope();
1561 pYZ = 1/TMath::Abs(trackParam->GetInverseBendingMomentum());
1562 pZ = pYZ/TMath::Sqrt(1+bendingSlope*bendingSlope);
1563 pX = pZ * nonBendingSlope;
1564 pY = pZ * bendingSlope;
1565 c = TMath::Sign(1.0, trackParam->GetInverseBendingMomentum());
1566 AliDebug(1, Form("Track parameters at Vertex z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1567 z, x, y, pX, pY, pZ, c));
1569 // track parameters at first hit
1570 trackParam1 = ((AliMUONTrackHit*)
1571 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1572 x = trackParam1->GetNonBendingCoor();
1573 y = trackParam1->GetBendingCoor();
1574 z = trackParam1->GetZ();
1575 bendingSlope = trackParam1->GetBendingSlope();
1576 nonBendingSlope = trackParam1->GetNonBendingSlope();
1577 pYZ = 1/TMath::Abs(trackParam1->GetInverseBendingMomentum());
1578 pZ = pYZ/TMath::Sqrt(1.0 + bendingSlope * bendingSlope);
1579 pX = pZ * nonBendingSlope;
1580 pY = pZ * bendingSlope;
1581 c = TMath::Sign(1.0, trackParam1->GetInverseBendingMomentum());
1582 AliDebug(1, Form("track parameters at z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1583 z, x, y, pX, pY, pZ, c));
1585 // informations about generated particles
1586 np = gAlice->GetMCApp()->GetNtrack();
1587 printf(" **** number of generated particles: %d \n", np);
1589 // for (Int_t iPart = 0; iPart < np; iPart++) {
1590 // p = gAlice->Particle(iPart);
1591 // printf(" particle %d: type= %d px= %f py= %f pz= %f pdg= %d\n",
1592 // iPart, p->GetPdgCode(), p->Px(), p->Py(), p->Pz(), p->GetPdgCode());
1598 //__________________________________________________________________________
1599 void AliMUONEventReconstructor::EventDumpTrigger(void)
1601 // Dump reconstructed trigger event
1602 // and the particle parameters
1604 AliMUONTriggerTrack *triggertrack ;
1605 Int_t nTriggerTracks = fMUONData->RecTriggerTracks()->GetEntriesFast();
1607 AliDebug(1, "****** enter EventDumpTrigger ******");
1608 AliDebug(1, Form("Number of Reconstructed tracks : %d ", nTriggerTracks));
1610 // Loop over reconstructed tracks
1611 for (Int_t trackIndex = 0; trackIndex < nTriggerTracks; trackIndex++) {
1612 triggertrack = (AliMUONTriggerTrack*)fMUONData->RecTriggerTracks()->At(trackIndex);
1613 printf(" trigger track number %i x11=%f y11=%f thetax=%f thetay=%f \n",
1615 triggertrack->GetX11(),triggertrack->GetY11(),
1616 triggertrack->GetThetax(),triggertrack->GetThetay());
1620 //__________________________________________________________________________
1621 void AliMUONEventReconstructor::FillEvent()
1623 // Create a new AliMUONRecoEvent, fill its track list, then add it as a
1624 // leaf in the Event branch of TreeRecoEvent tree
1625 cout << "Enter FillEvent() ...\n";
1628 fRecoEvent = new AliMUONRecoEvent();
1630 fRecoEvent->Clear();
1632 //save current directory
1633 TDirectory *current = gDirectory;
1634 if (!fTreeFile) fTreeFile = new TFile("tree_reco.root", "RECREATE");
1635 if (!fEventTree) fEventTree = new TTree("TreeRecoEvent", "MUON reconstructed events");
1636 //AZif (fRecoEvent->MakeDumpTracks(fRecTracksPtr)) {
1637 if (fRecoEvent->MakeDumpTracks(fMuons, fRecTracksPtr, this)) { //AZ
1638 if (AliLog::GetGlobalDebugLevel() > 1) fRecoEvent->EventInfo();
1639 TBranch *branch = fEventTree->GetBranch("Event");
1640 if (!branch) branch = fEventTree->Branch("Event", "AliMUONRecoEvent", &fRecoEvent, 64000);
1641 branch->SetAutoDelete();
1646 // restore directory
1650 //__________________________________________________________________________
1651 void AliMUONEventReconstructor::MakeTrackCandidatesK(void)
1653 // To make initial tracks for Kalman filter from the list of segments
1655 AliMUONSegment *segment;
1656 AliMUONTrackK *trackK;
1658 AliDebug(1,"Enter MakeTrackCandidatesK");
1659 // Reset the TClonesArray of reconstructed tracks
1660 if (fRecTracksPtr) fRecTracksPtr->Delete();
1661 // Delete in order that the Track destructors are called,
1662 // hence the space for the TClonesArray of pointers to TrackHit's is freed
1665 AliMUONTrackK a(this, fHitsForRecPtr); // bad idea ???
1666 // Loop over stations(1...) 5 and 4
1667 for (istat=4; istat>=3; istat--) {
1668 // Loop over segments in the station
1669 for (iseg=0; iseg<fNSegments[istat]; iseg++) {
1670 // Transform segments to tracks and evaluate covariance matrix
1671 segment = (AliMUONSegment*) ((*fSegmentsPtr[istat])[iseg]);
1672 trackK = new ((*fRecTracksPtr)[fNRecTracks]) AliMUONTrackK(segment);
1674 } // for (iseg=0;...)
1675 } // for (istat=4;...)
1679 //__________________________________________________________________________
1680 void AliMUONEventReconstructor::FollowTracksK(void)
1682 // Follow tracks using Kalman filter
1684 Int_t icand, ichamBeg = 0, ichamEnd, chamBits;
1685 Double_t zDipole1, zDipole2;
1686 AliMUONTrackK *trackK;
1687 AliMUONHitForRec *hit;
1688 AliMUONRawCluster *clus;
1689 TClonesArray *rawclusters;
1690 clus = 0; rawclusters = 0;
1692 //AZ(z->-z) zDipole1 = GetSimpleBPosition() - GetSimpleBLength()/2;
1693 //AZ(z->-z) zDipole2 = zDipole1 + GetSimpleBLength();
1694 zDipole1 = GetSimpleBPosition() + GetSimpleBLength()/2;
1695 zDipole2 = zDipole1 - GetSimpleBLength();
1698 trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[0]);
1699 if (trackK->DebugLevel() > 0) {
1700 for (Int_t i1=0; i1<fNHitsForRec; i1++) {
1701 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[i1]);
1702 //if (hit->GetTTRTrack() > 1 || hit->GetTrackRefSignal() == 0) continue;
1703 printf(" Hit # %d %10.4f %10.4f %10.4f",
1704 hit->GetChamberNumber(), hit->GetBendingCoor(),
1705 hit->GetNonBendingCoor(), hit->GetZ());
1706 if (fRecTrackRefHits) {
1707 // from track ref hits
1708 printf(" %3d %3d \n", hit->GetTrackRefSignal(), hit->GetTTRTrack());
1710 // from raw clusters
1711 rawclusters = fMUONData->RawClusters(hit->GetChamberNumber());
1712 clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(hit->
1714 printf("%3d", clus->GetTrack(1)-1);
1715 if (clus->GetTrack(2) != 0) printf("%3d \n", clus->GetTrack(2)-1);
1717 } // if (fRecTrackRefHits)
1719 } // if (trackK->DebugLevel() > 0)
1723 nSeeds = fNRecTracks; // starting number of seeds
1724 // Loop over track candidates
1725 while (icand < fNRecTracks-1) {
1727 if (trackK->DebugLevel()>0) cout << " *** Kalman track candidate No. " << icand << endl;
1728 trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]);
1729 if (trackK->GetRecover() < 0) continue; // failed track
1731 // Discard candidate which will produce the double track
1734 ok = CheckCandidateK(icand,nSeeds);
1736 trackK->SetRecover(-1); // mark candidate to be removed
1743 if (trackK->GetRecover() == 0) hit = (AliMUONHitForRec*)
1744 trackK->GetHitOnTrack()->Last(); // last hit
1745 //else hit = (AliMUONHitForRec*) (*trackK->GetHitOnTrack())[1]; // 2'nd hit
1746 else hit = trackK->GetHitLastOk(); // hit where track stopped
1747 if (hit) ichamBeg = hit->GetChamberNumber();
1749 // Check propagation direction
1750 if (hit == NULL) ichamBeg = ichamEnd;
1751 //AZ(z->-z) else if (trackK->GetTrackDir() > 0) {
1752 else if (trackK->GetTrackDir() < 0) {
1753 ichamEnd = 9; // forward propagation
1754 ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kFALSE,zDipole1,zDipole2);
1756 ichamBeg = ichamEnd;
1757 ichamEnd = 6; // backward propagation
1758 // Change weight matrix and zero fChi2 for backpropagation
1759 trackK->StartBack();
1760 //AZ trackK->SetTrackDir(-1);
1761 trackK->SetTrackDir(1);
1762 ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kTRUE,zDipole1,zDipole2);
1763 ichamBeg = ichamEnd;
1767 if (trackK->GetBPFlag()) {
1769 ichamEnd = 6; // backward propagation
1770 // Change weight matrix and zero fChi2 for backpropagation
1771 trackK->StartBack();
1772 ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kTRUE,zDipole1,zDipole2);
1773 ichamBeg = ichamEnd;
1779 //AZ trackK->SetTrackDir(-1);
1780 trackK->SetTrackDir(1);
1781 trackK->SetBPFlag(kFALSE);
1782 ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kFALSE,zDipole1,zDipole2);
1784 if (!ok) { trackK->SetRecover(-1); continue; } // mark candidate to be removed
1787 if (trackK->GetRecover() >= 0) {
1788 ok = trackK->Smooth();
1789 if (!ok) trackK->SetRecover(-1); // mark candidate to be removed
1792 // Majority 3 of 4 in first 2 stations
1795 Double_t chi2_max = 0;
1796 for (Int_t i=0; i<trackK->GetNTrackHits(); i++) {
1797 hit = (AliMUONHitForRec*) (*trackK->GetHitOnTrack())[i];
1798 chamBits |= BIT(hit->GetChamberNumber());
1799 if (trackK->GetChi2PerPoint(i) > chi2_max) chi2_max = trackK->GetChi2PerPoint(i);
1801 if (!((chamBits&3)==3 || (chamBits>>2&3)==3) && chi2_max > 25) {
1802 //trackK->Recover();
1803 trackK->SetRecover(-1); //mark candidate to be removed
1806 if (ok) trackK->SetTrackQuality(0); // compute "track quality"
1809 for (Int_t i=0; i<fNRecTracks; i++) {
1810 trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[i]);
1811 if (trackK->GetRecover() < 0) fRecTracksPtr->RemoveAt(i);
1814 // Compress TClonesArray
1815 fRecTracksPtr->Compress();
1816 fNRecTracks = fRecTracksPtr->GetEntriesFast();
1820 //__________________________________________________________________________
1821 Bool_t AliMUONEventReconstructor::CheckCandidateK(Int_t icand, Int_t nSeeds) const
1823 // Discards track candidate if it will produce the double track (having
1824 // the same seed segment hits as hits of a good track found before)
1825 AliMUONTrackK *track1, *track2;
1826 AliMUONHitForRec *hit1, *hit2, *hit;
1828 track1 = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]);
1829 hit1 = (AliMUONHitForRec*) (*track1->GetHitOnTrack())[0]; // 1'st hit
1830 hit2 = (AliMUONHitForRec*) (*track1->GetHitOnTrack())[1]; // 2'nd hit
1832 for (Int_t i=0; i<icand; i++) {
1833 track2 = (AliMUONTrackK*) ((*fRecTracksPtr)[i]);
1834 //if (track2->GetRecover() < 0) continue;
1835 if (track2->GetRecover() < 0 && icand >= nSeeds) continue;
1837 if (track1->GetStartSegment() == track2->GetStartSegment()) {
1841 for (Int_t j=0; j<track2->GetNTrackHits(); j++) {
1842 hit = (AliMUONHitForRec*) (*track2->GetHitOnTrack())[j];
1843 if (hit == hit1 || hit == hit2) {
1845 if (nSame == 2) return kFALSE;
1847 } // for (Int_t j=0;
1849 } // for (Int_t i=0;
1853 //__________________________________________________________________________
1854 void AliMUONEventReconstructor::RemoveDoubleTracksK(void)
1856 // Removes double tracks (sharing more than half of their hits). Keeps
1857 // the track with higher quality
1858 AliMUONTrackK *track1, *track2, *trackToKill;
1860 // Sort tracks according to their quality
1861 fRecTracksPtr->Sort();
1863 // Loop over first track of the pair
1864 track1 = (AliMUONTrackK*) fRecTracksPtr->First();
1865 Int_t debug = track1->DebugLevel();
1867 // Loop over second track of the pair
1868 track2 = (AliMUONTrackK*) fRecTracksPtr->After(track1);
1870 // Check whether or not to keep track2
1871 if (!track2->KeepTrack(track1)) {
1872 if (debug >= 0) cout << " Killed track: " << 1/(*track2->GetTrackParameters())(4,0) <<
1873 " " << track2->GetTrackQuality() << endl;
1874 trackToKill = track2;
1875 track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2);
1876 trackToKill->Kill();
1877 fRecTracksPtr->Compress();
1878 } else track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2);
1880 track1 = (AliMUONTrackK*) fRecTracksPtr->After(track1);
1883 fNRecTracks = fRecTracksPtr->GetEntriesFast();
1884 if (debug >= 0) cout << " Number of Kalman tracks: " << fNRecTracks << endl;
1887 //__________________________________________________________________________
1888 void AliMUONEventReconstructor::GoToVertex(void)
1890 // Propagates track to the vertex thru absorber
1891 // (using Branson correction for now)
1895 for (Int_t i=0; i<fNRecTracks; i++) {
1896 //((AliMUONTrackK*)(*fRecTracksPtr)[i])->Branson();
1897 ((AliMUONTrackK*)(*fRecTracksPtr)[i])->SetTrackQuality(1); // compute Chi2
1898 //((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToZ(zVertex); // w/out absorber
1899 ((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToVertex(1); // with absorber
1903 //__________________________________________________________________________
1904 void AliMUONEventReconstructor::SetTrackMethod(Int_t iTrackMethod)
1906 // Set track method and recreate track container if necessary
1908 fTrackMethod = iTrackMethod == 2 ? 2 : 1;
1909 if (fTrackMethod == 2) {
1910 if (fRecTracksPtr) delete fRecTracksPtr;
1911 fRecTracksPtr = new TClonesArray("AliMUONTrackK", 10);
1912 cout << " *** Tracking with the Kalman filter *** " << endl;
1913 } else cout << " *** Traditional tracking *** " << endl;