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