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 ////////////////////////////////////
35 #include <Riostream.h> // for cout
36 #include <stdlib.h> // for exit()
41 //#include "AliMUONChamber.h"
42 #include "AliMUONEventReconstructor.h"
43 #include "AliMUONHitForRec.h"
44 #include "AliMUONTriggerTrack.h"
45 //#include "AliMUONTriggerConstants.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 "AliConfig.h"
57 #include "AliRunLoader.h"
58 #include "AliLoader.h"
59 #include "AliMUONTrackK.h" //AZ
60 #include <TMatrixD.h> //AZ
63 //************* Defaults parameters for reconstruction
64 static const Double_t kDefaultMinBendingMomentum = 3.0;
65 static const Double_t kDefaultMaxBendingMomentum = 500.0;
66 static const Double_t kDefaultMaxChi2 = 100.0;
67 static const Double_t kDefaultMaxSigma2Distance = 16.0;
68 static const Double_t kDefaultBendingResolution = 0.01;
69 static const Double_t kDefaultNonBendingResolution = 0.144;
70 static const Double_t kDefaultChamberThicknessInX0 = 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 static const Double_t kDefaultSimpleBValue = 7.0;
78 static const Double_t kDefaultSimpleBLength = 428.0;
79 static const Double_t kDefaultSimpleBPosition = 1019.0;
80 static const Int_t kDefaultRecGeantHits = 0;
81 static const Double_t kDefaultEfficiency = 0.95;
83 static const Int_t kDefaultPrintLevel = -1;
85 ClassImp(AliMUONEventReconstructor) // Class implementation in ROOT context
87 //__________________________________________________________________________
88 AliMUONEventReconstructor::AliMUONEventReconstructor(AliLoader* loader)
90 // Constructor for class AliMUONEventReconstructor
91 SetReconstructionParametersToDefaults();
92 fTrackMethod = 1; //AZ - tracking method (1-default, 2-Kalman)
93 // Memory allocation for the TClonesArray of hits for reconstruction
94 // Is 10000 the right size ????
95 fHitsForRecPtr = new TClonesArray("AliMUONHitForRec", 10000);
96 fNHitsForRec = 0; // really needed or GetEntriesFast sufficient ????
97 // Memory allocation for the TClonesArray's of segments in stations
98 // Is 2000 the right size ????
99 for (Int_t st = 0; st < kMaxMuonTrackingStations; st++) {
100 fSegmentsPtr[st] = new TClonesArray("AliMUONSegment", 2000);
101 fNSegments[st] = 0; // really needed or GetEntriesFast sufficient ????
103 // Memory allocation for the TClonesArray of reconstructed tracks
104 // Is 10 the right size ????
105 fRecTracksPtr = new TClonesArray("AliMUONTrack", 10);
106 fNRecTracks = 0; // really needed or GetEntriesFast sufficient ????
108 fRecTriggerTracksPtr = new TClonesArray("AliMUONTriggerTrack", 10);
109 fNRecTriggerTracks = 0; // really needed or GetEntriesFast sufficient ????
110 // Memory allocation for the TClonesArray of hits on reconstructed tracks
111 // Is 100 the right size ????
112 fRecTrackHitsPtr = new TClonesArray("AliMUONTrack", 100);
113 fNRecTrackHits = 0; // really needed or GetEntriesFast sufficient ????
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 if (fPrintLevel >= 0) {
125 cout << "AliMUONEventReconstructor constructed with defaults" << endl; Dump();
126 cout << endl << "Magnetic field from root file:" << endl;
127 gAlice->Field()->Dump();
131 // Initializions for GEANT background events
134 fBkgGeantParticles = 0;
137 fBkgGeantEventNumber = -1;
139 // Initialize to 0 pointers to RecoEvent, tree and tree file
144 // initialize loader's
147 // initialize container
148 fMUONData = new AliMUONData(fLoader,"MUON","MUON");
150 // Loading AliRun master
151 AliRunLoader* runloader = fLoader->GetRunLoader();
152 if (runloader->GetAliRun() == 0x0) runloader->LoadgAlice();
153 gAlice = runloader->GetAliRun();
157 //__________________________________________________________________________
158 AliMUONEventReconstructor::AliMUONEventReconstructor (const AliMUONEventReconstructor& Reconstructor):TObject(Reconstructor)
160 // Dummy copy constructor
163 AliMUONEventReconstructor & AliMUONEventReconstructor::operator=(const AliMUONEventReconstructor& /*Reconstructor*/)
165 // Dummy assignment operator
169 //__________________________________________________________________________
170 AliMUONEventReconstructor::~AliMUONEventReconstructor(void)
172 // Destructor for class AliMUONEventReconstructor
177 // if (fEventTree) delete fEventTree;
178 if (fRecoEvent) delete fRecoEvent;
179 delete fHitsForRecPtr; // Correct destruction of everything ???? or delete [] ????
180 for (Int_t st = 0; st < kMaxMuonTrackingStations; st++)
181 delete fSegmentsPtr[st]; // Correct destruction of everything ????
185 //__________________________________________________________________________
186 void AliMUONEventReconstructor::SetReconstructionParametersToDefaults(void)
188 // Set reconstruction parameters to default values
189 // Would be much more convenient with a structure (or class) ????
190 fMinBendingMomentum = kDefaultMinBendingMomentum;
191 fMaxBendingMomentum = kDefaultMaxBendingMomentum;
192 fMaxChi2 = kDefaultMaxChi2;
193 fMaxSigma2Distance = kDefaultMaxSigma2Distance;
195 AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON");
196 // ******** Parameters for making HitsForRec
198 // like in TRACKF_STAT:
199 // 2 degrees for stations 1 and 2, or ch(0...) from 0 to 3;
200 // 30 cm for stations 3 to 5, or ch(0...) from 4 to 9
201 for (Int_t ch = 0; ch < kMaxMuonTrackingChambers; ch++) {
202 if (ch < 4) fRMin[ch] = TMath::Abs((&(pMUON->Chamber(ch)))->Z()) *
203 2.0 * TMath::Pi() / 180.0;
204 else fRMin[ch] = 30.0;
205 // maximum radius at 10 degrees and Z of chamber
206 fRMax[ch] = TMath::Abs((&(pMUON->Chamber(ch)))->Z()) *
207 10.0 * TMath::Pi() / 180.0;
210 // ******** Parameters for making segments
211 // should be parametrized ????
212 // according to interval between chambers in a station ????
213 // Maximum distance in non bending plane
214 // 5 * 0.22 just to remember the way it was made in TRACKF_STAT
216 for (Int_t st = 0; st < kMaxMuonTrackingStations; st++)
217 fSegmentMaxDistNonBending[st] = 5. * 0.22;
218 // Maximum distance in bending plane:
219 // values from TRACKF_STAT, corresponding to (J psi 20cm),
220 // scaled to the real distance between chambers in a station
221 fSegmentMaxDistBending[0] = TMath::Abs( 1.5 *
222 ((&(pMUON->Chamber(1)))->Z() - (&(pMUON->Chamber(0)))->Z()) / 20.0);
223 fSegmentMaxDistBending[1] = TMath::Abs( 1.5 *
224 ((&(pMUON->Chamber(3)))->Z() - (&(pMUON->Chamber(2)))->Z()) / 20.0);
225 fSegmentMaxDistBending[2] = TMath::Abs( 3.0 *
226 ((&(pMUON->Chamber(5)))->Z() - (&(pMUON->Chamber(4)))->Z()) / 20.0);
227 fSegmentMaxDistBending[3] = TMath::Abs( 6.0 *
228 ((&(pMUON->Chamber(7)))->Z() - (&(pMUON->Chamber(6)))->Z()) / 20.0);
229 fSegmentMaxDistBending[4] = TMath::Abs( 6.0 *
230 ((&(pMUON->Chamber(9)))->Z() - (&(pMUON->Chamber(8)))->Z()) / 20.0);
232 fBendingResolution = kDefaultBendingResolution;
233 fNonBendingResolution = kDefaultNonBendingResolution;
234 fChamberThicknessInX0 = kDefaultChamberThicknessInX0;
235 fSimpleBValue = kDefaultSimpleBValue;
236 fSimpleBLength = kDefaultSimpleBLength;
237 fSimpleBPosition = kDefaultSimpleBPosition;
238 fRecGeantHits = kDefaultRecGeantHits;
239 fEfficiency = kDefaultEfficiency;
240 fPrintLevel = kDefaultPrintLevel;
244 //__________________________________________________________________________
245 Double_t AliMUONEventReconstructor::GetImpactParamFromBendingMomentum(Double_t BendingMomentum) const
247 // Returns impact parameter at vertex in bending plane (cm),
248 // from the signed bending momentum "BendingMomentum" in bending plane (GeV/c),
249 // using simple values for dipole magnetic field.
250 // The sign of "BendingMomentum" is the sign of the charge.
251 return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
255 //__________________________________________________________________________
256 Double_t AliMUONEventReconstructor::GetBendingMomentumFromImpactParam(Double_t ImpactParam) const
258 // Returns signed bending momentum in bending plane (GeV/c),
259 // the sign being the sign of the charge for particles moving forward in Z,
260 // from the impact parameter "ImpactParam" at vertex in bending plane (cm),
261 // using simple values for dipole magnetic field.
262 return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
266 //__________________________________________________________________________
267 void AliMUONEventReconstructor::SetBkgGeantFile(Text_t *BkgGeantFileName)
269 // Set background file ... for GEANT hits
270 // Must be called after having loaded the firts signal event
271 if (fPrintLevel >= 0) {
272 cout << "Enter SetBkgGeantFile with BkgGeantFileName ``"
273 << BkgGeantFileName << "''" << endl;}
274 if (strlen(BkgGeantFileName)) {
275 // BkgGeantFileName not empty: try to open the file
276 if (fPrintLevel >= 2) {cout << "Before File(Bkg)" << endl; gDirectory->Dump();}
277 fBkgGeantFile = new TFile(BkgGeantFileName);
278 if (fPrintLevel >= 2) {cout << "After File(Bkg)" << endl; gDirectory->Dump();}
279 if (fBkgGeantFile-> IsOpen()) {
280 if (fPrintLevel >= 0) {
281 cout << "Background for GEANT hits in file: ``" << BkgGeantFileName
282 << "'' successfully opened" << endl;}
285 cout << "Background for GEANT hits in file: " << BkgGeantFileName << endl;
286 cout << "NOT FOUND: EXIT" << endl;
287 exit(0); // right instruction for exit ????
289 // Arrays for "particles" and "hits"
290 fBkgGeantParticles = new TClonesArray("TParticle", 200);
291 fBkgGeantHits = new TClonesArray("AliMUONHit", 2000);
292 // Event number to -1 for initialization
293 fBkgGeantEventNumber = -1;
294 // Back to the signal file:
295 // first signal event must have been loaded previously,
296 // otherwise, Segmentation violation at the next instruction
297 // How is it possible to do smething better ????
298 ((gAlice->TreeK())->GetCurrentFile())->cd();
299 if (fPrintLevel >= 2) {cout << "After cd(gAlice)" << endl; gDirectory->Dump();}
304 //__________________________________________________________________________
305 void AliMUONEventReconstructor::NextBkgGeantEvent(void)
307 // Get next event in background file for GEANT hits
308 // Goes back to event number 0 when end of file is reached
311 if (fPrintLevel >= 0) {
312 cout << "Enter NextBkgGeantEvent" << endl;}
313 // Clean previous event
314 if(fBkgGeantTK) delete fBkgGeantTK;
316 if(fBkgGeantParticles) fBkgGeantParticles->Clear();
317 if(fBkgGeantTH) delete fBkgGeantTH;
319 if(fBkgGeantHits) fBkgGeantHits->Clear();
320 // Increment event number
321 fBkgGeantEventNumber++;
322 // Get access to Particles and Hits for event from background file
323 if (fPrintLevel >= 2) {cout << "Before cd(Bkg)" << endl; gDirectory->Dump();}
325 if (fPrintLevel >= 2) {cout << "After cd(Bkg)" << endl; gDirectory->Dump();}
326 // Particles: TreeK for event and branch "Particles"
327 sprintf(treeName, "TreeK%d", fBkgGeantEventNumber);
328 fBkgGeantTK = (TTree*)gDirectory->Get(treeName);
330 if (fPrintLevel >= 0) {
331 cout << "Cannot find Kine Tree for background event: " <<
332 fBkgGeantEventNumber << endl;
333 cout << "Goes back to event 0" << endl;
335 fBkgGeantEventNumber = 0;
336 sprintf(treeName, "TreeK%d", fBkgGeantEventNumber);
337 fBkgGeantTK = (TTree*)gDirectory->Get(treeName);
339 cout << "ERROR: cannot find Kine Tree for background event: " <<
340 fBkgGeantEventNumber << endl;
345 fBkgGeantTK->SetBranchAddress("Particles", &fBkgGeantParticles);
346 fBkgGeantTK->GetEvent(0); // why event 0 ???? necessary ????
347 // Hits: TreeH for event and branch "MUON"
348 sprintf(treeName, "TreeH%d", fBkgGeantEventNumber);
349 fBkgGeantTH = (TTree*)gDirectory->Get(treeName);
351 cout << "ERROR: cannot find Hits Tree for background event: " <<
352 fBkgGeantEventNumber << endl;
355 if (fBkgGeantTH && fBkgGeantHits) {
356 branch = fBkgGeantTH->GetBranch("MUON");
357 if (branch) branch->SetAddress(&fBkgGeantHits);
359 fBkgGeantTH->GetEntries(); // necessary ????
360 // Back to the signal file
361 ((gAlice->TreeK())->GetCurrentFile())->cd();
362 if (fPrintLevel >= 2) {cout << "After cd(gAlice)" << endl; gDirectory->Dump();}
366 //__________________________________________________________________________
367 void AliMUONEventReconstructor::EventReconstruct(void)
369 // To reconstruct one event
370 if (fPrintLevel >= 1) cout << "enter EventReconstruct" << endl;
371 MakeEventToBeReconstructed();
374 if (fMUONData->IsTriggerTrackBranchesInTree())
375 ValidateTracksWithTrigger();
377 // Add tracks to MUON data container
378 for(Int_t i=0; i<GetNRecTracks(); i++) {
379 AliMUONTrack * track = (AliMUONTrack*) GetRecTracksPtr()->At(i);
380 fMUONData->AddRecTrack(*track);
386 //__________________________________________________________________________
387 void AliMUONEventReconstructor::EventReconstructTrigger(void)
389 // To reconstruct one event
390 if (fPrintLevel >= 1) cout << "enter EventReconstructTrigger" << endl;
395 //__________________________________________________________________________
396 void AliMUONEventReconstructor::ResetHitsForRec(void)
398 // To reset the array and the number of HitsForRec,
399 // and also the number of HitsForRec
400 // and the index of the first HitForRec per chamber
401 if (fHitsForRecPtr) fHitsForRecPtr->Clear();
403 for (Int_t ch = 0; ch < kMaxMuonTrackingChambers; ch++)
404 fNHitsForRecPerChamber[ch] = fIndexOfFirstHitForRecPerChamber[ch] = 0;
408 //__________________________________________________________________________
409 void AliMUONEventReconstructor::ResetSegments(void)
411 // To reset the TClonesArray of segments and the number of Segments
413 for (Int_t st = 0; st < kMaxMuonTrackingStations; st++) {
414 if (fSegmentsPtr[st]) fSegmentsPtr[st]->Clear();
420 //__________________________________________________________________________
421 void AliMUONEventReconstructor::ResetTracks(void)
423 // To reset the TClonesArray of reconstructed tracks
424 if (fRecTracksPtr) fRecTracksPtr->Delete();
425 // Delete in order that the Track destructors are called,
426 // hence the space for the TClonesArray of pointers to TrackHit's is freed
431 //__________________________________________________________________________
432 void AliMUONEventReconstructor::ResetTriggerTracks(void)
434 // To reset the TClonesArray of reconstructed trigger tracks
435 if (fRecTriggerTracksPtr) fRecTriggerTracksPtr->Delete();
436 // Delete in order that the Track destructors are called,
437 // hence the space for the TClonesArray of pointers to TrackHit's is freed
438 fNRecTriggerTracks = 0;
442 //__________________________________________________________________________
443 void AliMUONEventReconstructor::ResetTrackHits(void)
445 // To reset the TClonesArray of hits on reconstructed tracks
446 if (fRecTrackHitsPtr) fRecTrackHitsPtr->Clear();
451 //__________________________________________________________________________
452 void AliMUONEventReconstructor::MakeEventToBeReconstructed(void)
454 // To make the list of hits to be reconstructed,
455 // either from the GEANT hits or from the raw clusters
456 // according to the parameter set for the reconstructor
457 // TString evfoldname = AliConfig::fgkDefaultEventFolderName;//to be interfaced properly
459 // AliRunLoader* rl = AliRunLoader::GetRunLoader(evfoldname);
462 // Error("MakeEventToBeReconstructed",
463 // "Can not find Run Loader in Event Folder named %s.",
464 // evfoldname.Data());
467 // AliLoader* gime = rl->GetLoader("MUONLoader");
470 // Error("MakeEventToBeReconstructed","Can not get MUON Loader from Run Loader.");
474 if (fPrintLevel >= 1) cout << "enter MakeEventToBeReconstructed" << endl;
476 if (fRecGeantHits == 1) {
477 // Reconstruction from GEANT hits
478 // Back to the signal file
479 TTree* treeH = fLoader->TreeH();
482 Int_t retval = fLoader->LoadHits();
485 Error("MakeEventToBeReconstructed","Error occured while loading hits.");
488 treeH = fLoader->TreeH();
491 Error("MakeEventToBeReconstructed","Can not get TreeH");
495 fMUONData->SetTreeAddress("H");
496 AddHitsForRecFromGEANT(treeH);
499 AddHitsForRecFromBkgGEANT(fBkgGeantTH, fBkgGeantHits);
500 // Sort HitsForRec in increasing order with respect to chamber number
501 SortHitsForRecWithIncreasingChamber();
504 // Reconstruction from raw clusters
505 // AliMUON *MUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
506 // Security on MUON ????
507 // TreeR assumed to be be "prepared" in calling function
508 // by "MUON->GetTreeR(nev)" ????
509 TTree *treeR = fLoader->TreeR();
510 fMUONData->SetTreeAddress("RC");
511 AddHitsForRecFromRawClusters(treeR);
512 // No sorting: it is done automatically in the previous function
514 if (fPrintLevel >= 10) {
515 cout << "end of MakeEventToBeReconstructed" << endl;
516 cout << "NHitsForRec: " << fNHitsForRec << endl;
517 for (Int_t ch = 0; ch < kMaxMuonTrackingChambers; ch++) {
518 cout << "chamber(0...): " << ch
519 << " NHitsForRec: " << fNHitsForRecPerChamber[ch]
520 << " index(first HitForRec): " << fIndexOfFirstHitForRecPerChamber[ch]
522 for (Int_t hit = fIndexOfFirstHitForRecPerChamber[ch];
523 hit < fIndexOfFirstHitForRecPerChamber[ch] + fNHitsForRecPerChamber[ch];
525 cout << "HitForRec index(0...): " << hit << endl;
526 ((*fHitsForRecPtr)[hit])->Dump();
533 //__________________________________________________________________________
534 void AliMUONEventReconstructor::AddHitsForRecFromGEANT(TTree *TH)
536 // To add to the list of hits for reconstruction
537 // the GEANT signal hits from a hit tree TH.
538 Int_t hitBits, chamBits; //AZ
539 if (fPrintLevel >= 2)
540 cout << "enter AddHitsForRecFromGEANT with TH: " << TH << endl;
541 if (TH == NULL) return;
542 // AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
543 //AliMUONData * muondata = pMUON->GetMUONData();
544 // Security on MUON ????
545 // See whether it could be the same for signal and background ????
546 // Loop over tracks in tree
547 Int_t ntracks = (Int_t) TH->GetEntries();
548 if (fPrintLevel >= 2)
549 cout << "ntracks: " << ntracks << endl;
551 for (Int_t track = 0; track < ntracks; track++) {
552 fMUONData->ResetHits();
558 Int_t itrack = track; //AZ
561 nhits = (Int_t) fMUONData->Hits()->GetEntriesFast();
562 AliMUONHit* mHit=0x0;
564 for(ihit=0; ihit<nhits; ihit++) {
565 mHit = static_cast<AliMUONHit*>(fMUONData->Hits()->At(ihit));
566 Int_t ipart = TMath::Abs ((Int_t) mHit->Particle()); //AZ
567 if (NewHitForRecFromGEANT(mHit,track, hit, 1) && ipart == 13
568 //if (NewHitForRecFromGEANT(mHit,itrack-1, hit, 1) && ipart == 13
569 && itrack <= 2 && !BIT(mHit->Chamber()-1) ) chamBits |= BIT(mHit->Chamber()-1); //AZ - set bit
572 if (chamBits&3 && chamBits>>2&3 && chamBits>>4&3 && chamBits>>6&3 &&
573 chamBits>>8&3 && ((chamBits>>6&3)==3 || (chamBits>>8&3)==3))
575 //if (chamBits&3 && chamBits>>2&3 && chamBits>>4&3 && chamBits>>6&3 &&
576 // chamBits>>8&3 && ((chamBits>>6&3)==3 || (chamBits>>8&3)==3) &&
577 // ((chamBits&3)==3 || (chamBits>>2&3)==3)) fMuons += 1;
578 } // end of track loop
582 //__________________________________________________________________________
583 void AliMUONEventReconstructor::AddHitsForRecFromBkgGEANT(TTree *TH, TClonesArray *Hits)
585 // To add to the list of hits for reconstruction
586 // the GEANT background hits from a hit tree TH and a pointer Hits to a hit list.
587 // How to have only one function "AddHitsForRecFromGEANT" ????
588 if (fPrintLevel >= 2)
589 cout << "enter AddHitsForRecFromBkgGEANT with TH: " << TH << endl;
590 if (TH == NULL) return;
591 // Loop over tracks in tree
592 Int_t ntracks = (Int_t) TH->GetEntries();
593 if (fPrintLevel >= 2)
594 cout << "ntracks: " << ntracks << endl;
595 for (Int_t track = 0; track < ntracks; track++) {
596 if (Hits) Hits->Clear();
599 for (Int_t hit = 0; hit < Hits->GetEntriesFast(); hit++) {
600 NewHitForRecFromGEANT((AliMUONHit*) (*Hits)[hit], track, hit, 0);
602 } // end of track loop
606 //__________________________________________________________________________
607 AliMUONHitForRec* AliMUONEventReconstructor::NewHitForRecFromGEANT(AliMUONHit* Hit, Int_t TrackNumber, Int_t HitNumber, Int_t Signal)
609 // To make a new hit for reconstruction from a GEANT hit pointed to by "Hit",
610 // with hit number "HitNumber" in the track numbered "TrackNumber",
611 // either from signal ("Signal" = 1) or background ("Signal" = 0) event.
612 // Selects hits in tracking (not trigger) chambers.
613 // Takes into account the efficiency (fEfficiency)
614 // and the smearing from resolution (fBendingResolution and fNonBendingResolution).
615 // Adds a condition on the radius between RMin and RMax
616 // to better simulate the real chambers.
617 // Returns the pointer to the new hit for reconstruction,
618 // or NULL in case of inefficiency or non tracking chamber or bad radius.
619 // No condition on at most 20 cm from a muon from a resonance
620 // like in Fortran TRACKF_STAT.
621 AliMUONHitForRec* hitForRec;
622 Double_t bendCoor, nonBendCoor, radius;
623 Int_t chamber = Hit->Chamber() - 1; // chamber(0...)
624 // only in tracking chambers (fChamber starts at 1)
625 if (chamber >= kMaxMuonTrackingChambers) return NULL;
626 // only if hit is efficient (keep track for checking ????)
627 if (gRandom->Rndm() > fEfficiency) return NULL;
628 // only if radius between RMin and RMax
630 nonBendCoor = Hit->X();
631 radius = TMath::Sqrt((bendCoor * bendCoor) + (nonBendCoor * nonBendCoor));
632 // This cut is not needed with a realistic chamber geometry !!!!
633 // if ((radius < fRMin[chamber]) || (radius > fRMax[chamber])) return NULL;
634 // new AliMUONHitForRec from GEANT hit and increment number of AliMUONHitForRec's
635 hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(Hit);
637 // add smearing from resolution
638 hitForRec->SetBendingCoor(bendCoor + gRandom->Gaus(0., fBendingResolution));
639 hitForRec->SetNonBendingCoor(nonBendCoor
640 + gRandom->Gaus(0., fNonBendingResolution));
641 // // !!!! without smearing
642 // hitForRec->SetBendingCoor(bendCoor);
643 // hitForRec->SetNonBendingCoor(nonBendCoor);
644 // more information into HitForRec
645 // resolution: angular effect to be added here ????
646 hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
647 hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
649 hitForRec->SetHitNumber(HitNumber);
650 hitForRec->SetTHTrack(TrackNumber);
651 hitForRec->SetGeantSignal(Signal);
652 if (fPrintLevel >= 10) {
653 cout << "track: " << TrackNumber << " hit: " << HitNumber << endl;
655 cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
660 //__________________________________________________________________________
661 void AliMUONEventReconstructor::SortHitsForRecWithIncreasingChamber()
663 // Sort HitsForRec's in increasing order with respect to chamber number.
664 // Uses the function "Compare".
665 // Update the information for HitsForRec per chamber too.
666 Int_t ch, nhits, prevch;
667 fHitsForRecPtr->Sort();
668 for (ch = 0; ch < kMaxMuonTrackingChambers; ch++) {
669 fNHitsForRecPerChamber[ch] = 0;
670 fIndexOfFirstHitForRecPerChamber[ch] = 0;
672 prevch = 0; // previous chamber
673 nhits = 0; // number of hits in current chamber
674 // Loop over HitsForRec
675 for (Int_t hit = 0; hit < fNHitsForRec; hit++) {
676 // chamber number (0...)
677 ch = ((AliMUONHitForRec*) ((*fHitsForRecPtr)[hit]))->GetChamberNumber();
678 // increment number of hits in current chamber
679 (fNHitsForRecPerChamber[ch])++;
680 // update index of first HitForRec in current chamber
681 // if chamber number different from previous one
683 fIndexOfFirstHitForRecPerChamber[ch] = hit;
690 // //__________________________________________________________________________
691 // void AliMUONEventReconstructor::AddHitsForRecFromCathodeCorrelations(TTree* TC)
693 // // OLD VERSION WHEN ONE ONE WAS USING SO CALLED CATHODE CORRELATIONS
694 // // To add to the list of hits for reconstruction
695 // // the (cathode correlated) raw clusters
696 // // No condition added, like in Fortran TRACKF_STAT,
697 // // on the radius between RMin and RMax.
698 // AliMUONHitForRec *hitForRec;
699 // if (fPrintLevel >= 1) cout << "enter AddHitsForRecFromCathodeCorrelations" << endl;
700 // AliMUON *MUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
701 // // Security on MUON ????
702 // // Loop over tracking chambers
703 // for (Int_t ch = 0; ch < kMaxMuonTrackingChambers; ch++) {
704 // // number of HitsForRec to 0 for the chamber
705 // fNHitsForRecPerChamber[ch] = 0;
706 // // index of first HitForRec for the chamber
707 // if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0;
708 // else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec;
709 // TClonesArray *reconst_hits = MUON->ReconstHitsAddress(ch);
710 // MUON->ResetReconstHits();
712 // Int_t ncor = (Int_t)reconst_hits->GetEntries();
713 // // Loop over (cathode correlated) raw clusters
714 // for (Int_t cor = 0; cor < ncor; cor++) {
715 // AliMUONReconstHit * mCor =
716 // (AliMUONReconstHit*) reconst_hits->UncheckedAt(cor);
717 // // new AliMUONHitForRec from (cathode correlated) raw cluster
718 // // and increment number of AliMUONHitForRec's (total and in chamber)
719 // hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(mCor);
721 // (fNHitsForRecPerChamber[ch])++;
722 // // more information into HitForRec
723 // hitForRec->SetChamberNumber(ch);
724 // hitForRec->SetHitNumber(cor);
725 // // Z coordinate of the chamber (cm) with sign opposite to GEANT sign
726 // // could (should) be more exact from chamber geometry ????
727 // hitForRec->SetZ(-(&(MUON->Chamber(ch)))->Z());
728 // if (fPrintLevel >= 10) {
729 // cout << "chamber (0...): " << ch <<
730 // " cathcorrel (0...): " << cor << endl;
732 // cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
733 // hitForRec->Dump();}
734 // } // end of cluster loop
735 // } // end of chamber loop
739 //__________________________________________________________________________
740 void AliMUONEventReconstructor::AddHitsForRecFromRawClusters(TTree* TR)
742 // To add to the list of hits for reconstruction all the raw clusters
743 // No condition added, like in Fortran TRACKF_STAT,
744 // on the radius between RMin and RMax.
745 AliMUONHitForRec *hitForRec;
746 AliMUONRawCluster *clus;
747 Int_t iclus, nclus, nTRentries;
748 TClonesArray *rawclusters;
749 if (fPrintLevel >= 1) cout << "enter AddHitsForRecFromRawClusters" << endl;
751 // TString evfoldname = AliConfig::fgkDefaultEventFolderName;//to be interfaced properly
752 // AliRunLoader* rl = AliRunLoader::GetRunLoader(evfoldname);
755 // Error("MakeEventToBeReconstructed",
756 // "Can not find Run Loader in Event Folder named %s.",
757 // evfoldname.Data());
760 // AliLoader* gime = rl->GetLoader("MUONLoader");
763 // Error("MakeEventToBeReconstructed","Can not get MUON Loader from Run Loader.");
766 // // Loading AliRun master
768 // gAlice = rl->GetAliRun();
770 // Loading MUON subsystem
771 // AliMUON * pMUON = (AliMUON *) gAlice->GetDetector("MUON");
773 nTRentries = Int_t(TR->GetEntries());
774 if (nTRentries != 1) {
775 cout << "Error in AliMUONEventReconstructor::AddHitsForRecFromRawClusters"
777 cout << "nTRentries = " << nTRentries << " not equal to 1" << endl;
780 fLoader->TreeR()->GetEvent(0); // only one entry
782 // Loop over tracking chambers
783 for (Int_t ch = 0; ch < kMaxMuonTrackingChambers; ch++) {
784 // number of HitsForRec to 0 for the chamber
785 fNHitsForRecPerChamber[ch] = 0;
786 // index of first HitForRec for the chamber
787 if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0;
788 else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec;
789 rawclusters =fMUONData->RawClusters(ch);
790 // pMUON->ResetRawClusters();
791 // TR->GetEvent((Int_t) (TR->GetEntries()) - 1); // to be checked ????
792 nclus = (Int_t) (rawclusters->GetEntries());
793 // Loop over (cathode correlated) raw clusters
794 for (iclus = 0; iclus < nclus; iclus++) {
795 clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(iclus);
796 // new AliMUONHitForRec from raw cluster
797 // and increment number of AliMUONHitForRec's (total and in chamber)
798 hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(clus);
800 (fNHitsForRecPerChamber[ch])++;
801 // more information into HitForRec
802 // resolution: info should be already in raw cluster and taken from it ????
803 hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
804 hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
805 // original raw cluster
806 hitForRec->SetChamberNumber(ch);
807 hitForRec->SetHitNumber(iclus);
808 // Z coordinate of the raw cluster (cm)
809 hitForRec->SetZ(clus->fZ[0]);
810 if (fPrintLevel >= 10) {
811 cout << "chamber (0...): " << ch <<
812 " raw cluster (0...): " << iclus << endl;
814 cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
816 } // end of cluster loop
817 } // end of chamber loop
818 SortHitsForRecWithIncreasingChamber(); //AZ
822 //__________________________________________________________________________
823 void AliMUONEventReconstructor::MakeSegments(void)
825 // To make the list of segments in all stations,
826 // from the list of hits to be reconstructed
827 if (fPrintLevel >= 1) cout << "enter MakeSegments" << endl;
829 // Loop over stations
830 Int_t nb = (fTrackMethod == 2) ? 3 : 0; //AZ
831 //AZ for (Int_t st = 0; st < kMaxMuonTrackingStations; st++)
832 for (Int_t st = nb; st < kMaxMuonTrackingStations; st++) //AZ
833 MakeSegmentsPerStation(st);
834 if (fPrintLevel >= 10) {
835 cout << "end of MakeSegments" << endl;
836 for (Int_t st = 0; st < kMaxMuonTrackingStations; st++) {
837 cout << "station(0...): " << st
838 << " Segments: " << fNSegments[st]
841 seg < fNSegments[st];
843 cout << "Segment index(0...): " << seg << endl;
844 ((*fSegmentsPtr[st])[seg])->Dump();
851 //__________________________________________________________________________
852 void AliMUONEventReconstructor::MakeSegmentsPerStation(Int_t Station)
854 // To make the list of segments in station number "Station" (0...)
855 // from the list of hits to be reconstructed.
856 // Updates "fNSegments"[Station].
857 // Segments in stations 4 and 5 are sorted
858 // according to increasing absolute value of "impact parameter"
859 AliMUONHitForRec *hit1Ptr, *hit2Ptr;
860 AliMUONSegment *segment;
862 Double_t bendingSlope, distBend, distNonBend, extBendCoor, extNonBendCoor,
863 impactParam = 0., maxImpactParam = 0., minImpactParam = 0.; // =0 to avoid compilation warnings.
864 AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
865 if (fPrintLevel >= 1)
866 cout << "enter MakeSegmentsPerStation (0...) " << Station << endl;
867 // first and second chambers (0...) in the station
868 Int_t ch1 = 2 * Station;
870 // variable true for stations downstream of the dipole:
871 // Station(0..4) equal to 3 or 4
872 if ((Station == 3) || (Station == 4)) {
874 // maximum impact parameter (cm) according to fMinBendingMomentum...
876 TMath::Abs(GetImpactParamFromBendingMomentum(fMinBendingMomentum));
877 // minimum impact parameter (cm) according to fMaxBendingMomentum...
879 TMath::Abs(GetImpactParamFromBendingMomentum(fMaxBendingMomentum));
881 else last2st = kFALSE;
882 // extrapolation factor from Z of first chamber to Z of second chamber
883 // dZ to be changed to take into account fine structure of chambers ????
884 Double_t extrapFact =
885 (&(pMUON->Chamber(ch2)))->Z() / (&(pMUON->Chamber(ch1)))->Z();
886 // index for current segment
887 Int_t segmentIndex = 0;
888 // Loop over HitsForRec in the first chamber of the station
889 for (Int_t hit1 = fIndexOfFirstHitForRecPerChamber[ch1];
890 hit1 < fIndexOfFirstHitForRecPerChamber[ch1] + fNHitsForRecPerChamber[ch1];
892 // pointer to the HitForRec
893 hit1Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit1]);
895 // on the straight line joining the HitForRec to the vertex (0,0,0),
896 // to the Z of the second chamber of the station
897 extBendCoor = extrapFact * hit1Ptr->GetBendingCoor();
898 extNonBendCoor = extrapFact * hit1Ptr->GetNonBendingCoor();
899 // Loop over HitsForRec in the second chamber of the station
900 for (Int_t hit2 = fIndexOfFirstHitForRecPerChamber[ch2];
901 hit2 < fIndexOfFirstHitForRecPerChamber[ch2] + fNHitsForRecPerChamber[ch2];
903 // pointer to the HitForRec
904 hit2Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit2]);
905 // absolute values of distances, in bending and non bending planes,
906 // between the HitForRec in the second chamber
907 // and the previous extrapolation
908 distBend = TMath::Abs(hit2Ptr->GetBendingCoor() - extBendCoor);
909 distNonBend = TMath::Abs(hit2Ptr->GetNonBendingCoor() - extNonBendCoor);
912 bendingSlope = (hit1Ptr->GetBendingCoor() - hit2Ptr->GetBendingCoor()) /
913 (hit1Ptr->GetZ() - hit2Ptr->GetZ());
914 // absolute value of impact parameter
916 TMath::Abs(hit1Ptr->GetBendingCoor() - hit1Ptr->GetZ() * bendingSlope);
918 // check for distances not too large,
919 // and impact parameter not too big if stations downstream of the dipole.
920 // Conditions "distBend" and "impactParam" correlated for these stations ????
921 if ((distBend < fSegmentMaxDistBending[Station]) &&
922 (distNonBend < fSegmentMaxDistNonBending[Station]) &&
923 (!last2st || (impactParam < maxImpactParam)) &&
924 (!last2st || (impactParam > minImpactParam))) {
926 segment = new ((*fSegmentsPtr[Station])[segmentIndex])
927 AliMUONSegment(hit1Ptr, hit2Ptr);
928 // update "link" to this segment from the hit in the first chamber
929 if (hit1Ptr->GetNSegments() == 0)
930 hit1Ptr->SetIndexOfFirstSegment(segmentIndex);
931 hit1Ptr->SetNSegments(hit1Ptr->GetNSegments() + 1);
932 if (fPrintLevel >= 10) {
933 cout << "segmentIndex(0...): " << segmentIndex
934 << " distBend: " << distBend
935 << " distNonBend: " << distNonBend
938 cout << "HitForRec in first chamber" << endl;
940 cout << "HitForRec in second chamber" << endl;
943 // increment index for current segment
947 } // for (Int_t hit1...
948 fNSegments[Station] = segmentIndex;
949 // Sorting according to "impact parameter" if station(1..5) 4 or 5,
950 // i.e. Station(0..4) 3 or 4, using the function "Compare".
951 // After this sorting, it is impossible to use
952 // the "fNSegments" and "fIndexOfFirstSegment"
953 // of the HitForRec in the first chamber to explore all segments formed with it.
954 // Is this sorting really needed ????
955 if ((Station == 3) || (Station == 4)) (fSegmentsPtr[Station])->Sort();
956 if (fPrintLevel >= 1) cout << "Station: " << Station << " NSegments: "
957 << fNSegments[Station] << endl;
961 //__________________________________________________________________________
962 void AliMUONEventReconstructor::MakeTracks(void)
964 // To make the tracks,
965 // from the list of segments and points in all stations
966 if (fPrintLevel >= 1) cout << "enter MakeTracks" << endl;
967 // The order may be important for the following Reset's
970 if (fTrackMethod == 2) { //AZ - Kalman filter
971 MakeTrackCandidatesK();
972 // Follow tracks in stations(1..) 3, 2 and 1
974 // Remove double tracks
975 RemoveDoubleTracksK();
976 // Propagate tracks to the vertex thru absorber
979 // Look for candidates from at least 3 aligned points in stations(1..) 4 and 5
980 MakeTrackCandidates();
981 // Follow tracks in stations(1..) 3, 2 and 1
983 // Remove double tracks
984 RemoveDoubleTracks();
985 UpdateTrackParamAtHit();
990 //__________________________________________________________________________
991 void AliMUONEventReconstructor::ValidateTracksWithTrigger(void)
994 TClonesArray *recTriggerTracks = fMUONData->RecTriggerTracks();
995 if (recTriggerTracks == 0x0) {
996 fMUONData->SetTreeAddress("RL");
997 fMUONData->GetRecTriggerTracks();
998 recTriggerTracks = fMUONData->RecTriggerTracks();
1000 track = (AliMUONTrack*) fRecTracksPtr->First();
1002 track->MatchTriggerTrack(recTriggerTracks);
1003 track = (AliMUONTrack*) fRecTracksPtr->After(track);
1007 //__________________________________________________________________________
1008 Bool_t AliMUONEventReconstructor::MakeTriggerTracks(void)
1010 // To make the trigger tracks from Local Trigger
1011 if (fPrintLevel >= 1) cout << "enter MakeTriggerTracks" << endl;
1012 // ResetTriggerTracks();
1016 TClonesArray *localTrigger;
1017 TClonesArray *globalTrigger;
1018 AliMUONLocalTrigger *locTrg;
1019 AliMUONGlobalTrigger *gloTrg;
1020 AliMUONTriggerCircuit *circuit;
1021 AliMUONTriggerTrack *recTriggerTrack = 0;
1023 TTree* TR = fLoader->TreeR();
1025 // Loading MUON subsystem
1026 AliMUON * pMUON = (AliMUON *) gAlice->GetDetector("MUON");
1028 nTRentries = Int_t(TR->GetEntries());
1030 TR->GetEvent(0); // only one entry
1032 if (!(fMUONData->IsTriggerBranchesInTree())) {
1033 cout << "Warning in AliMUONEventReconstructor::MakeTriggerTracks"
1035 cout << "Trigger information is not avalaible, nTRentries = " << nTRentries << " not equal to 1" << endl;
1039 fMUONData->SetTreeAddress("GLT");
1040 fMUONData->GetTrigger();
1042 // global trigger for trigger pattern
1044 globalTrigger = fMUONData->GlobalTrigger();
1045 gloTrg = (AliMUONGlobalTrigger*)globalTrigger->UncheckedAt(0);
1046 if (gloTrg->SinglePlusLpt()) gloTrigPat|= 0x1;
1047 if (gloTrg->SinglePlusHpt()) gloTrigPat|= 0x2;
1048 if (gloTrg->SinglePlusApt()) gloTrigPat|= 0x4;
1050 if (gloTrg->SingleMinusLpt()) gloTrigPat|= 0x8;
1051 if (gloTrg->SingleMinusHpt()) gloTrigPat|= 0x10;
1052 if (gloTrg->SingleMinusApt()) gloTrigPat|= 0x20;
1054 if (gloTrg->SingleUndefLpt()) gloTrigPat|= 0x40;
1055 if (gloTrg->SingleUndefHpt()) gloTrigPat|= 0x80;
1056 if (gloTrg->SingleUndefApt()) gloTrigPat|= 0x100;
1058 if (gloTrg->PairUnlikeLpt()) gloTrigPat|= 0x200;
1059 if (gloTrg->PairUnlikeHpt()) gloTrigPat|= 0x400;
1060 if (gloTrg->PairUnlikeApt()) gloTrigPat|= 0x800;
1062 if (gloTrg->PairLikeLpt()) gloTrigPat|= 0x1000;
1063 if (gloTrg->PairLikeHpt()) gloTrigPat|= 0x2000;
1064 if (gloTrg->PairLikeApt()) gloTrigPat|= 0x4000;
1068 // local trigger for tracking
1069 localTrigger = fMUONData->LocalTrigger();
1070 Int_t nlocals = (Int_t) (localTrigger->GetEntries());
1071 Float_t z11 = ( &(pMUON->Chamber(10)) )->Z();
1072 Float_t z21 = ( &(pMUON->Chamber(12)) )->Z();
1074 for (Int_t i=0; i<nlocals; i++) { // loop on Local Trigger
1075 locTrg = (AliMUONLocalTrigger*)localTrigger->UncheckedAt(i);
1076 circuit = &(pMUON->TriggerCircuit(locTrg->LoCircuit()));
1077 Float_t y11 = circuit->GetY11Pos(locTrg->LoStripX());
1078 Int_t stripX21 = locTrg->LoStripX()+locTrg->LoDev()+1;
1079 Float_t y21 = circuit->GetY21Pos(stripX21);
1080 Float_t x11 = circuit->GetX11Pos(locTrg->LoStripY());
1081 Float_t thetax = TMath::ATan2( x11 , z11 );
1082 Float_t thetay = TMath::ATan2( (y21-y11) , (z21-z11) );
1084 recTriggerTrack = new AliMUONTriggerTrack(x11,y11,thetax,thetay,gloTrigPat,this);
1085 // since static statement does not work, set gloTrigPat for each track
1087 // fNRecTriggerTracks++;
1088 fMUONData->AddRecTriggerTrack(*recTriggerTrack);
1089 } // end of loop on Local Trigger
1093 //__________________________________________________________________________
1094 Int_t AliMUONEventReconstructor::MakeTrackCandidatesWithTwoSegments(AliMUONSegment *BegSegment)
1096 // To make track candidates with two segments in stations(1..) 4 and 5,
1097 // the first segment being pointed to by "BegSegment".
1098 // Returns the number of such track candidates.
1099 Int_t endStation, iEndSegment, nbCan2Seg;
1100 AliMUONSegment *endSegment, *extrapSegment;
1101 AliMUONTrack *recTrack;
1103 if (fPrintLevel >= 1) cout << "enter MakeTrackCandidatesWithTwoSegments" << endl;
1104 // Station for the end segment
1105 endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
1106 // multiple scattering factor corresponding to one chamber
1107 mcsFactor = 0.0136 /
1108 GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
1109 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1110 // linear extrapolation to end station
1112 BegSegment->CreateSegmentFromLinearExtrapToStation(endStation, mcsFactor);
1113 // number of candidates with 2 segments to 0
1115 // Loop over segments in the end station
1116 for (iEndSegment = 0; iEndSegment < fNSegments[endStation]; iEndSegment++) {
1117 // pointer to segment
1118 endSegment = (AliMUONSegment*) ((*fSegmentsPtr[endStation])[iEndSegment]);
1119 // test compatibility between current segment and "extrapSegment"
1120 // 4 because 4 quantities in chi2
1122 NormalizedChi2WithSegment(extrapSegment,
1123 fMaxSigma2Distance)) <= 4.0) {
1124 // both segments compatible:
1125 // make track candidate from "begSegment" and "endSegment"
1126 if (fPrintLevel >= 2)
1127 cout << "TrackCandidate with Segment " << iEndSegment <<
1128 " in Station(0..) " << endStation << endl;
1129 // flag for both segments in one track:
1130 // to be done in track constructor ????
1131 BegSegment->SetInTrack(kTRUE);
1132 endSegment->SetInTrack(kTRUE);
1133 recTrack = new ((*fRecTracksPtr)[fNRecTracks])
1134 AliMUONTrack(BegSegment, endSegment, this);
1136 if (fPrintLevel >= 10) recTrack->RecursiveDump();
1137 // increment number of track candidates with 2 segments
1140 } // for (iEndSegment = 0;...
1141 delete extrapSegment; // should not delete HitForRec's it points to !!!!
1145 //__________________________________________________________________________
1146 Int_t AliMUONEventReconstructor::MakeTrackCandidatesWithOneSegmentAndOnePoint(AliMUONSegment *BegSegment)
1148 // To make track candidates with one segment and one point
1149 // in stations(1..) 4 and 5,
1150 // the segment being pointed to by "BegSegment".
1151 Int_t ch, ch1, ch2, endStation, iHit, iHitMax, iHitMin, nbCan1Seg1Hit;
1152 AliMUONHitForRec *extrapHitForRec, *hit;
1153 AliMUONTrack *recTrack;
1155 if (fPrintLevel >= 1)
1156 cout << "enter MakeTrackCandidatesWithOneSegmentAndOnePoint" << endl;
1157 // station for the end point
1158 endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
1159 // multiple scattering factor corresponding to one chamber
1160 mcsFactor = 0.0136 /
1161 GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
1162 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1163 // first and second chambers(0..) in the end station
1164 ch1 = 2 * endStation;
1166 // number of candidates to 0
1168 // Loop over chambers of the end station
1169 for (ch = ch2; ch >= ch1; ch--) {
1170 // linear extrapolation to chamber
1172 BegSegment->CreateHitForRecFromLinearExtrapToChamber(ch, mcsFactor);
1173 // limits for the hit index in the loop
1174 iHitMin = fIndexOfFirstHitForRecPerChamber[ch];
1175 iHitMax = iHitMin + fNHitsForRecPerChamber[ch];
1176 // Loop over HitForRec's in the chamber
1177 for (iHit = iHitMin; iHit < iHitMax; iHit++) {
1178 // pointer to HitForRec
1179 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1180 // test compatibility between current HitForRec and "extrapHitForRec"
1181 // 2 because 2 quantities in chi2
1183 NormalizedChi2WithHitForRec(extrapHitForRec,
1184 fMaxSigma2Distance)) <= 2.0) {
1185 // both HitForRec's compatible:
1186 // make track candidate from begSegment and current HitForRec
1187 if (fPrintLevel >= 2)
1188 cout << "TrackCandidate with HitForRec " << iHit <<
1189 " in Chamber(0..) " << ch << endl;
1190 // flag for beginning segments in one track:
1191 // to be done in track constructor ????
1192 BegSegment->SetInTrack(kTRUE);
1193 recTrack = new ((*fRecTracksPtr)[fNRecTracks])
1194 AliMUONTrack(BegSegment, hit, this);
1195 // the right place to eliminate "double counting" ???? how ????
1197 if (fPrintLevel >= 10) recTrack->RecursiveDump();
1198 // increment number of track candidates
1201 } // for (iHit = iHitMin;...
1202 delete extrapHitForRec;
1203 } // for (ch = ch2;...
1204 return nbCan1Seg1Hit;
1207 //__________________________________________________________________________
1208 void AliMUONEventReconstructor::MakeTrackCandidates(void)
1210 // To make track candidates
1211 // with at least 3 aligned points in stations(1..) 4 and 5
1212 // (two Segment's or one Segment and one HitForRec)
1213 Int_t begStation, iBegSegment, nbCan1Seg1Hit, nbCan2Seg;
1214 AliMUONSegment *begSegment;
1215 if (fPrintLevel >= 1) cout << "enter MakeTrackCandidates" << endl;
1216 // Loop over stations(1..) 5 and 4 for the beginning segment
1217 for (begStation = 4; begStation > 2; begStation--) {
1218 // Loop over segments in the beginning station
1219 for (iBegSegment = 0; iBegSegment < fNSegments[begStation]; iBegSegment++) {
1220 // pointer to segment
1221 begSegment = (AliMUONSegment*) ((*fSegmentsPtr[begStation])[iBegSegment]);
1222 if (fPrintLevel >= 2)
1223 cout << "look for TrackCandidate's with Segment " << iBegSegment <<
1224 " in Station(0..) " << begStation << endl;
1225 // Look for track candidates with two segments,
1226 // "begSegment" and all compatible segments in other station.
1227 // Only for beginning station(1..) 5
1228 // because candidates with 2 segments have to looked for only once.
1229 if (begStation == 4)
1230 nbCan2Seg = MakeTrackCandidatesWithTwoSegments(begSegment);
1231 // Look for track candidates with one segment and one point,
1232 // "begSegment" and all compatible HitForRec's in other station.
1233 // Only if "begSegment" does not belong already to a track candidate.
1234 // Is that a too strong condition ????
1235 if (!(begSegment->GetInTrack()))
1236 nbCan1Seg1Hit = MakeTrackCandidatesWithOneSegmentAndOnePoint(begSegment);
1237 } // for (iBegSegment = 0;...
1238 } // for (begStation = 4;...
1242 //__________________________________________________________________________
1243 void AliMUONEventReconstructor::FollowTracks(void)
1245 // Follow tracks in stations(1..) 3, 2 and 1
1246 // too long: should be made more modular !!!!
1247 AliMUONHitForRec *bestHit, *extrapHit, *extrapCorrHit, *hit;
1248 AliMUONSegment *bestSegment, *extrapSegment, *extrapCorrSegment, *segment;
1249 AliMUONTrack *track, *nextTrack;
1250 AliMUONTrackParam *trackParam1, trackParam[2], trackParamVertex;
1251 // -1 to avoid compilation warnings
1252 Int_t ch = -1, chInStation, chBestHit = -1, iHit, iSegment, station, trackIndex;
1253 Double_t bestChi2, chi2, dZ1, dZ2, dZ3, maxSigma2Distance, mcsFactor;
1254 Double_t bendingMomentum, chi2Norm = 0.;
1255 AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
1256 // local maxSigma2Distance, for easy increase in testing
1257 maxSigma2Distance = fMaxSigma2Distance;
1258 if (fPrintLevel >= 2)
1259 cout << "enter FollowTracks" << endl;
1260 // Loop over track candidates
1261 track = (AliMUONTrack*) fRecTracksPtr->First();
1264 // Follow function for each track candidate ????
1266 nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track); // prepare next track
1267 if (fPrintLevel >= 2)
1268 cout << "FollowTracks: track candidate(0..): " << trackIndex << endl;
1269 // Fit track candidate
1270 track->SetFitMCS(0); // without multiple Coulomb scattering
1271 track->SetFitNParam(3); // with 3 parameters (X = Y = 0)
1272 track->SetFitStart(0); // from parameters at vertex
1274 if (fPrintLevel >= 10) {
1275 cout << "FollowTracks: track candidate(0..): " << trackIndex
1276 << " after fit in stations(0..) 3 and 4" << endl;
1277 track->RecursiveDump();
1279 // Loop over stations(1..) 3, 2 and 1
1280 // something SPECIAL for stations 2 and 1 for majority 3 coincidence ????
1281 // otherwise: majority coincidence 2 !!!!
1282 for (station = 2; station >= 0; station--) {
1283 // Track parameters at first track hit (smallest Z)
1284 trackParam1 = ((AliMUONTrackHit*)
1285 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1286 // extrapolation to station
1287 trackParam1->ExtrapToStation(station, trackParam);
1288 extrapSegment = new AliMUONSegment(); // empty segment
1289 extrapCorrSegment = new AliMUONSegment(); // empty corrected segment
1290 // multiple scattering factor corresponding to one chamber
1291 // and momentum in bending plane (not total)
1292 mcsFactor = 0.0136 * trackParam1->GetInverseBendingMomentum();
1293 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1294 // Z difference from previous station
1295 dZ1 = (&(pMUON->Chamber(2 * station)))->Z() -
1296 (&(pMUON->Chamber(2 * station + 2)))->Z();
1297 // Z difference between the two previous stations
1298 dZ2 = (&(pMUON->Chamber(2 * station + 2)))->Z() -
1299 (&(pMUON->Chamber(2 * station + 4)))->Z();
1300 // Z difference between the two chambers in the previous station
1301 dZ3 = (&(pMUON->Chamber(2 * station)))->Z() -
1302 (&(pMUON->Chamber(2 * station + 1)))->Z();
1303 extrapSegment->SetBendingCoorReso2(fBendingResolution * fBendingResolution);
1305 SetNonBendingCoorReso2(fNonBendingResolution * fNonBendingResolution);
1306 extrapSegment->UpdateFromStationTrackParam
1307 (trackParam, mcsFactor, dZ1, dZ2, dZ3, station,
1308 trackParam1->GetInverseBendingMomentum());
1309 // same thing for corrected segment
1310 // better to use copy constructor, after checking that it works properly !!!!
1311 extrapCorrSegment->SetBendingCoorReso2(fBendingResolution * fBendingResolution);
1313 SetNonBendingCoorReso2(fNonBendingResolution * fNonBendingResolution);
1314 extrapCorrSegment->UpdateFromStationTrackParam
1315 (trackParam, mcsFactor, dZ1, dZ2, dZ3, station,
1316 trackParam1->GetInverseBendingMomentum());
1319 if (fPrintLevel >= 10) {
1320 cout << "FollowTracks: track candidate(0..): " << trackIndex
1321 << " Look for segment in station(0..): " << station << endl;
1323 // Loop over segments in station
1324 for (iSegment = 0; iSegment < fNSegments[station]; iSegment++) {
1325 // Look for best compatible Segment in station
1326 // should consider all possibilities ????
1327 // multiple scattering ????
1328 // separation in 2 functions: Segment and HitForRec ????
1329 segment = (AliMUONSegment*) ((*fSegmentsPtr[station])[iSegment]);
1330 // correction of corrected segment (fBendingCoor and fNonBendingCoor)
1331 // according to real Z value of "segment" and slopes of "extrapSegment"
1333 SetBendingCoor(extrapSegment->GetBendingCoor() +
1334 extrapSegment->GetBendingSlope() *
1335 (segment->GetHitForRec1()->GetZ() -
1336 (&(pMUON->Chamber(2 * station)))->Z()));
1338 SetNonBendingCoor(extrapSegment->GetNonBendingCoor() +
1339 extrapSegment->GetNonBendingSlope() *
1340 (segment->GetHitForRec1()->GetZ() -
1341 (&(pMUON->Chamber(2 * station)))->Z()));
1343 NormalizedChi2WithSegment(extrapCorrSegment, maxSigma2Distance);
1344 if (chi2 < bestChi2) {
1345 // update best Chi2 and Segment if better found
1346 bestSegment = segment;
1351 // best segment found: add it to track candidate
1352 track->AddSegment(bestSegment);
1353 // set track parameters at these two TrakHit's
1354 track->SetTrackParamAtHit(track->GetNTrackHits() - 2, &(trackParam[0]));
1355 track->SetTrackParamAtHit(track->GetNTrackHits() - 1, &(trackParam[1]));
1356 if (fPrintLevel >= 10) {
1357 cout << "FollowTracks: track candidate(0..): " << trackIndex
1358 << " Added segment in station(0..): " << station << endl;
1359 track->RecursiveDump();
1363 // No best segment found:
1364 // Look for best compatible HitForRec in station:
1365 // should consider all possibilities ????
1366 // multiple scattering ???? do about like for extrapSegment !!!!
1367 extrapHit = new AliMUONHitForRec(); // empty hit
1368 extrapCorrHit = new AliMUONHitForRec(); // empty corrected hit
1371 if (fPrintLevel >= 10) {
1372 cout << "FollowTracks: track candidate(0..): " << trackIndex
1373 << " Segment not found, look for hit in station(0..): " << station
1376 // Loop over chambers of the station
1377 for (chInStation = 0; chInStation < 2; chInStation++) {
1378 // coordinates of extrapolated hit
1380 SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor());
1382 SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor());
1383 // resolutions from "extrapSegment"
1384 extrapHit->SetBendingReso2(extrapSegment->GetBendingCoorReso2());
1385 extrapHit->SetNonBendingReso2(extrapSegment->GetNonBendingCoorReso2());
1386 // same things for corrected hit
1387 // better to use copy constructor, after checking that it works properly !!!!
1389 SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor());
1391 SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor());
1392 extrapHit->SetBendingReso2(extrapSegment->GetBendingCoorReso2());
1393 extrapHit->SetNonBendingReso2(extrapSegment->GetNonBendingCoorReso2());
1394 // Loop over hits in the chamber
1395 ch = 2 * station + chInStation;
1396 for (iHit = fIndexOfFirstHitForRecPerChamber[ch];
1397 iHit < fIndexOfFirstHitForRecPerChamber[ch] +
1398 fNHitsForRecPerChamber[ch];
1400 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1401 // correction of corrected hit (fBendingCoor and fNonBendingCoor)
1402 // according to real Z value of "hit" and slopes of right "trackParam"
1404 SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor() +
1405 (&(trackParam[chInStation]))->GetBendingSlope() *
1407 (&(trackParam[chInStation]))->GetZ()));
1409 SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor() +
1410 (&(trackParam[chInStation]))->GetNonBendingSlope() *
1412 (&(trackParam[chInStation]))->GetZ()));
1413 // condition for hit not already in segment ????
1414 chi2 = hit->NormalizedChi2WithHitForRec(extrapHit, maxSigma2Distance);
1415 if (chi2 < bestChi2) {
1416 // update best Chi2 and HitForRec if better found
1419 chBestHit = chInStation;
1424 // best hit found: add it to track candidate
1425 track->AddHitForRec(bestHit);
1426 // set track parameters at this TrackHit
1427 track->SetTrackParamAtHit(track->GetNTrackHits() - 1,
1428 &(trackParam[chBestHit]));
1429 if (fPrintLevel >= 10) {
1430 cout << "FollowTracks: track candidate(0..): " << trackIndex
1431 << " Added hit in station(0..): " << station << endl;
1432 track->RecursiveDump();
1436 // Remove current track candidate
1437 // and corresponding TrackHit's, ...
1439 delete extrapSegment;
1440 delete extrapCorrSegment;
1442 delete extrapCorrHit;
1443 break; // stop the search for this candidate:
1444 // exit from the loop over station
1447 delete extrapCorrHit;
1449 delete extrapSegment;
1450 delete extrapCorrSegment;
1451 // Sort track hits according to increasing Z
1452 track->GetTrackHitsPtr()->Sort();
1453 // Update track parameters at first track hit (smallest Z)
1454 trackParam1 = ((AliMUONTrackHit*)
1455 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1456 bendingMomentum = 0.;
1457 if (TMath::Abs(trackParam1->GetInverseBendingMomentum()) > 0.)
1458 bendingMomentum = TMath::Abs(1/(trackParam1->GetInverseBendingMomentum()));
1459 // Track removed if bendingMomentum not in window [min, max]
1460 if ((bendingMomentum < fMinBendingMomentum) || (bendingMomentum > fMaxBendingMomentum)) {
1462 break; // stop the search for this candidate:
1463 // exit from the loop over station
1466 // with multiple Coulomb scattering if all stations
1467 if (station == 0) track->SetFitMCS(1);
1468 // without multiple Coulomb scattering if not all stations
1469 else track->SetFitMCS(0);
1470 track->SetFitNParam(5); // with 5 parameters (momentum and position)
1471 track->SetFitStart(1); // from parameters at first hit
1473 Double_t numberOfDegFree = (2.0 * track->GetNTrackHits() - 5);
1474 if (numberOfDegFree > 0) {
1475 chi2Norm = track->GetFitFMin() / numberOfDegFree;
1479 // Track removed if normalized chi2 too high
1480 if (chi2Norm > fMaxChi2) {
1482 break; // stop the search for this candidate:
1483 // exit from the loop over station
1485 if (fPrintLevel >= 10) {
1486 cout << "FollowTracks: track candidate(0..): " << trackIndex
1487 << " after fit from station(0..): " << station << " to 4" << endl;
1488 track->RecursiveDump();
1490 // Track extrapolation to the vertex through the absorber (Branson)
1491 // after going through the first station
1493 trackParamVertex = *trackParam1;
1494 (&trackParamVertex)->ExtrapToVertex();
1495 track->SetTrackParamAtVertex(&trackParamVertex);
1496 if (fPrintLevel >= 1) {
1497 cout << "FollowTracks: track candidate(0..): " << trackIndex
1498 << " after extrapolation to vertex" << endl;
1499 track->RecursiveDump();
1502 } // for (station = 2;...
1503 // go really to next track
1506 // Compression of track array (necessary after Remove ????)
1507 fRecTracksPtr->Compress();
1511 //__________________________________________________________________________
1512 void AliMUONEventReconstructor::RemoveDoubleTracks(void)
1514 // To remove double tracks.
1515 // Tracks are considered identical
1516 // if they have at least half of their hits in common.
1517 // Among two identical tracks, one keeps the track with the larger number of hits
1518 // or, if these numbers are equal, the track with the minimum Chi2.
1519 AliMUONTrack *track1, *track2, *trackToRemove;
1520 Bool_t identicalTracks;
1521 Int_t hitsInCommon, nHits1, nHits2;
1522 identicalTracks = kTRUE;
1523 while (identicalTracks) {
1524 identicalTracks = kFALSE;
1525 // Loop over first track of the pair
1526 track1 = (AliMUONTrack*) fRecTracksPtr->First();
1527 while (track1 && (!identicalTracks)) {
1528 nHits1 = track1->GetNTrackHits();
1529 // Loop over second track of the pair
1530 track2 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1531 while (track2 && (!identicalTracks)) {
1532 nHits2 = track2->GetNTrackHits();
1533 // number of hits in common between two tracks
1534 hitsInCommon = track1->HitsInCommon(track2);
1535 // check for identical tracks
1536 if ((4 * hitsInCommon) >= (nHits1 + nHits2)) {
1537 identicalTracks = kTRUE;
1538 // decide which track to remove
1539 if (nHits1 > nHits2) trackToRemove = track2;
1540 else if (nHits1 < nHits2) trackToRemove = track1;
1541 else if ((track1->GetFitFMin()) < (track2->GetFitFMin()))
1542 trackToRemove = track2;
1543 else trackToRemove = track1;
1545 trackToRemove->Remove();
1547 track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
1549 track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1555 //__________________________________________________________________________
1556 void AliMUONEventReconstructor::UpdateTrackParamAtHit()
1558 // Set track parameters after track fitting. Fill fTrackParamAtHit of AliMUONTrack's
1559 AliMUONTrack *track;
1560 AliMUONTrackHit *trackHit;
1561 AliMUONTrackParam *trackParam;
1562 track = (AliMUONTrack*) fRecTracksPtr->First();
1564 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->First();
1566 trackParam = trackHit->GetTrackParam();
1567 track->AddTrackParamAtHit(trackParam);
1568 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->After(trackHit);
1570 track = (AliMUONTrack*) fRecTracksPtr->After(track);
1575 //__________________________________________________________________________
1576 void AliMUONEventReconstructor::EventDump(void)
1578 // Dump reconstructed event (track parameters at vertex and at first hit),
1579 // and the particle parameters
1581 AliMUONTrack *track;
1582 AliMUONTrackParam *trackParam, *trackParam1;
1583 Double_t bendingSlope, nonBendingSlope, pYZ;
1584 Double_t pX, pY, pZ, x, y, z, c;
1585 Int_t np, trackIndex, nTrackHits;
1587 if (fPrintLevel >= 1) cout << "****** enter EventDump ******" << endl;
1588 if (fPrintLevel >= 1) {
1589 cout << " Number of Reconstructed tracks :" << fNRecTracks << endl;
1591 fRecTracksPtr->Compress(); // for simple loop without "Next" since no hole
1592 // Loop over reconstructed tracks
1593 for (trackIndex = 0; trackIndex < fNRecTracks; trackIndex++) {
1594 if (fTrackMethod != 1) continue; //AZ - skip the rest for now
1595 if (fPrintLevel >= 1)
1596 cout << " track number: " << trackIndex << endl;
1597 // function for each track for modularity ????
1598 track = (AliMUONTrack*) ((*fRecTracksPtr)[trackIndex]);
1599 nTrackHits = track->GetNTrackHits();
1600 if (fPrintLevel >= 1)
1601 cout << " number of track hits: " << nTrackHits << endl;
1602 // track parameters at Vertex
1603 trackParam = track->GetTrackParamAtVertex();
1604 x = trackParam->GetNonBendingCoor();
1605 y = trackParam->GetBendingCoor();
1606 z = trackParam->GetZ();
1607 bendingSlope = trackParam->GetBendingSlope();
1608 nonBendingSlope = trackParam->GetNonBendingSlope();
1609 pYZ = 1/TMath::Abs(trackParam->GetInverseBendingMomentum());
1610 pZ = pYZ/TMath::Sqrt(1+bendingSlope*bendingSlope);
1611 pX = pZ * nonBendingSlope;
1612 pY = pZ * bendingSlope;
1613 c = TMath::Sign(1.0, trackParam->GetInverseBendingMomentum());
1614 if (fPrintLevel >= 1)
1615 printf(" track parameters at Vertex z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1616 z, x, y, pX, pY, pZ, c);
1618 // track parameters at first hit
1619 trackParam1 = ((AliMUONTrackHit*)
1620 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1621 x = trackParam1->GetNonBendingCoor();
1622 y = trackParam1->GetBendingCoor();
1623 z = trackParam1->GetZ();
1624 bendingSlope = trackParam1->GetBendingSlope();
1625 nonBendingSlope = trackParam1->GetNonBendingSlope();
1626 pYZ = 1/TMath::Abs(trackParam1->GetInverseBendingMomentum());
1627 pZ = pYZ/TMath::Sqrt(1.0 + bendingSlope * bendingSlope);
1628 pX = pZ * nonBendingSlope;
1629 pY = pZ * bendingSlope;
1630 c = TMath::Sign(1.0, trackParam1->GetInverseBendingMomentum());
1631 if (fPrintLevel >= 1)
1632 printf(" track parameters at z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1633 z, x, y, pX, pY, pZ, c);
1635 // informations about generated particles
1636 np = gAlice->GetMCApp()->GetNtrack();
1637 printf(" **** number of generated particles: %d \n", np);
1639 // for (Int_t iPart = 0; iPart < np; iPart++) {
1640 // p = gAlice->Particle(iPart);
1641 // printf(" particle %d: type= %d px= %f py= %f pz= %f pdg= %d\n",
1642 // iPart, p->GetPdgCode(), p->Px(), p->Py(), p->Pz(), p->GetPdgCode());
1648 //__________________________________________________________________________
1649 void AliMUONEventReconstructor::EventDumpTrigger(void)
1651 // Dump reconstructed trigger event
1652 // and the particle parameters
1654 AliMUONTriggerTrack *triggertrack ;
1655 Int_t nTriggerTracks = fMUONData->RecTriggerTracks()->GetEntriesFast();
1657 if (fPrintLevel >= 1) {
1658 cout << "****** enter EventDumpTrigger ******" << endl;
1659 cout << " Number of Reconstructed tracks :" << nTriggerTracks << endl;
1661 // Loop over reconstructed tracks
1662 for (Int_t trackIndex = 0; trackIndex < nTriggerTracks; trackIndex++) {
1663 triggertrack = (AliMUONTriggerTrack*)fMUONData->RecTriggerTracks()->At(trackIndex);
1664 printf(" trigger track number %i x11=%f y11=%f thetax=%f thetay=%f \n",
1666 triggertrack->GetX11(),triggertrack->GetY11(),
1667 triggertrack->GetThetax(),triggertrack->GetThetay());
1671 //__________________________________________________________________________
1672 void AliMUONEventReconstructor::FillEvent()
1674 // Create a new AliMUONRecoEvent, fill its track list, then add it as a
1675 // leaf in the Event branch of TreeRecoEvent tree
1676 cout << "Enter FillEvent() ...\n";
1679 fRecoEvent = new AliMUONRecoEvent();
1681 fRecoEvent->Clear();
1683 //save current directory
1684 TDirectory *current = gDirectory;
1685 if (!fTreeFile) fTreeFile = new TFile("tree_reco.root", "RECREATE");
1686 if (!fEventTree) fEventTree = new TTree("TreeRecoEvent", "MUON reconstructed events");
1687 //AZif (fRecoEvent->MakeDumpTracks(fRecTracksPtr)) {
1688 if (fRecoEvent->MakeDumpTracks(fMuons, fRecTracksPtr, this)) { //AZ
1689 if (fPrintLevel > 1) fRecoEvent->EventInfo();
1690 TBranch *branch = fEventTree->GetBranch("Event");
1691 if (!branch) branch = fEventTree->Branch("Event", "AliMUONRecoEvent", &fRecoEvent, 64000);
1692 branch->SetAutoDelete();
1697 // restore directory
1701 //__________________________________________________________________________
1702 void AliMUONEventReconstructor::MakeTrackCandidatesK(void)
1704 // To make initial tracks for Kalman filter from the list of segments
1706 AliMUONSegment *segment;
1707 AliMUONTrackK *trackK;
1709 if (fPrintLevel >= 1) cout << "enter MakeTrackCandidatesK" << endl;
1710 // Reset the TClonesArray of reconstructed tracks
1711 if (fRecTracksPtr) fRecTracksPtr->Delete();
1712 // Delete in order that the Track destructors are called,
1713 // hence the space for the TClonesArray of pointers to TrackHit's is freed
1716 AliMUONTrackK a(this, fHitsForRecPtr); // bad idea ???
1717 // Loop over stations(1...) 5 and 4
1718 for (istat=4; istat>=3; istat--) {
1719 // Loop over segments in the station
1720 for (iseg=0; iseg<fNSegments[istat]; iseg++) {
1721 // Transform segments to tracks and evaluate covariance matrix
1722 segment = (AliMUONSegment*) ((*fSegmentsPtr[istat])[iseg]);
1723 trackK = new ((*fRecTracksPtr)[fNRecTracks]) AliMUONTrackK(segment);
1725 } // for (iseg=0;...)
1726 } // for (istat=4;...)
1730 //__________________________________________________________________________
1731 void AliMUONEventReconstructor::FollowTracksK(void)
1733 // Follow tracks using Kalman filter
1735 Int_t icand, ichamBeg, ichamEnd, chamBits;
1736 Double_t zDipole1, zDipole2;
1737 AliMUONTrackK *trackK;
1738 AliMUONHitForRec *hit;
1739 AliMUONRawCluster *clus;
1740 TClonesArray *rawclusters;
1742 clus = 0; rawclusters = 0;
1744 zDipole1 = GetSimpleBPosition() - GetSimpleBLength()/2;
1745 zDipole2 = zDipole1 + GetSimpleBLength();
1748 pMUON = (AliMUON*) gAlice->GetModule("MUON");
1749 for (Int_t i1=0; i1<fNHitsForRec; i1++) {
1750 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[i1]);
1751 //if (hit->GetTHTrack() > 1 || hit->GetGeantSignal() == 0) continue;
1753 cout << " Hit #" << hit->GetChamberNumber() << " ";
1754 cout << hit->GetBendingCoor() << " ";
1755 cout << hit->GetNonBendingCoor() << " ";
1756 cout << hit->GetZ() << " ";
1757 cout << hit->GetGeantSignal() << " ";
1758 cout << hit->GetTHTrack() << endl;
1761 printf(" Hit # %d %10.4f %10.4f %10.4f",
1762 hit->GetChamberNumber(), hit->GetBendingCoor(),
1763 hit->GetNonBendingCoor(), hit->GetZ());
1764 if (fRecGeantHits) {
1766 printf(" %3d %3d \n", hit->GetGeantSignal(), hit->GetTHTrack());
1768 // from raw clusters
1769 rawclusters = pMUON->RawClustAddress(hit->GetChamberNumber());
1770 clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(hit->
1772 printf("%3d", clus->fTracks[1]-1);
1773 if (clus->fTracks[2] != 0) printf("%3d \n", clus->fTracks[2]-1);
1780 Int_t nSeeds = fNRecTracks; // starting number of seeds
1781 // Loop over track candidates
1782 while (icand < fNRecTracks-1) {
1784 trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]);
1786 // Discard candidate which will produce the double track
1788 Ok = CheckCandidateK(icand,nSeeds);
1790 //trackK->SetRecover(-1); // mark candidate to be removed
1796 if (trackK->GetRecover() == 0) hit = (AliMUONHitForRec*)
1797 trackK->GetHitOnTrack()->Last(); // last hit
1798 else hit = (AliMUONHitForRec*) (*trackK->GetHitOnTrack())[1]; // 2'nd hit
1799 ichamBeg = hit->GetChamberNumber();
1801 // Check propagation direction
1802 if (trackK->GetTrackDir() > 0) {
1803 ichamEnd = 9; // forward propagation
1804 Ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kFALSE,zDipole1,zDipole2);
1806 ichamBeg = ichamEnd;
1807 ichamEnd = 6; // backward propagation
1808 // Change weight matrix and zero fChi2 for backpropagation
1809 trackK->StartBack();
1810 Ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kTRUE,zDipole1,zDipole2);
1811 ichamBeg = ichamEnd;
1815 if (trackK->GetBPFlag()) {
1817 ichamEnd = 6; // backward propagation
1818 // Change weight matrix and zero fChi2 for backpropagation
1819 trackK->StartBack();
1820 Ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kTRUE,zDipole1,zDipole2);
1821 ichamBeg = ichamEnd;
1827 trackK->SetTrackDir(-1);
1828 trackK->SetBPFlag(kFALSE);
1829 Ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kFALSE,zDipole1,zDipole2);
1831 if (Ok) trackK->SetTrackQuality(0); // compute "track quality"
1832 else trackK->SetRecover(-1); // mark candidate to be removed
1834 // Majority 3 of 4 in first 2 stations
1836 for (Int_t i=0; i<trackK->GetNTrackHits(); i++) {
1837 hit = (AliMUONHitForRec*) (*trackK->GetHitOnTrack())[i];
1838 chamBits |= BIT(hit->GetChamberNumber()-1);
1840 //if (!((chamBits&3)==3 || (chamBits>>2&3)==3)) trackK->SetRecover(-1);
1841 //mark candidate to be removed
1844 for (Int_t i=0; i<fNRecTracks; i++) {
1845 trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[i]);
1846 if (trackK->GetRecover() < 0) fRecTracksPtr->RemoveAt(i);
1849 // Compress TClonesArray
1850 fRecTracksPtr->Compress();
1851 fNRecTracks = fRecTracksPtr->GetEntriesFast();
1855 //__________________________________________________________________________
1856 Bool_t AliMUONEventReconstructor::CheckCandidateK(Int_t icand, Int_t nSeeds)
1858 // Discards track candidate if it will produce the double track (having
1859 // the same seed segment hits as hits of a good track found before)
1860 AliMUONTrackK *track1, *track2;
1861 AliMUONHitForRec *hit1, *hit2, *hit;
1863 track1 = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]);
1864 hit1 = (AliMUONHitForRec*) (*track1->GetHitOnTrack())[0]; // 1'st hit
1865 hit2 = (AliMUONHitForRec*) (*track1->GetHitOnTrack())[1]; // 2'nd hit
1867 for (Int_t i=0; i<icand; i++) {
1868 track2 = (AliMUONTrackK*) ((*fRecTracksPtr)[i]);
1869 //if (track2->GetRecover() < 0) continue;
1870 if (track2->GetRecover() < 0 && icand >= nSeeds) continue;
1872 if (track1->GetStartSegment() == track2->GetStartSegment()) {
1876 for (Int_t j=0; j<track2->GetNTrackHits(); j++) {
1877 hit = (AliMUONHitForRec*) (*track2->GetHitOnTrack())[j];
1878 if (hit == hit1 || hit == hit2) {
1880 if (nSame == 2) return kFALSE;
1882 } // for (Int_t j=0;
1884 } // for (Int_t i=0;
1888 //__________________________________________________________________________
1889 void AliMUONEventReconstructor::RemoveDoubleTracksK(void)
1891 // Removes double tracks (sharing more than half of their hits). Keeps
1892 // the track with higher quality
1893 AliMUONTrackK *track1, *track2, *trackToKill;
1895 // Sort tracks according to their quality
1896 fRecTracksPtr->Sort();
1898 // Loop over first track of the pair
1899 track1 = (AliMUONTrackK*) fRecTracksPtr->First();
1901 // Loop over second track of the pair
1902 track2 = (AliMUONTrackK*) fRecTracksPtr->After(track1);
1904 // Check whether or not to keep track2
1905 if (!track2->KeepTrack(track1)) {
1906 cout << " Killed track: " << 1/(*track2->GetTrackParameters())(4,0) <<
1907 " " << track2->GetTrackQuality() << endl;
1908 trackToKill = track2;
1909 track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2);
1910 trackToKill->Kill();
1911 fRecTracksPtr->Compress();
1912 } else track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2);
1914 track1 = (AliMUONTrackK*) fRecTracksPtr->After(track1);
1917 fNRecTracks = fRecTracksPtr->GetEntriesFast();
1918 cout << " Number of Kalman tracks: " << fNRecTracks << endl;
1921 //__________________________________________________________________________
1922 void AliMUONEventReconstructor::GoToVertex(void)
1924 // Propagates track to the vertex thru absorber
1925 // (using Branson correction for now)
1929 for (Int_t i=0; i<fNRecTracks; i++) {
1930 //((AliMUONTrackK*)(*fRecTracksPtr)[i])->Branson();
1931 //((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToZ(zVertex); // w/out absorber
1932 ((AliMUONTrackK*)(*fRecTracksPtr)[i])->SetTrackQuality(1); // compute Chi2
1933 ((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToVertex(); // with absorber