]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/AliMUONEventReconstructor.cxx
Corrections to violations of coding conventions
[u/mrichter/AliRoot.git] / MUON / AliMUONEventReconstructor.cxx
CommitLineData
a9e2aefa 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
16/*
17$Log$
3831f268 18Revision 1.24 2001/03/30 09:37:58 gosset
19Initialisations of pointers... for GEANT background events in the constructor
20
47a6b602 21Revision 1.23 2001/01/26 21:44:45 morsch
22Use access functions to AliMUONDigit, ... member data.
23
0a5f36d9 24Revision 1.22 2001/01/26 20:00:53 hristov
25Major upgrade of AliRoot code
3f5cf0b3 26Revision 1.20 2001/01/08 11:01:02 gosset
27Modifications used for addendum to Dimuon TDR (JP Cussonneau):
28*. MaxBendingMomentum to make both a segment and a track (default 500)
29*. MaxChi2 per degree of freedom to make a track (default 100)
30*. MinBendingMomentum used also to make a track
31 and not only a segment (default 3)
32*. wider roads for track search in stations 1 to 3
33*. extrapolation to actual Z instead of Z(chamber) in FollowTracks
34*. in track fit:
35 - limits on parameters X and Y (+/-500)
36 - covariance matrices in double precision
37 - normalization of covariance matrices before inversion
38 - suppression of Minuit printouts
39*. correction against memory leak (delete extrapHit) in FollowTracks
40*. RMax to 10 degrees with Z(chamber) instead of fixed values;
41 RMin and Rmax cuts suppressed in NewHitForRecFromGEANT,
42 because useless with realistic geometry
43
d0bfce8d 44Revision 1.19 2000/11/20 11:24:10 gosset
45New package for reconstructed tracks (A. Gheata):
46* tree output in the file "tree_reco.root"
47* display events and make histograms from this file
48
c7ba256d 49Revision 1.18 2000/10/26 12:47:03 gosset
50Real distance between chambers of each station taken into account
51for the reconstruction parameters "fSegmentMaxDistBending[5]"
52
860cef81 53Revision 1.17 2000/10/24 09:26:20 gosset
54Comments updated
55
ba5b68db 56Revision 1.16 2000/10/24 09:22:35 gosset
57Method AddHitsForRecFromRawClusters: real Z of raw cluster and not Z of chamber
58
cfe0107f 59Revision 1.15 2000/10/12 15:17:30 gosset
60Sign of fSimpleBValue corrected: sign ox Bx and not Bz (thanks to Galina)
61
065bd0e1 62Revision 1.14 2000/10/04 18:21:26 morsch
63Include stdlib.h
64
e7d455d5 65Revision 1.13 2000/10/02 21:28:09 fca
66Removal of useless dependecies via forward declarations
67
94de3818 68Revision 1.12 2000/10/02 16:58:29 egangler
69Cleaning of the code :
70-> coding conventions
71-> void Streamers
72-> some useless includes removed or replaced by "class" statement
73
ecfa008b 74Revision 1.11 2000/09/22 09:16:33 hristov
75Casting needed on DEC
76
8832c7b4 77Revision 1.10 2000/09/19 09:49:50 gosset
78AliMUONEventReconstructor package
79* track extrapolation independent from reco_muon.F, use of AliMagF...
80* possibility to use new magnetic field (automatic from generated root file)
81
a6f03ddb 82Revision 1.9 2000/07/20 12:45:27 gosset
83New "EventReconstructor..." structure,
84 hopefully more adapted to tree/streamer.
85"AliMUONEventReconstructor::RemoveDoubleTracks"
86 to keep only one track among similar ones.
87
8429a5e4 88Revision 1.8 2000/07/18 16:04:06 gosset
89AliMUONEventReconstructor package:
90* a few minor modifications and more comments
91* a few corrections
92 * right sign for Z of raw clusters
93 * right loop over chambers inside station
94 * symmetrized covariance matrix for measurements (TrackChi2MCS)
95 * right sign of charge in extrapolation (ExtrapToZ)
96 * right zEndAbsorber for Branson correction below 3 degrees
97* use of TVirtualFitter instead of TMinuit for AliMUONTrack::Fit
98* no parameter for AliMUONTrack::Fit() but more fit parameters in Track object
99
956019b6 100Revision 1.7 2000/07/03 12:28:06 gosset
101Printout at the right place after extrapolation to vertex
102
6ae236ee 103Revision 1.6 2000/06/30 12:01:06 gosset
104Correction for hit search in the right chamber (JPC)
105
d82671a0 106Revision 1.5 2000/06/30 10:15:48 gosset
107Changes to EventReconstructor...:
108precision fit with multiple Coulomb scattering;
109extrapolation to vertex with Branson correction in absorber (JPC)
110
04b5ea16 111Revision 1.4 2000/06/27 14:11:36 gosset
112Corrections against violations of coding conventions
113
9cbdf048 114Revision 1.3 2000/06/16 07:27:08 gosset
115To remove problem in running RuleChecker, like in MUON-dev
116
79e1e601 117Revision 1.1.2.5 2000/06/16 07:00:26 gosset
118To remove problem in running RuleChecker
119
a9e2aefa 120Revision 1.1.2.4 2000/06/12 08:00:07 morsch
121Dummy streamer to solve CINT compilation problem (to be investigated !)
122
123Revision 1.1.2.3 2000/06/09 20:59:57 morsch
124Make includes consistent with new file structure.
125
126Revision 1.1.2.2 2000/06/09 12:58:05 gosset
127Removed comment beginnings in Log sections of .cxx files
128Suppressed most violations of coding rules
129
130Revision 1.1.2.1 2000/06/07 14:44:53 gosset
131Addition of files for track reconstruction in C++
132*/
133
3831f268 134////////////////////////////////////
a9e2aefa 135//
136// MUON event reconstructor in ALICE
137//
138// This class contains as data:
139// * the parameters for the event reconstruction
140// * a pointer to the array of hits to be reconstructed (the event)
141// * a pointer to the array of segments made with these hits inside each station
142// * a pointer to the array of reconstructed tracks
143//
144// It contains as methods, among others:
145// * MakeEventToBeReconstructed to build the array of hits to be reconstructed
146// * MakeSegments to build the segments
147// * MakeTracks to build the tracks
3831f268 148//
149////////////////////////////////////
a9e2aefa 150
3831f268 151#include <iostream.h> // for cout
a9e2aefa 152
94de3818 153#include <TTree.h>
a9e2aefa 154
a9e2aefa 155#include "AliMUON.h"
3831f268 156#include "AliMUONChamber.h"
157#include "AliMUONEventReconstructor.h"
a9e2aefa 158#include "AliMUONHitForRec.h"
a9e2aefa 159#include "AliMUONRawCluster.h"
3831f268 160#include "AliMUONRecoEvent.h"
161#include "AliMUONSegment.h"
a9e2aefa 162#include "AliMUONTrack.h"
a9e2aefa 163#include "AliMUONTrackHit.h"
94de3818 164#include "AliMagF.h"
3831f268 165#include "AliRun.h" // for gAlice
a9e2aefa 166
a9e2aefa 167//************* Defaults parameters for reconstruction
9cbdf048 168static const Double_t kDefaultMinBendingMomentum = 3.0;
d0bfce8d 169static const Double_t kDefaultMaxBendingMomentum = 500.0;
170static const Double_t kDefaultMaxChi2 = 100.0;
9cbdf048 171static const Double_t kDefaultMaxSigma2Distance = 16.0;
172static const Double_t kDefaultBendingResolution = 0.01;
173static const Double_t kDefaultNonBendingResolution = 0.144;
174static const Double_t kDefaultChamberThicknessInX0 = 0.03;
a9e2aefa 175// Simple magnetic field:
176// Value taken from macro MUONtracking.C: 0.7 T, hence 7 kG
177// Length and Position from reco_muon.F, with opposite sign:
178// Length = ZMAGEND-ZCOIL
179// Position = (ZMAGEND+ZCOIL)/2
180// to be ajusted differently from real magnetic field ????
9cbdf048 181static const Double_t kDefaultSimpleBValue = 7.0;
182static const Double_t kDefaultSimpleBLength = 428.0;
183static const Double_t kDefaultSimpleBPosition = 1019.0;
184static const Int_t kDefaultRecGeantHits = 0;
185static const Double_t kDefaultEfficiency = 0.95;
a9e2aefa 186
9cbdf048 187static const Int_t kDefaultPrintLevel = 0;
a9e2aefa 188
189ClassImp(AliMUONEventReconstructor) // Class implementation in ROOT context
190
191 //__________________________________________________________________________
192AliMUONEventReconstructor::AliMUONEventReconstructor(void)
193{
194 // Constructor for class AliMUONEventReconstructor
195 SetReconstructionParametersToDefaults();
196 // Memory allocation for the TClonesArray of hits for reconstruction
197 // Is 10000 the right size ????
198 fHitsForRecPtr = new TClonesArray("AliMUONHitForRec", 10000);
199 fNHitsForRec = 0; // really needed or GetEntriesFast sufficient ????
200 // Memory allocation for the TClonesArray's of segments in stations
201 // Is 2000 the right size ????
9cbdf048 202 for (Int_t st = 0; st < kMaxMuonTrackingStations; st++) {
a9e2aefa 203 fSegmentsPtr[st] = new TClonesArray("AliMUONSegment", 2000);
204 fNSegments[st] = 0; // really needed or GetEntriesFast sufficient ????
205 }
206 // Memory allocation for the TClonesArray of reconstructed tracks
207 // Is 10 the right size ????
208 fRecTracksPtr = new TClonesArray("AliMUONTrack", 10);
209 fNRecTracks = 0; // really needed or GetEntriesFast sufficient ????
8429a5e4 210 // Memory allocation for the TClonesArray of hits on reconstructed tracks
211 // Is 100 the right size ????
212 fRecTrackHitsPtr = new TClonesArray("AliMUONTrack", 100);
213 fNRecTrackHits = 0; // really needed or GetEntriesFast sufficient ????
a9e2aefa 214
065bd0e1 215 // Sign of fSimpleBValue according to sign of Bx value at (50,50,950).
a6f03ddb 216 Float_t b[3], x[3];
217 x[0] = 50.; x[1] = 50.; x[2] = 950.;
218 gAlice->Field()->Field(x, b);
065bd0e1 219 fSimpleBValue = TMath::Sign(fSimpleBValue,(Double_t) b[0]);
a6f03ddb 220 // See how to get fSimple(BValue, BLength, BPosition)
221 // automatically calculated from the actual magnetic field ????
a9e2aefa 222
223 if (fPrintLevel >= 0) {
a6f03ddb 224 cout << "AliMUONEventReconstructor constructed with defaults" << endl; Dump();
225 cout << endl << "Magnetic field from root file:" << endl;
226 gAlice->Field()->Dump();
227 cout << endl;
228 }
c7ba256d 229
47a6b602 230 // Initializions for GEANT background events
231 fBkgGeantFile = 0;
232 fBkgGeantTK = 0;
233 fBkgGeantParticles = 0;
234 fBkgGeantTH = 0;
235 fBkgGeantHits = 0;
236 fBkgGeantEventNumber = -1;
237
c7ba256d 238 // Initialize to 0 pointers to RecoEvent, tree and tree file
239 fRecoEvent = 0;
240 fEventTree = 0;
241 fTreeFile = 0;
242
a9e2aefa 243 return;
244}
245
9cbdf048 246AliMUONEventReconstructor::AliMUONEventReconstructor (const AliMUONEventReconstructor& Reconstructor)
247{
248 // Dummy copy constructor
249}
250
251AliMUONEventReconstructor & AliMUONEventReconstructor::operator=(const AliMUONEventReconstructor& Reconstructor)
252{
253 // Dummy assignment operator
254 return *this;
255}
256
a9e2aefa 257 //__________________________________________________________________________
258AliMUONEventReconstructor::~AliMUONEventReconstructor(void)
259{
260 // Destructor for class AliMUONEventReconstructor
c7ba256d 261 if (fTreeFile) {
262 fTreeFile->Close();
263 delete fTreeFile;
264 }
265// if (fEventTree) delete fEventTree;
266 if (fRecoEvent) delete fRecoEvent;
a9e2aefa 267 delete fHitsForRecPtr; // Correct destruction of everything ???? or delete [] ????
9cbdf048 268 for (Int_t st = 0; st < kMaxMuonTrackingStations; st++)
a9e2aefa 269 delete fSegmentsPtr[st]; // Correct destruction of everything ????
270 return;
271}
272
273 //__________________________________________________________________________
274void AliMUONEventReconstructor::SetReconstructionParametersToDefaults(void)
275{
276 // Set reconstruction parameters to default values
277 // Would be much more convenient with a structure (or class) ????
9cbdf048 278 fMinBendingMomentum = kDefaultMinBendingMomentum;
d0bfce8d 279 fMaxBendingMomentum = kDefaultMaxBendingMomentum;
280 fMaxChi2 = kDefaultMaxChi2;
9cbdf048 281 fMaxSigma2Distance = kDefaultMaxSigma2Distance;
a9e2aefa 282
9cbdf048 283 AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON");
a9e2aefa 284 // ******** Parameters for making HitsForRec
285 // minimum radius,
286 // like in TRACKF_STAT:
287 // 2 degrees for stations 1 and 2, or ch(0...) from 0 to 3;
288 // 30 cm for stations 3 to 5, or ch(0...) from 4 to 9
9cbdf048 289 for (Int_t ch = 0; ch < kMaxMuonTrackingChambers; ch++) {
290 if (ch < 4) fRMin[ch] = TMath::Abs((&(pMUON->Chamber(ch)))->Z()) *
a9e2aefa 291 2.0 * TMath::Pi() / 180.0;
292 else fRMin[ch] = 30.0;
d0bfce8d 293 // maximum radius at 10 degrees and Z of chamber
294 fRMax[ch] = TMath::Abs((&(pMUON->Chamber(ch)))->Z()) *
295 10.0 * TMath::Pi() / 180.0;
a9e2aefa 296 }
a9e2aefa 297
298 // ******** Parameters for making segments
299 // should be parametrized ????
300 // according to interval between chambers in a station ????
301 // Maximum distance in non bending plane
302 // 5 * 0.22 just to remember the way it was made in TRACKF_STAT
303 // SIGCUT*DYMAX(IZ)
9cbdf048 304 for (Int_t st = 0; st < kMaxMuonTrackingStations; st++)
a9e2aefa 305 fSegmentMaxDistNonBending[st] = 5. * 0.22;
860cef81 306 // Maximum distance in bending plane:
307 // values from TRACKF_STAT, corresponding to (J psi 20cm),
308 // scaled to the real distance between chambers in a station
309 fSegmentMaxDistBending[0] = 1.5 *
310 ((&(pMUON->Chamber(1)))->Z() - (&(pMUON->Chamber(0)))->Z()) / 20.0;
311 fSegmentMaxDistBending[1] = 1.5 *
312 ((&(pMUON->Chamber(3)))->Z() - (&(pMUON->Chamber(2)))->Z()) / 20.0;
313 fSegmentMaxDistBending[2] = 3.0 *
314 ((&(pMUON->Chamber(5)))->Z() - (&(pMUON->Chamber(4)))->Z()) / 20.0;
315 fSegmentMaxDistBending[3] = 6.0 *
316 ((&(pMUON->Chamber(7)))->Z() - (&(pMUON->Chamber(6)))->Z()) / 20.0;
317 fSegmentMaxDistBending[4] = 6.0 *
318 ((&(pMUON->Chamber(9)))->Z() - (&(pMUON->Chamber(8)))->Z()) / 20.0;
a9e2aefa 319
9cbdf048 320 fBendingResolution = kDefaultBendingResolution;
321 fNonBendingResolution = kDefaultNonBendingResolution;
322 fChamberThicknessInX0 = kDefaultChamberThicknessInX0;
323 fSimpleBValue = kDefaultSimpleBValue;
324 fSimpleBLength = kDefaultSimpleBLength;
325 fSimpleBPosition = kDefaultSimpleBPosition;
326 fRecGeantHits = kDefaultRecGeantHits;
327 fEfficiency = kDefaultEfficiency;
328 fPrintLevel = kDefaultPrintLevel;
a9e2aefa 329 return;
330}
331
332//__________________________________________________________________________
3831f268 333Double_t AliMUONEventReconstructor::GetImpactParamFromBendingMomentum(Double_t BendingMomentum) const
a9e2aefa 334{
335 // Returns impact parameter at vertex in bending plane (cm),
336 // from the signed bending momentum "BendingMomentum" in bending plane (GeV/c),
337 // using simple values for dipole magnetic field.
065bd0e1 338 // The sign of "BendingMomentum" is the sign of the charge.
a9e2aefa 339 return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
340 BendingMomentum);
341}
342
343//__________________________________________________________________________
3831f268 344Double_t AliMUONEventReconstructor::GetBendingMomentumFromImpactParam(Double_t ImpactParam) const
a9e2aefa 345{
346 // Returns signed bending momentum in bending plane (GeV/c),
065bd0e1 347 // the sign being the sign of the charge for particles moving forward in Z,
a9e2aefa 348 // from the impact parameter "ImpactParam" at vertex in bending plane (cm),
349 // using simple values for dipole magnetic field.
a9e2aefa 350 return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition /
351 ImpactParam);
352}
353
354//__________________________________________________________________________
355void AliMUONEventReconstructor::SetBkgGeantFile(Text_t *BkgGeantFileName)
356{
357 // Set background file ... for GEANT hits
358 // Must be called after having loaded the firts signal event
359 if (fPrintLevel >= 0) {
360 cout << "Enter SetBkgGeantFile with BkgGeantFileName ``"
361 << BkgGeantFileName << "''" << endl;}
362 if (strlen(BkgGeantFileName)) {
363 // BkgGeantFileName not empty: try to open the file
364 if (fPrintLevel >= 2) {cout << "Before File(Bkg)" << endl; gDirectory->Dump();}
365 fBkgGeantFile = new TFile(BkgGeantFileName);
366 if (fPrintLevel >= 2) {cout << "After File(Bkg)" << endl; gDirectory->Dump();}
367 if (fBkgGeantFile-> IsOpen()) {
368 if (fPrintLevel >= 0) {
369 cout << "Background for GEANT hits in file: ``" << BkgGeantFileName
370 << "'' successfully opened" << endl;}
371 }
372 else {
373 cout << "Background for GEANT hits in file: " << BkgGeantFileName << endl;
374 cout << "NOT FOUND: EXIT" << endl;
375 exit(0); // right instruction for exit ????
376 }
377 // Arrays for "particles" and "hits"
378 fBkgGeantParticles = new TClonesArray("TParticle", 200);
379 fBkgGeantHits = new TClonesArray("AliMUONHit", 2000);
380 // Event number to -1 for initialization
381 fBkgGeantEventNumber = -1;
382 // Back to the signal file:
383 // first signal event must have been loaded previously,
384 // otherwise, Segmentation violation at the next instruction
385 // How is it possible to do smething better ????
386 ((gAlice->TreeK())->GetCurrentFile())->cd();
387 if (fPrintLevel >= 2) {cout << "After cd(gAlice)" << endl; gDirectory->Dump();}
388 }
389 return;
390}
391
392//__________________________________________________________________________
393void AliMUONEventReconstructor::NextBkgGeantEvent(void)
394{
395 // Get next event in background file for GEANT hits
396 // Goes back to event number 0 when end of file is reached
397 char treeName[20];
398 TBranch *branch;
399 if (fPrintLevel >= 0) {
400 cout << "Enter NextBkgGeantEvent" << endl;}
401 // Clean previous event
402 if(fBkgGeantTK) delete fBkgGeantTK;
403 fBkgGeantTK = NULL;
404 if(fBkgGeantParticles) fBkgGeantParticles->Clear();
405 if(fBkgGeantTH) delete fBkgGeantTH;
406 fBkgGeantTH = NULL;
407 if(fBkgGeantHits) fBkgGeantHits->Clear();
408 // Increment event number
409 fBkgGeantEventNumber++;
410 // Get access to Particles and Hits for event from background file
411 if (fPrintLevel >= 2) {cout << "Before cd(Bkg)" << endl; gDirectory->Dump();}
412 fBkgGeantFile->cd();
413 if (fPrintLevel >= 2) {cout << "After cd(Bkg)" << endl; gDirectory->Dump();}
414 // Particles: TreeK for event and branch "Particles"
415 sprintf(treeName, "TreeK%d", fBkgGeantEventNumber);
416 fBkgGeantTK = (TTree*)gDirectory->Get(treeName);
417 if (!fBkgGeantTK) {
418 if (fPrintLevel >= 0) {
419 cout << "Cannot find Kine Tree for background event: " <<
420 fBkgGeantEventNumber << endl;
421 cout << "Goes back to event 0" << endl;
422 }
423 fBkgGeantEventNumber = 0;
424 sprintf(treeName, "TreeK%d", fBkgGeantEventNumber);
425 fBkgGeantTK = (TTree*)gDirectory->Get(treeName);
426 if (!fBkgGeantTK) {
427 cout << "ERROR: cannot find Kine Tree for background event: " <<
428 fBkgGeantEventNumber << endl;
429 exit(0);
430 }
431 }
432 if (fBkgGeantTK)
433 fBkgGeantTK->SetBranchAddress("Particles", &fBkgGeantParticles);
434 fBkgGeantTK->GetEvent(0); // why event 0 ???? necessary ????
435 // Hits: TreeH for event and branch "MUON"
436 sprintf(treeName, "TreeH%d", fBkgGeantEventNumber);
437 fBkgGeantTH = (TTree*)gDirectory->Get(treeName);
438 if (!fBkgGeantTH) {
439 cout << "ERROR: cannot find Hits Tree for background event: " <<
440 fBkgGeantEventNumber << endl;
441 exit(0);
442 }
443 if (fBkgGeantTH && fBkgGeantHits) {
444 branch = fBkgGeantTH->GetBranch("MUON");
445 if (branch) branch->SetAddress(&fBkgGeantHits);
446 }
447 fBkgGeantTH->GetEntries(); // necessary ????
448 // Back to the signal file
449 ((gAlice->TreeK())->GetCurrentFile())->cd();
450 if (fPrintLevel >= 2) {cout << "After cd(gAlice)" << endl; gDirectory->Dump();}
451 return;
452}
453
454//__________________________________________________________________________
455void AliMUONEventReconstructor::EventReconstruct(void)
456{
457 // To reconstruct one event
458 if (fPrintLevel >= 1) cout << "enter EventReconstruct" << endl;
459 MakeEventToBeReconstructed();
460 MakeSegments();
461 MakeTracks();
462 return;
463}
464
465 //__________________________________________________________________________
466void AliMUONEventReconstructor::ResetHitsForRec(void)
467{
468 // To reset the array and the number of HitsForRec,
469 // and also the number of HitsForRec
470 // and the index of the first HitForRec per chamber
471 if (fHitsForRecPtr) fHitsForRecPtr->Clear();
472 fNHitsForRec = 0;
9cbdf048 473 for (Int_t ch = 0; ch < kMaxMuonTrackingChambers; ch++)
a9e2aefa 474 fNHitsForRecPerChamber[ch] = fIndexOfFirstHitForRecPerChamber[ch] = 0;
475 return;
476}
477
478 //__________________________________________________________________________
479void AliMUONEventReconstructor::ResetSegments(void)
480{
481 // To reset the TClonesArray of segments and the number of Segments
482 // for all stations
9cbdf048 483 for (Int_t st = 0; st < kMaxMuonTrackingStations; st++) {
a9e2aefa 484 if (fSegmentsPtr[st]) fSegmentsPtr[st]->Clear();
485 fNSegments[st] = 0;
486 }
487 return;
488}
489
490 //__________________________________________________________________________
491void AliMUONEventReconstructor::ResetTracks(void)
492{
493 // To reset the TClonesArray of reconstructed tracks
8429a5e4 494 if (fRecTracksPtr) fRecTracksPtr->Delete();
495 // Delete in order that the Track destructors are called,
496 // hence the space for the TClonesArray of pointers to TrackHit's is freed
a9e2aefa 497 fNRecTracks = 0;
498 return;
499}
500
8429a5e4 501 //__________________________________________________________________________
502void AliMUONEventReconstructor::ResetTrackHits(void)
503{
504 // To reset the TClonesArray of hits on reconstructed tracks
505 if (fRecTrackHitsPtr) fRecTrackHitsPtr->Clear();
506 fNRecTrackHits = 0;
507 return;
508}
509
a9e2aefa 510 //__________________________________________________________________________
511void AliMUONEventReconstructor::MakeEventToBeReconstructed(void)
512{
513 // To make the list of hits to be reconstructed,
514 // either from the GEANT hits or from the raw clusters
515 // according to the parameter set for the reconstructor
516 if (fPrintLevel >= 1) cout << "enter MakeEventToBeReconstructed" << endl;
517 ResetHitsForRec();
518 if (fRecGeantHits == 1) {
519 // Reconstruction from GEANT hits
520 // Back to the signal file
521 ((gAlice->TreeK())->GetCurrentFile())->cd();
522 // Signal hits
523 // AliMUON *MUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
524 // Security on MUON ????
525 AddHitsForRecFromGEANT(gAlice->TreeH());
526 // Background hits
527 AddHitsForRecFromBkgGEANT(fBkgGeantTH, fBkgGeantHits);
528 // Sort HitsForRec in increasing order with respect to chamber number
529 SortHitsForRecWithIncreasingChamber();
530 }
531 else {
532 // Reconstruction from raw clusters
533 // AliMUON *MUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
534 // Security on MUON ????
535 // TreeR assumed to be be "prepared" in calling function
536 // by "MUON->GetTreeR(nev)" ????
9cbdf048 537 TTree *treeR = gAlice->TreeR();
538 AddHitsForRecFromRawClusters(treeR);
a9e2aefa 539 // No sorting: it is done automatically in the previous function
540 }
541 if (fPrintLevel >= 10) {
542 cout << "end of MakeEventToBeReconstructed" << endl;
543 cout << "NHitsForRec: " << fNHitsForRec << endl;
9cbdf048 544 for (Int_t ch = 0; ch < kMaxMuonTrackingChambers; ch++) {
a9e2aefa 545 cout << "chamber(0...): " << ch
546 << " NHitsForRec: " << fNHitsForRecPerChamber[ch]
547 << " index(first HitForRec): " << fIndexOfFirstHitForRecPerChamber[ch]
548 << endl;
549 for (Int_t hit = fIndexOfFirstHitForRecPerChamber[ch];
550 hit < fIndexOfFirstHitForRecPerChamber[ch] + fNHitsForRecPerChamber[ch];
551 hit++) {
552 cout << "HitForRec index(0...): " << hit << endl;
553 ((*fHitsForRecPtr)[hit])->Dump();
554 }
555 }
556 }
557 return;
558}
559
560 //__________________________________________________________________________
561void AliMUONEventReconstructor::AddHitsForRecFromGEANT(TTree *TH)
562{
563 // To add to the list of hits for reconstruction
564 // the GEANT signal hits from a hit tree TH.
565 if (fPrintLevel >= 2)
566 cout << "enter AddHitsForRecFromGEANT with TH: " << TH << endl;
567 if (TH == NULL) return;
9cbdf048 568 AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
a9e2aefa 569 // Security on MUON ????
570 // See whether it could be the same for signal and background ????
571 // Loop over tracks in tree
572 Int_t ntracks = (Int_t) TH->GetEntries();
573 if (fPrintLevel >= 2)
574 cout << "ntracks: " << ntracks << endl;
575 for (Int_t track = 0; track < ntracks; track++) {
576 gAlice->ResetHits();
577 TH->GetEvent(track);
578 // Loop over hits
579 Int_t hit = 0;
9cbdf048 580 for (AliMUONHit* mHit = (AliMUONHit*) pMUON->FirstHit(-1);
a9e2aefa 581 mHit;
9cbdf048 582 mHit = (AliMUONHit*) pMUON->NextHit(), hit++) {
a9e2aefa 583 NewHitForRecFromGEANT(mHit,track, hit, 1);
584 } // end of hit loop
585 } // end of track loop
586 return;
587}
588
589 //__________________________________________________________________________
590void AliMUONEventReconstructor::AddHitsForRecFromBkgGEANT(TTree *TH, TClonesArray *Hits)
591{
592 // To add to the list of hits for reconstruction
593 // the GEANT background hits from a hit tree TH and a pointer Hits to a hit list.
594 // How to have only one function "AddHitsForRecFromGEANT" ????
595 if (fPrintLevel >= 2)
596 cout << "enter AddHitsForRecFromBkgGEANT with TH: " << TH << endl;
597 if (TH == NULL) return;
598 // Loop over tracks in tree
599 Int_t ntracks = (Int_t) TH->GetEntries();
600 if (fPrintLevel >= 2)
601 cout << "ntracks: " << ntracks << endl;
602 for (Int_t track = 0; track < ntracks; track++) {
603 if (Hits) Hits->Clear();
604 TH->GetEvent(track);
605 // Loop over hits
606 for (Int_t hit = 0; hit < Hits->GetEntriesFast(); hit++) {
607 NewHitForRecFromGEANT((AliMUONHit*) (*Hits)[hit], track, hit, 0);
608 } // end of hit loop
609 } // end of track loop
610 return;
611}
612
613 //__________________________________________________________________________
614AliMUONHitForRec* AliMUONEventReconstructor::NewHitForRecFromGEANT(AliMUONHit* Hit, Int_t TrackNumber, Int_t HitNumber, Int_t Signal)
615{
616 // To make a new hit for reconstruction from a GEANT hit pointed to by "Hit",
617 // with hit number "HitNumber" in the track numbered "TrackNumber",
618 // either from signal ("Signal" = 1) or background ("Signal" = 0) event.
619 // Selects hits in tracking (not trigger) chambers.
620 // Takes into account the efficiency (fEfficiency)
621 // and the smearing from resolution (fBendingResolution and fNonBendingResolution).
622 // Adds a condition on the radius between RMin and RMax
623 // to better simulate the real chambers.
624 // Returns the pointer to the new hit for reconstruction,
625 // or NULL in case of inefficiency or non tracking chamber or bad radius.
626 // No condition on at most 20 cm from a muon from a resonance
627 // like in Fortran TRACKF_STAT.
628 AliMUONHitForRec* hitForRec;
629 Double_t bendCoor, nonBendCoor, radius;
0a5f36d9 630 Int_t chamber = Hit->Chamber() - 1; // chamber(0...)
a9e2aefa 631 // only in tracking chambers (fChamber starts at 1)
9cbdf048 632 if (chamber >= kMaxMuonTrackingChambers) return NULL;
a9e2aefa 633 // only if hit is efficient (keep track for checking ????)
634 if (gRandom->Rndm() > fEfficiency) return NULL;
635 // only if radius between RMin and RMax
94de3818 636 bendCoor = Hit->Y();
637 nonBendCoor = Hit->X();
a9e2aefa 638 radius = TMath::Sqrt((bendCoor * bendCoor) + (nonBendCoor * nonBendCoor));
d0bfce8d 639 // This cut is not needed with a realistic chamber geometry !!!!
640// if ((radius < fRMin[chamber]) || (radius > fRMax[chamber])) return NULL;
a9e2aefa 641 // new AliMUONHitForRec from GEANT hit and increment number of AliMUONHitForRec's
642 hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(Hit);
643 fNHitsForRec++;
644 // add smearing from resolution
645 hitForRec->SetBendingCoor(bendCoor + gRandom->Gaus(0., fBendingResolution));
646 hitForRec->SetNonBendingCoor(nonBendCoor
647 + gRandom->Gaus(0., fNonBendingResolution));
d0bfce8d 648// // !!!! without smearing
649// hitForRec->SetBendingCoor(bendCoor);
650// hitForRec->SetNonBendingCoor(nonBendCoor);
a9e2aefa 651 // more information into HitForRec
652 // resolution: angular effect to be added here ????
653 hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
654 hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
655 // GEANT track info
656 hitForRec->SetHitNumber(HitNumber);
657 hitForRec->SetTHTrack(TrackNumber);
658 hitForRec->SetGeantSignal(Signal);
659 if (fPrintLevel >= 10) {
660 cout << "track: " << TrackNumber << " hit: " << HitNumber << endl;
661 Hit->Dump();
662 cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
663 hitForRec->Dump();}
664 return hitForRec;
665}
666
667 //__________________________________________________________________________
668void AliMUONEventReconstructor::SortHitsForRecWithIncreasingChamber()
669{
670 // Sort HitsForRec's in increasing order with respect to chamber number.
671 // Uses the function "Compare".
672 // Update the information for HitsForRec per chamber too.
673 Int_t ch, nhits, prevch;
674 fHitsForRecPtr->Sort();
9cbdf048 675 for (ch = 0; ch < kMaxMuonTrackingChambers; ch++) {
a9e2aefa 676 fNHitsForRecPerChamber[ch] = 0;
677 fIndexOfFirstHitForRecPerChamber[ch] = 0;
678 }
679 prevch = 0; // previous chamber
680 nhits = 0; // number of hits in current chamber
681 // Loop over HitsForRec
682 for (Int_t hit = 0; hit < fNHitsForRec; hit++) {
683 // chamber number (0...)
684 ch = ((AliMUONHitForRec*) ((*fHitsForRecPtr)[hit]))->GetChamberNumber();
685 // increment number of hits in current chamber
686 (fNHitsForRecPerChamber[ch])++;
687 // update index of first HitForRec in current chamber
688 // if chamber number different from previous one
689 if (ch != prevch) {
690 fIndexOfFirstHitForRecPerChamber[ch] = hit;
691 prevch = ch;
692 }
693 }
694 return;
695}
696
697// //__________________________________________________________________________
698// void AliMUONEventReconstructor::AddHitsForRecFromCathodeCorrelations(TTree* TC)
699// {
700// // OLD VERSION WHEN ONE ONE WAS USING SO CALLED CATHODE CORRELATIONS
701// // To add to the list of hits for reconstruction
702// // the (cathode correlated) raw clusters
703// // No condition added, like in Fortran TRACKF_STAT,
704// // on the radius between RMin and RMax.
705// AliMUONHitForRec *hitForRec;
706// if (fPrintLevel >= 1) cout << "enter AddHitsForRecFromCathodeCorrelations" << endl;
707// AliMUON *MUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
708// // Security on MUON ????
709// // Loop over tracking chambers
9cbdf048 710// for (Int_t ch = 0; ch < kMaxMuonTrackingChambers; ch++) {
a9e2aefa 711// // number of HitsForRec to 0 for the chamber
712// fNHitsForRecPerChamber[ch] = 0;
713// // index of first HitForRec for the chamber
714// if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0;
715// else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec;
716// TClonesArray *reconst_hits = MUON->ReconstHitsAddress(ch);
717// MUON->ResetReconstHits();
718// TC->GetEvent();
719// Int_t ncor = (Int_t)reconst_hits->GetEntries();
720// // Loop over (cathode correlated) raw clusters
721// for (Int_t cor = 0; cor < ncor; cor++) {
722// AliMUONReconstHit * mCor =
723// (AliMUONReconstHit*) reconst_hits->UncheckedAt(cor);
724// // new AliMUONHitForRec from (cathode correlated) raw cluster
725// // and increment number of AliMUONHitForRec's (total and in chamber)
726// hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(mCor);
727// fNHitsForRec++;
728// (fNHitsForRecPerChamber[ch])++;
729// // more information into HitForRec
730// hitForRec->SetChamberNumber(ch);
731// hitForRec->SetHitNumber(cor);
732// // Z coordinate of the chamber (cm) with sign opposite to GEANT sign
733// // could (should) be more exact from chamber geometry ????
734// hitForRec->SetZ(-(&(MUON->Chamber(ch)))->Z());
735// if (fPrintLevel >= 10) {
736// cout << "chamber (0...): " << ch <<
737// " cathcorrel (0...): " << cor << endl;
738// mCor->Dump();
739// cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
740// hitForRec->Dump();}
741// } // end of cluster loop
742// } // end of chamber loop
743// return;
744// }
745
746 //__________________________________________________________________________
747void AliMUONEventReconstructor::AddHitsForRecFromRawClusters(TTree* TR)
748{
749 // To add to the list of hits for reconstruction all the raw clusters
750 // No condition added, like in Fortran TRACKF_STAT,
751 // on the radius between RMin and RMax.
752 AliMUONHitForRec *hitForRec;
753 AliMUONRawCluster *clus;
754 Int_t iclus, nclus;
755 TClonesArray *rawclusters;
756 if (fPrintLevel >= 1) cout << "enter AddHitsForRecFromRawClusters" << endl;
9cbdf048 757 AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
a9e2aefa 758 // Security on MUON ????
759 // Loop over tracking chambers
9cbdf048 760 for (Int_t ch = 0; ch < kMaxMuonTrackingChambers; ch++) {
a9e2aefa 761 // number of HitsForRec to 0 for the chamber
762 fNHitsForRecPerChamber[ch] = 0;
763 // index of first HitForRec for the chamber
764 if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0;
765 else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec;
9cbdf048 766 rawclusters = pMUON->RawClustAddress(ch);
767 pMUON->ResetRawClusters();
a9e2aefa 768 TR->GetEvent((Int_t) (TR->GetEntries()) - 1); // to be checked ????
769 nclus = (Int_t) (rawclusters->GetEntries());
770 // Loop over (cathode correlated) raw clusters
771 for (iclus = 0; iclus < nclus; iclus++) {
772 clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(iclus);
773 // new AliMUONHitForRec from raw cluster
774 // and increment number of AliMUONHitForRec's (total and in chamber)
775 hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(clus);
776 fNHitsForRec++;
777 (fNHitsForRecPerChamber[ch])++;
778 // more information into HitForRec
779 // resolution: info should be already in raw cluster and taken from it ????
780 hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution);
781 hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution);
782 // original raw cluster
783 hitForRec->SetChamberNumber(ch);
784 hitForRec->SetHitNumber(iclus);
ba5b68db 785 // Z coordinate of the raw cluster (cm)
cfe0107f 786 hitForRec->SetZ(clus->fZ[0]);
a9e2aefa 787 if (fPrintLevel >= 10) {
788 cout << "chamber (0...): " << ch <<
789 " raw cluster (0...): " << iclus << endl;
790 clus->Dump();
791 cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
792 hitForRec->Dump();}
793 } // end of cluster loop
794 } // end of chamber loop
795 return;
796}
797
798 //__________________________________________________________________________
799void AliMUONEventReconstructor::MakeSegments(void)
800{
801 // To make the list of segments in all stations,
802 // from the list of hits to be reconstructed
803 if (fPrintLevel >= 1) cout << "enter MakeSegments" << endl;
804 ResetSegments();
805 // Loop over stations
9cbdf048 806 for (Int_t st = 0; st < kMaxMuonTrackingStations; st++)
a9e2aefa 807 MakeSegmentsPerStation(st);
808 if (fPrintLevel >= 10) {
809 cout << "end of MakeSegments" << endl;
9cbdf048 810 for (Int_t st = 0; st < kMaxMuonTrackingStations; st++) {
a9e2aefa 811 cout << "station(0...): " << st
812 << " Segments: " << fNSegments[st]
813 << endl;
814 for (Int_t seg = 0;
815 seg < fNSegments[st];
816 seg++) {
817 cout << "Segment index(0...): " << seg << endl;
818 ((*fSegmentsPtr[st])[seg])->Dump();
819 }
820 }
821 }
822 return;
823}
824
825 //__________________________________________________________________________
826void AliMUONEventReconstructor::MakeSegmentsPerStation(Int_t Station)
827{
828 // To make the list of segments in station number "Station" (0...)
829 // from the list of hits to be reconstructed.
830 // Updates "fNSegments"[Station].
831 // Segments in stations 4 and 5 are sorted
832 // according to increasing absolute value of "impact parameter"
833 AliMUONHitForRec *hit1Ptr, *hit2Ptr;
834 AliMUONSegment *segment;
835 Bool_t last2st;
836 Double_t bendingSlope, distBend, distNonBend, extBendCoor, extNonBendCoor,
d0bfce8d 837 impactParam = 0., maxImpactParam = 0., minImpactParam = 0.; // =0 to avoid compilation warnings.
9cbdf048 838 AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
a9e2aefa 839 if (fPrintLevel >= 1)
840 cout << "enter MakeSegmentsPerStation (0...) " << Station << endl;
841 // first and second chambers (0...) in the station
842 Int_t ch1 = 2 * Station;
843 Int_t ch2 = ch1 + 1;
844 // variable true for stations downstream of the dipole:
845 // Station(0..4) equal to 3 or 4
846 if ((Station == 3) || (Station == 4)) {
847 last2st = kTRUE;
848 // maximum impact parameter (cm) according to fMinBendingMomentum...
849 maxImpactParam =
850 TMath::Abs(GetImpactParamFromBendingMomentum(fMinBendingMomentum));
d0bfce8d 851 // minimum impact parameter (cm) according to fMaxBendingMomentum...
852 minImpactParam =
853 TMath::Abs(GetImpactParamFromBendingMomentum(fMaxBendingMomentum));
a9e2aefa 854 }
855 else last2st = kFALSE;
856 // extrapolation factor from Z of first chamber to Z of second chamber
857 // dZ to be changed to take into account fine structure of chambers ????
858 Double_t extrapFact =
9cbdf048 859 (&(pMUON->Chamber(ch2)))->Z() / (&(pMUON->Chamber(ch1)))->Z();
a9e2aefa 860 // index for current segment
861 Int_t segmentIndex = 0;
862 // Loop over HitsForRec in the first chamber of the station
863 for (Int_t hit1 = fIndexOfFirstHitForRecPerChamber[ch1];
864 hit1 < fIndexOfFirstHitForRecPerChamber[ch1] + fNHitsForRecPerChamber[ch1];
865 hit1++) {
866 // pointer to the HitForRec
867 hit1Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit1]);
868 // extrapolation,
869 // on the straight line joining the HitForRec to the vertex (0,0,0),
870 // to the Z of the second chamber of the station
871 extBendCoor = extrapFact * hit1Ptr->GetBendingCoor();
872 extNonBendCoor = extrapFact * hit1Ptr->GetNonBendingCoor();
873 // Loop over HitsForRec in the second chamber of the station
874 for (Int_t hit2 = fIndexOfFirstHitForRecPerChamber[ch2];
875 hit2 < fIndexOfFirstHitForRecPerChamber[ch2] + fNHitsForRecPerChamber[ch2];
876 hit2++) {
877 // pointer to the HitForRec
878 hit2Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit2]);
879 // absolute values of distances, in bending and non bending planes,
880 // between the HitForRec in the second chamber
881 // and the previous extrapolation
882 distBend = TMath::Abs(hit2Ptr->GetBendingCoor() - extBendCoor);
883 distNonBend = TMath::Abs(hit2Ptr->GetNonBendingCoor() - extNonBendCoor);
884 if (last2st) {
885 // bending slope
886 bendingSlope = (hit1Ptr->GetBendingCoor() - hit2Ptr->GetBendingCoor()) /
887 (hit1Ptr->GetZ() - hit2Ptr->GetZ());
888 // absolute value of impact parameter
889 impactParam =
890 TMath::Abs(hit1Ptr->GetBendingCoor() - hit2Ptr->GetZ() * bendingSlope);
891 }
892 // check for distances not too large,
893 // and impact parameter not too big if stations downstream of the dipole.
894 // Conditions "distBend" and "impactParam" correlated for these stations ????
895 if ((distBend < fSegmentMaxDistBending[Station]) &&
896 (distNonBend < fSegmentMaxDistNonBending[Station]) &&
d0bfce8d 897 (!last2st || (impactParam < maxImpactParam)) &&
898 (!last2st || (impactParam > minImpactParam))) {
a9e2aefa 899 // make new segment
900 segment = new ((*fSegmentsPtr[Station])[segmentIndex])
901 AliMUONSegment(hit1Ptr, hit2Ptr);
902 // update "link" to this segment from the hit in the first chamber
903 if (hit1Ptr->GetNSegments() == 0)
904 hit1Ptr->SetIndexOfFirstSegment(segmentIndex);
905 hit1Ptr->SetNSegments(hit1Ptr->GetNSegments() + 1);
906 if (fPrintLevel >= 10) {
907 cout << "segmentIndex(0...): " << segmentIndex
908 << " distBend: " << distBend
909 << " distNonBend: " << distNonBend
910 << endl;
911 segment->Dump();
912 cout << "HitForRec in first chamber" << endl;
913 hit1Ptr->Dump();
914 cout << "HitForRec in second chamber" << endl;
915 hit2Ptr->Dump();
916 };
917 // increment index for current segment
918 segmentIndex++;
919 }
920 } //for (Int_t hit2
921 } // for (Int_t hit1...
922 fNSegments[Station] = segmentIndex;
923 // Sorting according to "impact parameter" if station(1..5) 4 or 5,
924 // i.e. Station(0..4) 3 or 4, using the function "Compare".
925 // After this sorting, it is impossible to use
926 // the "fNSegments" and "fIndexOfFirstSegment"
927 // of the HitForRec in the first chamber to explore all segments formed with it.
928 // Is this sorting really needed ????
929 if ((Station == 3) || (Station == 4)) (fSegmentsPtr[Station])->Sort();
930 if (fPrintLevel >= 1) cout << "Station: " << Station << " NSegments: "
931 << fNSegments[Station] << endl;
932 return;
933}
934
935 //__________________________________________________________________________
936void AliMUONEventReconstructor::MakeTracks(void)
937{
938 // To make the tracks,
939 // from the list of segments and points in all stations
940 if (fPrintLevel >= 1) cout << "enter MakeTracks" << endl;
8429a5e4 941 // The order may be important for the following Reset's
a9e2aefa 942 ResetTracks();
8429a5e4 943 ResetTrackHits();
a9e2aefa 944 // Look for candidates from at least 3 aligned points in stations(1..) 4 and 5
945 MakeTrackCandidates();
946 // Follow tracks in stations(1..) 3, 2 and 1
947 FollowTracks();
8429a5e4 948 // Remove double tracks
949 RemoveDoubleTracks();
a9e2aefa 950 return;
951}
952
953 //__________________________________________________________________________
954Int_t AliMUONEventReconstructor::MakeTrackCandidatesWithTwoSegments(AliMUONSegment *BegSegment)
955{
956 // To make track candidates with two segments in stations(1..) 4 and 5,
957 // the first segment being pointed to by "BegSegment".
958 // Returns the number of such track candidates.
959 Int_t endStation, iEndSegment, nbCan2Seg;
960 AliMUONSegment *endSegment, *extrapSegment;
961 AliMUONTrack *recTrack;
962 Double_t mcsFactor;
963 if (fPrintLevel >= 1) cout << "enter MakeTrackCandidatesWithTwoSegments" << endl;
964 // Station for the end segment
965 endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
966 // multiple scattering factor corresponding to one chamber
967 mcsFactor = 0.0136 /
968 GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
969 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
970 // linear extrapolation to end station
971 extrapSegment =
972 BegSegment->CreateSegmentFromLinearExtrapToStation(endStation, mcsFactor);
973 // number of candidates with 2 segments to 0
974 nbCan2Seg = 0;
975 // Loop over segments in the end station
976 for (iEndSegment = 0; iEndSegment < fNSegments[endStation]; iEndSegment++) {
977 // pointer to segment
978 endSegment = (AliMUONSegment*) ((*fSegmentsPtr[endStation])[iEndSegment]);
979 // test compatibility between current segment and "extrapSegment"
04b5ea16 980 // 4 because 4 quantities in chi2
a9e2aefa 981 if ((endSegment->
982 NormalizedChi2WithSegment(extrapSegment,
983 fMaxSigma2Distance)) <= 4.0) {
984 // both segments compatible:
985 // make track candidate from "begSegment" and "endSegment"
986 if (fPrintLevel >= 2)
987 cout << "TrackCandidate with Segment " << iEndSegment <<
988 " in Station(0..) " << endStation << endl;
989 // flag for both segments in one track:
990 // to be done in track constructor ????
991 BegSegment->SetInTrack(kTRUE);
992 endSegment->SetInTrack(kTRUE);
993 recTrack = new ((*fRecTracksPtr)[fNRecTracks])
994 AliMUONTrack(BegSegment, endSegment, this);
995 fNRecTracks++;
996 if (fPrintLevel >= 10) recTrack->RecursiveDump();
997 // increment number of track candidates with 2 segments
998 nbCan2Seg++;
999 }
1000 } // for (iEndSegment = 0;...
1001 delete extrapSegment; // should not delete HitForRec's it points to !!!!
1002 return nbCan2Seg;
1003}
1004
1005 //__________________________________________________________________________
1006Int_t AliMUONEventReconstructor::MakeTrackCandidatesWithOneSegmentAndOnePoint(AliMUONSegment *BegSegment)
1007{
1008 // To make track candidates with one segment and one point
1009 // in stations(1..) 4 and 5,
1010 // the segment being pointed to by "BegSegment".
1011 Int_t ch, ch1, ch2, endStation, iHit, iHitMax, iHitMin, nbCan1Seg1Hit;
1012 AliMUONHitForRec *extrapHitForRec, *hit;
1013 AliMUONTrack *recTrack;
1014 Double_t mcsFactor;
1015 if (fPrintLevel >= 1)
1016 cout << "enter MakeTrackCandidatesWithOneSegmentAndOnePoint" << endl;
1017 // station for the end point
1018 endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
1019 // multiple scattering factor corresponding to one chamber
1020 mcsFactor = 0.0136 /
1021 GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
1022 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1023 // first and second chambers(0..) in the end station
1024 ch1 = 2 * endStation;
1025 ch2 = ch1 + 1;
1026 // number of candidates to 0
1027 nbCan1Seg1Hit = 0;
1028 // Loop over chambers of the end station
1029 for (ch = ch2; ch >= ch1; ch--) {
1030 // linear extrapolation to chamber
1031 extrapHitForRec =
1032 BegSegment->CreateHitForRecFromLinearExtrapToChamber(ch, mcsFactor);
1033 // limits for the hit index in the loop
1034 iHitMin = fIndexOfFirstHitForRecPerChamber[ch];
1035 iHitMax = iHitMin + fNHitsForRecPerChamber[ch];
1036 // Loop over HitForRec's in the chamber
1037 for (iHit = iHitMin; iHit < iHitMax; iHit++) {
1038 // pointer to HitForRec
1039 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
1040 // test compatibility between current HitForRec and "extrapHitForRec"
04b5ea16 1041 // 2 because 2 quantities in chi2
a9e2aefa 1042 if ((hit->
1043 NormalizedChi2WithHitForRec(extrapHitForRec,
1044 fMaxSigma2Distance)) <= 2.0) {
1045 // both HitForRec's compatible:
1046 // make track candidate from begSegment and current HitForRec
1047 if (fPrintLevel >= 2)
1048 cout << "TrackCandidate with HitForRec " << iHit <<
1049 " in Chamber(0..) " << ch << endl;
1050 // flag for beginning segments in one track:
1051 // to be done in track constructor ????
1052 BegSegment->SetInTrack(kTRUE);
1053 recTrack = new ((*fRecTracksPtr)[fNRecTracks])
1054 AliMUONTrack(BegSegment, hit, this);
1055 // the right place to eliminate "double counting" ???? how ????
1056 fNRecTracks++;
1057 if (fPrintLevel >= 10) recTrack->RecursiveDump();
1058 // increment number of track candidates
1059 nbCan1Seg1Hit++;
1060 }
1061 } // for (iHit = iHitMin;...
1062 delete extrapHitForRec;
1063 } // for (ch = ch2;...
1064 return nbCan1Seg1Hit;
1065}
1066
1067 //__________________________________________________________________________
1068void AliMUONEventReconstructor::MakeTrackCandidates(void)
1069{
1070 // To make track candidates
1071 // with at least 3 aligned points in stations(1..) 4 and 5
1072 // (two Segment's or one Segment and one HitForRec)
1073 Int_t begStation, iBegSegment, nbCan1Seg1Hit, nbCan2Seg;
1074 AliMUONSegment *begSegment;
1075 if (fPrintLevel >= 1) cout << "enter MakeTrackCandidates" << endl;
1076 // Loop over stations(1..) 5 and 4 for the beginning segment
1077 for (begStation = 4; begStation > 2; begStation--) {
1078 // Loop over segments in the beginning station
1079 for (iBegSegment = 0; iBegSegment < fNSegments[begStation]; iBegSegment++) {
1080 // pointer to segment
1081 begSegment = (AliMUONSegment*) ((*fSegmentsPtr[begStation])[iBegSegment]);
1082 if (fPrintLevel >= 2)
1083 cout << "look for TrackCandidate's with Segment " << iBegSegment <<
1084 " in Station(0..) " << begStation << endl;
1085 // Look for track candidates with two segments,
1086 // "begSegment" and all compatible segments in other station.
1087 // Only for beginning station(1..) 5
1088 // because candidates with 2 segments have to looked for only once.
1089 if (begStation == 4)
1090 nbCan2Seg = MakeTrackCandidatesWithTwoSegments(begSegment);
956019b6 1091 // Look for track candidates with one segment and one point,
a9e2aefa 1092 // "begSegment" and all compatible HitForRec's in other station.
1093 // Only if "begSegment" does not belong already to a track candidate.
1094 // Is that a too strong condition ????
1095 if (!(begSegment->GetInTrack()))
1096 nbCan1Seg1Hit = MakeTrackCandidatesWithOneSegmentAndOnePoint(begSegment);
1097 } // for (iBegSegment = 0;...
1098 } // for (begStation = 4;...
1099 return;
1100}
1101
1102 //__________________________________________________________________________
1103void AliMUONEventReconstructor::FollowTracks(void)
1104{
1105 // Follow tracks in stations(1..) 3, 2 and 1
04b5ea16 1106 // too long: should be made more modular !!!!
d0bfce8d 1107 AliMUONHitForRec *bestHit, *extrapHit, *extrapCorrHit, *hit;
1108 AliMUONSegment *bestSegment, *extrapSegment, *extrapCorrSegment, *segment;
04b5ea16 1109 AliMUONTrack *track, *nextTrack;
1110 AliMUONTrackParam *trackParam1, trackParam[2], trackParamVertex;
ecfa008b 1111 // -1 to avoid compilation warnings
1112 Int_t ch = -1, chInStation, chBestHit = -1, iHit, iSegment, station, trackIndex;
04b5ea16 1113 Double_t bestChi2, chi2, dZ1, dZ2, dZ3, maxSigma2Distance, mcsFactor;
d0bfce8d 1114 Double_t bendingMomentum, chi2Norm = 0.;
9cbdf048 1115 AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
04b5ea16 1116 // local maxSigma2Distance, for easy increase in testing
1117 maxSigma2Distance = fMaxSigma2Distance;
a9e2aefa 1118 if (fPrintLevel >= 2)
1119 cout << "enter FollowTracks" << endl;
1120 // Loop over track candidates
04b5ea16 1121 track = (AliMUONTrack*) fRecTracksPtr->First();
1122 trackIndex = -1;
1123 while (track) {
1124 // Follow function for each track candidate ????
1125 trackIndex++;
1126 nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track); // prepare next track
a9e2aefa 1127 if (fPrintLevel >= 2)
1128 cout << "FollowTracks: track candidate(0..): " << trackIndex << endl;
956019b6 1129 // Fit track candidate
1130 track->SetFitMCS(0); // without multiple Coulomb scattering
1131 track->SetFitNParam(3); // with 3 parameters (X = Y = 0)
1132 track->SetFitStart(0); // from parameters at vertex
1133 track->Fit();
a9e2aefa 1134 if (fPrintLevel >= 10) {
1135 cout << "FollowTracks: track candidate(0..): " << trackIndex
956019b6 1136 << " after fit in stations(0..) 3 and 4" << endl;
a9e2aefa 1137 track->RecursiveDump();
1138 }
1139 // Loop over stations(1..) 3, 2 and 1
04b5ea16 1140 // something SPECIAL for stations 2 and 1 for majority 3 coincidence ????
1141 // otherwise: majority coincidence 2 !!!!
a9e2aefa 1142 for (station = 2; station >= 0; station--) {
1143 // Track parameters at first track hit (smallest Z)
1144 trackParam1 = ((AliMUONTrackHit*)
1145 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1146 // extrapolation to station
1147 trackParam1->ExtrapToStation(station, trackParam);
79e1e601 1148 extrapSegment = new AliMUONSegment(); // empty segment
d0bfce8d 1149 extrapCorrSegment = new AliMUONSegment(); // empty corrected segment
a9e2aefa 1150 // multiple scattering factor corresponding to one chamber
1151 // and momentum in bending plane (not total)
1152 mcsFactor = 0.0136 * trackParam1->GetInverseBendingMomentum();
1153 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
1154 // Z difference from previous station
9cbdf048 1155 dZ1 = (&(pMUON->Chamber(2 * station)))->Z() -
1156 (&(pMUON->Chamber(2 * station + 2)))->Z();
a9e2aefa 1157 // Z difference between the two previous stations
9cbdf048 1158 dZ2 = (&(pMUON->Chamber(2 * station + 2)))->Z() -
1159 (&(pMUON->Chamber(2 * station + 4)))->Z();
04b5ea16 1160 // Z difference between the two chambers in the previous station
1161 dZ3 = (&(pMUON->Chamber(2 * station)))->Z() -
1162 (&(pMUON->Chamber(2 * station + 1)))->Z();
1163 extrapSegment->SetBendingCoorReso2(fBendingResolution * fBendingResolution);
1164 extrapSegment->
1165 SetNonBendingCoorReso2(fNonBendingResolution * fNonBendingResolution);
1166 extrapSegment->UpdateFromStationTrackParam
1167 (trackParam, mcsFactor, dZ1, dZ2, dZ3, station,
1168 trackParam1->GetInverseBendingMomentum());
d0bfce8d 1169 // same thing for corrected segment
1170 // better to use copy constructor, after checking that it works properly !!!!
1171 extrapCorrSegment->SetBendingCoorReso2(fBendingResolution * fBendingResolution);
1172 extrapCorrSegment->
1173 SetNonBendingCoorReso2(fNonBendingResolution * fNonBendingResolution);
1174 extrapCorrSegment->UpdateFromStationTrackParam
1175 (trackParam, mcsFactor, dZ1, dZ2, dZ3, station,
1176 trackParam1->GetInverseBendingMomentum());
a9e2aefa 1177 bestChi2 = 5.0;
1178 bestSegment = NULL;
1179 if (fPrintLevel >= 10) {
1180 cout << "FollowTracks: track candidate(0..): " << trackIndex
1181 << " Look for segment in station(0..): " << station << endl;
1182 }
1183 // Loop over segments in station
1184 for (iSegment = 0; iSegment < fNSegments[station]; iSegment++) {
1185 // Look for best compatible Segment in station
1186 // should consider all possibilities ????
1187 // multiple scattering ????
1188 // separation in 2 functions: Segment and HitForRec ????
1189 segment = (AliMUONSegment*) ((*fSegmentsPtr[station])[iSegment]);
d0bfce8d 1190 // correction of corrected segment (fBendingCoor and fNonBendingCoor)
1191 // according to real Z value of "segment" and slopes of "extrapSegment"
1192 extrapCorrSegment->
1193 SetBendingCoor(extrapSegment->GetBendingCoor() +
1194 extrapSegment->GetBendingSlope() *
1195 (segment->GetHitForRec1()->GetZ() -
1196 (&(pMUON->Chamber(2 * station)))->Z()));
1197 extrapCorrSegment->
1198 SetNonBendingCoor(extrapSegment->GetNonBendingCoor() +
1199 extrapSegment->GetNonBendingSlope() *
1200 (segment->GetHitForRec1()->GetZ() -
1201 (&(pMUON->Chamber(2 * station)))->Z()));
1202 chi2 = segment->
1203 NormalizedChi2WithSegment(extrapCorrSegment, maxSigma2Distance);
a9e2aefa 1204 if (chi2 < bestChi2) {
1205 // update best Chi2 and Segment if better found
1206 bestSegment = segment;
1207 bestChi2 = chi2;
1208 }
1209 }
1210 if (bestSegment) {
1211 // best segment found: add it to track candidate
1212 track->AddSegment(bestSegment);
1213 // set track parameters at these two TrakHit's
1214 track->SetTrackParamAtHit(track->GetNTrackHits() - 2, &(trackParam[0]));
1215 track->SetTrackParamAtHit(track->GetNTrackHits() - 1, &(trackParam[1]));
1216 if (fPrintLevel >= 10) {
1217 cout << "FollowTracks: track candidate(0..): " << trackIndex
1218 << " Added segment in station(0..): " << station << endl;
1219 track->RecursiveDump();
1220 }
1221 }
1222 else {
1223 // No best segment found:
1224 // Look for best compatible HitForRec in station:
1225 // should consider all possibilities ????
1226 // multiple scattering ???? do about like for extrapSegment !!!!
79e1e601 1227 extrapHit = new AliMUONHitForRec(); // empty hit
d0bfce8d 1228 extrapCorrHit = new AliMUONHitForRec(); // empty corrected hit
a9e2aefa 1229 bestChi2 = 3.0;
1230 bestHit = NULL;
1231 if (fPrintLevel >= 10) {
1232 cout << "FollowTracks: track candidate(0..): " << trackIndex
1233 << " Segment not found, look for hit in station(0..): " << station
1234 << endl;
1235 }
1236 // Loop over chambers of the station
d82671a0 1237 for (chInStation = 0; chInStation < 2; chInStation++) {
a9e2aefa 1238 // coordinates of extrapolated hit
d82671a0 1239 extrapHit->
1240 SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor());
1241 extrapHit->
1242 SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor());
a9e2aefa 1243 // resolutions from "extrapSegment"
1244 extrapHit->SetBendingReso2(extrapSegment->GetBendingCoorReso2());
1245 extrapHit->SetNonBendingReso2(extrapSegment->GetNonBendingCoorReso2());
d0bfce8d 1246 // same things for corrected hit
1247 // better to use copy constructor, after checking that it works properly !!!!
1248 extrapCorrHit->
1249 SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor());
1250 extrapCorrHit->
1251 SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor());
1252 extrapHit->SetBendingReso2(extrapSegment->GetBendingCoorReso2());
1253 extrapHit->SetNonBendingReso2(extrapSegment->GetNonBendingCoorReso2());
a9e2aefa 1254 // Loop over hits in the chamber
956019b6 1255 ch = 2 * station + chInStation;
a9e2aefa 1256 for (iHit = fIndexOfFirstHitForRecPerChamber[ch];
1257 iHit < fIndexOfFirstHitForRecPerChamber[ch] +
1258 fNHitsForRecPerChamber[ch];
1259 iHit++) {
1260 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
d0bfce8d 1261 // correction of corrected hit (fBendingCoor and fNonBendingCoor)
1262 // according to real Z value of "hit" and slopes of right "trackParam"
1263 extrapCorrHit->
1264 SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor() +
1265 (&(trackParam[chInStation]))->GetBendingSlope() *
1266 (hit->GetZ() -
1267 (&(trackParam[chInStation]))->GetZ()));
1268 extrapCorrHit->
1269 SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor() +
1270 (&(trackParam[chInStation]))->GetNonBendingSlope() *
1271 (hit->GetZ() -
1272 (&(trackParam[chInStation]))->GetZ()));
a9e2aefa 1273 // condition for hit not already in segment ????
1274 chi2 = hit->NormalizedChi2WithHitForRec(extrapHit, maxSigma2Distance);
1275 if (chi2 < bestChi2) {
1276 // update best Chi2 and HitForRec if better found
1277 bestHit = hit;
1278 bestChi2 = chi2;
d82671a0 1279 chBestHit = chInStation;
a9e2aefa 1280 }
1281 }
1282 }
1283 if (bestHit) {
1284 // best hit found: add it to track candidate
1285 track->AddHitForRec(bestHit);
956019b6 1286 // set track parameters at this TrackHit
a9e2aefa 1287 track->SetTrackParamAtHit(track->GetNTrackHits() - 1,
1288 &(trackParam[chBestHit]));
1289 if (fPrintLevel >= 10) {
1290 cout << "FollowTracks: track candidate(0..): " << trackIndex
1291 << " Added hit in station(0..): " << station << endl;
1292 track->RecursiveDump();
1293 }
1294 }
1295 else {
8429a5e4 1296 // Remove current track candidate
1297 // and corresponding TrackHit's, ...
1298 track->Remove();
04b5ea16 1299 delete extrapSegment;
d0bfce8d 1300 delete extrapCorrSegment;
1301 delete extrapHit;
1302 delete extrapCorrHit;
a9e2aefa 1303 break; // stop the search for this candidate:
1304 // exit from the loop over station
1305 }
d0bfce8d 1306 delete extrapHit;
1307 delete extrapCorrHit;
a9e2aefa 1308 }
04b5ea16 1309 delete extrapSegment;
d0bfce8d 1310 delete extrapCorrSegment;
a9e2aefa 1311 // Sort track hits according to increasing Z
1312 track->GetTrackHitsPtr()->Sort();
1313 // Update track parameters at first track hit (smallest Z)
1314 trackParam1 = ((AliMUONTrackHit*)
1315 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
d0bfce8d 1316 bendingMomentum = 0.;
1317 if (TMath::Abs(trackParam1->GetInverseBendingMomentum()) > 0.)
1318 bendingMomentum = TMath::Abs(1/(trackParam1->GetInverseBendingMomentum()));
1319 // Track removed if bendingMomentum not in window [min, max]
1320 if ((bendingMomentum < fMinBendingMomentum) || (bendingMomentum > fMaxBendingMomentum)) {
1321 track->Remove();
1322 break; // stop the search for this candidate:
1323 // exit from the loop over station
1324 }
956019b6 1325 // Track fit
04b5ea16 1326 // with multiple Coulomb scattering if all stations
1327 if (station == 0) track->SetFitMCS(1);
956019b6 1328 // without multiple Coulomb scattering if not all stations
04b5ea16 1329 else track->SetFitMCS(0);
956019b6 1330 track->SetFitNParam(5); // with 5 parameters (momentum and position)
1331 track->SetFitStart(1); // from parameters at first hit
1332 track->Fit();
d0bfce8d 1333 Double_t numberOfDegFree = (2.0 * track->GetNTrackHits() - 5);
1334 if (numberOfDegFree > 0) {
1335 chi2Norm = track->GetFitFMin() / numberOfDegFree;
1336 } else {
1337 chi2Norm = 1.e10;
1338 }
1339 // Track removed if normalized chi2 too high
1340 if (chi2Norm > fMaxChi2) {
1341 track->Remove();
1342 break; // stop the search for this candidate:
1343 // exit from the loop over station
1344 }
a9e2aefa 1345 if (fPrintLevel >= 10) {
1346 cout << "FollowTracks: track candidate(0..): " << trackIndex
1347 << " after fit from station(0..): " << station << " to 4" << endl;
1348 track->RecursiveDump();
1349 }
04b5ea16 1350 // Track extrapolation to the vertex through the absorber (Branson)
956019b6 1351 // after going through the first station
04b5ea16 1352 if (station == 0) {
1353 trackParamVertex = *trackParam1;
1354 (&trackParamVertex)->ExtrapToVertex();
1355 track->SetTrackParamAtVertex(&trackParamVertex);
6ae236ee 1356 if (fPrintLevel >= 1) {
1357 cout << "FollowTracks: track candidate(0..): " << trackIndex
1358 << " after extrapolation to vertex" << endl;
1359 track->RecursiveDump();
1360 }
04b5ea16 1361 }
a9e2aefa 1362 } // for (station = 2;...
04b5ea16 1363 // go really to next track
1364 track = nextTrack;
1365 } // while (track)
1366 // Compression of track array (necessary after Remove ????)
1367 fRecTracksPtr->Compress();
1368 return;
1369}
1370
1371 //__________________________________________________________________________
8429a5e4 1372void AliMUONEventReconstructor::RemoveDoubleTracks(void)
1373{
1374 // To remove double tracks.
1375 // Tracks are considered identical
1376 // if they have at least half of their hits in common.
1377 // Among two identical tracks, one keeps the track with the larger number of hits
1378 // or, if these numbers are equal, the track with the minimum Chi2.
1379 AliMUONTrack *track1, *track2, *trackToRemove;
1380 Bool_t identicalTracks;
1381 Int_t hitsInCommon, nHits1, nHits2;
1382 identicalTracks = kTRUE;
1383 while (identicalTracks) {
1384 identicalTracks = kFALSE;
1385 // Loop over first track of the pair
1386 track1 = (AliMUONTrack*) fRecTracksPtr->First();
1387 while (track1 && (!identicalTracks)) {
1388 nHits1 = track1->GetNTrackHits();
1389 // Loop over second track of the pair
1390 track2 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1391 while (track2 && (!identicalTracks)) {
1392 nHits2 = track2->GetNTrackHits();
1393 // number of hits in common between two tracks
1394 hitsInCommon = track1->HitsInCommon(track2);
1395 // check for identical tracks
1396 if ((4 * hitsInCommon) >= (nHits1 + nHits2)) {
1397 identicalTracks = kTRUE;
1398 // decide which track to remove
1399 if (nHits1 > nHits2) trackToRemove = track2;
1400 else if (nHits1 < nHits2) trackToRemove = track1;
1401 else if ((track1->GetFitFMin()) < (track2->GetFitFMin()))
1402 trackToRemove = track2;
1403 else trackToRemove = track1;
1404 // remove it
1405 trackToRemove->Remove();
1406 }
1407 track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
1408 } // track2
1409 track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1410 } // track1
1411 }
1412 return;
1413}
1414
1415 //__________________________________________________________________________
04b5ea16 1416void AliMUONEventReconstructor::EventDump(void)
1417{
1418 // Dump reconstructed event (track parameters at vertex and at first hit),
1419 // and the particle parameters
1420
1421 AliMUONTrack *track;
1422 AliMUONTrackParam *trackParam, *trackParam1;
04b5ea16 1423 TParticle *p;
1424 Double_t bendingSlope, nonBendingSlope, pYZ;
1425 Double_t pX, pY, pZ, x, y, z, c;
1426 Int_t np, trackIndex, nTrackHits;
1427
1428 if (fPrintLevel >= 1) cout << "****** enter EventDump ******" << endl;
1429 if (fPrintLevel >= 1) {
1430 cout << " Number of Reconstructed tracks :" << fNRecTracks << endl;
1431 }
1432 fRecTracksPtr->Compress(); // for simple loop without "Next" since no hole
1433 // Loop over reconstructed tracks
1434 for (trackIndex = 0; trackIndex < fNRecTracks; trackIndex++) {
1435 if (fPrintLevel >= 1)
1436 cout << " track number: " << trackIndex << endl;
1437 // function for each track for modularity ????
1438 track = (AliMUONTrack*) ((*fRecTracksPtr)[trackIndex]);
1439 nTrackHits = track->GetNTrackHits();
1440 if (fPrintLevel >= 1)
1441 cout << " number of track hits: " << nTrackHits << endl;
1442 // track parameters at Vertex
1443 trackParam = track->GetTrackParamAtVertex();
1444 x = trackParam->GetNonBendingCoor();
1445 y = trackParam->GetBendingCoor();
1446 z = trackParam->GetZ();
1447 bendingSlope = trackParam->GetBendingSlope();
1448 nonBendingSlope = trackParam->GetNonBendingSlope();
1449 pYZ = 1/TMath::Abs(trackParam->GetInverseBendingMomentum());
1450 pZ = pYZ/TMath::Sqrt(1+bendingSlope*bendingSlope);
1451 pX = pZ * nonBendingSlope;
1452 pY = pZ * bendingSlope;
1453 c = TMath::Sign(1.0, trackParam->GetInverseBendingMomentum());
1454 if (fPrintLevel >= 1)
1455 printf(" track parameters at Vertex z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1456 z, x, y, pX, pY, pZ, c);
1457
1458 // track parameters at first hit
1459 trackParam1 = ((AliMUONTrackHit*)
1460 (track->GetTrackHitsPtr()->First()))->GetTrackParam();
1461 x = trackParam1->GetNonBendingCoor();
1462 y = trackParam1->GetBendingCoor();
1463 z = trackParam1->GetZ();
1464 bendingSlope = trackParam1->GetBendingSlope();
1465 nonBendingSlope = trackParam1->GetNonBendingSlope();
1466 pYZ = 1/TMath::Abs(trackParam1->GetInverseBendingMomentum());
1467 pZ = pYZ/TMath::Sqrt(1.0 + bendingSlope * bendingSlope);
1468 pX = pZ * nonBendingSlope;
1469 pY = pZ * bendingSlope;
1470 c = TMath::Sign(1.0, trackParam1->GetInverseBendingMomentum());
1471 if (fPrintLevel >= 1)
1472 printf(" track parameters at z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1473 z, x, y, pX, pY, pZ, c);
1474 }
1475 // informations about generated particles
2ab0c725 1476 np = gAlice->GetNtrack();
04b5ea16 1477 printf(" **** number of generated particles: %d \n", np);
1478
1479 for (Int_t iPart = 0; iPart < np; iPart++) {
2ab0c725 1480 p = gAlice->Particle(iPart);
04b5ea16 1481 printf(" particle %d: type= %d px= %f py= %f pz= %f pdg= %d\n",
1482 iPart, p->GetPdgCode(), p->Px(), p->Py(), p->Pz(), p->GetPdgCode());
1483 }
a9e2aefa 1484 return;
1485}
1486
c7ba256d 1487void AliMUONEventReconstructor::FillEvent()
1488{
1489// Create a new AliMUONRecoEvent, fill its track list, then add it as a
1490// leaf in the Event branch of TreeRecoEvent tree
1491 cout << "Enter FillEvent() ...\n";
1492
1493 if (!fRecoEvent) {
1494 fRecoEvent = new AliMUONRecoEvent();
1495 } else {
1496 fRecoEvent->Clear();
1497 }
1498 //save current directory
1499 TDirectory *current = gDirectory;
1500 if (!fTreeFile) fTreeFile = new TFile("tree_reco.root", "RECREATE");
1501 if (!fEventTree) fEventTree = new TTree("TreeRecoEvent", "MUON reconstructed events");
1502 if (fRecoEvent->MakeDumpTracks(fRecTracksPtr)) {
1503 if (fPrintLevel > 1) fRecoEvent->EventInfo();
1504 TBranch *branch = fEventTree->GetBranch("Event");
1505 if (!branch) branch = fEventTree->Branch("Event", "AliMUONRecoEvent", &fRecoEvent, 64000,1);
1506 branch->SetAutoDelete();
1507 fTreeFile->cd();
1508 fEventTree->Fill();
1509 fTreeFile->Write();
1510 }
1511 // restore directory
1512 current->cd();
1513}