]>
Commit | Line | Data |
---|---|---|
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 | ||
cc87ebcd | 16 | /* $Id$ */ |
17 | ||
3831f268 | 18 | //////////////////////////////////// |
a9e2aefa | 19 | // |
29f1b13a | 20 | // MUON track reconstructor in ALICE (class renamed from AliMUONEventReconstructor) |
a9e2aefa | 21 | // |
22 | // This class contains as data: | |
29f1b13a | 23 | // * the parameters for the track reconstruction |
a9e2aefa | 24 | // * a pointer to the array of hits to be reconstructed (the event) |
25 | // * a pointer to the array of segments made with these hits inside each station | |
26 | // * a pointer to the array of reconstructed tracks | |
27 | // | |
28 | // It contains as methods, among others: | |
29 | // * MakeEventToBeReconstructed to build the array of hits to be reconstructed | |
30 | // * MakeSegments to build the segments | |
31 | // * MakeTracks to build the tracks | |
3831f268 | 32 | // |
33 | //////////////////////////////////// | |
a9e2aefa | 34 | |
ae17f568 | 35 | #include <stdlib.h> |
30178c30 | 36 | #include <Riostream.h> |
37 | #include <TDirectory.h> | |
38 | #include <TFile.h> | |
cc87ebcd | 39 | #include <TMatrixD.h> |
a9e2aefa | 40 | |
29f1b13a | 41 | #include "AliMUONTrackReconstructor.h" |
30178c30 | 42 | #include "AliMUON.h" |
29fc2c86 | 43 | #include "AliMUONConstants.h" |
a9e2aefa | 44 | #include "AliMUONHitForRec.h" |
276c44b7 | 45 | #include "AliMUONTriggerTrack.h" |
276c44b7 | 46 | #include "AliMUONTriggerCircuit.h" |
bc4a7ff4 | 47 | #include "AliMUONTriggerCircuitNew.h" |
a9e2aefa | 48 | #include "AliMUONRawCluster.h" |
276c44b7 | 49 | #include "AliMUONLocalTrigger.h" |
52c9bc11 | 50 | #include "AliMUONGlobalTrigger.h" |
3831f268 | 51 | #include "AliMUONSegment.h" |
a9e2aefa | 52 | #include "AliMUONTrack.h" |
a9e2aefa | 53 | #include "AliMUONTrackHit.h" |
94de3818 | 54 | #include "AliMagF.h" |
3831f268 | 55 | #include "AliRun.h" // for gAlice |
88cb7938 | 56 | #include "AliRunLoader.h" |
57 | #include "AliLoader.h" | |
cc87ebcd | 58 | #include "AliMUONTrackK.h" |
5d12ce38 | 59 | #include "AliMC.h" |
8c343c7c | 60 | #include "AliLog.h" |
29fc2c86 | 61 | #include "AliTrackReference.h" |
a9e2aefa | 62 | |
a9e2aefa | 63 | //************* Defaults parameters for reconstruction |
29f1b13a | 64 | const Double_t AliMUONTrackReconstructor::fgkDefaultMinBendingMomentum = 3.0; |
7b96283c | 65 | const Double_t AliMUONTrackReconstructor::fgkDefaultMaxBendingMomentum = 3000.0; |
29f1b13a | 66 | const Double_t AliMUONTrackReconstructor::fgkDefaultMaxChi2 = 100.0; |
67 | const Double_t AliMUONTrackReconstructor::fgkDefaultMaxSigma2Distance = 16.0; | |
68 | const Double_t AliMUONTrackReconstructor::fgkDefaultBendingResolution = 0.01; | |
69 | const Double_t AliMUONTrackReconstructor::fgkDefaultNonBendingResolution = 0.144; | |
70 | const Double_t AliMUONTrackReconstructor::fgkDefaultChamberThicknessInX0 = 0.03; | |
a9e2aefa | 71 | // Simple magnetic field: |
72 | // Value taken from macro MUONtracking.C: 0.7 T, hence 7 kG | |
73 | // Length and Position from reco_muon.F, with opposite sign: | |
74 | // Length = ZMAGEND-ZCOIL | |
75 | // Position = (ZMAGEND+ZCOIL)/2 | |
76 | // to be ajusted differently from real magnetic field ???? | |
29f1b13a | 77 | const Double_t AliMUONTrackReconstructor::fgkDefaultSimpleBValue = 7.0; |
78 | const Double_t AliMUONTrackReconstructor::fgkDefaultSimpleBLength = 428.0; | |
79 | const Double_t AliMUONTrackReconstructor::fgkDefaultSimpleBPosition = 1019.0; | |
80 | const Int_t AliMUONTrackReconstructor::fgkDefaultRecTrackRefHits = 0; | |
81 | const Double_t AliMUONTrackReconstructor::fgkDefaultEfficiency = 0.95; | |
a9e2aefa | 82 | |
29f1b13a | 83 | ClassImp(AliMUONTrackReconstructor) // Class implementation in ROOT context |
a9e2aefa | 84 | |
3bc8b580 | 85 | //__________________________________________________________________________ |
86 | AliMUONTrackReconstructor::AliMUONTrackReconstructor(AliLoader* loader, AliMUONData* data) | |
54d7ba50 | 87 | : TObject(), |
88 | fTrackMethod(1), //AZ - tracking method (1-default, 2-Kalman) | |
89 | fMinBendingMomentum(fgkDefaultMinBendingMomentum), | |
90 | fMaxBendingMomentum(fgkDefaultMaxBendingMomentum), | |
91 | fMaxChi2(fgkDefaultMaxChi2), | |
92 | fMaxSigma2Distance(fgkDefaultMaxSigma2Distance), | |
93 | fBendingResolution(fgkDefaultBendingResolution), | |
94 | fNonBendingResolution(fgkDefaultNonBendingResolution), | |
95 | fChamberThicknessInX0(fgkDefaultChamberThicknessInX0), | |
96 | fSimpleBValue(fgkDefaultSimpleBValue), | |
97 | fSimpleBLength(fgkDefaultSimpleBLength), | |
98 | fSimpleBPosition(fgkDefaultSimpleBPosition), | |
99 | fRecTrackRefHits(fgkDefaultRecTrackRefHits), | |
100 | fEfficiency(fgkDefaultEfficiency), | |
101 | // Initializions for track ref. background events | |
102 | fBkgTrackRefFile(0x0), | |
103 | fBkgTrackRefTK(0x0), | |
104 | fBkgTrackRefParticles(0x0), | |
105 | fBkgTrackRefTTR(0x0), | |
106 | fBkgTrackRefEventNumber(-1), | |
107 | fHitsForRecPtr(0x0), | |
108 | fNHitsForRec(0), | |
109 | fRecTracksPtr(0x0), | |
110 | fNRecTracks(0), | |
111 | fRecTrackHitsPtr(0x0), | |
112 | fNRecTrackHits(0), | |
113 | fMUONData(data), | |
114 | fLoader(loader), | |
115 | fMuons(0) | |
a9e2aefa | 116 | { |
29f1b13a | 117 | // Constructor for class AliMUONTrackReconstructor |
a9e2aefa | 118 | SetReconstructionParametersToDefaults(); |
54d7ba50 | 119 | |
a9e2aefa | 120 | // Memory allocation for the TClonesArray of hits for reconstruction |
121 | // Is 10000 the right size ???? | |
122 | fHitsForRecPtr = new TClonesArray("AliMUONHitForRec", 10000); | |
54d7ba50 | 123 | |
a9e2aefa | 124 | // Memory allocation for the TClonesArray's of segments in stations |
125 | // Is 2000 the right size ???? | |
190f97ea | 126 | for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) { |
a9e2aefa | 127 | fSegmentsPtr[st] = new TClonesArray("AliMUONSegment", 2000); |
128 | fNSegments[st] = 0; // really needed or GetEntriesFast sufficient ???? | |
129 | } | |
130 | // Memory allocation for the TClonesArray of reconstructed tracks | |
131 | // Is 10 the right size ???? | |
132 | fRecTracksPtr = new TClonesArray("AliMUONTrack", 10); | |
54d7ba50 | 133 | |
8429a5e4 | 134 | // Memory allocation for the TClonesArray of hits on reconstructed tracks |
135 | // Is 100 the right size ???? | |
b035a148 | 136 | fRecTrackHitsPtr = new TClonesArray("AliMUONTrackHit", 100); |
a9e2aefa | 137 | |
5b64e914 | 138 | // Sign of fSimpleBValue according to sign of Bx value at (50,50,-950). |
a6f03ddb | 139 | Float_t b[3], x[3]; |
5b64e914 | 140 | x[0] = 50.; x[1] = 50.; x[2] = -950.; |
a6f03ddb | 141 | gAlice->Field()->Field(x, b); |
065bd0e1 | 142 | fSimpleBValue = TMath::Sign(fSimpleBValue,(Double_t) b[0]); |
5b64e914 | 143 | fSimpleBPosition = TMath::Sign(fSimpleBPosition,(Double_t) x[2]); |
a6f03ddb | 144 | // See how to get fSimple(BValue, BLength, BPosition) |
145 | // automatically calculated from the actual magnetic field ???? | |
a9e2aefa | 146 | |
29f1b13a | 147 | AliDebug(1,"AliMUONTrackReconstructor constructed with defaults"); |
b840996b | 148 | if ( AliLog::GetGlobalDebugLevel()>0) Dump(); |
149 | AliDebug(1,"Magnetic field from root file:"); | |
150 | if ( AliLog::GetGlobalDebugLevel()>0) gAlice->Field()->Dump(); | |
151 | ||
52c9bc11 | 152 | // initialize container |
3bc8b580 | 153 | // fMUONData = new AliMUONData(fLoader,"MUON","MUON"); |
54d7ba50 | 154 | //fMUONData = data; |
52c9bc11 | 155 | |
a9e2aefa | 156 | return; |
157 | } | |
9cbdf048 | 158 | |
a9e2aefa | 159 | //__________________________________________________________________________ |
29f1b13a | 160 | AliMUONTrackReconstructor::~AliMUONTrackReconstructor(void) |
a9e2aefa | 161 | { |
29f1b13a | 162 | // Destructor for class AliMUONTrackReconstructor |
a9e2aefa | 163 | delete fHitsForRecPtr; // Correct destruction of everything ???? or delete [] ???? |
190f97ea | 164 | for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) |
a9e2aefa | 165 | delete fSegmentsPtr[st]; // Correct destruction of everything ???? |
166 | return; | |
167 | } | |
168 | ||
169 | //__________________________________________________________________________ | |
29f1b13a | 170 | void AliMUONTrackReconstructor::SetReconstructionParametersToDefaults(void) |
a9e2aefa | 171 | { |
172 | // Set reconstruction parameters to default values | |
173 | // Would be much more convenient with a structure (or class) ???? | |
a9e2aefa | 174 | |
a9e2aefa | 175 | // ******** Parameters for making HitsForRec |
176 | // minimum radius, | |
177 | // like in TRACKF_STAT: | |
178 | // 2 degrees for stations 1 and 2, or ch(0...) from 0 to 3; | |
179 | // 30 cm for stations 3 to 5, or ch(0...) from 4 to 9 | |
190f97ea | 180 | for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) { |
06c11448 | 181 | if (ch < 4) fRMin[ch] = TMath::Abs(AliMUONConstants::DefaultChamberZ(ch)) * |
a9e2aefa | 182 | 2.0 * TMath::Pi() / 180.0; |
183 | else fRMin[ch] = 30.0; | |
d0bfce8d | 184 | // maximum radius at 10 degrees and Z of chamber |
06c11448 | 185 | fRMax[ch] = TMath::Abs(AliMUONConstants::DefaultChamberZ(ch)) * |
d0bfce8d | 186 | 10.0 * TMath::Pi() / 180.0; |
a9e2aefa | 187 | } |
a9e2aefa | 188 | |
189 | // ******** Parameters for making segments | |
190 | // should be parametrized ???? | |
191 | // according to interval between chambers in a station ???? | |
192 | // Maximum distance in non bending plane | |
193 | // 5 * 0.22 just to remember the way it was made in TRACKF_STAT | |
194 | // SIGCUT*DYMAX(IZ) | |
190f97ea | 195 | for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) |
a9e2aefa | 196 | fSegmentMaxDistNonBending[st] = 5. * 0.22; |
860cef81 | 197 | // Maximum distance in bending plane: |
198 | // values from TRACKF_STAT, corresponding to (J psi 20cm), | |
199 | // scaled to the real distance between chambers in a station | |
ef21c79e | 200 | fSegmentMaxDistBending[0] = TMath::Abs( 1.5 * |
06c11448 | 201 | (AliMUONConstants::DefaultChamberZ(1) - AliMUONConstants::DefaultChamberZ(0)) / 20.0); |
e516b01d | 202 | fSegmentMaxDistBending[1] = TMath::Abs( 1.5 * |
06c11448 | 203 | (AliMUONConstants::DefaultChamberZ(3) - AliMUONConstants::DefaultChamberZ(2)) / 20.0); |
e516b01d | 204 | fSegmentMaxDistBending[2] = TMath::Abs( 3.0 * |
06c11448 | 205 | (AliMUONConstants::DefaultChamberZ(5) - AliMUONConstants::DefaultChamberZ(4)) / 20.0); |
e516b01d | 206 | fSegmentMaxDistBending[3] = TMath::Abs( 6.0 * |
06c11448 | 207 | (AliMUONConstants::DefaultChamberZ(7) - AliMUONConstants::DefaultChamberZ(6)) / 20.0); |
e516b01d | 208 | fSegmentMaxDistBending[4] = TMath::Abs( 6.0 * |
06c11448 | 209 | (AliMUONConstants::DefaultChamberZ(9) - AliMUONConstants::DefaultChamberZ(8)) / 20.0); |
210 | ||
a9e2aefa | 211 | return; |
212 | } | |
213 | ||
214 | //__________________________________________________________________________ | |
29f1b13a | 215 | Double_t AliMUONTrackReconstructor::GetImpactParamFromBendingMomentum(Double_t BendingMomentum) const |
a9e2aefa | 216 | { |
217 | // Returns impact parameter at vertex in bending plane (cm), | |
218 | // from the signed bending momentum "BendingMomentum" in bending plane (GeV/c), | |
219 | // using simple values for dipole magnetic field. | |
065bd0e1 | 220 | // The sign of "BendingMomentum" is the sign of the charge. |
a9e2aefa | 221 | return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition / |
222 | BendingMomentum); | |
223 | } | |
224 | ||
225 | //__________________________________________________________________________ | |
29f1b13a | 226 | Double_t AliMUONTrackReconstructor::GetBendingMomentumFromImpactParam(Double_t ImpactParam) const |
a9e2aefa | 227 | { |
228 | // Returns signed bending momentum in bending plane (GeV/c), | |
065bd0e1 | 229 | // the sign being the sign of the charge for particles moving forward in Z, |
a9e2aefa | 230 | // from the impact parameter "ImpactParam" at vertex in bending plane (cm), |
231 | // using simple values for dipole magnetic field. | |
a9e2aefa | 232 | return (-0.0003 * fSimpleBValue * fSimpleBLength * fSimpleBPosition / |
233 | ImpactParam); | |
234 | } | |
235 | ||
236 | //__________________________________________________________________________ | |
29f1b13a | 237 | void AliMUONTrackReconstructor::SetBkgTrackRefFile(Text_t *BkgTrackRefFileName) |
a9e2aefa | 238 | { |
29fc2c86 | 239 | // Set background file ... for track ref. hits |
a9e2aefa | 240 | // Must be called after having loaded the firts signal event |
b840996b | 241 | AliDebug(1,Form("Enter SetBkgTrackRefFile with BkgTrackRefFileName %s",BkgTrackRefFileName)); |
242 | ||
29fc2c86 | 243 | if (strlen(BkgTrackRefFileName)) { |
244 | // BkgTrackRefFileName not empty: try to open the file | |
b840996b | 245 | |
246 | if(AliLog::GetGlobalDebugLevel()>1) { | |
247 | cout << "Before File(Bkg)" << endl; gDirectory->Dump(); | |
248 | } | |
29fc2c86 | 249 | fBkgTrackRefFile = new TFile(BkgTrackRefFileName); |
b840996b | 250 | if(AliLog::GetGlobalDebugLevel()>1) { |
251 | cout << "After File(Bkg)" << endl; gDirectory->Dump(); | |
252 | } | |
29fc2c86 | 253 | if (fBkgTrackRefFile-> IsOpen()) { |
b840996b | 254 | if(AliLog::GetGlobalDebugLevel()>0) { |
29fc2c86 | 255 | cout << "Background for Track ref. hits in file: ``" << BkgTrackRefFileName |
a9e2aefa | 256 | << "'' successfully opened" << endl;} |
257 | } | |
258 | else { | |
29fc2c86 | 259 | cout << "Background for Track Ref. hits in file: " << BkgTrackRefFileName << endl; |
a9e2aefa | 260 | cout << "NOT FOUND: EXIT" << endl; |
261 | exit(0); // right instruction for exit ???? | |
262 | } | |
263 | // Arrays for "particles" and "hits" | |
29fc2c86 | 264 | fBkgTrackRefParticles = new TClonesArray("TParticle", 200); |
a9e2aefa | 265 | // Event number to -1 for initialization |
29fc2c86 | 266 | fBkgTrackRefEventNumber = -1; |
a9e2aefa | 267 | // Back to the signal file: |
268 | // first signal event must have been loaded previously, | |
269 | // otherwise, Segmentation violation at the next instruction | |
270 | // How is it possible to do smething better ???? | |
271 | ((gAlice->TreeK())->GetCurrentFile())->cd(); | |
b840996b | 272 | if(AliLog::GetGlobalDebugLevel()>1) cout << "After cd(gAlice)" << endl; gDirectory->Dump(); |
a9e2aefa | 273 | } |
274 | return; | |
275 | } | |
276 | ||
277 | //__________________________________________________________________________ | |
29f1b13a | 278 | void AliMUONTrackReconstructor::NextBkgTrackRefEvent(void) |
a9e2aefa | 279 | { |
29fc2c86 | 280 | // Get next event in background file for track ref. hits |
a9e2aefa | 281 | // Goes back to event number 0 when end of file is reached |
282 | char treeName[20]; | |
b840996b | 283 | |
284 | AliDebug(1,"Enter NextBkgTrackRefEvent"); | |
a9e2aefa | 285 | // Clean previous event |
29fc2c86 | 286 | if(fBkgTrackRefTK) delete fBkgTrackRefTK; |
287 | fBkgTrackRefTK = NULL; | |
288 | if(fBkgTrackRefParticles) fBkgTrackRefParticles->Clear(); | |
289 | if(fBkgTrackRefTTR) delete fBkgTrackRefTTR; | |
290 | fBkgTrackRefTTR = NULL; | |
a9e2aefa | 291 | // Increment event number |
29fc2c86 | 292 | fBkgTrackRefEventNumber++; |
a9e2aefa | 293 | // Get access to Particles and Hits for event from background file |
b840996b | 294 | if (AliLog::GetGlobalDebugLevel()>1) cout << "Before cd(Bkg)" << endl; gDirectory->Dump(); |
29fc2c86 | 295 | fBkgTrackRefFile->cd(); |
b840996b | 296 | if (AliLog::GetGlobalDebugLevel()>1) cout << "After cd(Bkg)" << endl; gDirectory->Dump(); |
a9e2aefa | 297 | // Particles: TreeK for event and branch "Particles" |
29fc2c86 | 298 | sprintf(treeName, "TreeK%d", fBkgTrackRefEventNumber); |
299 | fBkgTrackRefTK = (TTree*)gDirectory->Get(treeName); | |
300 | if (!fBkgTrackRefTK) { | |
b840996b | 301 | |
302 | AliDebug(1,Form("Cannot find Kine Tree for background event: %d",fBkgTrackRefEventNumber)); | |
303 | AliDebug(1,"Goes back to event 0"); | |
304 | ||
29fc2c86 | 305 | fBkgTrackRefEventNumber = 0; |
306 | sprintf(treeName, "TreeK%d", fBkgTrackRefEventNumber); | |
307 | fBkgTrackRefTK = (TTree*)gDirectory->Get(treeName); | |
308 | if (!fBkgTrackRefTK) { | |
309 | AliError(Form("cannot find Kine Tree for background event: %d",fBkgTrackRefEventNumber)); | |
8c343c7c | 310 | exit(0); |
a9e2aefa | 311 | } |
312 | } | |
29fc2c86 | 313 | if (fBkgTrackRefTK) |
314 | fBkgTrackRefTK->SetBranchAddress("Particles", &fBkgTrackRefParticles); | |
315 | fBkgTrackRefTK->GetEvent(0); // why event 0 ???? necessary ???? | |
a9e2aefa | 316 | // Hits: TreeH for event and branch "MUON" |
29fc2c86 | 317 | sprintf(treeName, "TreeTR%d", fBkgTrackRefEventNumber); |
318 | fBkgTrackRefTTR = (TTree*)gDirectory->Get(treeName); | |
319 | if (!fBkgTrackRefTTR) { | |
320 | AliError(Form("cannot find Hits Tree for background event: %d",fBkgTrackRefEventNumber)); | |
a9e2aefa | 321 | exit(0); |
322 | } | |
29fc2c86 | 323 | fBkgTrackRefTTR->GetEntries(); // necessary ???? |
a9e2aefa | 324 | // Back to the signal file |
325 | ((gAlice->TreeK())->GetCurrentFile())->cd(); | |
b840996b | 326 | if (AliLog::GetGlobalDebugLevel()>1) |
327 | cout << "After cd(gAlice)" << endl; gDirectory->Dump(); | |
a9e2aefa | 328 | return; |
329 | } | |
330 | ||
331 | //__________________________________________________________________________ | |
29f1b13a | 332 | void AliMUONTrackReconstructor::EventReconstruct(void) |
a9e2aefa | 333 | { |
334 | // To reconstruct one event | |
b840996b | 335 | AliDebug(1,"Enter EventReconstruct"); |
fff097e0 | 336 | ResetTracks(); //AZ |
337 | ResetTrackHits(); //AZ | |
338 | ResetSegments(); //AZ | |
339 | ResetHitsForRec(); //AZ | |
a9e2aefa | 340 | MakeEventToBeReconstructed(); |
341 | MakeSegments(); | |
342 | MakeTracks(); | |
d837040f | 343 | if (fMUONData->IsTriggerTrackBranchesInTree()) |
344 | ValidateTracksWithTrigger(); | |
345 | ||
346 | // Add tracks to MUON data container | |
347 | for(Int_t i=0; i<GetNRecTracks(); i++) { | |
348 | AliMUONTrack * track = (AliMUONTrack*) GetRecTracksPtr()->At(i); | |
349 | fMUONData->AddRecTrack(*track); | |
350 | } | |
351 | ||
a9e2aefa | 352 | return; |
353 | } | |
354 | ||
276c44b7 | 355 | //__________________________________________________________________________ |
29f1b13a | 356 | void AliMUONTrackReconstructor::EventReconstructTrigger(void) |
276c44b7 | 357 | { |
358 | // To reconstruct one event | |
b840996b | 359 | AliDebug(1,"Enter EventReconstructTrigger"); |
276c44b7 | 360 | MakeTriggerTracks(); |
361 | return; | |
362 | } | |
363 | ||
a9e2aefa | 364 | //__________________________________________________________________________ |
29f1b13a | 365 | void AliMUONTrackReconstructor::ResetHitsForRec(void) |
a9e2aefa | 366 | { |
367 | // To reset the array and the number of HitsForRec, | |
368 | // and also the number of HitsForRec | |
369 | // and the index of the first HitForRec per chamber | |
0ce28091 | 370 | if (fHitsForRecPtr) fHitsForRecPtr->Delete(); |
a9e2aefa | 371 | fNHitsForRec = 0; |
190f97ea | 372 | for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) |
a9e2aefa | 373 | fNHitsForRecPerChamber[ch] = fIndexOfFirstHitForRecPerChamber[ch] = 0; |
374 | return; | |
375 | } | |
376 | ||
377 | //__________________________________________________________________________ | |
29f1b13a | 378 | void AliMUONTrackReconstructor::ResetSegments(void) |
a9e2aefa | 379 | { |
380 | // To reset the TClonesArray of segments and the number of Segments | |
381 | // for all stations | |
190f97ea | 382 | for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) { |
a9e2aefa | 383 | if (fSegmentsPtr[st]) fSegmentsPtr[st]->Clear(); |
384 | fNSegments[st] = 0; | |
385 | } | |
386 | return; | |
387 | } | |
388 | ||
389 | //__________________________________________________________________________ | |
29f1b13a | 390 | void AliMUONTrackReconstructor::ResetTracks(void) |
a9e2aefa | 391 | { |
392 | // To reset the TClonesArray of reconstructed tracks | |
8429a5e4 | 393 | if (fRecTracksPtr) fRecTracksPtr->Delete(); |
394 | // Delete in order that the Track destructors are called, | |
395 | // hence the space for the TClonesArray of pointers to TrackHit's is freed | |
a9e2aefa | 396 | fNRecTracks = 0; |
397 | return; | |
398 | } | |
399 | ||
276c44b7 | 400 | |
8429a5e4 | 401 | //__________________________________________________________________________ |
29f1b13a | 402 | void AliMUONTrackReconstructor::ResetTrackHits(void) |
8429a5e4 | 403 | { |
404 | // To reset the TClonesArray of hits on reconstructed tracks | |
405 | if (fRecTrackHitsPtr) fRecTrackHitsPtr->Clear(); | |
406 | fNRecTrackHits = 0; | |
407 | return; | |
408 | } | |
409 | ||
a9e2aefa | 410 | //__________________________________________________________________________ |
29f1b13a | 411 | void AliMUONTrackReconstructor::MakeEventToBeReconstructed(void) |
a9e2aefa | 412 | { |
413 | // To make the list of hits to be reconstructed, | |
29fc2c86 | 414 | // either from the track ref. hits or from the raw clusters |
a9e2aefa | 415 | // according to the parameter set for the reconstructor |
e191bb57 | 416 | // TString evfoldname = AliConfig::GetDefaultEventFolderName();//to be interfaced properly |
88cb7938 | 417 | |
52c9bc11 | 418 | // AliRunLoader* rl = AliRunLoader::GetRunLoader(evfoldname); |
419 | // if (rl == 0x0) | |
420 | // { | |
421 | // Error("MakeEventToBeReconstructed", | |
422 | // "Can not find Run Loader in Event Folder named %s.", | |
423 | // evfoldname.Data()); | |
424 | // return; | |
425 | // } | |
426 | // AliLoader* gime = rl->GetLoader("MUONLoader"); | |
427 | // if (gime == 0x0) | |
428 | // { | |
429 | // Error("MakeEventToBeReconstructed","Can not get MUON Loader from Run Loader."); | |
430 | // return; | |
431 | // } | |
29fc2c86 | 432 | AliRunLoader *runLoader = fLoader->GetRunLoader(); |
433 | ||
b840996b | 434 | AliDebug(1,"Enter MakeEventToBeReconstructed"); |
fff097e0 | 435 | //AZ ResetHitsForRec(); |
29fc2c86 | 436 | if (fRecTrackRefHits == 1) { |
437 | // Reconstruction from track ref. hits | |
a9e2aefa | 438 | // Back to the signal file |
29fc2c86 | 439 | TTree* treeTR = runLoader->TreeTR(); |
440 | if (treeTR == 0x0) | |
441 | { | |
442 | Int_t retval = runLoader->LoadTrackRefs("READ"); | |
443 | if ( retval) | |
88cb7938 | 444 | { |
8c343c7c | 445 | AliError("Error occured while loading hits."); |
88cb7938 | 446 | return; |
447 | } | |
29fc2c86 | 448 | treeTR = runLoader->TreeTR(); |
449 | if (treeTR == 0x0) | |
88cb7938 | 450 | { |
29fc2c86 | 451 | AliError("Can not get TreeTR"); |
452 | return; | |
88cb7938 | 453 | } |
29fc2c86 | 454 | } |
455 | ||
456 | AddHitsForRecFromTrackRef(treeTR,1); | |
88cb7938 | 457 | |
a9e2aefa | 458 | // Background hits |
29fc2c86 | 459 | AddHitsForRecFromTrackRef(fBkgTrackRefTTR,0); |
a9e2aefa | 460 | // Sort HitsForRec in increasing order with respect to chamber number |
461 | SortHitsForRecWithIncreasingChamber(); | |
462 | } | |
463 | else { | |
464 | // Reconstruction from raw clusters | |
465 | // AliMUON *MUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ???? | |
466 | // Security on MUON ???? | |
467 | // TreeR assumed to be be "prepared" in calling function | |
468 | // by "MUON->GetTreeR(nev)" ???? | |
52c9bc11 | 469 | TTree *treeR = fLoader->TreeR(); |
cc87ebcd | 470 | //AZ? fMUONData->SetTreeAddress("RC"); |
9cbdf048 | 471 | AddHitsForRecFromRawClusters(treeR); |
a9e2aefa | 472 | // No sorting: it is done automatically in the previous function |
473 | } | |
b840996b | 474 | |
475 | AliDebug(1,"End of MakeEventToBeReconstructed"); | |
476 | if (AliLog::GetGlobalDebugLevel() > 0) { | |
477 | AliDebug(1, Form("NHitsForRec: %d",fNHitsForRec)); | |
190f97ea | 478 | for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) { |
b840996b | 479 | AliDebug(1, Form("Chamber(0...): %d",ch)); |
480 | AliDebug(1, Form("NHitsForRec: %d", fNHitsForRecPerChamber[ch])); | |
481 | AliDebug(1, Form("Index(first HitForRec): %d", fIndexOfFirstHitForRecPerChamber[ch])); | |
482 | for (Int_t hit = fIndexOfFirstHitForRecPerChamber[ch]; | |
483 | hit < fIndexOfFirstHitForRecPerChamber[ch] + fNHitsForRecPerChamber[ch]; | |
484 | hit++) { | |
485 | AliDebug(1, Form("HitForRec index(0...): %d",hit)); | |
486 | ((*fHitsForRecPtr)[hit])->Dump(); | |
a9e2aefa | 487 | } |
488 | } | |
489 | } | |
490 | return; | |
491 | } | |
492 | ||
493 | //__________________________________________________________________________ | |
29f1b13a | 494 | void AliMUONTrackReconstructor::AddHitsForRecFromTrackRef(TTree *TTR, Int_t signal) |
a9e2aefa | 495 | { |
496 | // To add to the list of hits for reconstruction | |
29fc2c86 | 497 | // the signal hits from a track reference tree TreeTR. |
498 | TClonesArray *listOfTrackRefs = NULL; | |
499 | AliTrackReference *trackRef; | |
500 | ||
501 | Int_t track = 0; | |
b840996b | 502 | AliDebug(2,Form("Enter AddHitsForRecFromTrackRef with TTR: %d", TTR)); |
29fc2c86 | 503 | if (TTR == NULL) return; |
504 | ||
505 | listOfTrackRefs = CleanTrackRefs(TTR); | |
1391e633 | 506 | |
29fc2c86 | 507 | Int_t ntracks = listOfTrackRefs->GetEntriesFast(); |
a9e2aefa | 508 | |
b840996b | 509 | AliDebug(2,Form("ntracks: %d", ntracks)); |
29fc2c86 | 510 | |
511 | for (Int_t index = 0; index < ntracks; index++) { | |
512 | trackRef = (AliTrackReference*) listOfTrackRefs->At(index); | |
513 | track = trackRef->GetTrack(); | |
514 | ||
515 | NewHitForRecFromTrackRef(trackRef,track,signal); | |
516 | } // end of track ref. | |
517 | ||
518 | listOfTrackRefs->Delete(); | |
519 | delete listOfTrackRefs; | |
a9e2aefa | 520 | return; |
521 | } | |
522 | ||
29fc2c86 | 523 | |
a9e2aefa | 524 | //__________________________________________________________________________ |
29f1b13a | 525 | AliMUONHitForRec* AliMUONTrackReconstructor::NewHitForRecFromTrackRef(AliTrackReference* Hit, Int_t TrackNumber, Int_t Signal) |
a9e2aefa | 526 | { |
29fc2c86 | 527 | // To make a new hit for reconstruction from a track ref. hit pointed to by "Hit", |
528 | // with the track numbered "TrackNumber", | |
a9e2aefa | 529 | // either from signal ("Signal" = 1) or background ("Signal" = 0) event. |
530 | // Selects hits in tracking (not trigger) chambers. | |
531 | // Takes into account the efficiency (fEfficiency) | |
532 | // and the smearing from resolution (fBendingResolution and fNonBendingResolution). | |
533 | // Adds a condition on the radius between RMin and RMax | |
534 | // to better simulate the real chambers. | |
535 | // Returns the pointer to the new hit for reconstruction, | |
536 | // or NULL in case of inefficiency or non tracking chamber or bad radius. | |
537 | // No condition on at most 20 cm from a muon from a resonance | |
538 | // like in Fortran TRACKF_STAT. | |
539 | AliMUONHitForRec* hitForRec; | |
540 | Double_t bendCoor, nonBendCoor, radius; | |
30d3ea8a | 541 | Int_t detElemId = Hit->UserId(), chamber = 0; |
542 | if (detElemId) chamber = detElemId / 100 - 1; | |
543 | else chamber = AliMUONConstants::ChamberNumber(Hit->Z()); // chamber(0...) | |
29fc2c86 | 544 | if (chamber < 0) return NULL; |
a9e2aefa | 545 | // only in tracking chambers (fChamber starts at 1) |
190f97ea | 546 | if (chamber >= AliMUONConstants::NTrackingCh()) return NULL; |
a9e2aefa | 547 | // only if hit is efficient (keep track for checking ????) |
548 | if (gRandom->Rndm() > fEfficiency) return NULL; | |
549 | // only if radius between RMin and RMax | |
94de3818 | 550 | bendCoor = Hit->Y(); |
551 | nonBendCoor = Hit->X(); | |
a9e2aefa | 552 | radius = TMath::Sqrt((bendCoor * bendCoor) + (nonBendCoor * nonBendCoor)); |
d0bfce8d | 553 | // This cut is not needed with a realistic chamber geometry !!!! |
554 | // if ((radius < fRMin[chamber]) || (radius > fRMax[chamber])) return NULL; | |
29fc2c86 | 555 | // new AliMUONHitForRec from track ref. hit and increment number of AliMUONHitForRec's |
a9e2aefa | 556 | hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(Hit); |
557 | fNHitsForRec++; | |
558 | // add smearing from resolution | |
559 | hitForRec->SetBendingCoor(bendCoor + gRandom->Gaus(0., fBendingResolution)); | |
560 | hitForRec->SetNonBendingCoor(nonBendCoor | |
561 | + gRandom->Gaus(0., fNonBendingResolution)); | |
d0bfce8d | 562 | // // !!!! without smearing |
563 | // hitForRec->SetBendingCoor(bendCoor); | |
564 | // hitForRec->SetNonBendingCoor(nonBendCoor); | |
a9e2aefa | 565 | // more information into HitForRec |
566 | // resolution: angular effect to be added here ???? | |
567 | hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution); | |
568 | hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution); | |
29fc2c86 | 569 | // track ref. info |
570 | hitForRec->SetTTRTrack(TrackNumber); | |
571 | hitForRec->SetTrackRefSignal(Signal); | |
b840996b | 572 | if (AliLog::GetGlobalDebugLevel()> 1) { |
573 | AliDebug(2,Form("Track: %d", TrackNumber)); | |
a9e2aefa | 574 | Hit->Dump(); |
575 | cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl; | |
b840996b | 576 | hitForRec->Dump(); |
577 | } | |
a9e2aefa | 578 | return hitForRec; |
579 | } | |
29fc2c86 | 580 | //__________________________________________________________________________ |
29f1b13a | 581 | TClonesArray* AliMUONTrackReconstructor::CleanTrackRefs(TTree *treeTR) |
29fc2c86 | 582 | { |
583 | // Make hits from track ref.. | |
584 | // Re-calculate hits parameters because two AliTrackReferences are recorded for | |
585 | // each chamber (one when particle is entering + one when particle is leaving | |
586 | // the sensitive volume) | |
587 | ||
588 | AliTrackReference *trackReference; | |
589 | Float_t x1, y1, z1, pX1, pY1, pZ1; | |
590 | Float_t x2, y2, z2, pX2, pY2, pZ2; | |
591 | Int_t track1, track2; | |
592 | Int_t nRec = 0; | |
593 | Float_t maxGasGap = 1.; // cm | |
594 | Int_t iHit1 = 0; | |
595 | Int_t iHitMin = 0; | |
596 | ||
597 | AliTrackReference *trackReferenceNew = new AliTrackReference(); | |
598 | ||
599 | TClonesArray* trackRefs = new TClonesArray("AliTrackReference", 10); | |
600 | TClonesArray* cleanTrackRefs = new TClonesArray("AliTrackReference", 10); | |
601 | ||
602 | if (treeTR == NULL) return NULL; | |
603 | TBranch* branch = treeTR->GetBranch("MUON"); | |
604 | if (branch == NULL) return NULL; | |
605 | branch->SetAddress(&trackRefs); | |
606 | ||
607 | Int_t nTrack = (Int_t)treeTR->GetEntries(); | |
608 | for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) { | |
609 | treeTR->GetEntry(iTrack); | |
610 | iHitMin = 0; | |
611 | iHit1 = 0; | |
612 | while (iHit1 < trackRefs->GetEntries()) { | |
613 | trackReference = (AliTrackReference*)trackRefs->At(iHit1); | |
614 | x1 = trackReference->X(); | |
615 | y1 = trackReference->Y(); | |
616 | z1 = trackReference->Z(); | |
617 | pX1 = trackReference->Px(); | |
618 | pY1 = trackReference->Py(); | |
619 | pZ1 = trackReference->Pz(); | |
620 | track1 = trackReference->GetTrack(); | |
621 | nRec = 1; | |
622 | iHitMin = iHit1+1; | |
623 | for (Int_t iHit2 = iHit1+1; iHit2 < trackRefs->GetEntries(); iHit2++) { | |
624 | trackReference = (AliTrackReference*)trackRefs->At(iHit2); | |
625 | x2 = trackReference->X(); | |
626 | y2 = trackReference->Y(); | |
627 | z2 = trackReference->Z(); | |
628 | pX2 = trackReference->Px(); | |
629 | pY2 = trackReference->Py(); | |
630 | pZ2 = trackReference->Pz(); | |
631 | track2 = trackReference->GetTrack(); | |
632 | if (track2 == track1 && TMath::Abs(z2-z1) < maxGasGap ) { | |
633 | nRec++; | |
634 | x1 += x2; | |
635 | y1 += y2; | |
636 | z1 += z2; | |
637 | pX1 += pX2; | |
638 | pY1 += pY2; | |
639 | pZ1 += pZ2; | |
640 | iHitMin = iHit2+1; | |
641 | } | |
642 | ||
643 | } // for iHit2 | |
644 | x1 /= (Float_t)nRec; | |
645 | y1 /= (Float_t)nRec; | |
646 | z1 /= (Float_t)nRec; | |
647 | pX1 /= (Float_t)nRec; | |
648 | pY1 /= (Float_t)nRec; | |
649 | pZ1 /= (Float_t)nRec; | |
650 | trackReferenceNew->SetPosition(x1,y1,z1); | |
651 | trackReferenceNew->SetMomentum(pX1,pY1,pZ1); | |
652 | trackReferenceNew->SetTrack(track1); | |
653 | {new ((*cleanTrackRefs)[cleanTrackRefs->GetEntriesFast()]) AliTrackReference(*trackReferenceNew);} | |
654 | iHit1 = iHitMin; | |
655 | } // while iHit1 | |
656 | } // for track | |
657 | ||
658 | trackRefs->Delete(); | |
659 | delete trackRefs; | |
660 | delete trackReferenceNew; | |
661 | return cleanTrackRefs; | |
662 | ||
663 | } | |
a9e2aefa | 664 | //__________________________________________________________________________ |
29f1b13a | 665 | void AliMUONTrackReconstructor::SortHitsForRecWithIncreasingChamber() |
a9e2aefa | 666 | { |
667 | // Sort HitsForRec's in increasing order with respect to chamber number. | |
668 | // Uses the function "Compare". | |
669 | // Update the information for HitsForRec per chamber too. | |
670 | Int_t ch, nhits, prevch; | |
671 | fHitsForRecPtr->Sort(); | |
190f97ea | 672 | for (ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) { |
a9e2aefa | 673 | fNHitsForRecPerChamber[ch] = 0; |
674 | fIndexOfFirstHitForRecPerChamber[ch] = 0; | |
675 | } | |
676 | prevch = 0; // previous chamber | |
677 | nhits = 0; // number of hits in current chamber | |
678 | // Loop over HitsForRec | |
679 | for (Int_t hit = 0; hit < fNHitsForRec; hit++) { | |
680 | // chamber number (0...) | |
681 | ch = ((AliMUONHitForRec*) ((*fHitsForRecPtr)[hit]))->GetChamberNumber(); | |
682 | // increment number of hits in current chamber | |
683 | (fNHitsForRecPerChamber[ch])++; | |
684 | // update index of first HitForRec in current chamber | |
685 | // if chamber number different from previous one | |
686 | if (ch != prevch) { | |
687 | fIndexOfFirstHitForRecPerChamber[ch] = hit; | |
688 | prevch = ch; | |
689 | } | |
690 | } | |
691 | return; | |
692 | } | |
693 | ||
a9e2aefa | 694 | //__________________________________________________________________________ |
29f1b13a | 695 | void AliMUONTrackReconstructor::AddHitsForRecFromRawClusters(TTree* TR) |
a9e2aefa | 696 | { |
697 | // To add to the list of hits for reconstruction all the raw clusters | |
698 | // No condition added, like in Fortran TRACKF_STAT, | |
699 | // on the radius between RMin and RMax. | |
700 | AliMUONHitForRec *hitForRec; | |
701 | AliMUONRawCluster *clus; | |
f69d51b5 | 702 | Int_t iclus, nclus, nTRentries; |
a9e2aefa | 703 | TClonesArray *rawclusters; |
b840996b | 704 | AliDebug(1,"Enter AddHitsForRecFromRawClusters"); |
88cb7938 | 705 | |
cc87ebcd | 706 | if (fTrackMethod != 3) { //AZ |
707 | fMUONData->SetTreeAddress("RC"); //AZ | |
708 | nTRentries = Int_t(TR->GetEntries()); | |
709 | if (nTRentries != 1) { | |
710 | AliError(Form("nTRentries = %d not equal to 1 ",nTRentries)); | |
711 | exit(0); | |
712 | } | |
713 | fLoader->TreeR()->GetEvent(0); // only one entry | |
f69d51b5 | 714 | } |
52c9bc11 | 715 | |
a9e2aefa | 716 | // Loop over tracking chambers |
190f97ea | 717 | for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) { |
a9e2aefa | 718 | // number of HitsForRec to 0 for the chamber |
719 | fNHitsForRecPerChamber[ch] = 0; | |
720 | // index of first HitForRec for the chamber | |
721 | if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0; | |
722 | else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec; | |
52c9bc11 | 723 | rawclusters =fMUONData->RawClusters(ch); |
f69d51b5 | 724 | // pMUON->ResetRawClusters(); |
725 | // TR->GetEvent((Int_t) (TR->GetEntries()) - 1); // to be checked ???? | |
a9e2aefa | 726 | nclus = (Int_t) (rawclusters->GetEntries()); |
727 | // Loop over (cathode correlated) raw clusters | |
728 | for (iclus = 0; iclus < nclus; iclus++) { | |
729 | clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(iclus); | |
730 | // new AliMUONHitForRec from raw cluster | |
731 | // and increment number of AliMUONHitForRec's (total and in chamber) | |
732 | hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(clus); | |
733 | fNHitsForRec++; | |
734 | (fNHitsForRecPerChamber[ch])++; | |
735 | // more information into HitForRec | |
736 | // resolution: info should be already in raw cluster and taken from it ???? | |
fff097e0 | 737 | //hitForRec->SetBendingReso2(fBendingResolution * fBendingResolution); |
738 | //hitForRec->SetNonBendingReso2(fNonBendingResolution * fNonBendingResolution); | |
739 | hitForRec->SetBendingReso2(clus->GetErrY() * clus->GetErrY()); | |
740 | hitForRec->SetNonBendingReso2(clus->GetErrX() * clus->GetErrX()); | |
a9e2aefa | 741 | // original raw cluster |
742 | hitForRec->SetChamberNumber(ch); | |
743 | hitForRec->SetHitNumber(iclus); | |
ba5b68db | 744 | // Z coordinate of the raw cluster (cm) |
ba12c242 | 745 | hitForRec->SetZ(clus->GetZ(0)); |
b840996b | 746 | if (AliLog::GetGlobalDebugLevel() > 1) { |
a9e2aefa | 747 | cout << "chamber (0...): " << ch << |
748 | " raw cluster (0...): " << iclus << endl; | |
749 | clus->Dump(); | |
750 | cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl; | |
751 | hitForRec->Dump();} | |
752 | } // end of cluster loop | |
753 | } // end of chamber loop | |
cc87ebcd | 754 | SortHitsForRecWithIncreasingChamber(); |
a9e2aefa | 755 | return; |
756 | } | |
757 | ||
758 | //__________________________________________________________________________ | |
29f1b13a | 759 | void AliMUONTrackReconstructor::MakeSegments(void) |
a9e2aefa | 760 | { |
761 | // To make the list of segments in all stations, | |
762 | // from the list of hits to be reconstructed | |
b840996b | 763 | AliDebug(1,"Enter MakeSegments"); |
fff097e0 | 764 | //AZ ResetSegments(); |
a9e2aefa | 765 | // Loop over stations |
cc87ebcd | 766 | Int_t nb = (fTrackMethod != 1) ? 3 : 0; //AZ |
767 | for (Int_t st = nb; st < AliMUONConstants::NTrackingCh()/2; st++) | |
a9e2aefa | 768 | MakeSegmentsPerStation(st); |
b840996b | 769 | if (AliLog::GetGlobalDebugLevel() > 1) { |
a9e2aefa | 770 | cout << "end of MakeSegments" << endl; |
190f97ea | 771 | for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) { |
a9e2aefa | 772 | cout << "station(0...): " << st |
773 | << " Segments: " << fNSegments[st] | |
774 | << endl; | |
775 | for (Int_t seg = 0; | |
776 | seg < fNSegments[st]; | |
777 | seg++) { | |
778 | cout << "Segment index(0...): " << seg << endl; | |
779 | ((*fSegmentsPtr[st])[seg])->Dump(); | |
780 | } | |
781 | } | |
782 | } | |
783 | return; | |
784 | } | |
785 | ||
786 | //__________________________________________________________________________ | |
29f1b13a | 787 | void AliMUONTrackReconstructor::MakeSegmentsPerStation(Int_t Station) |
a9e2aefa | 788 | { |
789 | // To make the list of segments in station number "Station" (0...) | |
790 | // from the list of hits to be reconstructed. | |
791 | // Updates "fNSegments"[Station]. | |
792 | // Segments in stations 4 and 5 are sorted | |
793 | // according to increasing absolute value of "impact parameter" | |
794 | AliMUONHitForRec *hit1Ptr, *hit2Ptr; | |
795 | AliMUONSegment *segment; | |
796 | Bool_t last2st; | |
797 | Double_t bendingSlope, distBend, distNonBend, extBendCoor, extNonBendCoor, | |
d0bfce8d | 798 | impactParam = 0., maxImpactParam = 0., minImpactParam = 0.; // =0 to avoid compilation warnings. |
b840996b | 799 | AliDebug(1,Form("Enter MakeSegmentsPerStation (0...) %d",Station)); |
a9e2aefa | 800 | // first and second chambers (0...) in the station |
801 | Int_t ch1 = 2 * Station; | |
802 | Int_t ch2 = ch1 + 1; | |
803 | // variable true for stations downstream of the dipole: | |
804 | // Station(0..4) equal to 3 or 4 | |
805 | if ((Station == 3) || (Station == 4)) { | |
806 | last2st = kTRUE; | |
807 | // maximum impact parameter (cm) according to fMinBendingMomentum... | |
808 | maxImpactParam = | |
809 | TMath::Abs(GetImpactParamFromBendingMomentum(fMinBendingMomentum)); | |
d0bfce8d | 810 | // minimum impact parameter (cm) according to fMaxBendingMomentum... |
811 | minImpactParam = | |
812 | TMath::Abs(GetImpactParamFromBendingMomentum(fMaxBendingMomentum)); | |
a9e2aefa | 813 | } |
814 | else last2st = kFALSE; | |
815 | // extrapolation factor from Z of first chamber to Z of second chamber | |
816 | // dZ to be changed to take into account fine structure of chambers ???? | |
e516b01d | 817 | Double_t extrapFact; |
a9e2aefa | 818 | // index for current segment |
819 | Int_t segmentIndex = 0; | |
820 | // Loop over HitsForRec in the first chamber of the station | |
821 | for (Int_t hit1 = fIndexOfFirstHitForRecPerChamber[ch1]; | |
822 | hit1 < fIndexOfFirstHitForRecPerChamber[ch1] + fNHitsForRecPerChamber[ch1]; | |
823 | hit1++) { | |
824 | // pointer to the HitForRec | |
825 | hit1Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit1]); | |
826 | // extrapolation, | |
827 | // on the straight line joining the HitForRec to the vertex (0,0,0), | |
828 | // to the Z of the second chamber of the station | |
a9e2aefa | 829 | // Loop over HitsForRec in the second chamber of the station |
830 | for (Int_t hit2 = fIndexOfFirstHitForRecPerChamber[ch2]; | |
831 | hit2 < fIndexOfFirstHitForRecPerChamber[ch2] + fNHitsForRecPerChamber[ch2]; | |
832 | hit2++) { | |
833 | // pointer to the HitForRec | |
834 | hit2Ptr = (AliMUONHitForRec*) ((*fHitsForRecPtr)[hit2]); | |
835 | // absolute values of distances, in bending and non bending planes, | |
836 | // between the HitForRec in the second chamber | |
837 | // and the previous extrapolation | |
e516b01d | 838 | extrapFact = hit2Ptr->GetZ()/ hit1Ptr->GetZ(); |
839 | extBendCoor = extrapFact * hit1Ptr->GetBendingCoor(); | |
840 | extNonBendCoor = extrapFact * hit1Ptr->GetNonBendingCoor(); | |
a9e2aefa | 841 | distBend = TMath::Abs(hit2Ptr->GetBendingCoor() - extBendCoor); |
842 | distNonBend = TMath::Abs(hit2Ptr->GetNonBendingCoor() - extNonBendCoor); | |
843 | if (last2st) { | |
844 | // bending slope | |
8999d47d | 845 | if ( hit1Ptr->GetZ() - hit2Ptr->GetZ() != 0.0 ) { |
846 | bendingSlope = (hit1Ptr->GetBendingCoor() - hit2Ptr->GetBendingCoor()) / | |
847 | (hit1Ptr->GetZ() - hit2Ptr->GetZ()); | |
848 | // absolute value of impact parameter | |
849 | impactParam = | |
850 | TMath::Abs(hit1Ptr->GetBendingCoor() - hit1Ptr->GetZ() * bendingSlope); | |
851 | } | |
852 | else { | |
853 | AliWarning("hit1Ptr->GetZ() = hit2Ptr->GetZ(): impactParam set to maxImpactParam"); | |
854 | impactParam = maxImpactParam; | |
855 | } | |
a9e2aefa | 856 | } |
857 | // check for distances not too large, | |
858 | // and impact parameter not too big if stations downstream of the dipole. | |
859 | // Conditions "distBend" and "impactParam" correlated for these stations ???? | |
860 | if ((distBend < fSegmentMaxDistBending[Station]) && | |
861 | (distNonBend < fSegmentMaxDistNonBending[Station]) && | |
d0bfce8d | 862 | (!last2st || (impactParam < maxImpactParam)) && |
863 | (!last2st || (impactParam > minImpactParam))) { | |
a9e2aefa | 864 | // make new segment |
865 | segment = new ((*fSegmentsPtr[Station])[segmentIndex]) | |
866 | AliMUONSegment(hit1Ptr, hit2Ptr); | |
867 | // update "link" to this segment from the hit in the first chamber | |
868 | if (hit1Ptr->GetNSegments() == 0) | |
869 | hit1Ptr->SetIndexOfFirstSegment(segmentIndex); | |
870 | hit1Ptr->SetNSegments(hit1Ptr->GetNSegments() + 1); | |
b840996b | 871 | if (AliLog::GetGlobalDebugLevel() > 1) { |
a9e2aefa | 872 | cout << "segmentIndex(0...): " << segmentIndex |
873 | << " distBend: " << distBend | |
874 | << " distNonBend: " << distNonBend | |
875 | << endl; | |
876 | segment->Dump(); | |
877 | cout << "HitForRec in first chamber" << endl; | |
878 | hit1Ptr->Dump(); | |
879 | cout << "HitForRec in second chamber" << endl; | |
880 | hit2Ptr->Dump(); | |
881 | }; | |
882 | // increment index for current segment | |
883 | segmentIndex++; | |
884 | } | |
885 | } //for (Int_t hit2 | |
886 | } // for (Int_t hit1... | |
887 | fNSegments[Station] = segmentIndex; | |
888 | // Sorting according to "impact parameter" if station(1..5) 4 or 5, | |
889 | // i.e. Station(0..4) 3 or 4, using the function "Compare". | |
890 | // After this sorting, it is impossible to use | |
891 | // the "fNSegments" and "fIndexOfFirstSegment" | |
892 | // of the HitForRec in the first chamber to explore all segments formed with it. | |
893 | // Is this sorting really needed ???? | |
894 | if ((Station == 3) || (Station == 4)) (fSegmentsPtr[Station])->Sort(); | |
b840996b | 895 | AliDebug(1,Form("Station: %d NSegments: %d ", Station, fNSegments[Station])); |
a9e2aefa | 896 | return; |
897 | } | |
898 | ||
899 | //__________________________________________________________________________ | |
29f1b13a | 900 | void AliMUONTrackReconstructor::MakeTracks(void) |
a9e2aefa | 901 | { |
902 | // To make the tracks, | |
903 | // from the list of segments and points in all stations | |
b840996b | 904 | AliDebug(1,"Enter MakeTracks"); |
8429a5e4 | 905 | // The order may be important for the following Reset's |
fff097e0 | 906 | //AZ ResetTracks(); |
907 | //AZ ResetTrackHits(); | |
cc87ebcd | 908 | if (fTrackMethod != 1) { //AZ - Kalman filter |
83dbc640 | 909 | MakeTrackCandidatesK(); |
b035a148 | 910 | if (fRecTracksPtr->GetEntriesFast() == 0) return; |
83dbc640 | 911 | // Follow tracks in stations(1..) 3, 2 and 1 |
912 | FollowTracksK(); | |
913 | // Remove double tracks | |
914 | RemoveDoubleTracksK(); | |
915 | // Propagate tracks to the vertex thru absorber | |
916 | GoToVertex(); | |
b035a148 | 917 | // Fill AliMUONTrack data members |
918 | FillMUONTrack(); | |
d837040f | 919 | } else { |
83dbc640 | 920 | // Look for candidates from at least 3 aligned points in stations(1..) 4 and 5 |
921 | MakeTrackCandidates(); | |
922 | // Follow tracks in stations(1..) 3, 2 and 1 | |
923 | FollowTracks(); | |
924 | // Remove double tracks | |
925 | RemoveDoubleTracks(); | |
d837040f | 926 | UpdateTrackParamAtHit(); |
b8dc484b | 927 | UpdateHitForRecAtHit(); |
52c9bc11 | 928 | } |
a9e2aefa | 929 | return; |
930 | } | |
931 | ||
276c44b7 | 932 | //__________________________________________________________________________ |
29f1b13a | 933 | void AliMUONTrackReconstructor::ValidateTracksWithTrigger(void) |
d837040f | 934 | { |
7e6d7e07 | 935 | // Try to match track from tracking system with trigger track |
d837040f | 936 | AliMUONTrack *track; |
24ab3324 | 937 | TClonesArray *recTriggerTracks; |
938 | ||
24ab3324 | 939 | fMUONData->SetTreeAddress("RL"); |
940 | fMUONData->GetRecTriggerTracks(); | |
941 | recTriggerTracks = fMUONData->RecTriggerTracks(); | |
e8765288 | 942 | |
d837040f | 943 | track = (AliMUONTrack*) fRecTracksPtr->First(); |
944 | while (track) { | |
945 | track->MatchTriggerTrack(recTriggerTracks); | |
946 | track = (AliMUONTrack*) fRecTracksPtr->After(track); | |
e8765288 | 947 | } |
24ab3324 | 948 | |
d837040f | 949 | } |
950 | ||
951 | //__________________________________________________________________________ | |
29f1b13a | 952 | Bool_t AliMUONTrackReconstructor::MakeTriggerTracks(void) |
276c44b7 | 953 | { |
954 | // To make the trigger tracks from Local Trigger | |
b840996b | 955 | AliDebug(1, "Enter MakeTriggerTracks"); |
276c44b7 | 956 | |
957 | Int_t nTRentries; | |
52c9bc11 | 958 | Long_t gloTrigPat; |
276c44b7 | 959 | TClonesArray *localTrigger; |
9131b4fe | 960 | TClonesArray *globalTrigger; |
276c44b7 | 961 | AliMUONLocalTrigger *locTrg; |
9131b4fe | 962 | AliMUONGlobalTrigger *gloTrg; |
bc4a7ff4 | 963 | |
52c9bc11 | 964 | AliMUONTriggerTrack *recTriggerTrack = 0; |
d837040f | 965 | |
30178c30 | 966 | TTree* treeR = fLoader->TreeR(); |
276c44b7 | 967 | |
276c44b7 | 968 | // Loading MUON subsystem |
969 | AliMUON * pMUON = (AliMUON *) gAlice->GetDetector("MUON"); | |
bc4a7ff4 | 970 | // Not really clean, but for the moment we must check whether the |
971 | // trigger uses new or old TriggerCircuit | |
972 | Bool_t newTrigger=kFALSE; | |
973 | if ( pMUON->DigitizerType().Contains("NewTrigger") ) newTrigger = kTRUE; | |
974 | ||
30178c30 | 975 | nTRentries = Int_t(treeR->GetEntries()); |
d837040f | 976 | |
30178c30 | 977 | treeR->GetEvent(0); // only one entry |
d837040f | 978 | |
979 | if (!(fMUONData->IsTriggerBranchesInTree())) { | |
8c343c7c | 980 | AliWarning(Form("Trigger information is not avalaible, nTRentries = %d not equal to 1",nTRentries)); |
d837040f | 981 | return kFALSE; |
276c44b7 | 982 | } |
276c44b7 | 983 | |
ce3e25a8 | 984 | fMUONData->SetTreeAddress("TC"); |
52c9bc11 | 985 | fMUONData->GetTrigger(); |
9131b4fe | 986 | |
52c9bc11 | 987 | // global trigger for trigger pattern |
9131b4fe | 988 | gloTrigPat = 0; |
52c9bc11 | 989 | globalTrigger = fMUONData->GlobalTrigger(); |
9131b4fe | 990 | gloTrg = (AliMUONGlobalTrigger*)globalTrigger->UncheckedAt(0); |
8efb0686 | 991 | if (gloTrg) { |
992 | if (gloTrg->SinglePlusLpt()) gloTrigPat|= 0x1; | |
993 | if (gloTrg->SinglePlusHpt()) gloTrigPat|= 0x2; | |
994 | if (gloTrg->SinglePlusApt()) gloTrigPat|= 0x4; | |
995 | ||
996 | if (gloTrg->SingleMinusLpt()) gloTrigPat|= 0x8; | |
997 | if (gloTrg->SingleMinusHpt()) gloTrigPat|= 0x10; | |
998 | if (gloTrg->SingleMinusApt()) gloTrigPat|= 0x20; | |
999 | ||
1000 | if (gloTrg->SingleUndefLpt()) gloTrigPat|= 0x40; | |
1001 | if (gloTrg->SingleUndefHpt()) gloTrigPat|= 0x80; | |
1002 | if (gloTrg->SingleUndefApt()) gloTrigPat|= 0x100; | |
1003 | ||
1004 | if (gloTrg->PairUnlikeLpt()) gloTrigPat|= 0x200; | |
1005 | if (gloTrg->PairUnlikeHpt()) gloTrigPat|= 0x400; | |
1006 | if (gloTrg->PairUnlikeApt()) gloTrigPat|= 0x800; | |
1007 | ||
1008 | if (gloTrg->PairLikeLpt()) gloTrigPat|= 0x1000; | |
1009 | if (gloTrg->PairLikeHpt()) gloTrigPat|= 0x2000; | |
1010 | if (gloTrg->PairLikeApt()) gloTrigPat|= 0x4000; | |
1011 | } | |
9131b4fe | 1012 | |
52c9bc11 | 1013 | // local trigger for tracking |
1014 | localTrigger = fMUONData->LocalTrigger(); | |
276c44b7 | 1015 | Int_t nlocals = (Int_t) (localTrigger->GetEntries()); |
06c11448 | 1016 | |
1017 | Float_t z11 = AliMUONConstants::DefaultChamberZ(10); | |
1018 | Float_t z21 = AliMUONConstants::DefaultChamberZ(12); | |
276c44b7 | 1019 | |
bc4a7ff4 | 1020 | Float_t y11 = 0.; |
1021 | Int_t stripX21 = 0; | |
1022 | Float_t y21 = 0.; | |
1023 | Float_t x11 = 0.; | |
1024 | ||
276c44b7 | 1025 | for (Int_t i=0; i<nlocals; i++) { // loop on Local Trigger |
52c9bc11 | 1026 | locTrg = (AliMUONLocalTrigger*)localTrigger->UncheckedAt(i); |
bc4a7ff4 | 1027 | |
1028 | if (!newTrigger) { // old trigger | |
1029 | // printf("AliMUONTrackReconstructor::MakeTriggerTrack using OLD trigger \n"); | |
1030 | AliMUONTriggerCircuit * circuit = &(pMUON->TriggerCircuit(locTrg->LoCircuit())); | |
1031 | y11 = circuit->GetY11Pos(locTrg->LoStripX()); | |
1032 | stripX21 = locTrg->LoStripX()+locTrg->LoDev()+1; | |
1033 | y21 = circuit->GetY21Pos(stripX21); | |
1034 | x11 = circuit->GetX11Pos(locTrg->LoStripY()); | |
1035 | } else { // new trigger | |
1036 | // printf("AliMUONTrackReconstructor::MakeTriggerTrack using NEW trigger \n"); | |
1037 | AliMUONTriggerCircuitNew * circuit = | |
1038 | &(pMUON->TriggerCircuitNew(locTrg->LoCircuit()-1)); // -1 !!! | |
1039 | y11 = circuit->GetY11Pos(locTrg->LoStripX()); | |
1040 | stripX21 = locTrg->LoStripX()+locTrg->LoDev()+1; | |
1041 | y21 = circuit->GetY21Pos(stripX21); | |
1042 | x11 = circuit->GetX11Pos(locTrg->LoStripY()); | |
1043 | } | |
1044 | // printf(" MakeTriggerTrack %d %d %d %d %d %f %f %f \n",i,locTrg->LoCircuit(),locTrg->LoStripX(),locTrg->LoStripX()+locTrg->LoDev()+1,locTrg->LoStripY(),y11, y21, x11); | |
1045 | ||
52c9bc11 | 1046 | Float_t thetax = TMath::ATan2( x11 , z11 ); |
1047 | Float_t thetay = TMath::ATan2( (y21-y11) , (z21-z11) ); | |
bc4a7ff4 | 1048 | |
7fe0032c | 1049 | recTriggerTrack = new AliMUONTriggerTrack(x11,y11,thetax,thetay,gloTrigPat); |
34f1bfa0 | 1050 | |
52c9bc11 | 1051 | // since static statement does not work, set gloTrigPat for each track |
bc4a7ff4 | 1052 | |
52c9bc11 | 1053 | fMUONData->AddRecTriggerTrack(*recTriggerTrack); |
34f1bfa0 | 1054 | delete recTriggerTrack; |
276c44b7 | 1055 | } // end of loop on Local Trigger |
d837040f | 1056 | return kTRUE; |
276c44b7 | 1057 | } |
1058 | ||
a9e2aefa | 1059 | //__________________________________________________________________________ |
29f1b13a | 1060 | Int_t AliMUONTrackReconstructor::MakeTrackCandidatesWithTwoSegments(AliMUONSegment *BegSegment) |
a9e2aefa | 1061 | { |
1062 | // To make track candidates with two segments in stations(1..) 4 and 5, | |
1063 | // the first segment being pointed to by "BegSegment". | |
1064 | // Returns the number of such track candidates. | |
1065 | Int_t endStation, iEndSegment, nbCan2Seg; | |
e516b01d | 1066 | AliMUONSegment *endSegment; |
1067 | AliMUONSegment *extrapSegment = NULL; | |
a9e2aefa | 1068 | AliMUONTrack *recTrack; |
1069 | Double_t mcsFactor; | |
b840996b | 1070 | AliDebug(1,"Enter MakeTrackCandidatesWithTwoSegments"); |
a9e2aefa | 1071 | // Station for the end segment |
1072 | endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2; | |
1073 | // multiple scattering factor corresponding to one chamber | |
1074 | mcsFactor = 0.0136 / | |
1075 | GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact()); | |
1076 | mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor; | |
1077 | // linear extrapolation to end station | |
a9e2aefa | 1078 | // number of candidates with 2 segments to 0 |
1079 | nbCan2Seg = 0; | |
1080 | // Loop over segments in the end station | |
1081 | for (iEndSegment = 0; iEndSegment < fNSegments[endStation]; iEndSegment++) { | |
1082 | // pointer to segment | |
1083 | endSegment = (AliMUONSegment*) ((*fSegmentsPtr[endStation])[iEndSegment]); | |
1084 | // test compatibility between current segment and "extrapSegment" | |
04b5ea16 | 1085 | // 4 because 4 quantities in chi2 |
e516b01d | 1086 | extrapSegment = |
1087 | BegSegment->CreateSegmentFromLinearExtrapToStation(endSegment->GetZ(), mcsFactor); | |
a9e2aefa | 1088 | if ((endSegment-> |
1089 | NormalizedChi2WithSegment(extrapSegment, | |
1090 | fMaxSigma2Distance)) <= 4.0) { | |
1091 | // both segments compatible: | |
1092 | // make track candidate from "begSegment" and "endSegment" | |
b840996b | 1093 | AliDebug(2,Form("TrackCandidate with Segment %d in Station(0..) %d", iEndSegment, endStation)); |
a9e2aefa | 1094 | // flag for both segments in one track: |
1095 | // to be done in track constructor ???? | |
1096 | BegSegment->SetInTrack(kTRUE); | |
1097 | endSegment->SetInTrack(kTRUE); | |
1098 | recTrack = new ((*fRecTracksPtr)[fNRecTracks]) | |
1099 | AliMUONTrack(BegSegment, endSegment, this); | |
1100 | fNRecTracks++; | |
b840996b | 1101 | if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump(); |
a9e2aefa | 1102 | // increment number of track candidates with 2 segments |
1103 | nbCan2Seg++; | |
1104 | } | |
ddf8948c | 1105 | delete extrapSegment; // should not delete HitForRec's it points to !!!! |
a9e2aefa | 1106 | } // for (iEndSegment = 0;... |
a9e2aefa | 1107 | return nbCan2Seg; |
1108 | } | |
1109 | ||
1110 | //__________________________________________________________________________ | |
29f1b13a | 1111 | Int_t AliMUONTrackReconstructor::MakeTrackCandidatesWithOneSegmentAndOnePoint(AliMUONSegment *BegSegment) |
a9e2aefa | 1112 | { |
1113 | // To make track candidates with one segment and one point | |
1114 | // in stations(1..) 4 and 5, | |
1115 | // the segment being pointed to by "BegSegment". | |
1116 | Int_t ch, ch1, ch2, endStation, iHit, iHitMax, iHitMin, nbCan1Seg1Hit; | |
e516b01d | 1117 | AliMUONHitForRec *extrapHitForRec= NULL; |
1118 | AliMUONHitForRec *hit; | |
a9e2aefa | 1119 | AliMUONTrack *recTrack; |
1120 | Double_t mcsFactor; | |
b840996b | 1121 | AliDebug(1,"Enter MakeTrackCandidatesWithOneSegmentAndOnePoint"); |
a9e2aefa | 1122 | // station for the end point |
1123 | endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2; | |
1124 | // multiple scattering factor corresponding to one chamber | |
1125 | mcsFactor = 0.0136 / | |
1126 | GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact()); | |
1127 | mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor; | |
1128 | // first and second chambers(0..) in the end station | |
1129 | ch1 = 2 * endStation; | |
1130 | ch2 = ch1 + 1; | |
1131 | // number of candidates to 0 | |
1132 | nbCan1Seg1Hit = 0; | |
1133 | // Loop over chambers of the end station | |
1134 | for (ch = ch2; ch >= ch1; ch--) { | |
a9e2aefa | 1135 | // limits for the hit index in the loop |
1136 | iHitMin = fIndexOfFirstHitForRecPerChamber[ch]; | |
1137 | iHitMax = iHitMin + fNHitsForRecPerChamber[ch]; | |
1138 | // Loop over HitForRec's in the chamber | |
1139 | for (iHit = iHitMin; iHit < iHitMax; iHit++) { | |
1140 | // pointer to HitForRec | |
1141 | hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]); | |
1142 | // test compatibility between current HitForRec and "extrapHitForRec" | |
04b5ea16 | 1143 | // 2 because 2 quantities in chi2 |
e516b01d | 1144 | // linear extrapolation to chamber |
1145 | extrapHitForRec = | |
1146 | BegSegment->CreateHitForRecFromLinearExtrapToChamber( hit->GetZ(), mcsFactor); | |
a9e2aefa | 1147 | if ((hit-> |
1148 | NormalizedChi2WithHitForRec(extrapHitForRec, | |
1149 | fMaxSigma2Distance)) <= 2.0) { | |
1150 | // both HitForRec's compatible: | |
1151 | // make track candidate from begSegment and current HitForRec | |
b840996b | 1152 | AliDebug(2, Form("TrackCandidate with HitForRec %d in Chamber(0..) %d", iHit, ch)); |
a9e2aefa | 1153 | // flag for beginning segments in one track: |
1154 | // to be done in track constructor ???? | |
1155 | BegSegment->SetInTrack(kTRUE); | |
1156 | recTrack = new ((*fRecTracksPtr)[fNRecTracks]) | |
1157 | AliMUONTrack(BegSegment, hit, this); | |
1158 | // the right place to eliminate "double counting" ???? how ???? | |
1159 | fNRecTracks++; | |
b840996b | 1160 | if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump(); |
a9e2aefa | 1161 | // increment number of track candidates |
1162 | nbCan1Seg1Hit++; | |
1163 | } | |
ddf8948c | 1164 | delete extrapHitForRec; |
a9e2aefa | 1165 | } // for (iHit = iHitMin;... |
a9e2aefa | 1166 | } // for (ch = ch2;... |
1167 | return nbCan1Seg1Hit; | |
1168 | } | |
1169 | ||
1170 | //__________________________________________________________________________ | |
29f1b13a | 1171 | void AliMUONTrackReconstructor::MakeTrackCandidates(void) |
a9e2aefa | 1172 | { |
1173 | // To make track candidates | |
1174 | // with at least 3 aligned points in stations(1..) 4 and 5 | |
1175 | // (two Segment's or one Segment and one HitForRec) | |
1176 | Int_t begStation, iBegSegment, nbCan1Seg1Hit, nbCan2Seg; | |
1177 | AliMUONSegment *begSegment; | |
b840996b | 1178 | AliDebug(1,"Enter MakeTrackCandidates"); |
a9e2aefa | 1179 | // Loop over stations(1..) 5 and 4 for the beginning segment |
1180 | for (begStation = 4; begStation > 2; begStation--) { | |
1181 | // Loop over segments in the beginning station | |
1182 | for (iBegSegment = 0; iBegSegment < fNSegments[begStation]; iBegSegment++) { | |
1183 | // pointer to segment | |
1184 | begSegment = (AliMUONSegment*) ((*fSegmentsPtr[begStation])[iBegSegment]); | |
b840996b | 1185 | AliDebug(2,Form("Look for TrackCandidate's with Segment %d in Station(0..) %d", iBegSegment, begStation)); |
a9e2aefa | 1186 | // Look for track candidates with two segments, |
1187 | // "begSegment" and all compatible segments in other station. | |
1188 | // Only for beginning station(1..) 5 | |
1189 | // because candidates with 2 segments have to looked for only once. | |
1190 | if (begStation == 4) | |
1191 | nbCan2Seg = MakeTrackCandidatesWithTwoSegments(begSegment); | |
956019b6 | 1192 | // Look for track candidates with one segment and one point, |
a9e2aefa | 1193 | // "begSegment" and all compatible HitForRec's in other station. |
1194 | // Only if "begSegment" does not belong already to a track candidate. | |
1195 | // Is that a too strong condition ???? | |
1196 | if (!(begSegment->GetInTrack())) | |
1197 | nbCan1Seg1Hit = MakeTrackCandidatesWithOneSegmentAndOnePoint(begSegment); | |
1198 | } // for (iBegSegment = 0;... | |
1199 | } // for (begStation = 4;... | |
1200 | return; | |
1201 | } | |
1202 | ||
1203 | //__________________________________________________________________________ | |
29f1b13a | 1204 | void AliMUONTrackReconstructor::FollowTracks(void) |
a9e2aefa | 1205 | { |
1206 | // Follow tracks in stations(1..) 3, 2 and 1 | |
04b5ea16 | 1207 | // too long: should be made more modular !!!! |
e516b01d | 1208 | AliMUONHitForRec *bestHit, *extrapHit, *hit; |
1209 | AliMUONSegment *bestSegment, *extrapSegment, *segment; | |
04b5ea16 | 1210 | AliMUONTrack *track, *nextTrack; |
1211 | AliMUONTrackParam *trackParam1, trackParam[2], trackParamVertex; | |
ecfa008b | 1212 | // -1 to avoid compilation warnings |
1213 | Int_t ch = -1, chInStation, chBestHit = -1, iHit, iSegment, station, trackIndex; | |
04b5ea16 | 1214 | Double_t bestChi2, chi2, dZ1, dZ2, dZ3, maxSigma2Distance, mcsFactor; |
d0bfce8d | 1215 | Double_t bendingMomentum, chi2Norm = 0.; |
cc87ebcd | 1216 | |
04b5ea16 | 1217 | // local maxSigma2Distance, for easy increase in testing |
1218 | maxSigma2Distance = fMaxSigma2Distance; | |
b840996b | 1219 | AliDebug(2,"Enter FollowTracks"); |
a9e2aefa | 1220 | // Loop over track candidates |
04b5ea16 | 1221 | track = (AliMUONTrack*) fRecTracksPtr->First(); |
1222 | trackIndex = -1; | |
1223 | while (track) { | |
1224 | // Follow function for each track candidate ???? | |
1225 | trackIndex++; | |
1226 | nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track); // prepare next track | |
b840996b | 1227 | AliDebug(2,Form("FollowTracks: track candidate(0..): %d", trackIndex)); |
956019b6 | 1228 | // Fit track candidate |
1229 | track->SetFitMCS(0); // without multiple Coulomb scattering | |
1230 | track->SetFitNParam(3); // with 3 parameters (X = Y = 0) | |
1231 | track->SetFitStart(0); // from parameters at vertex | |
1232 | track->Fit(); | |
b840996b | 1233 | if (AliLog::GetGlobalDebugLevel()> 2) { |
a9e2aefa | 1234 | cout << "FollowTracks: track candidate(0..): " << trackIndex |
956019b6 | 1235 | << " after fit in stations(0..) 3 and 4" << endl; |
a9e2aefa | 1236 | track->RecursiveDump(); |
1237 | } | |
1238 | // Loop over stations(1..) 3, 2 and 1 | |
04b5ea16 | 1239 | // something SPECIAL for stations 2 and 1 for majority 3 coincidence ???? |
1240 | // otherwise: majority coincidence 2 !!!! | |
a9e2aefa | 1241 | for (station = 2; station >= 0; station--) { |
1242 | // Track parameters at first track hit (smallest Z) | |
1243 | trackParam1 = ((AliMUONTrackHit*) | |
1244 | (track->GetTrackHitsPtr()->First()))->GetTrackParam(); | |
1245 | // extrapolation to station | |
1246 | trackParam1->ExtrapToStation(station, trackParam); | |
79e1e601 | 1247 | extrapSegment = new AliMUONSegment(); // empty segment |
a9e2aefa | 1248 | // multiple scattering factor corresponding to one chamber |
1249 | // and momentum in bending plane (not total) | |
1250 | mcsFactor = 0.0136 * trackParam1->GetInverseBendingMomentum(); | |
1251 | mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor; | |
1252 | // Z difference from previous station | |
f93b9bd8 | 1253 | dZ1 = AliMUONConstants::DefaultChamberZ(2 * station) - |
1254 | AliMUONConstants::DefaultChamberZ(2 * station + 2); | |
a9e2aefa | 1255 | // Z difference between the two previous stations |
f93b9bd8 | 1256 | dZ2 = AliMUONConstants::DefaultChamberZ(2 * station + 2) - |
1257 | AliMUONConstants::DefaultChamberZ(2 * station + 4); | |
04b5ea16 | 1258 | // Z difference between the two chambers in the previous station |
f93b9bd8 | 1259 | dZ3 = AliMUONConstants::DefaultChamberZ(2 * station) - |
1260 | AliMUONConstants::DefaultChamberZ(2 * station + 1); | |
04b5ea16 | 1261 | extrapSegment->SetBendingCoorReso2(fBendingResolution * fBendingResolution); |
1262 | extrapSegment-> | |
1263 | SetNonBendingCoorReso2(fNonBendingResolution * fNonBendingResolution); | |
1264 | extrapSegment->UpdateFromStationTrackParam | |
1265 | (trackParam, mcsFactor, dZ1, dZ2, dZ3, station, | |
1266 | trackParam1->GetInverseBendingMomentum()); | |
a9e2aefa | 1267 | bestChi2 = 5.0; |
1268 | bestSegment = NULL; | |
b840996b | 1269 | if (AliLog::GetGlobalDebugLevel() > 2) { |
a9e2aefa | 1270 | cout << "FollowTracks: track candidate(0..): " << trackIndex |
1271 | << " Look for segment in station(0..): " << station << endl; | |
1272 | } | |
1273 | // Loop over segments in station | |
1274 | for (iSegment = 0; iSegment < fNSegments[station]; iSegment++) { | |
1275 | // Look for best compatible Segment in station | |
1276 | // should consider all possibilities ???? | |
1277 | // multiple scattering ???? | |
1278 | // separation in 2 functions: Segment and HitForRec ???? | |
1279 | segment = (AliMUONSegment*) ((*fSegmentsPtr[station])[iSegment]); | |
d0bfce8d | 1280 | // correction of corrected segment (fBendingCoor and fNonBendingCoor) |
1281 | // according to real Z value of "segment" and slopes of "extrapSegment" | |
e516b01d | 1282 | (&(trackParam[0]))->ExtrapToZ(segment->GetZ()); |
1283 | (&(trackParam[1]))->ExtrapToZ(segment->GetZ()); | |
1284 | extrapSegment->SetBendingCoor((&(trackParam[0]))->GetBendingCoor()); | |
1285 | extrapSegment->SetNonBendingCoor((&(trackParam[0]))->GetNonBendingCoor()); | |
1286 | extrapSegment->SetBendingSlope((&(trackParam[0]))->GetBendingSlope()); | |
1287 | extrapSegment->SetNonBendingSlope((&(trackParam[0]))->GetNonBendingSlope()); | |
d0bfce8d | 1288 | chi2 = segment-> |
e516b01d | 1289 | NormalizedChi2WithSegment(extrapSegment, maxSigma2Distance); |
a9e2aefa | 1290 | if (chi2 < bestChi2) { |
1291 | // update best Chi2 and Segment if better found | |
1292 | bestSegment = segment; | |
1293 | bestChi2 = chi2; | |
1294 | } | |
1295 | } | |
1296 | if (bestSegment) { | |
1297 | // best segment found: add it to track candidate | |
e516b01d | 1298 | (&(trackParam[0]))->ExtrapToZ(bestSegment->GetZ()); |
1299 | (&(trackParam[1]))->ExtrapToZ(bestSegment->GetZ()); | |
a9e2aefa | 1300 | track->AddSegment(bestSegment); |
1301 | // set track parameters at these two TrakHit's | |
1302 | track->SetTrackParamAtHit(track->GetNTrackHits() - 2, &(trackParam[0])); | |
1303 | track->SetTrackParamAtHit(track->GetNTrackHits() - 1, &(trackParam[1])); | |
b840996b | 1304 | AliDebug(3, Form("FollowTracks: track candidate(0..): %d Added segment in station(0..): %d", trackIndex, station)); |
1305 | if (AliLog::GetGlobalDebugLevel()>2) track->RecursiveDump(); | |
a9e2aefa | 1306 | } |
1307 | else { | |
1308 | // No best segment found: | |
1309 | // Look for best compatible HitForRec in station: | |
1310 | // should consider all possibilities ???? | |
1311 | // multiple scattering ???? do about like for extrapSegment !!!! | |
79e1e601 | 1312 | extrapHit = new AliMUONHitForRec(); // empty hit |
a9e2aefa | 1313 | bestChi2 = 3.0; |
1314 | bestHit = NULL; | |
b840996b | 1315 | AliDebug(3, Form("FollowTracks: track candidate(0..): %d Segment not found, look for hit in station(0..): %d ", |
1316 | trackIndex, station)); | |
1317 | ||
a9e2aefa | 1318 | // Loop over chambers of the station |
d82671a0 | 1319 | for (chInStation = 0; chInStation < 2; chInStation++) { |
956019b6 | 1320 | ch = 2 * station + chInStation; |
a9e2aefa | 1321 | for (iHit = fIndexOfFirstHitForRecPerChamber[ch]; |
1322 | iHit < fIndexOfFirstHitForRecPerChamber[ch] + | |
1323 | fNHitsForRecPerChamber[ch]; | |
1324 | iHit++) { | |
1325 | hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]); | |
e516b01d | 1326 | // coordinates of extrapolated hit |
1327 | (&(trackParam[chInStation]))->ExtrapToZ(hit->GetZ()); | |
1328 | extrapHit-> | |
1329 | SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor()); | |
1330 | extrapHit-> | |
1331 | SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor()); | |
1332 | // resolutions from "extrapSegment" | |
1333 | extrapHit->SetBendingReso2(extrapSegment->GetBendingCoorReso2()); | |
1334 | extrapHit->SetNonBendingReso2(extrapSegment->GetNonBendingCoorReso2()); | |
1335 | // Loop over hits in the chamber | |
a9e2aefa | 1336 | // condition for hit not already in segment ???? |
1337 | chi2 = hit->NormalizedChi2WithHitForRec(extrapHit, maxSigma2Distance); | |
1338 | if (chi2 < bestChi2) { | |
1339 | // update best Chi2 and HitForRec if better found | |
1340 | bestHit = hit; | |
1341 | bestChi2 = chi2; | |
d82671a0 | 1342 | chBestHit = chInStation; |
a9e2aefa | 1343 | } |
1344 | } | |
1345 | } | |
1346 | if (bestHit) { | |
1347 | // best hit found: add it to track candidate | |
e516b01d | 1348 | (&(trackParam[chBestHit]))->ExtrapToZ(bestHit->GetZ()); |
a9e2aefa | 1349 | track->AddHitForRec(bestHit); |
956019b6 | 1350 | // set track parameters at this TrackHit |
a9e2aefa | 1351 | track->SetTrackParamAtHit(track->GetNTrackHits() - 1, |
1352 | &(trackParam[chBestHit])); | |
b840996b | 1353 | if (AliLog::GetGlobalDebugLevel() > 2) { |
a9e2aefa | 1354 | cout << "FollowTracks: track candidate(0..): " << trackIndex |
1355 | << " Added hit in station(0..): " << station << endl; | |
1356 | track->RecursiveDump(); | |
1357 | } | |
1358 | } | |
1359 | else { | |
8429a5e4 | 1360 | // Remove current track candidate |
1361 | // and corresponding TrackHit's, ... | |
1362 | track->Remove(); | |
04b5ea16 | 1363 | delete extrapSegment; |
d0bfce8d | 1364 | delete extrapHit; |
a9e2aefa | 1365 | break; // stop the search for this candidate: |
1366 | // exit from the loop over station | |
1367 | } | |
d0bfce8d | 1368 | delete extrapHit; |
a9e2aefa | 1369 | } |
04b5ea16 | 1370 | delete extrapSegment; |
a9e2aefa | 1371 | // Sort track hits according to increasing Z |
1372 | track->GetTrackHitsPtr()->Sort(); | |
1373 | // Update track parameters at first track hit (smallest Z) | |
1374 | trackParam1 = ((AliMUONTrackHit*) | |
1375 | (track->GetTrackHitsPtr()->First()))->GetTrackParam(); | |
d0bfce8d | 1376 | bendingMomentum = 0.; |
1377 | if (TMath::Abs(trackParam1->GetInverseBendingMomentum()) > 0.) | |
1378 | bendingMomentum = TMath::Abs(1/(trackParam1->GetInverseBendingMomentum())); | |
1379 | // Track removed if bendingMomentum not in window [min, max] | |
1380 | if ((bendingMomentum < fMinBendingMomentum) || (bendingMomentum > fMaxBendingMomentum)) { | |
1381 | track->Remove(); | |
1382 | break; // stop the search for this candidate: | |
1383 | // exit from the loop over station | |
1384 | } | |
956019b6 | 1385 | // Track fit |
04b5ea16 | 1386 | // with multiple Coulomb scattering if all stations |
1387 | if (station == 0) track->SetFitMCS(1); | |
956019b6 | 1388 | // without multiple Coulomb scattering if not all stations |
04b5ea16 | 1389 | else track->SetFitMCS(0); |
956019b6 | 1390 | track->SetFitNParam(5); // with 5 parameters (momentum and position) |
1391 | track->SetFitStart(1); // from parameters at first hit | |
1392 | track->Fit(); | |
d0bfce8d | 1393 | Double_t numberOfDegFree = (2.0 * track->GetNTrackHits() - 5); |
1394 | if (numberOfDegFree > 0) { | |
1395 | chi2Norm = track->GetFitFMin() / numberOfDegFree; | |
1396 | } else { | |
1397 | chi2Norm = 1.e10; | |
1398 | } | |
1399 | // Track removed if normalized chi2 too high | |
1400 | if (chi2Norm > fMaxChi2) { | |
1401 | track->Remove(); | |
1402 | break; // stop the search for this candidate: | |
1403 | // exit from the loop over station | |
1404 | } | |
b840996b | 1405 | if (AliLog::GetGlobalDebugLevel() > 2) { |
a9e2aefa | 1406 | cout << "FollowTracks: track candidate(0..): " << trackIndex |
1407 | << " after fit from station(0..): " << station << " to 4" << endl; | |
1408 | track->RecursiveDump(); | |
1409 | } | |
04b5ea16 | 1410 | // Track extrapolation to the vertex through the absorber (Branson) |
956019b6 | 1411 | // after going through the first station |
04b5ea16 | 1412 | if (station == 0) { |
1413 | trackParamVertex = *trackParam1; | |
889a0215 | 1414 | (&trackParamVertex)->ExtrapToVertex(0.,0.,0.); |
04b5ea16 | 1415 | track->SetTrackParamAtVertex(&trackParamVertex); |
b840996b | 1416 | if (AliLog::GetGlobalDebugLevel() > 0) { |
6ae236ee | 1417 | cout << "FollowTracks: track candidate(0..): " << trackIndex |
1418 | << " after extrapolation to vertex" << endl; | |
1419 | track->RecursiveDump(); | |
1420 | } | |
04b5ea16 | 1421 | } |
a9e2aefa | 1422 | } // for (station = 2;... |
04b5ea16 | 1423 | // go really to next track |
1424 | track = nextTrack; | |
1425 | } // while (track) | |
1426 | // Compression of track array (necessary after Remove ????) | |
1427 | fRecTracksPtr->Compress(); | |
1428 | return; | |
1429 | } | |
1430 | ||
1431 | //__________________________________________________________________________ | |
29f1b13a | 1432 | void AliMUONTrackReconstructor::RemoveDoubleTracks(void) |
8429a5e4 | 1433 | { |
1434 | // To remove double tracks. | |
1435 | // Tracks are considered identical | |
1436 | // if they have at least half of their hits in common. | |
1437 | // Among two identical tracks, one keeps the track with the larger number of hits | |
1438 | // or, if these numbers are equal, the track with the minimum Chi2. | |
1439 | AliMUONTrack *track1, *track2, *trackToRemove; | |
1440 | Bool_t identicalTracks; | |
1441 | Int_t hitsInCommon, nHits1, nHits2; | |
1442 | identicalTracks = kTRUE; | |
1443 | while (identicalTracks) { | |
1444 | identicalTracks = kFALSE; | |
1445 | // Loop over first track of the pair | |
1446 | track1 = (AliMUONTrack*) fRecTracksPtr->First(); | |
1447 | while (track1 && (!identicalTracks)) { | |
1448 | nHits1 = track1->GetNTrackHits(); | |
1449 | // Loop over second track of the pair | |
1450 | track2 = (AliMUONTrack*) fRecTracksPtr->After(track1); | |
1451 | while (track2 && (!identicalTracks)) { | |
1452 | nHits2 = track2->GetNTrackHits(); | |
1453 | // number of hits in common between two tracks | |
1454 | hitsInCommon = track1->HitsInCommon(track2); | |
1455 | // check for identical tracks | |
1456 | if ((4 * hitsInCommon) >= (nHits1 + nHits2)) { | |
1457 | identicalTracks = kTRUE; | |
1458 | // decide which track to remove | |
1459 | if (nHits1 > nHits2) trackToRemove = track2; | |
1460 | else if (nHits1 < nHits2) trackToRemove = track1; | |
1461 | else if ((track1->GetFitFMin()) < (track2->GetFitFMin())) | |
1462 | trackToRemove = track2; | |
1463 | else trackToRemove = track1; | |
1464 | // remove it | |
1465 | trackToRemove->Remove(); | |
1466 | } | |
1467 | track2 = (AliMUONTrack*) fRecTracksPtr->After(track2); | |
1468 | } // track2 | |
1469 | track1 = (AliMUONTrack*) fRecTracksPtr->After(track1); | |
1470 | } // track1 | |
1471 | } | |
1472 | return; | |
1473 | } | |
1474 | ||
d837040f | 1475 | //__________________________________________________________________________ |
29f1b13a | 1476 | void AliMUONTrackReconstructor::UpdateTrackParamAtHit() |
d837040f | 1477 | { |
1478 | // Set track parameters after track fitting. Fill fTrackParamAtHit of AliMUONTrack's | |
1479 | AliMUONTrack *track; | |
1480 | AliMUONTrackHit *trackHit; | |
1481 | AliMUONTrackParam *trackParam; | |
1482 | track = (AliMUONTrack*) fRecTracksPtr->First(); | |
1483 | while (track) { | |
1484 | trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->First(); | |
1485 | while (trackHit) { | |
1486 | trackParam = trackHit->GetTrackParam(); | |
1487 | track->AddTrackParamAtHit(trackParam); | |
1488 | trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->After(trackHit); | |
1489 | } // trackHit | |
1490 | track = (AliMUONTrack*) fRecTracksPtr->After(track); | |
1491 | } // track | |
1492 | return; | |
1493 | } | |
1494 | ||
b8dc484b | 1495 | //__________________________________________________________________________ |
29f1b13a | 1496 | void AliMUONTrackReconstructor::UpdateHitForRecAtHit() |
b8dc484b | 1497 | { |
1498 | // Set cluster parameterss after track fitting. Fill fHitForRecAtHit of AliMUONTrack's | |
1499 | AliMUONTrack *track; | |
1500 | AliMUONTrackHit *trackHit; | |
1501 | AliMUONHitForRec *hitForRec; | |
1502 | track = (AliMUONTrack*) fRecTracksPtr->First(); | |
1503 | while (track) { | |
1504 | trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->First(); | |
1505 | while (trackHit) { | |
1506 | hitForRec = trackHit->GetHitForRecPtr(); | |
1507 | track->AddHitForRecAtHit(hitForRec); | |
1508 | trackHit = (AliMUONTrackHit*) (track->GetTrackHitsPtr())->After(trackHit); | |
1509 | } // trackHit | |
1510 | track = (AliMUONTrack*) fRecTracksPtr->After(track); | |
1511 | } // track | |
1512 | return; | |
1513 | } | |
1514 | ||
b035a148 | 1515 | //__________________________________________________________________________ |
29f1b13a | 1516 | void AliMUONTrackReconstructor::FillMUONTrack() |
b035a148 | 1517 | { |
1518 | // Set track parameters at hits for Kalman track. Fill fTrackParamAtHit of AliMUONTrack's | |
1519 | AliMUONTrackK *track; | |
1520 | track = (AliMUONTrackK*) fRecTracksPtr->First(); | |
1521 | while (track) { | |
1522 | track->FillMUONTrack(); | |
1523 | track = (AliMUONTrackK*) fRecTracksPtr->After(track); | |
1524 | } | |
1525 | return; | |
1526 | } | |
1527 | ||
8429a5e4 | 1528 | //__________________________________________________________________________ |
29f1b13a | 1529 | void AliMUONTrackReconstructor::EventDump(void) |
04b5ea16 | 1530 | { |
1531 | // Dump reconstructed event (track parameters at vertex and at first hit), | |
1532 | // and the particle parameters | |
1533 | ||
1534 | AliMUONTrack *track; | |
1535 | AliMUONTrackParam *trackParam, *trackParam1; | |
04b5ea16 | 1536 | Double_t bendingSlope, nonBendingSlope, pYZ; |
1537 | Double_t pX, pY, pZ, x, y, z, c; | |
1538 | Int_t np, trackIndex, nTrackHits; | |
1539 | ||
b840996b | 1540 | AliDebug(1,"****** enter EventDump ******"); |
1541 | AliDebug(1, Form("Number of Reconstructed tracks : %d", fNRecTracks)); | |
1542 | ||
04b5ea16 | 1543 | fRecTracksPtr->Compress(); // for simple loop without "Next" since no hole |
1544 | // Loop over reconstructed tracks | |
1545 | for (trackIndex = 0; trackIndex < fNRecTracks; trackIndex++) { | |
83dbc640 | 1546 | if (fTrackMethod != 1) continue; //AZ - skip the rest for now |
b840996b | 1547 | AliDebug(1, Form("track number: %d", trackIndex)); |
04b5ea16 | 1548 | // function for each track for modularity ???? |
1549 | track = (AliMUONTrack*) ((*fRecTracksPtr)[trackIndex]); | |
1550 | nTrackHits = track->GetNTrackHits(); | |
b840996b | 1551 | AliDebug(1, Form("Number of track hits: %d ", nTrackHits)); |
04b5ea16 | 1552 | // track parameters at Vertex |
1553 | trackParam = track->GetTrackParamAtVertex(); | |
1554 | x = trackParam->GetNonBendingCoor(); | |
1555 | y = trackParam->GetBendingCoor(); | |
1556 | z = trackParam->GetZ(); | |
1557 | bendingSlope = trackParam->GetBendingSlope(); | |
1558 | nonBendingSlope = trackParam->GetNonBendingSlope(); | |
1559 | pYZ = 1/TMath::Abs(trackParam->GetInverseBendingMomentum()); | |
1560 | pZ = pYZ/TMath::Sqrt(1+bendingSlope*bendingSlope); | |
1561 | pX = pZ * nonBendingSlope; | |
1562 | pY = pZ * bendingSlope; | |
1563 | c = TMath::Sign(1.0, trackParam->GetInverseBendingMomentum()); | |
b840996b | 1564 | AliDebug(1, Form("Track parameters at Vertex z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n", |
1565 | z, x, y, pX, pY, pZ, c)); | |
04b5ea16 | 1566 | |
1567 | // track parameters at first hit | |
1568 | trackParam1 = ((AliMUONTrackHit*) | |
1569 | (track->GetTrackHitsPtr()->First()))->GetTrackParam(); | |
1570 | x = trackParam1->GetNonBendingCoor(); | |
1571 | y = trackParam1->GetBendingCoor(); | |
1572 | z = trackParam1->GetZ(); | |
1573 | bendingSlope = trackParam1->GetBendingSlope(); | |
1574 | nonBendingSlope = trackParam1->GetNonBendingSlope(); | |
1575 | pYZ = 1/TMath::Abs(trackParam1->GetInverseBendingMomentum()); | |
1576 | pZ = pYZ/TMath::Sqrt(1.0 + bendingSlope * bendingSlope); | |
1577 | pX = pZ * nonBendingSlope; | |
1578 | pY = pZ * bendingSlope; | |
1579 | c = TMath::Sign(1.0, trackParam1->GetInverseBendingMomentum()); | |
b840996b | 1580 | AliDebug(1, Form("track parameters at z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n", |
1581 | z, x, y, pX, pY, pZ, c)); | |
04b5ea16 | 1582 | } |
1583 | // informations about generated particles | |
5d12ce38 | 1584 | np = gAlice->GetMCApp()->GetNtrack(); |
04b5ea16 | 1585 | printf(" **** number of generated particles: %d \n", np); |
1586 | ||
88cb7938 | 1587 | // for (Int_t iPart = 0; iPart < np; iPart++) { |
1588 | // p = gAlice->Particle(iPart); | |
1589 | // printf(" particle %d: type= %d px= %f py= %f pz= %f pdg= %d\n", | |
1590 | // iPart, p->GetPdgCode(), p->Px(), p->Py(), p->Pz(), p->GetPdgCode()); | |
1591 | // } | |
a9e2aefa | 1592 | return; |
1593 | } | |
1594 | ||
276c44b7 | 1595 | |
1596 | //__________________________________________________________________________ | |
29f1b13a | 1597 | void AliMUONTrackReconstructor::EventDumpTrigger(void) |
276c44b7 | 1598 | { |
1599 | // Dump reconstructed trigger event | |
1600 | // and the particle parameters | |
1601 | ||
52c9bc11 | 1602 | AliMUONTriggerTrack *triggertrack ; |
1603 | Int_t nTriggerTracks = fMUONData->RecTriggerTracks()->GetEntriesFast(); | |
276c44b7 | 1604 | |
b840996b | 1605 | AliDebug(1, "****** enter EventDumpTrigger ******"); |
1606 | AliDebug(1, Form("Number of Reconstructed tracks : %d ", nTriggerTracks)); | |
1607 | ||
276c44b7 | 1608 | // Loop over reconstructed tracks |
52c9bc11 | 1609 | for (Int_t trackIndex = 0; trackIndex < nTriggerTracks; trackIndex++) { |
1610 | triggertrack = (AliMUONTriggerTrack*)fMUONData->RecTriggerTracks()->At(trackIndex); | |
276c44b7 | 1611 | printf(" trigger track number %i x11=%f y11=%f thetax=%f thetay=%f \n", |
1612 | trackIndex, | |
1613 | triggertrack->GetX11(),triggertrack->GetY11(), | |
1614 | triggertrack->GetThetax(),triggertrack->GetThetay()); | |
1615 | } | |
1616 | } | |
1617 | ||
83dbc640 | 1618 | //__________________________________________________________________________ |
29f1b13a | 1619 | void AliMUONTrackReconstructor::MakeTrackCandidatesK(void) |
83dbc640 | 1620 | { |
1621 | // To make initial tracks for Kalman filter from the list of segments | |
1622 | Int_t istat, iseg; | |
1623 | AliMUONSegment *segment; | |
1624 | AliMUONTrackK *trackK; | |
1625 | ||
b840996b | 1626 | AliDebug(1,"Enter MakeTrackCandidatesK"); |
83dbc640 | 1627 | |
fff097e0 | 1628 | AliMUONTrackK a(this, fHitsForRecPtr); |
83dbc640 | 1629 | // Loop over stations(1...) 5 and 4 |
1630 | for (istat=4; istat>=3; istat--) { | |
1631 | // Loop over segments in the station | |
1632 | for (iseg=0; iseg<fNSegments[istat]; iseg++) { | |
1633 | // Transform segments to tracks and evaluate covariance matrix | |
1634 | segment = (AliMUONSegment*) ((*fSegmentsPtr[istat])[iseg]); | |
cc87ebcd | 1635 | trackK = new ((*fRecTracksPtr)[fNRecTracks++]) AliMUONTrackK(segment); |
83dbc640 | 1636 | } // for (iseg=0;...) |
1637 | } // for (istat=4;...) | |
1638 | return; | |
1639 | } | |
1640 | ||
1641 | //__________________________________________________________________________ | |
29f1b13a | 1642 | void AliMUONTrackReconstructor::FollowTracksK(void) |
83dbc640 | 1643 | { |
1644 | // Follow tracks using Kalman filter | |
30178c30 | 1645 | Bool_t ok; |
b035a148 | 1646 | Int_t icand, ichamBeg = 0, ichamEnd, chamBits; |
83dbc640 | 1647 | Double_t zDipole1, zDipole2; |
1648 | AliMUONTrackK *trackK; | |
1649 | AliMUONHitForRec *hit; | |
1650 | AliMUONRawCluster *clus; | |
1651 | TClonesArray *rawclusters; | |
83dbc640 | 1652 | clus = 0; rawclusters = 0; |
1653 | ||
b035a148 | 1654 | zDipole1 = GetSimpleBPosition() + GetSimpleBLength()/2; |
1655 | zDipole2 = zDipole1 - GetSimpleBLength(); | |
83dbc640 | 1656 | |
1657 | // Print hits | |
b035a148 | 1658 | trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[0]); |
1659 | if (trackK->DebugLevel() > 0) { | |
1660 | for (Int_t i1=0; i1<fNHitsForRec; i1++) { | |
1661 | hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[i1]); | |
29fc2c86 | 1662 | //if (hit->GetTTRTrack() > 1 || hit->GetTrackRefSignal() == 0) continue; |
b035a148 | 1663 | printf(" Hit # %d %10.4f %10.4f %10.4f", |
1664 | hit->GetChamberNumber(), hit->GetBendingCoor(), | |
1665 | hit->GetNonBendingCoor(), hit->GetZ()); | |
29fc2c86 | 1666 | if (fRecTrackRefHits) { |
1667 | // from track ref hits | |
1668 | printf(" %3d %3d \n", hit->GetTrackRefSignal(), hit->GetTTRTrack()); | |
b035a148 | 1669 | } else { |
1670 | // from raw clusters | |
1671 | rawclusters = fMUONData->RawClusters(hit->GetChamberNumber()); | |
1672 | clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(hit-> | |
1673 | GetHitNumber()); | |
cc87ebcd | 1674 | printf(" %d", clus->GetTrack(1)); |
1675 | if (clus->GetTrack(2) != -1) printf(" %d \n", clus->GetTrack(2)); | |
b035a148 | 1676 | else printf("\n"); |
29fc2c86 | 1677 | } // if (fRecTrackRefHits) |
83dbc640 | 1678 | } |
b035a148 | 1679 | } // if (trackK->DebugLevel() > 0) |
83dbc640 | 1680 | |
1681 | icand = -1; | |
b035a148 | 1682 | Int_t nSeeds; |
1683 | nSeeds = fNRecTracks; // starting number of seeds | |
83dbc640 | 1684 | // Loop over track candidates |
1685 | while (icand < fNRecTracks-1) { | |
1686 | icand ++; | |
b035a148 | 1687 | if (trackK->DebugLevel()>0) cout << " *** Kalman track candidate No. " << icand << endl; |
83dbc640 | 1688 | trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]); |
b035a148 | 1689 | if (trackK->GetRecover() < 0) continue; // failed track |
83dbc640 | 1690 | |
1691 | // Discard candidate which will produce the double track | |
b035a148 | 1692 | /* |
83dbc640 | 1693 | if (icand > 0) { |
30178c30 | 1694 | ok = CheckCandidateK(icand,nSeeds); |
1695 | if (!ok) { | |
b035a148 | 1696 | trackK->SetRecover(-1); // mark candidate to be removed |
1697 | continue; | |
83dbc640 | 1698 | } |
1699 | } | |
b035a148 | 1700 | */ |
83dbc640 | 1701 | |
30178c30 | 1702 | ok = kTRUE; |
83dbc640 | 1703 | if (trackK->GetRecover() == 0) hit = (AliMUONHitForRec*) |
fff097e0 | 1704 | trackK->GetTrackHits()->Last(); // last hit |
b035a148 | 1705 | else hit = trackK->GetHitLastOk(); // hit where track stopped |
1706 | if (hit) ichamBeg = hit->GetChamberNumber(); | |
83dbc640 | 1707 | ichamEnd = 0; |
1708 | // Check propagation direction | |
cc87ebcd | 1709 | if (!hit) { ichamBeg = ichamEnd; AliFatal(" ??? "); } |
b035a148 | 1710 | else if (trackK->GetTrackDir() < 0) { |
83dbc640 | 1711 | ichamEnd = 9; // forward propagation |
30178c30 | 1712 | ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kFALSE,zDipole1,zDipole2); |
1713 | if (ok) { | |
83dbc640 | 1714 | ichamBeg = ichamEnd; |
1715 | ichamEnd = 6; // backward propagation | |
1716 | // Change weight matrix and zero fChi2 for backpropagation | |
1717 | trackK->StartBack(); | |
b035a148 | 1718 | trackK->SetTrackDir(1); |
30178c30 | 1719 | ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kTRUE,zDipole1,zDipole2); |
83dbc640 | 1720 | ichamBeg = ichamEnd; |
1721 | ichamEnd = 0; | |
1722 | } | |
1723 | } else { | |
1724 | if (trackK->GetBPFlag()) { | |
1725 | // backpropagation | |
1726 | ichamEnd = 6; // backward propagation | |
1727 | // Change weight matrix and zero fChi2 for backpropagation | |
1728 | trackK->StartBack(); | |
30178c30 | 1729 | ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kTRUE,zDipole1,zDipole2); |
83dbc640 | 1730 | ichamBeg = ichamEnd; |
1731 | ichamEnd = 0; | |
1732 | } | |
1733 | } | |
1734 | ||
30178c30 | 1735 | if (ok) { |
b035a148 | 1736 | trackK->SetTrackDir(1); |
83dbc640 | 1737 | trackK->SetBPFlag(kFALSE); |
30178c30 | 1738 | ok = trackK->KalmanFilter(ichamBeg,ichamEnd,kFALSE,zDipole1,zDipole2); |
83dbc640 | 1739 | } |
b035a148 | 1740 | if (!ok) { trackK->SetRecover(-1); continue; } // mark candidate to be removed |
1741 | ||
1742 | // Apply smoother | |
1743 | if (trackK->GetRecover() >= 0) { | |
1744 | ok = trackK->Smooth(); | |
1745 | if (!ok) trackK->SetRecover(-1); // mark candidate to be removed | |
1746 | } | |
83dbc640 | 1747 | |
1748 | // Majority 3 of 4 in first 2 stations | |
b035a148 | 1749 | if (!ok) continue; |
83dbc640 | 1750 | chamBits = 0; |
f29ba3e1 | 1751 | Double_t chi2max = 0; |
83dbc640 | 1752 | for (Int_t i=0; i<trackK->GetNTrackHits(); i++) { |
fff097e0 | 1753 | hit = (AliMUONHitForRec*) (*trackK->GetTrackHits())[i]; |
b035a148 | 1754 | chamBits |= BIT(hit->GetChamberNumber()); |
f29ba3e1 | 1755 | if (trackK->GetChi2PerPoint(i) > chi2max) chi2max = trackK->GetChi2PerPoint(i); |
83dbc640 | 1756 | } |
f29ba3e1 | 1757 | if (!((chamBits&3)==3 || (chamBits>>2&3)==3) && chi2max > 25) { |
b035a148 | 1758 | //trackK->Recover(); |
1759 | trackK->SetRecover(-1); //mark candidate to be removed | |
1760 | continue; | |
1761 | } | |
1762 | if (ok) trackK->SetTrackQuality(0); // compute "track quality" | |
83dbc640 | 1763 | } // while |
1764 | ||
1765 | for (Int_t i=0; i<fNRecTracks; i++) { | |
1766 | trackK = (AliMUONTrackK*) ((*fRecTracksPtr)[i]); | |
1767 | if (trackK->GetRecover() < 0) fRecTracksPtr->RemoveAt(i); | |
1768 | } | |
1769 | ||
1770 | // Compress TClonesArray | |
1771 | fRecTracksPtr->Compress(); | |
1772 | fNRecTracks = fRecTracksPtr->GetEntriesFast(); | |
1773 | return; | |
1774 | } | |
1775 | ||
1776 | //__________________________________________________________________________ | |
29f1b13a | 1777 | Bool_t AliMUONTrackReconstructor::CheckCandidateK(Int_t icand, Int_t nSeeds) const |
83dbc640 | 1778 | { |
1779 | // Discards track candidate if it will produce the double track (having | |
1780 | // the same seed segment hits as hits of a good track found before) | |
1781 | AliMUONTrackK *track1, *track2; | |
1782 | AliMUONHitForRec *hit1, *hit2, *hit; | |
1783 | ||
1784 | track1 = (AliMUONTrackK*) ((*fRecTracksPtr)[icand]); | |
fff097e0 | 1785 | hit1 = (AliMUONHitForRec*) (*track1->GetTrackHits())[0]; // 1'st hit |
1786 | hit2 = (AliMUONHitForRec*) (*track1->GetTrackHits())[1]; // 2'nd hit | |
83dbc640 | 1787 | |
1788 | for (Int_t i=0; i<icand; i++) { | |
1789 | track2 = (AliMUONTrackK*) ((*fRecTracksPtr)[i]); | |
1790 | //if (track2->GetRecover() < 0) continue; | |
1791 | if (track2->GetRecover() < 0 && icand >= nSeeds) continue; | |
1792 | ||
1793 | if (track1->GetStartSegment() == track2->GetStartSegment()) { | |
1794 | return kFALSE; | |
1795 | } else { | |
1796 | Int_t nSame = 0; | |
1797 | for (Int_t j=0; j<track2->GetNTrackHits(); j++) { | |
fff097e0 | 1798 | hit = (AliMUONHitForRec*) (*track2->GetTrackHits())[j]; |
83dbc640 | 1799 | if (hit == hit1 || hit == hit2) { |
1800 | nSame++; | |
1801 | if (nSame == 2) return kFALSE; | |
1802 | } | |
1803 | } // for (Int_t j=0; | |
1804 | } | |
1805 | } // for (Int_t i=0; | |
1806 | return kTRUE; | |
1807 | } | |
1808 | ||
1809 | //__________________________________________________________________________ | |
29f1b13a | 1810 | void AliMUONTrackReconstructor::RemoveDoubleTracksK(void) |
83dbc640 | 1811 | { |
1812 | // Removes double tracks (sharing more than half of their hits). Keeps | |
1813 | // the track with higher quality | |
1814 | AliMUONTrackK *track1, *track2, *trackToKill; | |
1815 | ||
1816 | // Sort tracks according to their quality | |
1817 | fRecTracksPtr->Sort(); | |
1818 | ||
1819 | // Loop over first track of the pair | |
1820 | track1 = (AliMUONTrackK*) fRecTracksPtr->First(); | |
b035a148 | 1821 | Int_t debug = track1->DebugLevel(); |
83dbc640 | 1822 | while (track1) { |
1823 | // Loop over second track of the pair | |
1824 | track2 = (AliMUONTrackK*) fRecTracksPtr->After(track1); | |
1825 | while (track2) { | |
1826 | // Check whether or not to keep track2 | |
1827 | if (!track2->KeepTrack(track1)) { | |
b035a148 | 1828 | if (debug >= 0) cout << " Killed track: " << 1/(*track2->GetTrackParameters())(4,0) << |
83dbc640 | 1829 | " " << track2->GetTrackQuality() << endl; |
1830 | trackToKill = track2; | |
1831 | track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2); | |
1832 | trackToKill->Kill(); | |
1833 | fRecTracksPtr->Compress(); | |
1834 | } else track2 = (AliMUONTrackK*) fRecTracksPtr->After(track2); | |
1835 | } // track2 | |
1836 | track1 = (AliMUONTrackK*) fRecTracksPtr->After(track1); | |
1837 | } // track1 | |
1838 | ||
1839 | fNRecTracks = fRecTracksPtr->GetEntriesFast(); | |
b035a148 | 1840 | if (debug >= 0) cout << " Number of Kalman tracks: " << fNRecTracks << endl; |
83dbc640 | 1841 | } |
1842 | ||
1843 | //__________________________________________________________________________ | |
29f1b13a | 1844 | void AliMUONTrackReconstructor::GoToVertex(void) |
83dbc640 | 1845 | { |
1846 | // Propagates track to the vertex thru absorber | |
1847 | // (using Branson correction for now) | |
1848 | ||
1849 | Double_t zVertex; | |
1850 | zVertex = 0; | |
1851 | for (Int_t i=0; i<fNRecTracks; i++) { | |
1852 | //((AliMUONTrackK*)(*fRecTracksPtr)[i])->Branson(); | |
83dbc640 | 1853 | ((AliMUONTrackK*)(*fRecTracksPtr)[i])->SetTrackQuality(1); // compute Chi2 |
b035a148 | 1854 | //((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToZ(zVertex); // w/out absorber |
1855 | ((AliMUONTrackK*)(*fRecTracksPtr)[i])->GoToVertex(1); // with absorber | |
83dbc640 | 1856 | } |
1857 | } | |
b035a148 | 1858 | |
1859 | //__________________________________________________________________________ | |
29f1b13a | 1860 | void AliMUONTrackReconstructor::SetTrackMethod(Int_t iTrackMethod) |
b035a148 | 1861 | { |
1862 | // Set track method and recreate track container if necessary | |
1863 | ||
cc87ebcd | 1864 | fTrackMethod = TMath::Min (iTrackMethod, 3); |
1865 | fTrackMethod = TMath::Max (fTrackMethod, 1); | |
1866 | if (fTrackMethod != 1) { | |
b035a148 | 1867 | if (fRecTracksPtr) delete fRecTracksPtr; |
1868 | fRecTracksPtr = new TClonesArray("AliMUONTrackK", 10); | |
cc87ebcd | 1869 | if (fTrackMethod == 2) cout << " *** Tracking with the Kalman filter *** " << endl; |
1870 | else cout << " *** Combined cluster / track finder ***" << endl; | |
fff097e0 | 1871 | } else cout << " *** Original tracking *** " << endl; |
b035a148 | 1872 | |
1873 | } |