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 "AliMUONHit.h"
44 #include "AliMUONHitForRec.h"
45 #include "AliMUONTriggerTrack.h"
46 #include "AliMUONTriggerCircuit.h"
47 #include "AliMUONRawCluster.h"
48 #include "AliMUONLocalTrigger.h"
49 #include "AliMUONGlobalTrigger.h"
50 #include "AliMUONRecoEvent.h"
51 #include "AliMUONSegment.h"
52 #include "AliMUONTrack.h"
53 #include "AliMUONTrackHit.h"
55 #include "AliRun.h" // for gAlice
56 #include "AliRunLoader.h"
57 #include "AliLoader.h"
58 #include "AliMUONTrackK.h" //AZ
61 //************* Defaults parameters for reconstruction
62 const Double_t AliMUONEventReconstructor::fgkDefaultMinBendingMomentum = 3.0;
63 const Double_t AliMUONEventReconstructor::fgkDefaultMaxBendingMomentum = 500.0;
64 const Double_t AliMUONEventReconstructor::fgkDefaultMaxChi2 = 100.0;
65 const Double_t AliMUONEventReconstructor::fgkDefaultMaxSigma2Distance = 16.0;
66 const Double_t AliMUONEventReconstructor::fgkDefaultBendingResolution = 0.01;
67 const Double_t AliMUONEventReconstructor::fgkDefaultNonBendingResolution = 0.144;
68 const Double_t AliMUONEventReconstructor::fgkDefaultChamberThicknessInX0 = 0.03;
69 // Simple magnetic field:
70 // Value taken from macro MUONtracking.C: 0.7 T, hence 7 kG
71 // Length and Position from reco_muon.F, with opposite sign:
72 // Length = ZMAGEND-ZCOIL
73 // Position = (ZMAGEND+ZCOIL)/2
74 // to be ajusted differently from real magnetic field ????
75 const Double_t AliMUONEventReconstructor::fgkDefaultSimpleBValue = 7.0;
76 const Double_t AliMUONEventReconstructor::fgkDefaultSimpleBLength = 428.0;
77 const Double_t AliMUONEventReconstructor::fgkDefaultSimpleBPosition = 1019.0;
78 const Int_t AliMUONEventReconstructor::fgkDefaultRecGeantHits = 0;
79 const Double_t AliMUONEventReconstructor::fgkDefaultEfficiency = 0.95;
81 const Int_t AliMUONEventReconstructor::fgkDefaultPrintLevel = -1;
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 < fgkMaxMuonTrackingStations; 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 fRecTrackHitsPtr = new TClonesArray("AliMUONTrack", 100);
112 fNRecTrackHits = 0; // really needed or GetEntriesFast sufficient ????
114 // Sign of fSimpleBValue according to sign of Bx value at (50,50,-950).
116 x[0] = 50.; x[1] = 50.; x[2] = -950.;
117 gAlice->Field()->Field(x, b);
118 fSimpleBValue = TMath::Sign(fSimpleBValue,(Double_t) b[0]);
119 fSimpleBPosition = TMath::Sign(fSimpleBPosition,(Double_t) x[2]);
120 // See how to get fSimple(BValue, BLength, BPosition)
121 // automatically calculated from the actual magnetic field ????
123 if (fPrintLevel >= 0) {
124 cout << "AliMUONEventReconstructor constructed with defaults" << endl; Dump();
125 cout << endl << "Magnetic field from root file:" << endl;
126 gAlice->Field()->Dump();
130 // Initializions for GEANT background events
133 fBkgGeantParticles = 0;
136 fBkgGeantEventNumber = -1;
138 // Initialize to 0 pointers to RecoEvent, tree and tree file
143 // initialize loader's
146 // initialize container
147 fMUONData = new AliMUONData(fLoader,"MUON","MUON");
149 // Loading AliRun master
150 AliRunLoader* runloader = fLoader->GetRunLoader();
151 if (runloader->GetAliRun() == 0x0) runloader->LoadgAlice();
152 gAlice = runloader->GetAliRun();
156 //__________________________________________________________________________
157 AliMUONEventReconstructor::AliMUONEventReconstructor (const AliMUONEventReconstructor& rhs)
160 // Protected copy constructor
162 Fatal("AliMUONEventReconstructor", "Not implemented.");
165 AliMUONEventReconstructor &
166 AliMUONEventReconstructor::operator=(const AliMUONEventReconstructor& rhs)
168 // Protected assignement operator
170 if (this == &rhs) return *this;
172 Fatal("operator=", "Not implemented.");
177 //__________________________________________________________________________
178 AliMUONEventReconstructor::~AliMUONEventReconstructor(void)
180 // Destructor for class AliMUONEventReconstructor
185 // if (fEventTree) delete fEventTree;
186 if (fRecoEvent) delete fRecoEvent;
187 delete fHitsForRecPtr; // Correct destruction of everything ???? or delete [] ????
188 for (Int_t st = 0; st < fgkMaxMuonTrackingStations; st++)
189 delete fSegmentsPtr[st]; // Correct destruction of everything ????
193 //__________________________________________________________________________
194 void AliMUONEventReconstructor::SetReconstructionParametersToDefaults(void)
196 // Set reconstruction parameters to default values
197 // Would be much more convenient with a structure (or class) ????
198 fMinBendingMomentum = fgkDefaultMinBendingMomentum;
199 fMaxBendingMomentum = fgkDefaultMaxBendingMomentum;
200 fMaxChi2 = fgkDefaultMaxChi2;
201 fMaxSigma2Distance = fgkDefaultMaxSigma2Distance;
203 AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON");
204 // ******** Parameters for making HitsForRec
206 // like in TRACKF_STAT:
207 // 2 degrees for stations 1 and 2, or ch(0...) from 0 to 3;
208 // 30 cm for stations 3 to 5, or ch(0...) from 4 to 9
209 for (Int_t ch = 0; ch < fgkMaxMuonTrackingChambers; ch++) {
210 if (ch < 4) fRMin[ch] = TMath::Abs((&(pMUON->Chamber(ch)))->Z()) *
211 2.0 * TMath::Pi() / 180.0;
212 else fRMin[ch] = 30.0;
213 // maximum radius at 10 degrees and Z of chamber
214 fRMax[ch] = TMath::Abs((&(pMUON->Chamber(ch)))->Z()) *
215 10.0 * TMath::Pi() / 180.0;
218 // ******** Parameters for making segments
219 // should be parametrized ????
220 // according to interval between chambers in a station ????
221 // Maximum distance in non bending plane
222 // 5 * 0.22 just to remember the way it was made in TRACKF_STAT
224 for (Int_t st = 0; st < fgkMaxMuonTrackingStations; st++)
225 fSegmentMaxDistNonBending[st] = 5. * 0.22;
226 // Maximum distance in bending plane:
227 // values from TRACKF_STAT, corresponding to (J psi 20cm),
228 // scaled to the real distance between chambers in a station
229 fSegmentMaxDistBending[0] = TMath::Abs( 1.5 *
230 ((&(pMUON->Chamber(1)))->Z() - (&(pMUON->Chamber(0)))->Z()) / 20.0);
231 fSegmentMaxDistBending[1] = TMath::Abs( 1.5 *
232 ((&(pMUON->Chamber(3)))->Z() - (&(pMUON->Chamber(2)))->Z()) / 20.0);
233 fSegmentMaxDistBending[2] = TMath::Abs( 3.0 *
234 ((&(pMUON->Chamber(5)))->Z() - (&(pMUON->Chamber(4)))->Z()) / 20.0);
235 fSegmentMaxDistBending[3] = TMath::Abs( 6.0 *
236 ((&(pMUON->Chamber(7)))->Z() - (&(pMUON->Chamber(6)))->Z()) / 20.0);
237 fSegmentMaxDistBending[4] = TMath::Abs( 6.0 *
238 ((&(pMUON->Chamber(9)))->Z() - (&(pMUON->Chamber(8)))->Z()) / 20.0);
240 fBendingResolution = fgkDefaultBendingResolution;
241 fNonBendingResolution = fgkDefaultNonBendingResolution;
242 fChamberThicknessInX0 = fgkDefaultChamberThicknessInX0;
243 fSimpleBValue = fgkDefaultSimpleBValue;
244 fSimpleBLength = fgkDefaultSimpleBLength;
245 fSimpleBPosition = fgkDefaultSimpleBPosition;
246 fRecGeantHits = fgkDefaultRecGeantHits;
247 fEfficiency = fgkDefaultEfficiency;
248 fPrintLevel = fgkDefaultPrintLevel;
252 //__________________________________________________________________________
253 Double_t AliMUONEventReconstructor::GetImpactParamFromBendingMomentum(Double_t BendingMomentum) const
255 // Returns impact parameter at vertex in bending plane (cm),
256 // from the signed bending momentum "BendingMomentum" in bending plane (GeV/c),
257 // using simple values for dipole magnetic field.
258 // The sign of "BendingMomentum" is the sign of the charge.
259 return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
263 //__________________________________________________________________________
264 Double_t AliMUONEventReconstructor::GetBendingMomentumFromImpactParam(Double_t ImpactParam) const
266 // Returns signed bending momentum in bending plane (GeV/c),
267 // the sign being the sign of the charge for particles moving forward in Z,
268 // from the impact parameter "ImpactParam" at vertex in bending plane (cm),
269 // using simple values for dipole magnetic field.
270 return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
274 //__________________________________________________________________________
275 void AliMUONEventReconstructor::SetBkgGeantFile(Text_t *BkgGeantFileName)
277 // Set background file ... for GEANT hits
278 // Must be called after having loaded the firts signal event
279 if (fPrintLevel >= 0) {
280 cout << "Enter SetBkgGeantFile with BkgGeantFileName ``"
281 << BkgGeantFileName << "''" << endl;}
282 if (strlen(BkgGeantFileName)) {
283 // BkgGeantFileName not empty: try to open the file
284 if (fPrintLevel >= 2) {cout << "Before File(Bkg)" << endl; gDirectory->Dump();}
285 fBkgGeantFile = new TFile(BkgGeantFileName);
286 if (fPrintLevel >= 2) {cout << "After File(Bkg)" << endl; gDirectory->Dump();}
287 if (fBkgGeantFile-> IsOpen()) {
288 if (fPrintLevel >= 0) {
289 cout << "Background for GEANT hits in file: ``" << BkgGeantFileName
290 << "'' successfully opened" << endl;}
293 cout << "Background for GEANT hits in file: " << BkgGeantFileName << endl;
294 cout << "NOT FOUND: EXIT" << endl;
295 exit(0); // right instruction for exit ????
297 // Arrays for "particles" and "hits"
298 fBkgGeantParticles = new TClonesArray("TParticle", 200);
299 fBkgGeantHits = new TClonesArray("AliMUONHit", 2000);
300 // Event number to -1 for initialization
301 fBkgGeantEventNumber = -1;
302 // Back to the signal file:
303 // first signal event must have been loaded previously,
304 // otherwise, Segmentation violation at the next instruction
305 // How is it possible to do smething better ????
306 ((gAlice->TreeK())->GetCurrentFile())->cd();
307 if (fPrintLevel >= 2) {cout << "After cd(gAlice)" << endl; gDirectory->Dump();}
312 //__________________________________________________________________________
313 void AliMUONEventReconstructor::NextBkgGeantEvent(void)
315 // Get next event in background file for GEANT hits
316 // Goes back to event number 0 when end of file is reached
319 if (fPrintLevel >= 0) {
320 cout << "Enter NextBkgGeantEvent" << endl;}
321 // Clean previous event
322 if(fBkgGeantTK) delete fBkgGeantTK;
324 if(fBkgGeantParticles) fBkgGeantParticles->Clear();
325 if(fBkgGeantTH) delete fBkgGeantTH;
327 if(fBkgGeantHits) fBkgGeantHits->Clear();
328 // Increment event number
329 fBkgGeantEventNumber++;
330 // Get access to Particles and Hits for event from background file
331 if (fPrintLevel >= 2) {cout << "Before cd(Bkg)" << endl; gDirectory->Dump();}
333 if (fPrintLevel >= 2) {cout << "After cd(Bkg)" << endl; gDirectory->Dump();}
334 // Particles: TreeK for event and branch "Particles"
335 sprintf(treeName, "TreeK%d", fBkgGeantEventNumber);
336 fBkgGeantTK = (TTree*)gDirectory->Get(treeName);
338 if (fPrintLevel >= 0) {
339 cout << "Cannot find Kine Tree for background event: " <<
340 fBkgGeantEventNumber << endl;
341 cout << "Goes back to event 0" << endl;
343 fBkgGeantEventNumber = 0;
344 sprintf(treeName, "TreeK%d", fBkgGeantEventNumber);
345 fBkgGeantTK = (TTree*)gDirectory->Get(treeName);
347 cout << "ERROR: cannot find Kine Tree for background event: " <<
348 fBkgGeantEventNumber << endl;
353 fBkgGeantTK->SetBranchAddress("Particles", &fBkgGeantParticles);
354 fBkgGeantTK->GetEvent(0); // why event 0 ???? necessary ????
355 // Hits: TreeH for event and branch "MUON"
356 sprintf(treeName, "TreeH%d", fBkgGeantEventNumber);
357 fBkgGeantTH = (TTree*)gDirectory->Get(treeName);
359 cout << "ERROR: cannot find Hits Tree for background event: " <<
360 fBkgGeantEventNumber << endl;
363 if (fBkgGeantTH && fBkgGeantHits) {
364 branch = fBkgGeantTH->GetBranch("MUON");
365 if (branch) branch->SetAddress(&fBkgGeantHits);
367 fBkgGeantTH->GetEntries(); // necessary ????
368 // Back to the signal file
369 ((gAlice->TreeK())->GetCurrentFile())->cd();
370 if (fPrintLevel >= 2) {cout << "After cd(gAlice)" << endl; gDirectory->Dump();}
374 //__________________________________________________________________________
375 void AliMUONEventReconstructor::EventReconstruct(void)
377 // To reconstruct one event
378 if (fPrintLevel >= 1) cout << "enter EventReconstruct" << endl;
379 MakeEventToBeReconstructed();
382 if (fMUONData->IsTriggerTrackBranchesInTree())
383 ValidateTracksWithTrigger();
385 // Add tracks to MUON data container
386 for(Int_t i=0; i<GetNRecTracks(); i++) {
387 AliMUONTrack * track = (AliMUONTrack*) GetRecTracksPtr()->At(i);
388 fMUONData->AddRecTrack(*track);
394 //__________________________________________________________________________
395 void AliMUONEventReconstructor::EventReconstructTrigger(void)
397 // To reconstruct one event
398 if (fPrintLevel >= 1) cout << "enter EventReconstructTrigger" << endl;
403 //__________________________________________________________________________
404 void AliMUONEventReconstructor::ResetHitsForRec(void)
406 // To reset the array and the number of HitsForRec,
407 // and also the number of HitsForRec
408 // and the index of the first HitForRec per chamber
409 if (fHitsForRecPtr) fHitsForRecPtr->Clear();
411 for (Int_t ch = 0; ch < fgkMaxMuonTrackingChambers; ch++)
412 fNHitsForRecPerChamber[ch] = fIndexOfFirstHitForRecPerChamber[ch] = 0;
416 //__________________________________________________________________________
417 void AliMUONEventReconstructor::ResetSegments(void)
419 // To reset the TClonesArray of segments and the number of Segments
421 for (Int_t st = 0; st < fgkMaxMuonTrackingStations; st++) {
422 if (fSegmentsPtr[st]) fSegmentsPtr[st]->Clear();
428 //__________________________________________________________________________
429 void AliMUONEventReconstructor::ResetTracks(void)
431 // To reset the TClonesArray of reconstructed tracks
432 if (fRecTracksPtr) fRecTracksPtr->Delete();
433 // Delete in order that the Track destructors are called,
434 // hence the space for the TClonesArray of pointers to TrackHit's is freed
439 //__________________________________________________________________________
440 void AliMUONEventReconstructor::ResetTriggerTracks(void)
442 // To reset the TClonesArray of reconstructed trigger tracks
443 if (fRecTriggerTracksPtr) fRecTriggerTracksPtr->Delete();
444 // Delete in order that the Track destructors are called,
445 // hence the space for the TClonesArray of pointers to TrackHit's is freed
446 fNRecTriggerTracks = 0;
450 //__________________________________________________________________________
451 void AliMUONEventReconstructor::ResetTrackHits(void)
453 // To reset the TClonesArray of hits on reconstructed tracks
454 if (fRecTrackHitsPtr) fRecTrackHitsPtr->Clear();
459 //__________________________________________________________________________
460 void AliMUONEventReconstructor::MakeEventToBeReconstructed(void)
462 // To make the list of hits to be reconstructed,
463 // either from the GEANT hits or from the raw clusters
464 // according to the parameter set for the reconstructor
465 // TString evfoldname = AliConfig::GetDefaultEventFolderName();//to be interfaced properly
467 // AliRunLoader* rl = AliRunLoader::GetRunLoader(evfoldname);
470 // Error("MakeEventToBeReconstructed",
471 // "Can not find Run Loader in Event Folder named %s.",
472 // evfoldname.Data());
475 // AliLoader* gime = rl->GetLoader("MUONLoader");
478 // Error("MakeEventToBeReconstructed","Can not get MUON Loader from Run Loader.");
482 if (fPrintLevel >= 1) cout << "enter MakeEventToBeReconstructed" << endl;
484 if (fRecGeantHits == 1) {
485 // Reconstruction from GEANT hits
486 // Back to the signal file
487 TTree* treeH = fLoader->TreeH();
490 Int_t retval = fLoader->LoadHits();
493 Error("MakeEventToBeReconstructed","Error occured while loading hits.");
496 treeH = fLoader->TreeH();
499 Error("MakeEventToBeReconstructed","Can not get TreeH");
503 fMUONData->SetTreeAddress("H");
504 AddHitsForRecFromGEANT(treeH);
507 AddHitsForRecFromBkgGEANT(fBkgGeantTH, fBkgGeantHits);
508 // Sort HitsForRec in increasing order with respect to chamber number
509 SortHitsForRecWithIncreasingChamber();
512 // Reconstruction from raw clusters
513 // AliMUON *MUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
514 // Security on MUON ????
515 // TreeR assumed to be be "prepared" in calling function
516 // by "MUON->GetTreeR(nev)" ????
517 TTree *treeR = fLoader->TreeR();
518 fMUONData->SetTreeAddress("RC");
519 AddHitsForRecFromRawClusters(treeR);
520 // No sorting: it is done automatically in the previous function
522 if (fPrintLevel >= 10) {
523 cout << "end of MakeEventToBeReconstructed" << endl;
524 cout << "NHitsForRec: " << fNHitsForRec << endl;
525 for (Int_t ch = 0; ch < fgkMaxMuonTrackingChambers; ch++) {
526 cout << "chamber(0...): " << ch
527 << " NHitsForRec: " << fNHitsForRecPerChamber[ch]
528 << " index(first HitForRec): " << fIndexOfFirstHitForRecPerChamber[ch]
530 for (Int_t hit = fIndexOfFirstHitForRecPerChamber[ch];
531 hit < fIndexOfFirstHitForRecPerChamber[ch] + fNHitsForRecPerChamber[ch];
533 cout << "HitForRec index(0...): " << hit << endl;
534 ((*fHitsForRecPtr)[hit])->Dump();
541 //__________________________________________________________________________
542 void AliMUONEventReconstructor::AddHitsForRecFromGEANT(TTree *TH)
544 // To add to the list of hits for reconstruction
545 // the GEANT signal hits from a hit tree TH.
546 Int_t hitBits, chamBits; //AZ
547 if (fPrintLevel >= 2)
548 cout << "enter AddHitsForRecFromGEANT with TH: " << TH << endl;
549 if (TH == NULL) return;
550 // AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
551 //AliMUONData * muondata = pMUON->GetMUONData();
552 // Security on MUON ????
553 // See whether it could be the same for signal and background ????
554 // Loop over tracks in tree
555 Int_t ntracks = (Int_t) TH->GetEntries();
556 if (fPrintLevel >= 2)
557 cout << "ntracks: " << ntracks << endl;
559 for (Int_t track = 0; track < ntracks; track++) {
560 fMUONData->ResetHits();
566 Int_t itrack = track; //AZ
569 nhits = (Int_t) fMUONData->Hits()->GetEntriesFast();
570 AliMUONHit* mHit=0x0;
572 for(ihit=0; ihit<nhits; ihit++) {
573 mHit = static_cast<AliMUONHit*>(fMUONData->Hits()->At(ihit));
574 Int_t ipart = TMath::Abs ((Int_t) mHit->Particle()); //AZ
575 if (NewHitForRecFromGEANT(mHit,track, hit, 1) && ipart == 13
576 //if (NewHitForRecFromGEANT(mHit,itrack-1, hit, 1) && ipart == 13
577 && itrack <= 2 && !BIT(mHit->Chamber()-1) ) chamBits |= BIT(mHit->Chamber()-1); //AZ - set bit
580 if (chamBits&3 && chamBits>>2&3 && chamBits>>4&3 && chamBits>>6&3 &&
581 chamBits>>8&3 && ((chamBits>>6&3)==3 || (chamBits>>8&3)==3))
583 //if (chamBits&3 && chamBits>>2&3 && chamBits>>4&3 && chamBits>>6&3 &&
584 // chamBits>>8&3 && ((chamBits>>6&3)==3 || (chamBits>>8&3)==3) &&
585 // ((chamBits&3)==3 || (chamBits>>2&3)==3)) fMuons += 1;
586 } // end of track loop
590 //__________________________________________________________________________
591 void AliMUONEventReconstructor::AddHitsForRecFromBkgGEANT(TTree *TH, TClonesArray *Hits)
593 // To add to the list of hits for reconstruction
594 // the GEANT background hits from a hit tree TH and a pointer Hits to a hit list.
595 // How to have only one function "AddHitsForRecFromGEANT" ????
596 if (fPrintLevel >= 2)
597 cout << "enter AddHitsForRecFromBkgGEANT with TH: " << TH << endl;
598 if (TH == NULL) return;
599 // Loop over tracks in tree
600 Int_t ntracks = (Int_t) TH->GetEntries();
601 if (fPrintLevel >= 2)
602 cout << "ntracks: " << ntracks << endl;
603 for (Int_t track = 0; track < ntracks; track++) {
604 if (Hits) Hits->Clear();
607 for (Int_t hit = 0; hit < Hits->GetEntriesFast(); hit++) {
608 NewHitForRecFromGEANT((AliMUONHit*) (*Hits)[hit], track, hit, 0);
610 } // end of track loop
614 //__________________________________________________________________________
615 AliMUONHitForRec* AliMUONEventReconstructor::NewHitForRecFromGEANT(AliMUONHit* Hit, Int_t TrackNumber, Int_t HitNumber, Int_t Signal)
617 // To make a new hit for reconstruction from a GEANT hit pointed to by "Hit",
618 // with hit number "HitNumber" in the track numbered "TrackNumber",
619 // either from signal ("Signal" = 1) or background ("Signal" = 0) event.
620 // Selects hits in tracking (not trigger) chambers.
621 // Takes into account the efficiency (fEfficiency)
622 // and the smearing from resolution (fBendingResolution and fNonBendingResolution).
623 // Adds a condition on the radius between RMin and RMax
624 // to better simulate the real chambers.
625 // Returns the pointer to the new hit for reconstruction,
626 // or NULL in case of inefficiency or non tracking chamber or bad radius.
627 // No condition on at most 20 cm from a muon from a resonance
628 // like in Fortran TRACKF_STAT.
629 AliMUONHitForRec* hitForRec;
630 Double_t bendCoor, nonBendCoor, radius;
631 Int_t chamber = Hit->Chamber() - 1; // chamber(0...)
632 // only in tracking chambers (fChamber starts at 1)
633 if (chamber >= fgkMaxMuonTrackingChambers) return NULL;
634 // only if hit is efficient (keep track for checking ????)
635 if (gRandom->Rndm() > fEfficiency) return NULL;
636 // only if radius between RMin and RMax
638 nonBendCoor = Hit->X();
639 radius = TMath::Sqrt((bendCoor * bendCoor) + (nonBendCoor * nonBendCoor));
640 // This cut is not needed with a realistic chamber geometry !!!!
641 // if ((radius < fRMin[chamber]) || (radius > fRMax[chamber])) return NULL;
642 // new AliMUONHitForRec from GEANT hit and increment number of AliMUONHitForRec's
643 hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(Hit);
645 // add smearing from resolution
646 hitForRec->SetBendingCoor(bendCoor + gRandom->Gaus(0., fBendingResolution));
647 hitForRec->SetNonBendingCoor(nonBendCoor
648 + gRandom->Gaus(0., fNonBendingResolution));
649 // // !!!! without smearing
650 // hitForRec->SetBendingCoor(bendCoor);
651 // hitForRec->SetNonBendingCoor(nonBendCoor);
652 // more information into HitForRec
653 // resolution: angular effect to be added here ????
654 hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
655 hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
657 hitForRec->SetHitNumber(HitNumber);
658 hitForRec->SetTHTrack(TrackNumber);
659 hitForRec->SetGeantSignal(Signal);
660 if (fPrintLevel >= 10) {
661 cout << "track: " << TrackNumber << " hit: " << HitNumber << endl;
663 cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
668 //__________________________________________________________________________
669 void AliMUONEventReconstructor::SortHitsForRecWithIncreasingChamber()
671 // Sort HitsForRec's in increasing order with respect to chamber number.
672 // Uses the function "Compare".
673 // Update the information for HitsForRec per chamber too.
674 Int_t ch, nhits, prevch;
675 fHitsForRecPtr->Sort();
676 for (ch = 0; ch < fgkMaxMuonTrackingChambers; ch++) {
677 fNHitsForRecPerChamber[ch] = 0;
678 fIndexOfFirstHitForRecPerChamber[ch] = 0;
680 prevch = 0; // previous chamber
681 nhits = 0; // number of hits in current chamber
682 // Loop over HitsForRec
683 for (Int_t hit = 0; hit < fNHitsForRec; hit++) {
684 // chamber number (0...)
685 ch = ((AliMUONHitForRec*) ((*fHitsForRecPtr)[hit]))->GetChamberNumber();
686 // increment number of hits in current chamber
687 (fNHitsForRecPerChamber[ch])++;
688 // update index of first HitForRec in current chamber
689 // if chamber number different from previous one
691 fIndexOfFirstHitForRecPerChamber[ch] = hit;
698 // //__________________________________________________________________________
699 // void AliMUONEventReconstructor::AddHitsForRecFromCathodeCorrelations(TTree* TC)
701 // // OLD VERSION WHEN ONE ONE WAS USING SO CALLED CATHODE CORRELATIONS
702 // // To add to the list of hits for reconstruction
703 // // the (cathode correlated) raw clusters
704 // // No condition added, like in Fortran TRACKF_STAT,
705 // // on the radius between RMin and RMax.
706 // AliMUONHitForRec *hitForRec;
707 // if (fPrintLevel >= 1) cout << "enter AddHitsForRecFromCathodeCorrelations" << endl;
708 // AliMUON *MUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
709 // // Security on MUON ????
710 // // Loop over tracking chambers
711 // for (Int_t ch = 0; ch < fgkMaxMuonTrackingChambers; ch++) {
712 // // number of HitsForRec to 0 for the chamber
713 // fNHitsForRecPerChamber[ch] = 0;
714 // // index of first HitForRec for the chamber
715 // if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0;
716 // else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec;
717 // TClonesArray *reconst_hits = MUON->ReconstHitsAddress(ch);
718 // MUON->ResetReconstHits();
720 // Int_t ncor = (Int_t)reconst_hits->GetEntries();
721 // // Loop over (cathode correlated) raw clusters
722 // for (Int_t cor = 0; cor < ncor; cor++) {
723 // AliMUONReconstHit * mCor =
724 // (AliMUONReconstHit*) reconst_hits->UncheckedAt(cor);
725 // // new AliMUONHitForRec from (cathode correlated) raw cluster
726 // // and increment number of AliMUONHitForRec's (total and in chamber)
727 // hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(mCor);
729 // (fNHitsForRecPerChamber[ch])++;
730 // // more information into HitForRec
731 // hitForRec->SetChamberNumber(ch);
732 // hitForRec->SetHitNumber(cor);
733 // // Z coordinate of the chamber (cm) with sign opposite to GEANT sign
734 // // could (should) be more exact from chamber geometry ????
735 // hitForRec->SetZ(-(&(MUON->Chamber(ch)))->Z());
736 // if (fPrintLevel >= 10) {
737 // cout << "chamber (0...): " << ch <<
738 // " cathcorrel (0...): " << cor << endl;
740 // cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
741 // hitForRec->Dump();}
742 // } // end of cluster loop
743 // } // end of chamber loop
747 //__________________________________________________________________________
748 void AliMUONEventReconstructor::AddHitsForRecFromRawClusters(TTree* TR)
750 // To add to the list of hits for reconstruction all the raw clusters
751 // No condition added, like in Fortran TRACKF_STAT,
752 // on the radius between RMin and RMax.
753 AliMUONHitForRec *hitForRec;
754 AliMUONRawCluster *clus;
755 Int_t iclus, nclus, nTRentries;
756 TClonesArray *rawclusters;
757 if (fPrintLevel >= 1) cout << "enter AddHitsForRecFromRawClusters" << endl;
759 // TString evfoldname = AliConfig::GetDefaultEventFolderName();//to be interfaced properly
760 // AliRunLoader* rl = AliRunLoader::GetRunLoader(evfoldname);
763 // Error("MakeEventToBeReconstructed",
764 // "Can not find Run Loader in Event Folder named %s.",
765 // evfoldname.Data());
768 // AliLoader* gime = rl->GetLoader("MUONLoader");
771 // Error("MakeEventToBeReconstructed","Can not get MUON Loader from Run Loader.");
774 // // Loading AliRun master
776 // gAlice = rl->GetAliRun();
778 // Loading MUON subsystem
779 // AliMUON * pMUON = (AliMUON *) gAlice->GetDetector("MUON");
781 nTRentries = Int_t(TR->GetEntries());
782 if (nTRentries != 1) {
783 cout << "Error in AliMUONEventReconstructor::AddHitsForRecFromRawClusters"
785 cout << "nTRentries = " << nTRentries << " not equal to 1" << endl;
788 fLoader->TreeR()->GetEvent(0); // only one entry
790 // Loop over tracking chambers
791 for (Int_t ch = 0; ch < fgkMaxMuonTrackingChambers; ch++) {
792 // number of HitsForRec to 0 for the chamber
793 fNHitsForRecPerChamber[ch] = 0;
794 // index of first HitForRec for the chamber
795 if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0;
796 else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec;
797 rawclusters =fMUONData->RawClusters(ch);
798 // pMUON->ResetRawClusters();
799 // TR->GetEvent((Int_t) (TR->GetEntries()) - 1); // to be checked ????
800 nclus = (Int_t) (rawclusters->GetEntries());
801 // Loop over (cathode correlated) raw clusters
802 for (iclus = 0; iclus < nclus; iclus++) {
803 clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(iclus);
804 // new AliMUONHitForRec from raw cluster
805 // and increment number of AliMUONHitForRec's (total and in chamber)
806 hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(clus);
808 (fNHitsForRecPerChamber[ch])++;
809 // more information into HitForRec
810 // resolution: info should be already in raw cluster and taken from it ????
811 hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
812 hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
813 // original raw cluster
814 hitForRec->SetChamberNumber(ch);
815 hitForRec->SetHitNumber(iclus);
816 // Z coordinate of the raw cluster (cm)
817 hitForRec->SetZ(clus->GetZ(0));
818 if (fPrintLevel >= 10) {
819 cout << "chamber (0...): " << ch <<
820 " raw cluster (0...): " << iclus << endl;
822 cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
824 } // end of cluster loop
825 } // end of chamber loop
826 SortHitsForRecWithIncreasingChamber(); //AZ
830 //__________________________________________________________________________
831 void AliMUONEventReconstructor::MakeSegments(void)
833 // To make the list of segments in all stations,
834 // from the list of hits to be reconstructed
835 if (fPrintLevel >= 1) cout << "enter MakeSegments" << endl;
837 // Loop over stations
838 Int_t nb = (fTrackMethod == 2) ? 3 : 0; //AZ
839 //AZ for (Int_t st = 0; st < fgkMaxMuonTrackingStations; st++)
840 for (Int_t st = nb; st < fgkMaxMuonTrackingStations; st++) //AZ
841 MakeSegmentsPerStation(st);
842 if (fPrintLevel >= 10) {
843 cout << "end of MakeSegments" << endl;
844 for (Int_t st = 0; st < fgkMaxMuonTrackingStations; st++) {
845 cout << "station(0...): " << st
846 << " Segments: " << fNSegments[st]
849 seg < fNSegments[st];
851 cout << "Segment index(0...): " << seg << endl;
852 ((*fSegmentsPtr[st])[seg])->Dump();
859 //__________________________________________________________________________
860 void AliMUONEventReconstructor::MakeSegmentsPerStation(Int_t Station)
862 // To make the list of segments in station number "Station" (0...)
863 // from the list of hits to be reconstructed.
864 // Updates "fNSegments"[Station].
865 // Segments in stations 4 and 5 are sorted
866 // according to increasing absolute value of "impact parameter"
867 AliMUONHitForRec *hit1Ptr, *hit2Ptr;
868 AliMUONSegment *segment;
870 Double_t bendingSlope, distBend, distNonBend, extBendCoor, extNonBendCoor,
871 impactParam = 0., maxImpactParam = 0., minImpactParam = 0.; // =0 to avoid compilation warnings.
872 AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
873 if (fPrintLevel >= 1)
874 cout << "enter MakeSegmentsPerStation (0...) " << Station << endl;
875 // first and second chambers (0...) in the station
876 Int_t ch1 = 2 * Station;
878 // variable true for stations downstream of the dipole:
879 // Station(0..4) equal to 3 or 4
880 if ((Station == 3) || (Station == 4)) {
882 // maximum impact parameter (cm) according to fMinBendingMomentum...
884 TMath::Abs(GetImpactParamFromBendingMomentum(fMinBendingMomentum));
885 // minimum impact parameter (cm) according to fMaxBendingMomentum...
887 TMath::Abs(GetImpactParamFromBendingMomentum(fMaxBendingMomentum));
889 else last2st = kFALSE;
890 // extrapolation factor from Z of first chamber to Z of second chamber
891 // dZ to be changed to take into account fine structure of chambers ????
892 Double_t extrapFact =
893 (&(pMUON->Chamber(ch2)))->Z() / (&(pMUON->Chamber(ch1)))->Z();
894 // index for current segment
895 Int_t segmentIndex = 0;
896 // Loop over HitsForRec in the first chamber of the station
897 for (Int_t hit1 = fIndexOfFirstHitForRecPerChamber[ch1];
898 hit1 < fIndexOfFirstHitForRecPerChamber[ch1] + fNHitsForRecPerChamber[ch1];
900 // pointer to the HitForRec
901 hit1Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit1]);
903 // on the straight line joining the HitForRec to the vertex (0,0,0),
904 // to the Z of the second chamber of the station
905 extBendCoor = extrapFact * hit1Ptr->GetBendingCoor();
906 extNonBendCoor = extrapFact * hit1Ptr->GetNonBendingCoor();
907 // Loop over HitsForRec in the second chamber of the station
908 for (Int_t hit2 = fIndexOfFirstHitForRecPerChamber[ch2];
909 hit2 < fIndexOfFirstHitForRecPerChamber[ch2] + fNHitsForRecPerChamber[ch2];
911 // pointer to the HitForRec
912 hit2Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit2]);
913 // absolute values of distances, in bending and non bending planes,
914 // between the HitForRec in the second chamber
915 // and the previous extrapolation
916 distBend = TMath::Abs(hit2Ptr->GetBendingCoor() - extBendCoor);
917 distNonBend = TMath::Abs(hit2Ptr->GetNonBendingCoor() - extNonBendCoor);
920 bendingSlope = (hit1Ptr->GetBendingCoor() - hit2Ptr->GetBendingCoor()) /
921 (hit1Ptr->GetZ() - hit2Ptr->GetZ());
922 // absolute value of impact parameter
924 TMath::Abs(hit1Ptr->GetBendingCoor() - hit1Ptr->GetZ() * bendingSlope);
926 // check for distances not too large,
927 // and impact parameter not too big if stations downstream of the dipole.
928 // Conditions "distBend" and "impactParam" correlated for these stations ????
929 if ((distBend < fSegmentMaxDistBending[Station]) &&
930 (distNonBend < fSegmentMaxDistNonBending[Station]) &&
931 (!last2st || (impactParam < maxImpactParam)) &&
932 (!last2st || (impactParam > minImpactParam))) {
934 segment = new ((*fSegmentsPtr[Station])[segmentIndex])
935 AliMUONSegment(hit1Ptr, hit2Ptr);
936 // update "link" to this segment from the hit in the first chamber
937 if (hit1Ptr->GetNSegments() == 0)
938 hit1Ptr->SetIndexOfFirstSegment(segmentIndex);
939 hit1Ptr->SetNSegments(hit1Ptr->GetNSegments() + 1);
940 if (fPrintLevel >= 10) {
941 cout << "segmentIndex(0...): " << segmentIndex
942 << " distBend: " << distBend
943 << " distNonBend: " << distNonBend
946 cout << "HitForRec in first chamber" << endl;
948 cout << "HitForRec in second chamber" << endl;
951 // increment index for current segment
955 } // for (Int_t hit1...
956 fNSegments[Station] = segmentIndex;
957 // Sorting according to "impact parameter" if station(1..5) 4 or 5,
958 // i.e. Station(0..4) 3 or 4, using the function "Compare".
959 // After this sorting, it is impossible to use
960 // the "fNSegments" and "fIndexOfFirstSegment"
961 // of the HitForRec in the first chamber to explore all segments formed with it.
962 // Is this sorting really needed ????
963 if ((Station == 3) || (Station == 4)) (fSegmentsPtr[Station])->Sort();
964 if (fPrintLevel >= 1) cout << "Station: " << Station << " NSegments: "
965 << fNSegments[Station] << endl;
969 //__________________________________________________________________________
970 void AliMUONEventReconstructor::MakeTracks(void)
972 // To make the tracks,
973 // from the list of segments and points in all stations
974 if (fPrintLevel >= 1) cout << "enter MakeTracks" << endl;
975 // The order may be important for the following Reset's
978 if (fTrackMethod == 2) { //AZ - Kalman filter
979 MakeTrackCandidatesK();
980 // Follow tracks in stations(1..) 3, 2 and 1
982 // Remove double tracks
983 RemoveDoubleTracksK();
984 // Propagate tracks to the vertex thru absorber
987 // Look for candidates from at least 3 aligned points in stations(1..) 4 and 5
988 MakeTrackCandidates();
989 // Follow tracks in stations(1..) 3, 2 and 1
991 // Remove double tracks
992 RemoveDoubleTracks();
993 UpdateTrackParamAtHit();
998 //__________________________________________________________________________
999 void AliMUONEventReconstructor::ValidateTracksWithTrigger(void)
1001 AliMUONTrack *track;
1002 TClonesArray *recTriggerTracks;
1004 fMUONData->ResetRecTriggerTracks();
1005 fMUONData->SetTreeAddress("RL");
1006 fMUONData->GetRecTriggerTracks();
1007 recTriggerTracks = fMUONData->RecTriggerTracks();
1009 track = (AliMUONTrack*) fRecTracksPtr->First();
1011 track->MatchTriggerTrack(recTriggerTracks);
1012 track = (AliMUONTrack*) fRecTracksPtr->After(track);
1017 //__________________________________________________________________________
1018 Bool_t AliMUONEventReconstructor::MakeTriggerTracks(void)
1020 // To make the trigger tracks from Local Trigger
1021 if (fPrintLevel >= 1) cout << "enter MakeTriggerTracks" << endl;
1022 // ResetTriggerTracks();
1026 TClonesArray *localTrigger;
1027 TClonesArray *globalTrigger;
1028 AliMUONLocalTrigger *locTrg;
1029 AliMUONGlobalTrigger *gloTrg;
1030 AliMUONTriggerCircuit *circuit;
1031 AliMUONTriggerTrack *recTriggerTrack = 0;
1033 TTree* treeR = fLoader->TreeR();
1035 // Loading MUON subsystem
1036 AliMUON * pMUON = (AliMUON *) gAlice->GetDetector("MUON");
1038 nTRentries = Int_t(treeR->GetEntries());
1040 treeR->GetEvent(0); // only one entry
1042 if (!(fMUONData->IsTriggerBranchesInTree())) {
1043 cout << "Warning in AliMUONEventReconstructor::MakeTriggerTracks"
1045 cout << "Trigger information is not avalaible, nTRentries = " << nTRentries << " not equal to 1" << endl;
1049 fMUONData->SetTreeAddress("GLT");
1050 fMUONData->GetTrigger();
1052 // global trigger for trigger pattern
1054 globalTrigger = fMUONData->GlobalTrigger();
1055 gloTrg = (AliMUONGlobalTrigger*)globalTrigger->UncheckedAt(0);
1056 if (gloTrg->SinglePlusLpt()) gloTrigPat|= 0x1;
1057 if (gloTrg->SinglePlusHpt()) gloTrigPat|= 0x2;
1058 if (gloTrg->SinglePlusApt()) gloTrigPat|= 0x4;
1060 if (gloTrg->SingleMinusLpt()) gloTrigPat|= 0x8;
1061 if (gloTrg->SingleMinusHpt()) gloTrigPat|= 0x10;
1062 if (gloTrg->SingleMinusApt()) gloTrigPat|= 0x20;
1064 if (gloTrg->SingleUndefLpt()) gloTrigPat|= 0x40;
1065 if (gloTrg->SingleUndefHpt()) gloTrigPat|= 0x80;
1066 if (gloTrg->SingleUndefApt()) gloTrigPat|= 0x100;
1068 if (gloTrg->PairUnlikeLpt()) gloTrigPat|= 0x200;
1069 if (gloTrg->PairUnlikeHpt()) gloTrigPat|= 0x400;
1070 if (gloTrg->PairUnlikeApt()) gloTrigPat|= 0x800;
1072 if (gloTrg->PairLikeLpt()) gloTrigPat|= 0x1000;
1073 if (gloTrg->PairLikeHpt()) gloTrigPat|= 0x2000;
1074 if (gloTrg->PairLikeApt()) gloTrigPat|= 0x4000;
1078 // local trigger for tracking
1079 localTrigger = fMUONData->LocalTrigger();
1080 Int_t nlocals = (Int_t) (localTrigger->GetEntries());
1081 Float_t z11 = ( &(pMUON->Chamber(10)) )->Z();
1082 Float_t z21 = ( &(pMUON->Chamber(12)) )->Z();
1084 for (Int_t i=0; i<nlocals; i++) { // loop on Local Trigger
1085 locTrg = (AliMUONLocalTrigger*)localTrigger->UncheckedAt(i);
1086 circuit = &(pMUON->TriggerCircuit(locTrg->LoCircuit()));
1087 Float_t y11 = circuit->GetY11Pos(locTrg->LoStripX());
1088 Int_t stripX21 = locTrg->LoStripX()+locTrg->LoDev()+1;
1089 Float_t y21 = circuit->GetY21Pos(stripX21);
1090 Float_t x11 = circuit->GetX11Pos(locTrg->LoStripY());
1091 Float_t thetax = TMath::ATan2( x11 , z11 );
1092 Float_t thetay = TMath::ATan2( (y21-y11) , (z21-z11) );
1094 recTriggerTrack = new AliMUONTriggerTrack(x11,y11,thetax,thetay,gloTrigPat,this);
1095 // since static statement does not work, set gloTrigPat for each track
1097 // fNRecTriggerTracks++;
1098 fMUONData->AddRecTriggerTrack(*recTriggerTrack);
1099 } // end of loop on Local Trigger
1103 //__________________________________________________________________________
1104 Int_t AliMUONEventReconstructor::MakeTrackCandidatesWithTwoSegments(AliMUONSegment *BegSegment)
1106 // To make track candidates with two segments in stations(1..) 4 and 5,
1107 // the first segment being pointed to by "BegSegment".
1108 // Returns the number of such track candidates.
1109 Int_t endStation, iEndSegment, nbCan2Seg;
1110 AliMUONSegment *endSegment, *extrapSegment;
1111 AliMUONTrack *recTrack;
1113 if (fPrintLevel >= 1) cout << "enter MakeTrackCandidatesWithTwoSegments" << endl;
1114 // Station for the end segment
1115 endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
1116 // multiple scattering factor corresponding to one chamber
1117 mcsFactor = 0.0136 /
1118 GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
1119 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1120 // linear extrapolation to end station
1122 BegSegment->CreateSegmentFromLinearExtrapToStation(endStation, mcsFactor);
1123 // number of candidates with 2 segments to 0
1125 // Loop over segments in the end station
1126 for (iEndSegment = 0; iEndSegment < fNSegments[endStation]; iEndSegment++) {
1127 // pointer to segment
1128 endSegment = (AliMUONSegment*) ((*fSegmentsPtr[endStation])[iEndSegment]);
1129 // test compatibility between current segment and "extrapSegment"
1130 // 4 because 4 quantities in chi2
1132 NormalizedChi2WithSegment(extrapSegment,
1133 fMaxSigma2Distance)) <= 4.0) {
1134 // both segments compatible:
1135 // make track candidate from "begSegment" and "endSegment"
1136 if (fPrintLevel >= 2)
1137 cout << "TrackCandidate with Segment " << iEndSegment <<
1138 " in Station(0..) " << endStation << endl;
1139 // flag for both segments in one track:
1140 // to be done in track constructor ????
1141 BegSegment->SetInTrack(kTRUE);
1142 endSegment->SetInTrack(kTRUE);
1143 recTrack = new ((*fRecTracksPtr)[fNRecTracks])
1144 AliMUONTrack(BegSegment, endSegment, this);
1146 if (fPrintLevel >= 10) recTrack->RecursiveDump();
1147 // increment number of track candidates with 2 segments
1150 } // for (iEndSegment = 0;...
1151 delete extrapSegment; // should not delete HitForRec's it points to !!!!
1155 //__________________________________________________________________________
1156 Int_t AliMUONEventReconstructor::MakeTrackCandidatesWithOneSegmentAndOnePoint(AliMUONSegment *BegSegment)
1158 // To make track candidates with one segment and one point
1159 // in stations(1..) 4 and 5,
1160 // the segment being pointed to by "BegSegment".
1161 Int_t ch, ch1, ch2, endStation, iHit, iHitMax, iHitMin, nbCan1Seg1Hit;
1162 AliMUONHitForRec *extrapHitForRec, *hit;
1163 AliMUONTrack *recTrack;
1165 if (fPrintLevel >= 1)
1166 cout << "enter MakeTrackCandidatesWithOneSegmentAndOnePoint" << endl;
1167 // station for the end point
1168 endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
1169 // multiple scattering factor corresponding to one chamber
1170 mcsFactor = 0.0136 /
1171 GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
1172 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1173 // first and second chambers(0..) in the end station
1174 ch1 = 2 * endStation;
1176 // number of candidates to 0
1178 // Loop over chambers of the end station
1179 for (ch = ch2; ch >= ch1; ch--) {
1180 // linear extrapolation to chamber
1182 BegSegment->CreateHitForRecFromLinearExtrapToChamber(ch, mcsFactor);
1183 // limits for the hit index in the loop
1184 iHitMin = fIndexOfFirstHitForRecPerChamber[ch];
1185 iHitMax = iHitMin + fNHitsForRecPerChamber[ch];
1186 // Loop over HitForRec's in the chamber
1187 for (iHit = iHitMin; iHit < iHitMax; iHit++) {
1188 // pointer to HitForRec
1189 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1190 // test compatibility between current HitForRec and "extrapHitForRec"
1191 // 2 because 2 quantities in chi2
1193 NormalizedChi2WithHitForRec(extrapHitForRec,
1194 fMaxSigma2Distance)) <= 2.0) {
1195 // both HitForRec's compatible:
1196 // make track candidate from begSegment and current HitForRec
1197 if (fPrintLevel >= 2)
1198 cout << "TrackCandidate with HitForRec " << iHit <<
1199 " in Chamber(0..) " << ch << endl;
1200 // flag for beginning segments in one track:
1201 // to be done in track constructor ????
1202 BegSegment->SetInTrack(kTRUE);
1203 recTrack = new ((*fRecTracksPtr)[fNRecTracks])
1204 AliMUONTrack(BegSegment, hit, this);
1205 // the right place to eliminate "double counting" ???? how ????
1207 if (fPrintLevel >= 10) recTrack->RecursiveDump();
1208 // increment number of track candidates
1211 } // for (iHit = iHitMin;...
1212 delete extrapHitForRec;
1213 } // for (ch = ch2;...
1214 return nbCan1Seg1Hit;
1217 //__________________________________________________________________________
1218 void AliMUONEventReconstructor::MakeTrackCandidates(void)
1220 // To make track candidates
1221 // with at least 3 aligned points in stations(1..) 4 and 5
1222 // (two Segment's or one Segment and one HitForRec)
1223 Int_t begStation, iBegSegment, nbCan1Seg1Hit, nbCan2Seg;
1224 AliMUONSegment *begSegment;
1225 if (fPrintLevel >= 1) cout << "enter MakeTrackCandidates" << endl;
1226 // Loop over stations(1..) 5 and 4 for the beginning segment
1227 for (begStation = 4; begStation > 2; begStation--) {
1228 // Loop over segments in the beginning station
1229 for (iBegSegment = 0; iBegSegment < fNSegments[begStation]; iBegSegment++) {
1230 // pointer to segment
1231 begSegment = (AliMUONSegment*) ((*fSegmentsPtr[begStation])[iBegSegment]);
1232 if (fPrintLevel >= 2)
1233 cout << "look for TrackCandidate's with Segment " << iBegSegment <<
1234 " in Station(0..) " << begStation << endl;
1235 // Look for track candidates with two segments,
1236 // "begSegment" and all compatible segments in other station.
1237 // Only for beginning station(1..) 5
1238 // because candidates with 2 segments have to looked for only once.
1239 if (begStation == 4)
1240 nbCan2Seg = MakeTrackCandidatesWithTwoSegments(begSegment);
1241 // Look for track candidates with one segment and one point,
1242 // "begSegment" and all compatible HitForRec's in other station.
1243 // Only if "begSegment" does not belong already to a track candidate.
1244 // Is that a too strong condition ????
1245 if (!(begSegment->GetInTrack()))
1246 nbCan1Seg1Hit = MakeTrackCandidatesWithOneSegmentAndOnePoint(begSegment);
1247 } // for (iBegSegment = 0;...
1248 } // for (begStation = 4;...
1252 //__________________________________________________________________________
1253 void AliMUONEventReconstructor::FollowTracks(void)
1255 // Follow tracks in stations(1..) 3, 2 and 1
1256 // too long: should be made more modular !!!!
1257 AliMUONHitForRec *bestHit, *extrapHit, *extrapCorrHit, *hit;
1258 AliMUONSegment *bestSegment, *extrapSegment, *extrapCorrSegment, *segment;
1259 AliMUONTrack *track, *nextTrack;
1260 AliMUONTrackParam *trackParam1, trackParam[2], trackParamVertex;
1261 // -1 to avoid compilation warnings
1262 Int_t ch = -1, chInStation, chBestHit = -1, iHit, iSegment, station, trackIndex;
1263 Double_t bestChi2, chi2, dZ1, dZ2, dZ3, maxSigma2Distance, mcsFactor;
1264 Double_t bendingMomentum, chi2Norm = 0.;
1265 AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
1266 // local maxSigma2Distance, for easy increase in testing
1267 maxSigma2Distance = fMaxSigma2Distance;
1268 if (fPrintLevel >= 2)
1269 cout << "enter FollowTracks" << endl;
1270 // Loop over track candidates
1271 track = (AliMUONTrack*) fRecTracksPtr->First();
1274 // Follow function for each track candidate ????
1276 nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track); // prepare next track
1277 if (fPrintLevel >= 2)
1278 cout << "FollowTracks: track candidate(0..): " << trackIndex << endl;
1279 // Fit track candidate
1280 track->SetFitMCS(0); // without multiple Coulomb scattering
1281 track->SetFitNParam(3); // with 3 parameters (X = Y = 0)
1282 track->SetFitStart(0); // from parameters at vertex
1284 if (fPrintLevel >= 10) {
1285 cout << "FollowTracks: track candidate(0..): " << trackIndex
1286 << " after fit in stations(0..) 3 and 4" << endl;
1287 track->RecursiveDump();
1289 // Loop over stations(1..) 3, 2 and 1
1290 // something SPECIAL for stations 2 and 1 for majority 3 coincidence ????
1291 // otherwise: majority coincidence 2 !!!!
1292 for (station = 2; station >= 0; station--) {
1293 // Track parameters at first track hit (smallest Z)
1294 trackParam1 = ((AliMUONTrackHit*)
1295 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1296 // extrapolation to station
1297 trackParam1->ExtrapToStation(station, trackParam);
1298 extrapSegment = new AliMUONSegment(); // empty segment
1299 extrapCorrSegment = new AliMUONSegment(); // empty corrected segment
1300 // multiple scattering factor corresponding to one chamber
1301 // and momentum in bending plane (not total)
1302 mcsFactor = 0.0136 * trackParam1->GetInverseBendingMomentum();
1303 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1304 // Z difference from previous station
1305 dZ1 = (&(pMUON->Chamber(2 * station)))->Z() -
1306 (&(pMUON->Chamber(2 * station + 2)))->Z();
1307 // Z difference between the two previous stations
1308 dZ2 = (&(pMUON->Chamber(2 * station + 2)))->Z() -
1309 (&(pMUON->Chamber(2 * station + 4)))->Z();
1310 // Z difference between the two chambers in the previous station
1311 dZ3 = (&(pMUON->Chamber(2 * station)))->Z() -
1312 (&(pMUON->Chamber(2 * station + 1)))->Z();
1313 extrapSegment->SetBendingCoorReso2(fBendingResolution * fBendingResolution);
1315 SetNonBendingCoorReso2(fNonBendingResolution * fNonBendingResolution);
1316 extrapSegment->UpdateFromStationTrackParam
1317 (trackParam, mcsFactor, dZ1, dZ2, dZ3, station,
1318 trackParam1->GetInverseBendingMomentum());
1319 // same thing for corrected segment
1320 // better to use copy constructor, after checking that it works properly !!!!
1321 extrapCorrSegment->SetBendingCoorReso2(fBendingResolution * fBendingResolution);
1323 SetNonBendingCoorReso2(fNonBendingResolution * fNonBendingResolution);
1324 extrapCorrSegment->UpdateFromStationTrackParam
1325 (trackParam, mcsFactor, dZ1, dZ2, dZ3, station,
1326 trackParam1->GetInverseBendingMomentum());
1329 if (fPrintLevel >= 10) {
1330 cout << "FollowTracks: track candidate(0..): " << trackIndex
1331 << " Look for segment in station(0..): " << station << endl;
1333 // Loop over segments in station
1334 for (iSegment = 0; iSegment < fNSegments[station]; iSegment++) {
1335 // Look for best compatible Segment in station
1336 // should consider all possibilities ????
1337 // multiple scattering ????
1338 // separation in 2 functions: Segment and HitForRec ????
1339 segment = (AliMUONSegment*) ((*fSegmentsPtr[station])[iSegment]);
1340 // correction of corrected segment (fBendingCoor and fNonBendingCoor)
1341 // according to real Z value of "segment" and slopes of "extrapSegment"
1343 SetBendingCoor(extrapSegment->GetBendingCoor() +
1344 extrapSegment->GetBendingSlope() *
1345 (segment->GetHitForRec1()->GetZ() -
1346 (&(pMUON->Chamber(2 * station)))->Z()));
1348 SetNonBendingCoor(extrapSegment->GetNonBendingCoor() +
1349 extrapSegment->GetNonBendingSlope() *
1350 (segment->GetHitForRec1()->GetZ() -
1351 (&(pMUON->Chamber(2 * station)))->Z()));
1353 NormalizedChi2WithSegment(extrapCorrSegment, maxSigma2Distance);
1354 if (chi2 < bestChi2) {
1355 // update best Chi2 and Segment if better found
1356 bestSegment = segment;
1361 // best segment found: add it to track candidate
1362 track->AddSegment(bestSegment);
1363 // set track parameters at these two TrakHit's
1364 track->SetTrackParamAtHit(track->GetNTrackHits() - 2, &(trackParam[0]));
1365 track->SetTrackParamAtHit(track->GetNTrackHits() - 1, &(trackParam[1]));
1366 if (fPrintLevel >= 10) {
1367 cout << "FollowTracks: track candidate(0..): " << trackIndex
1368 << " Added segment in station(0..): " << station << endl;
1369 track->RecursiveDump();
1373 // No best segment found:
1374 // Look for best compatible HitForRec in station:
1375 // should consider all possibilities ????
1376 // multiple scattering ???? do about like for extrapSegment !!!!
1377 extrapHit = new AliMUONHitForRec(); // empty hit
1378 extrapCorrHit = new AliMUONHitForRec(); // empty corrected hit
1381 if (fPrintLevel >= 10) {
1382 cout << "FollowTracks: track candidate(0..): " << trackIndex
1383 << " Segment not found, look for hit in station(0..): " << station
1386 // Loop over chambers of the station
1387 for (chInStation = 0; chInStation < 2; chInStation++) {
1388 // coordinates of extrapolated hit
1390 SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor());
1392 SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor());
1393 // resolutions from "extrapSegment"
1394 extrapHit->SetBendingReso2(extrapSegment->GetBendingCoorReso2());
1395 extrapHit->SetNonBendingReso2(extrapSegment->GetNonBendingCoorReso2());
1396 // same things for corrected hit
1397 // better to use copy constructor, after checking that it works properly !!!!
1399 SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor());
1401 SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor());
1402 extrapHit->SetBendingReso2(extrapSegment->GetBendingCoorReso2());
1403 extrapHit->SetNonBendingReso2(extrapSegment->GetNonBendingCoorReso2());
1404 // Loop over hits in the chamber
1405 ch = 2 * station + chInStation;
1406 for (iHit = fIndexOfFirstHitForRecPerChamber[ch];
1407 iHit < fIndexOfFirstHitForRecPerChamber[ch] +
1408 fNHitsForRecPerChamber[ch];
1410 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1411 // correction of corrected hit (fBendingCoor and fNonBendingCoor)
1412 // according to real Z value of "hit" and slopes of right "trackParam"
1414 SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor() +
1415 (&(trackParam[chInStation]))->GetBendingSlope() *
1417 (&(trackParam[chInStation]))->GetZ()));
1419 SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor() +
1420 (&(trackParam[chInStation]))->GetNonBendingSlope() *
1422 (&(trackParam[chInStation]))->GetZ()));
1423 // condition for hit not already in segment ????
1424 chi2 = hit->NormalizedChi2WithHitForRec(extrapHit, maxSigma2Distance);
1425 if (chi2 < bestChi2) {
1426 // update best Chi2 and HitForRec if better found
1429 chBestHit = chInStation;
1434 // best hit found: add it to track candidate
1435 track->AddHitForRec(bestHit);
1436 // set track parameters at this TrackHit
1437 track->SetTrackParamAtHit(track->GetNTrackHits() - 1,
1438 &(trackParam[chBestHit]));
1439 if (fPrintLevel >= 10) {
1440 cout << "FollowTracks: track candidate(0..): " << trackIndex
1441 << " Added hit in station(0..): " << station << endl;
1442 track->RecursiveDump();
1446 // Remove current track candidate
1447 // and corresponding TrackHit's, ...
1449 delete extrapSegment;
1450 delete extrapCorrSegment;
1452 delete extrapCorrHit;
1453 break; // stop the search for this candidate:
1454 // exit from the loop over station
1457 delete extrapCorrHit;
1459 delete extrapSegment;
1460 delete extrapCorrSegment;
1461 // Sort track hits according to increasing Z
1462 track->GetTrackHitsPtr()->Sort();
1463 // Update track parameters at first track hit (smallest Z)
1464 trackParam1 = ((AliMUONTrackHit*)
1465 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1466 bendingMomentum = 0.;
1467 if (TMath::Abs(trackParam1->GetInverseBendingMomentum()) > 0.)
1468 bendingMomentum = TMath::Abs(1/(trackParam1->GetInverseBendingMomentum()));
1469 // Track removed if bendingMomentum not in window [min, max]
1470 if ((bendingMomentum < fMinBendingMomentum) || (bendingMomentum > fMaxBendingMomentum)) {
1472 break; // stop the search for this candidate:
1473 // exit from the loop over station
1476 // with multiple Coulomb scattering if all stations
1477 if (station == 0) track->SetFitMCS(1);
1478 // without multiple Coulomb scattering if not all stations
1479 else track->SetFitMCS(0);
1480 track->SetFitNParam(5); // with 5 parameters (momentum and position)
1481 track->SetFitStart(1); // from parameters at first hit
1483 Double_t numberOfDegFree = (2.0 * track->GetNTrackHits() - 5);
1484 if (numberOfDegFree > 0) {
1485 chi2Norm = track->GetFitFMin() / numberOfDegFree;
1489 // Track removed if normalized chi2 too high
1490 if (chi2Norm > fMaxChi2) {
1492 break; // stop the search for this candidate:
1493 // exit from the loop over station
1495 if (fPrintLevel >= 10) {
1496 cout << "FollowTracks: track candidate(0..): " << trackIndex
1497 << " after fit from station(0..): " << station << " to 4" << endl;
1498 track->RecursiveDump();
1500 // Track extrapolation to the vertex through the absorber (Branson)
1501 // after going through the first station
1503 trackParamVertex = *trackParam1;
1504 (&trackParamVertex)->ExtrapToVertex();
1505 track->SetTrackParamAtVertex(&trackParamVertex);
1506 if (fPrintLevel >= 1) {
1507 cout << "FollowTracks: track candidate(0..): " << trackIndex
1508 << " after extrapolation to vertex" << endl;
1509 track->RecursiveDump();
1512 } // for (station = 2;...
1513 // go really to next track
1516 // Compression of track array (necessary after Remove ????)
1517 fRecTracksPtr->Compress();
1521 //__________________________________________________________________________
1522 void AliMUONEventReconstructor::RemoveDoubleTracks(void)
1524 // To remove double tracks.
1525 // Tracks are considered identical
1526 // if they have at least half of their hits in common.
1527 // Among two identical tracks, one keeps the track with the larger number of hits
1528 // or, if these numbers are equal, the track with the minimum Chi2.
1529 AliMUONTrack *track1, *track2, *trackToRemove;
1530 Bool_t identicalTracks;
1531 Int_t hitsInCommon, nHits1, nHits2;
1532 identicalTracks = kTRUE;
1533 while (identicalTracks) {
1534 identicalTracks = kFALSE;
1535 // Loop over first track of the pair
1536 track1 = (AliMUONTrack*) fRecTracksPtr->First();
1537 while (track1 && (!identicalTracks)) {
1538 nHits1 = track1->GetNTrackHits();
1539 // Loop over second track of the pair
1540 track2 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1541 while (track2 && (!identicalTracks)) {
1542 nHits2 = track2->GetNTrackHits();
1543 // number of hits in common between two tracks
1544 hitsInCommon = track1->HitsInCommon(track2);
1545 // check for identical tracks
1546 if ((4 * hitsInCommon) >= (nHits1 + nHits2)) {
1547 identicalTracks = kTRUE;
1548 // decide which track to remove
1549 if (nHits1 > nHits2) trackToRemove = track2;
1550 else if (nHits1 < nHits2) trackToRemove = track1;
1551 else if ((track1->GetFitFMin()) < (track2->GetFitFMin()))
1552 trackToRemove = track2;
1553 else trackToRemove = track1;
1555 trackToRemove->Remove();
1557 track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
1559 track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1565 //__________________________________________________________________________
1566 void AliMUONEventReconstructor::UpdateTrackParamAtHit()
1568 // Set track parameters after track fitting. Fill fTrackParamAtHit of AliMUONTrack's
1569 AliMUONTrack *track;
1570 AliMUONTrackHit *trackHit;
1571 AliMUONTrackParam *trackParam;
1572 track = (AliMUONTrack*) fRecTracksPtr->First();
1574 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->First();
1576 trackParam = trackHit->GetTrackParam();
1577 track->AddTrackParamAtHit(trackParam);
1578 trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->After(trackHit);
1580 track = (AliMUONTrack*) fRecTracksPtr->After(track);
1585 //__________________________________________________________________________
1586 void AliMUONEventReconstructor::EventDump(void)
1588 // Dump reconstructed event (track parameters at vertex and at first hit),
1589 // and the particle parameters
1591 AliMUONTrack *track;
1592 AliMUONTrackParam *trackParam, *trackParam1;
1593 Double_t bendingSlope, nonBendingSlope, pYZ;
1594 Double_t pX, pY, pZ, x, y, z, c;
1595 Int_t np, trackIndex, nTrackHits;
1597 if (fPrintLevel >= 1) cout << "****** enter EventDump ******" << endl;
1598 if (fPrintLevel >= 1) {
1599 cout << " Number of Reconstructed tracks :" << fNRecTracks << endl;
1601 fRecTracksPtr->Compress(); // for simple loop without "Next" since no hole
1602 // Loop over reconstructed tracks
1603 for (trackIndex = 0; trackIndex < fNRecTracks; trackIndex++) {
1604 if (fTrackMethod != 1) continue; //AZ - skip the rest for now
1605 if (fPrintLevel >= 1)
1606 cout << " track number: " << trackIndex << endl;
1607 // function for each track for modularity ????
1608 track = (AliMUONTrack*) ((*fRecTracksPtr)[trackIndex]);
1609 nTrackHits = track->GetNTrackHits();
1610 if (fPrintLevel >= 1)
1611 cout << " number of track hits: " << nTrackHits << endl;
1612 // track parameters at Vertex
1613 trackParam = track->GetTrackParamAtVertex();
1614 x = trackParam->GetNonBendingCoor();
1615 y = trackParam->GetBendingCoor();
1616 z = trackParam->GetZ();
1617 bendingSlope = trackParam->GetBendingSlope();
1618 nonBendingSlope = trackParam->GetNonBendingSlope();
1619 pYZ = 1/TMath::Abs(trackParam->GetInverseBendingMomentum());
1620 pZ = pYZ/TMath::Sqrt(1+bendingSlope*bendingSlope);
1621 pX = pZ * nonBendingSlope;
1622 pY = pZ * bendingSlope;
1623 c = TMath::Sign(1.0, trackParam->GetInverseBendingMomentum());
1624 if (fPrintLevel >= 1)
1625 printf(" track parameters at Vertex z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1626 z, x, y, pX, pY, pZ, c);
1628 // track parameters at first hit
1629 trackParam1 = ((AliMUONTrackHit*)
1630 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1631 x = trackParam1->GetNonBendingCoor();
1632 y = trackParam1->GetBendingCoor();
1633 z = trackParam1->GetZ();
1634 bendingSlope = trackParam1->GetBendingSlope();
1635 nonBendingSlope = trackParam1->GetNonBendingSlope();
1636 pYZ = 1/TMath::Abs(trackParam1->GetInverseBendingMomentum());
1637 pZ = pYZ/TMath::Sqrt(1.0 + bendingSlope * bendingSlope);
1638 pX = pZ * nonBendingSlope;
1639 pY = pZ * bendingSlope;
1640 c = TMath::Sign(1.0, trackParam1->GetInverseBendingMomentum());
1641 if (fPrintLevel >= 1)
1642 printf(" track parameters at z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1643 z, x, y, pX, pY, pZ, c);
1645 // informations about generated particles
1646 np = gAlice->GetMCApp()->GetNtrack();
1647 printf(" **** number of generated particles: %d \n", np);
1649 // for (Int_t iPart = 0; iPart < np; iPart++) {
1650 // p = gAlice->Particle(iPart);
1651 // printf(" particle %d: type= %d px= %f py= %f pz= %f pdg= %d\n",
1652 // iPart, p->GetPdgCode(), p->Px(), p->Py(), p->Pz(), p->GetPdgCode());
1658 //__________________________________________________________________________
1659 void AliMUONEventReconstructor::EventDumpTrigger(void)
1661 // Dump reconstructed trigger event
1662 // and the particle parameters
1664 AliMUONTriggerTrack *triggertrack ;
1665 Int_t nTriggerTracks = fMUONData->RecTriggerTracks()->GetEntriesFast();
1667 if (fPrintLevel >= 1) {
1668 cout << "****** enter EventDumpTrigger ******" << endl;
1669 cout << " Number of Reconstructed tracks :" << nTriggerTracks << endl;
1671 // Loop over reconstructed tracks
1672 for (Int_t trackIndex = 0; trackIndex < nTriggerTracks; trackIndex++) {
1673 triggertrack = (AliMUONTriggerTrack*)fMUONData->RecTriggerTracks()->At(trackIndex);
1674 printf(" trigger track number %i x11=%f y11=%f thetax=%f thetay=%f \n",
1676 triggertrack->GetX11(),triggertrack->GetY11(),
1677 triggertrack->GetThetax(),triggertrack->GetThetay());
1681 //__________________________________________________________________________
1682 void AliMUONEventReconstructor::FillEvent()
1684 // Create a new AliMUONRecoEvent, fill its track list, then add it as a
1685 // leaf in the Event branch of TreeRecoEvent tree
1686 cout << "Enter FillEvent() ...\n";
1689 fRecoEvent = new AliMUONRecoEvent();
1691 fRecoEvent->Clear();
1693 //save current directory
1694 TDirectory *current = gDirectory;
1695 if (!fTreeFile) fTreeFile = new TFile("tree_reco.root", "RECREATE");
1696 if (!fEventTree) fEventTree = new TTree("TreeRecoEvent", "MUON reconstructed events");
1697 //AZif (fRecoEvent->MakeDumpTracks(fRecTracksPtr)) {
1698 if (fRecoEvent->MakeDumpTracks(fMuons, fRecTracksPtr, this)) { //AZ
1699 if (fPrintLevel > 1) fRecoEvent->EventInfo();
1700 TBranch *branch = fEventTree->GetBranch("Event");
1701 if (!branch) branch = fEventTree->Branch("Event", "AliMUONRecoEvent", &fRecoEvent, 64000);
1702 branch->SetAutoDelete();
1707 // restore directory
1711 //__________________________________________________________________________
1712 void AliMUONEventReconstructor::MakeTrackCandidatesK(void)
1714 // To make initial tracks for Kalman filter from the list of segments
1716 AliMUONSegment *segment;
1717 AliMUONTrackK *trackK;
1719 if (fPrintLevel >= 1) cout << "enter MakeTrackCandidatesK" << endl;
1720 // Reset the TClonesArray of reconstructed tracks
1721 if (fRecTracksPtr) fRecTracksPtr->Delete();
1722 // Delete in order that the Track destructors are called,
1723 // hence the space for the TClonesArray of pointers to TrackHit's is freed
1726 AliMUONTrackK a(this, fHitsForRecPtr); // bad idea ???
1727 // Loop over stations(1...) 5 and 4
1728 for (istat=4; istat>=3; istat--) {
1729 // Loop over segments in the station
1730 for (iseg=0; iseg<fNSegments[istat]; iseg++) {
1731 // Transform segments to tracks and evaluate covariance matrix
1732 segment = (AliMUONSegment*) ((*fSegmentsPtr[istat])[iseg]);
1733 trackK = new ((*fRecTracksPtr)[fNRecTracks]) AliMUONTrackK(segment);
1735 } // for (iseg=0;...)
1736 } // for (istat=4;...)
1740 //__________________________________________________________________________
1741 void AliMUONEventReconstructor::FollowTracksK(void)
1743 // Follow tracks using Kalman filter
1745 Int_t icand, ichamBeg, ichamEnd, chamBits;
1746 Double_t zDipole1, zDipole2;
1747 AliMUONTrackK *trackK;
1748 AliMUONHitForRec *hit;
1749 AliMUONRawCluster *clus;
1750 TClonesArray *rawclusters;
1752 clus = 0; rawclusters = 0;
1754 zDipole1 = GetSimpleBPosition() - GetSimpleBLength()/2;
1755 zDipole2 = zDipole1 + GetSimpleBLength();
1758 pMUON = (AliMUON*) gAlice->GetModule("MUON");
1759 for (Int_t i1=0; i1<fNHitsForRec; i1++) {
1760 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[i1]);
1761 //if (hit->GetTHTrack() > 1 || hit->GetGeantSignal() == 0) continue;
1763 cout << " Hit #" << hit->GetChamberNumber() << " ";
1764 cout << hit->GetBendingCoor() << " ";
1765 cout << hit->GetNonBendingCoor() << " ";
1766 cout << hit->GetZ() << " ";
1767 cout << hit->GetGeantSignal() << " ";
1768 cout << hit->GetTHTrack() << endl;
1771 printf(" Hit # %d %10.4f %10.4f %10.4f",
1772 hit->GetChamberNumber(), hit->GetBendingCoor(),
1773 hit->GetNonBendingCoor(), hit->GetZ());
1774 if (fRecGeantHits) {
1776 printf(" %3d %3d \n", hit->GetGeantSignal(), hit->GetTHTrack());
1778 // from raw clusters
1779 rawclusters = pMUON->RawClustAddress(hit->GetChamberNumber());
1780 clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(hit->
1782 printf("%3d", clus->fTracks[1]-1);
1783 if (clus->fTracks[2] != 0) printf("%3d \n", clus->fTracks[2]-1);
1790 Int_t nSeeds = fNRecTracks; // starting number of seeds
1791 // Loop over track candidates
1792 while (icand < fNRecTracks-1) {
1794 trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]);
1796 // Discard candidate which will produce the double track
1798 ok = CheckCandidateK(icand,nSeeds);
1800 //trackK->SetRecover(-1); // mark candidate to be removed
1806 if (trackK->GetRecover() == 0) hit = (AliMUONHitForRec*)
1807 trackK->GetHitOnTrack()->Last(); // last hit
1808 else hit = (AliMUONHitForRec*) (*trackK->GetHitOnTrack())[1]; // 2'nd hit
1809 ichamBeg = hit->GetChamberNumber();
1811 // Check propagation direction
1812 if (trackK->GetTrackDir() > 0) {
1813 ichamEnd = 9; // forward propagation
1814 ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kFALSE,zDipole1,zDipole2);
1816 ichamBeg = ichamEnd;
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;
1825 if (trackK->GetBPFlag()) {
1827 ichamEnd = 6; // backward propagation
1828 // Change weight matrix and zero fChi2 for backpropagation
1829 trackK->StartBack();
1830 ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kTRUE,zDipole1,zDipole2);
1831 ichamBeg = ichamEnd;
1837 trackK->SetTrackDir(-1);
1838 trackK->SetBPFlag(kFALSE);
1839 ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kFALSE,zDipole1,zDipole2);
1841 if (ok) trackK->SetTrackQuality(0); // compute "track quality"
1842 else trackK->SetRecover(-1); // mark candidate to be removed
1844 // Majority 3 of 4 in first 2 stations
1846 for (Int_t i=0; i<trackK->GetNTrackHits(); i++) {
1847 hit = (AliMUONHitForRec*) (*trackK->GetHitOnTrack())[i];
1848 chamBits |= BIT(hit->GetChamberNumber()-1);
1850 //if (!((chamBits&3)==3 || (chamBits>>2&3)==3)) trackK->SetRecover(-1);
1851 //mark candidate to be removed
1854 for (Int_t i=0; i<fNRecTracks; i++) {
1855 trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[i]);
1856 if (trackK->GetRecover() < 0) fRecTracksPtr->RemoveAt(i);
1859 // Compress TClonesArray
1860 fRecTracksPtr->Compress();
1861 fNRecTracks = fRecTracksPtr->GetEntriesFast();
1865 //__________________________________________________________________________
1866 Bool_t AliMUONEventReconstructor::CheckCandidateK(Int_t icand, Int_t nSeeds) const
1868 // Discards track candidate if it will produce the double track (having
1869 // the same seed segment hits as hits of a good track found before)
1870 AliMUONTrackK *track1, *track2;
1871 AliMUONHitForRec *hit1, *hit2, *hit;
1873 track1 = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]);
1874 hit1 = (AliMUONHitForRec*) (*track1->GetHitOnTrack())[0]; // 1'st hit
1875 hit2 = (AliMUONHitForRec*) (*track1->GetHitOnTrack())[1]; // 2'nd hit
1877 for (Int_t i=0; i<icand; i++) {
1878 track2 = (AliMUONTrackK*) ((*fRecTracksPtr)[i]);
1879 //if (track2->GetRecover() < 0) continue;
1880 if (track2->GetRecover() < 0 && icand >= nSeeds) continue;
1882 if (track1->GetStartSegment() == track2->GetStartSegment()) {
1886 for (Int_t j=0; j<track2->GetNTrackHits(); j++) {
1887 hit = (AliMUONHitForRec*) (*track2->GetHitOnTrack())[j];
1888 if (hit == hit1 || hit == hit2) {
1890 if (nSame == 2) return kFALSE;
1892 } // for (Int_t j=0;
1894 } // for (Int_t i=0;
1898 //__________________________________________________________________________
1899 void AliMUONEventReconstructor::RemoveDoubleTracksK(void)
1901 // Removes double tracks (sharing more than half of their hits). Keeps
1902 // the track with higher quality
1903 AliMUONTrackK *track1, *track2, *trackToKill;
1905 // Sort tracks according to their quality
1906 fRecTracksPtr->Sort();
1908 // Loop over first track of the pair
1909 track1 = (AliMUONTrackK*) fRecTracksPtr->First();
1911 // Loop over second track of the pair
1912 track2 = (AliMUONTrackK*) fRecTracksPtr->After(track1);
1914 // Check whether or not to keep track2
1915 if (!track2->KeepTrack(track1)) {
1916 cout << " Killed track: " << 1/(*track2->GetTrackParameters())(4,0) <<
1917 " " << track2->GetTrackQuality() << endl;
1918 trackToKill = track2;
1919 track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2);
1920 trackToKill->Kill();
1921 fRecTracksPtr->Compress();
1922 } else track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2);
1924 track1 = (AliMUONTrackK*) fRecTracksPtr->After(track1);
1927 fNRecTracks = fRecTracksPtr->GetEntriesFast();
1928 cout << " Number of Kalman tracks: " << fNRecTracks << endl;
1931 //__________________________________________________________________________
1932 void AliMUONEventReconstructor::GoToVertex(void)
1934 // Propagates track to the vertex thru absorber
1935 // (using Branson correction for now)
1939 for (Int_t i=0; i<fNRecTracks; i++) {
1940 //((AliMUONTrackK*)(*fRecTracksPtr)[i])->Branson();
1941 //((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToZ(zVertex); // w/out absorber
1942 ((AliMUONTrackK*)(*fRecTracksPtr)[i])->SetTrackQuality(1); // compute Chi2
1943 ((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToVertex(); // with absorber