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