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