]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/AliMUONTrackReconstructor.cxx
removing unneeded include
[u/mrichter/AliRoot.git] / MUON / AliMUONTrackReconstructor.cxx
CommitLineData
a9e2aefa 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
cc87ebcd 16/* $Id$ */
17
7945aae7 18/// \class AliMUONTrackReconstructor
19/// MUON track reconstructor using the original method
20///
21/// This class contains as data:
22/// - the parameters for the track reconstruction
23///
24/// It contains as methods, among others:
25/// - MakeTracks to build the tracks
26///
a9e2aefa 27
ae17f568 28#include <stdlib.h>
30178c30 29#include <Riostream.h>
cc87ebcd 30#include <TMatrixD.h>
a9e2aefa 31
2457f726 32#include "AliMUONVTrackReconstructor.h"
29f1b13a 33#include "AliMUONTrackReconstructor.h"
5e671e06 34#include "AliMUONData.h"
29fc2c86 35#include "AliMUONConstants.h"
a9e2aefa 36#include "AliMUONRawCluster.h"
2457f726 37#include "AliMUONHitForRec.h"
3831f268 38#include "AliMUONSegment.h"
a9e2aefa 39#include "AliMUONTrack.h"
37827b29 40#include "AliMUONTrackParam.h"
41#include "AliMUONTrackExtrap.h"
8c343c7c 42#include "AliLog.h"
de2cd600 43#include <TVirtualFitter.h>
a9e2aefa 44
de2cd600 45// Functions to be minimized with Minuit
46void TrackChi2(Int_t &NParam, Double_t *Gradient, Double_t &Chi2, Double_t *Param, Int_t Flag);
47void TrackChi2MCS(Int_t &NParam, Double_t *Gradient, Double_t &Chi2, Double_t *Param, Int_t Flag);
48
49void mnvertLocal(Double_t* a, Int_t l, Int_t m, Int_t n, Int_t& ifail);
50
51Double_t MultipleScatteringAngle2(AliMUONTrackParam *param);
52
7945aae7 53/// \cond CLASSIMP
29f1b13a 54ClassImp(AliMUONTrackReconstructor) // Class implementation in ROOT context
7945aae7 55/// \endcond
a9e2aefa 56
2457f726 57//************* Defaults parameters for reconstruction
58const Double_t AliMUONTrackReconstructor::fgkDefaultMaxChi2 = 100.0;
59
60TVirtualFitter* AliMUONTrackReconstructor::fgFitter = 0x0;
61
3bc8b580 62//__________________________________________________________________________
cc9e7528 63AliMUONTrackReconstructor::AliMUONTrackReconstructor(AliMUONData* data)
2457f726 64 : AliMUONVTrackReconstructor(data),
65 fMaxChi2(fgkDefaultMaxChi2)
a9e2aefa 66{
2457f726 67 /// Constructor for class AliMUONTrackReconstructor
de2cd600 68
a9e2aefa 69 // Memory allocation for the TClonesArray of reconstructed tracks
a9e2aefa 70 fRecTracksPtr = new TClonesArray("AliMUONTrack", 10);
a9e2aefa 71}
9cbdf048 72
a9e2aefa 73 //__________________________________________________________________________
29f1b13a 74AliMUONTrackReconstructor::~AliMUONTrackReconstructor(void)
a9e2aefa 75{
2457f726 76 /// Destructor for class AliMUONTrackReconstructor
de2cd600 77 delete fRecTracksPtr;
a9e2aefa 78}
a9e2aefa 79
80 //__________________________________________________________________________
2457f726 81void AliMUONTrackReconstructor::AddHitsForRecFromRawClusters()
a9e2aefa 82{
2457f726 83 /// To add to the list of hits for reconstruction all the raw clusters
e6b25a6e 84 TTree *treeR;
a9e2aefa 85 AliMUONHitForRec *hitForRec;
86 AliMUONRawCluster *clus;
f69d51b5 87 Int_t iclus, nclus, nTRentries;
a9e2aefa 88 TClonesArray *rawclusters;
b840996b 89 AliDebug(1,"Enter AddHitsForRecFromRawClusters");
2457f726 90
e6b25a6e 91 treeR = fMUONData->TreeR();
92 if (!treeR) {
93 AliError("TreeR must be loaded");
2457f726 94 exit(0);
f69d51b5 95 }
e6b25a6e 96
97 nTRentries = Int_t(treeR->GetEntries());
98
99 if (!(fMUONData->IsRawClusterBranchesInTree())) {
100 AliError(Form("RawCluster information is not avalaible, nTRentries = %d not equal to 1",nTRentries));
101 exit(0);
102 }
103
104 fMUONData->SetTreeAddress("RC");
2457f726 105 fMUONData->GetRawClusters(); // only one entry
106
a9e2aefa 107 // Loop over tracking chambers
190f97ea 108 for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
a9e2aefa 109 // number of HitsForRec to 0 for the chamber
110 fNHitsForRecPerChamber[ch] = 0;
111 // index of first HitForRec for the chamber
112 if (ch == 0) fIndexOfFirstHitForRecPerChamber[ch] = 0;
113 else fIndexOfFirstHitForRecPerChamber[ch] = fNHitsForRec;
52c9bc11 114 rawclusters =fMUONData->RawClusters(ch);
a9e2aefa 115 nclus = (Int_t) (rawclusters->GetEntries());
116 // Loop over (cathode correlated) raw clusters
117 for (iclus = 0; iclus < nclus; iclus++) {
118 clus = (AliMUONRawCluster*) rawclusters->UncheckedAt(iclus);
119 // new AliMUONHitForRec from raw cluster
120 // and increment number of AliMUONHitForRec's (total and in chamber)
121 hitForRec = new ((*fHitsForRecPtr)[fNHitsForRec]) AliMUONHitForRec(clus);
122 fNHitsForRec++;
123 (fNHitsForRecPerChamber[ch])++;
124 // more information into HitForRec
fff097e0 125 hitForRec->SetBendingReso2(clus->GetErrY() * clus->GetErrY());
126 hitForRec->SetNonBendingReso2(clus->GetErrX() * clus->GetErrX());
a9e2aefa 127 // original raw cluster
128 hitForRec->SetChamberNumber(ch);
129 hitForRec->SetHitNumber(iclus);
ba5b68db 130 // Z coordinate of the raw cluster (cm)
ba12c242 131 hitForRec->SetZ(clus->GetZ(0));
80e33fc5 132
133 StdoutToAliDebug(3,
134 cout << "Chamber " << ch <<
135 " raw cluster " << iclus << " : " << endl;
136 clus->Print("full");
137 cout << "AliMUONHitForRec number (1...): " << fNHitsForRec << endl;
138 hitForRec->Print("full");
139 );
a9e2aefa 140 } // end of cluster loop
141 } // end of chamber loop
cc87ebcd 142 SortHitsForRecWithIncreasingChamber();
2457f726 143
144 AliDebug(1,"End of AddHitsForRecFromRawClusters");
145 if (AliLog::GetGlobalDebugLevel() > 0) {
146 AliDebug(1, Form("NHitsForRec: %d",fNHitsForRec));
147 for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
148 AliDebug(1, Form("Chamber(0...): %d",ch));
149 AliDebug(1, Form("NHitsForRec: %d", fNHitsForRecPerChamber[ch]));
150 AliDebug(1, Form("Index(first HitForRec): %d", fIndexOfFirstHitForRecPerChamber[ch]));
151 for (Int_t hit = fIndexOfFirstHitForRecPerChamber[ch];
152 hit < fIndexOfFirstHitForRecPerChamber[ch] + fNHitsForRecPerChamber[ch];
153 hit++) {
154 AliDebug(1, Form("HitForRec index(0...): %d",hit));
155 ((*fHitsForRecPtr)[hit])->Dump();
156 }
157 }
158 }
159
a9e2aefa 160 return;
161}
162
163 //__________________________________________________________________________
29f1b13a 164void AliMUONTrackReconstructor::MakeSegments(void)
a9e2aefa 165{
2457f726 166 /// To make the list of segments in all stations,
167 /// from the list of hits to be reconstructed
b840996b 168 AliDebug(1,"Enter MakeSegments");
a9e2aefa 169 // Loop over stations
2457f726 170 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++) MakeSegmentsPerStation(st);
80e33fc5 171
172 StdoutToAliDebug(3,
a9e2aefa 173 cout << "end of MakeSegments" << endl;
80e33fc5 174 for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++)
175 {
176 cout << "station " << st
177 << " has " << fNSegments[st] << " segments:"
178 << endl;
179 for (Int_t seg = 0; seg < fNSegments[st]; seg++)
180 {
181 ((*fSegmentsPtr[st])[seg])->Print();
a9e2aefa 182 }
183 }
80e33fc5 184 );
a9e2aefa 185 return;
186}
187
188 //__________________________________________________________________________
29f1b13a 189void AliMUONTrackReconstructor::MakeTracks(void)
a9e2aefa 190{
2457f726 191 /// To make the tracks from the list of segments and points in all stations
b840996b 192 AliDebug(1,"Enter MakeTracks");
2457f726 193 // Look for candidates from at least 3 aligned points in stations(1..) 4 and 5
194 MakeTrackCandidates();
195 // Follow tracks in stations(1..) 3, 2 and 1
196 FollowTracks();
197 // Remove double tracks
198 RemoveDoubleTracks();
199 UpdateHitForRecAtHit();
276c44b7 200}
201
de2cd600 202 //__________________________________________________________________________
203void AliMUONTrackReconstructor::MakeTrackCandidates(void)
204{
2457f726 205 /// To make track candidates
206 /// with at least 3 aligned points in stations(1..) 4 and 5
207 /// (two Segment's or one Segment and one HitForRec)
de2cd600 208 Int_t begStation, iBegSegment, nbCan1Seg1Hit, nbCan2Seg;
209 AliMUONSegment *begSegment;
210 AliDebug(1,"Enter MakeTrackCandidates");
211 // Loop over stations(1..) 5 and 4 for the beginning segment
212 for (begStation = 4; begStation > 2; begStation--) {
213 // Loop over segments in the beginning station
214 for (iBegSegment = 0; iBegSegment < fNSegments[begStation]; iBegSegment++) {
215 // pointer to segment
216 begSegment = (AliMUONSegment*) ((*fSegmentsPtr[begStation])[iBegSegment]);
217 AliDebug(2,Form("Look for TrackCandidate's with Segment %d in Station(0..) %d", iBegSegment, begStation));
218 // Look for track candidates with two segments,
219 // "begSegment" and all compatible segments in other station.
220 // Only for beginning station(1..) 5
221 // because candidates with 2 segments have to looked for only once.
222 if (begStation == 4)
223 nbCan2Seg = MakeTrackCandidatesWithTwoSegments(begSegment);
224 // Look for track candidates with one segment and one point,
225 // "begSegment" and all compatible HitForRec's in other station.
226 // Only if "begSegment" does not belong already to a track candidate.
227 // Is that a too strong condition ????
228 if (!(begSegment->GetInTrack()))
229 nbCan1Seg1Hit = MakeTrackCandidatesWithOneSegmentAndOnePoint(begSegment);
230 } // for (iBegSegment = 0;...
231 } // for (begStation = 4;...
232 return;
233}
234
a9e2aefa 235 //__________________________________________________________________________
29f1b13a 236Int_t AliMUONTrackReconstructor::MakeTrackCandidatesWithTwoSegments(AliMUONSegment *BegSegment)
a9e2aefa 237{
2457f726 238 /// To make track candidates with two segments in stations(1..) 4 and 5,
239 /// the first segment being pointed to by "BegSegment".
240 /// Returns the number of such track candidates.
a9e2aefa 241 Int_t endStation, iEndSegment, nbCan2Seg;
e516b01d 242 AliMUONSegment *endSegment;
243 AliMUONSegment *extrapSegment = NULL;
a9e2aefa 244 AliMUONTrack *recTrack;
245 Double_t mcsFactor;
b840996b 246 AliDebug(1,"Enter MakeTrackCandidatesWithTwoSegments");
a9e2aefa 247 // Station for the end segment
248 endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
249 // multiple scattering factor corresponding to one chamber
250 mcsFactor = 0.0136 /
251 GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
252 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
253 // linear extrapolation to end station
a9e2aefa 254 // number of candidates with 2 segments to 0
255 nbCan2Seg = 0;
256 // Loop over segments in the end station
257 for (iEndSegment = 0; iEndSegment < fNSegments[endStation]; iEndSegment++) {
258 // pointer to segment
259 endSegment = (AliMUONSegment*) ((*fSegmentsPtr[endStation])[iEndSegment]);
260 // test compatibility between current segment and "extrapSegment"
04b5ea16 261 // 4 because 4 quantities in chi2
e516b01d 262 extrapSegment =
263 BegSegment->CreateSegmentFromLinearExtrapToStation(endSegment->GetZ(), mcsFactor);
a9e2aefa 264 if ((endSegment->
265 NormalizedChi2WithSegment(extrapSegment,
266 fMaxSigma2Distance)) <= 4.0) {
267 // both segments compatible:
268 // make track candidate from "begSegment" and "endSegment"
b840996b 269 AliDebug(2,Form("TrackCandidate with Segment %d in Station(0..) %d", iEndSegment, endStation));
a9e2aefa 270 // flag for both segments in one track:
271 // to be done in track constructor ????
272 BegSegment->SetInTrack(kTRUE);
273 endSegment->SetInTrack(kTRUE);
de2cd600 274 recTrack = new ((*fRecTracksPtr)[fNRecTracks]) AliMUONTrack(BegSegment, endSegment);
275 // Set track parameters at vertex from last stations 4 & 5
276 CalcTrackParamAtVertex(recTrack);
a9e2aefa 277 fNRecTracks++;
b840996b 278 if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump();
a9e2aefa 279 // increment number of track candidates with 2 segments
280 nbCan2Seg++;
281 }
ddf8948c 282 delete extrapSegment; // should not delete HitForRec's it points to !!!!
a9e2aefa 283 } // for (iEndSegment = 0;...
a9e2aefa 284 return nbCan2Seg;
285}
286
287 //__________________________________________________________________________
29f1b13a 288Int_t AliMUONTrackReconstructor::MakeTrackCandidatesWithOneSegmentAndOnePoint(AliMUONSegment *BegSegment)
a9e2aefa 289{
2457f726 290 /// To make track candidates with one segment and one point
291 /// in stations(1..) 4 and 5,
292 /// the segment being pointed to by "BegSegment".
a9e2aefa 293 Int_t ch, ch1, ch2, endStation, iHit, iHitMax, iHitMin, nbCan1Seg1Hit;
e516b01d 294 AliMUONHitForRec *extrapHitForRec= NULL;
295 AliMUONHitForRec *hit;
a9e2aefa 296 AliMUONTrack *recTrack;
297 Double_t mcsFactor;
b840996b 298 AliDebug(1,"Enter MakeTrackCandidatesWithOneSegmentAndOnePoint");
a9e2aefa 299 // station for the end point
300 endStation = 7 - (BegSegment->GetHitForRec1())->GetChamberNumber() / 2;
301 // multiple scattering factor corresponding to one chamber
302 mcsFactor = 0.0136 /
303 GetBendingMomentumFromImpactParam(BegSegment->GetBendingImpact());
304 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
305 // first and second chambers(0..) in the end station
306 ch1 = 2 * endStation;
307 ch2 = ch1 + 1;
308 // number of candidates to 0
309 nbCan1Seg1Hit = 0;
310 // Loop over chambers of the end station
311 for (ch = ch2; ch >= ch1; ch--) {
a9e2aefa 312 // limits for the hit index in the loop
313 iHitMin = fIndexOfFirstHitForRecPerChamber[ch];
314 iHitMax = iHitMin + fNHitsForRecPerChamber[ch];
315 // Loop over HitForRec's in the chamber
316 for (iHit = iHitMin; iHit < iHitMax; iHit++) {
317 // pointer to HitForRec
318 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
319 // test compatibility between current HitForRec and "extrapHitForRec"
04b5ea16 320 // 2 because 2 quantities in chi2
e516b01d 321 // linear extrapolation to chamber
322 extrapHitForRec =
323 BegSegment->CreateHitForRecFromLinearExtrapToChamber( hit->GetZ(), mcsFactor);
a9e2aefa 324 if ((hit->
325 NormalizedChi2WithHitForRec(extrapHitForRec,
326 fMaxSigma2Distance)) <= 2.0) {
327 // both HitForRec's compatible:
328 // make track candidate from begSegment and current HitForRec
b840996b 329 AliDebug(2, Form("TrackCandidate with HitForRec %d in Chamber(0..) %d", iHit, ch));
a9e2aefa 330 // flag for beginning segments in one track:
331 // to be done in track constructor ????
332 BegSegment->SetInTrack(kTRUE);
de2cd600 333 recTrack = new ((*fRecTracksPtr)[fNRecTracks]) AliMUONTrack(BegSegment, hit);
334 // Set track parameters at vertex from last stations 4 & 5
335 CalcTrackParamAtVertex(recTrack);
a9e2aefa 336 // the right place to eliminate "double counting" ???? how ????
337 fNRecTracks++;
b840996b 338 if (AliLog::GetGlobalDebugLevel() > 1) recTrack->RecursiveDump();
a9e2aefa 339 // increment number of track candidates
340 nbCan1Seg1Hit++;
341 }
ddf8948c 342 delete extrapHitForRec;
a9e2aefa 343 } // for (iHit = iHitMin;...
a9e2aefa 344 } // for (ch = ch2;...
345 return nbCan1Seg1Hit;
346}
347
348 //__________________________________________________________________________
2457f726 349void AliMUONTrackReconstructor::CalcTrackParamAtVertex(AliMUONTrack *Track) const
a9e2aefa 350{
2457f726 351 /// Set track parameters at vertex.
352 /// TrackHit's are assumed to be only in stations(1..) 4 and 5,
353 /// and sorted according to increasing Z.
354 /// Parameters are calculated from information in HitForRec's
355 /// of first and last TrackHit's.
de2cd600 356 AliMUONTrackParam *trackParamVertex = Track->GetTrackParamAtVertex(); // pointer to track parameters at vertex
357 // Pointer to HitForRec attached to first TrackParamAtHit
358 AliMUONHitForRec *firstHit = ((AliMUONTrackParam*) (Track->GetTrackParamAtHit()->First()))->GetHitForRecPtr();
359 // Pointer to HitForRec attached to last TrackParamAtHit
360 AliMUONHitForRec *lastHit = ((AliMUONTrackParam*) (Track->GetTrackParamAtHit()->Last()))->GetHitForRecPtr();
361 // Z difference between first and last hits
362 Double_t deltaZ = firstHit->GetZ() - lastHit->GetZ();
363 // bending slope in stations(1..) 4 and 5
364 Double_t bendingSlope = (firstHit->GetBendingCoor() - lastHit->GetBendingCoor()) / deltaZ;
365 trackParamVertex->SetBendingSlope(bendingSlope);
366 // impact parameter
367 Double_t impactParam = firstHit->GetBendingCoor() - bendingSlope * firstHit->GetZ();
368 // signed bending momentum
369 trackParamVertex->SetInverseBendingMomentum(1.0 / GetBendingMomentumFromImpactParam(impactParam));
370 // bending slope at vertex
2457f726 371 trackParamVertex->SetBendingSlope(bendingSlope + impactParam / fSimpleBPosition);
de2cd600 372 // non bending slope
373 trackParamVertex->SetNonBendingSlope((firstHit->GetNonBendingCoor() - lastHit->GetNonBendingCoor()) / deltaZ);
374 // vertex coordinates at (0,0,0)
375 trackParamVertex->SetZ(0.0);
376 trackParamVertex->SetBendingCoor(0.0);
377 trackParamVertex->SetNonBendingCoor(0.0);
a9e2aefa 378}
379
380 //__________________________________________________________________________
29f1b13a 381void AliMUONTrackReconstructor::FollowTracks(void)
a9e2aefa 382{
2457f726 383 /// Follow tracks in stations(1..) 3, 2 and 1
04b5ea16 384 // too long: should be made more modular !!!!
e516b01d 385 AliMUONHitForRec *bestHit, *extrapHit, *hit;
386 AliMUONSegment *bestSegment, *extrapSegment, *segment;
04b5ea16 387 AliMUONTrack *track, *nextTrack;
388 AliMUONTrackParam *trackParam1, trackParam[2], trackParamVertex;
ecfa008b 389 // -1 to avoid compilation warnings
390 Int_t ch = -1, chInStation, chBestHit = -1, iHit, iSegment, station, trackIndex;
04b5ea16 391 Double_t bestChi2, chi2, dZ1, dZ2, dZ3, maxSigma2Distance, mcsFactor;
d0bfce8d 392 Double_t bendingMomentum, chi2Norm = 0.;
cc87ebcd 393
1a38e749 394
04b5ea16 395 // local maxSigma2Distance, for easy increase in testing
396 maxSigma2Distance = fMaxSigma2Distance;
b840996b 397 AliDebug(2,"Enter FollowTracks");
a9e2aefa 398 // Loop over track candidates
04b5ea16 399 track = (AliMUONTrack*) fRecTracksPtr->First();
400 trackIndex = -1;
401 while (track) {
04b5ea16 402 trackIndex++;
403 nextTrack = (AliMUONTrack*) fRecTracksPtr->After(track); // prepare next track
b840996b 404 AliDebug(2,Form("FollowTracks: track candidate(0..): %d", trackIndex));
de2cd600 405 // Fit track candidate from parameters at vertex
406 // -> with 3 parameters (X_vertex and Y_vertex are fixed)
407 // without multiple Coulomb scattering
408 Fit(track,0,0);
b840996b 409 if (AliLog::GetGlobalDebugLevel()> 2) {
a9e2aefa 410 cout << "FollowTracks: track candidate(0..): " << trackIndex
956019b6 411 << " after fit in stations(0..) 3 and 4" << endl;
a9e2aefa 412 track->RecursiveDump();
413 }
414 // Loop over stations(1..) 3, 2 and 1
04b5ea16 415 // something SPECIAL for stations 2 and 1 for majority 3 coincidence ????
416 // otherwise: majority coincidence 2 !!!!
a9e2aefa 417 for (station = 2; station >= 0; station--) {
418 // Track parameters at first track hit (smallest Z)
de2cd600 419 trackParam1 = (AliMUONTrackParam*) (track->GetTrackParamAtHit()->First());
a9e2aefa 420 // extrapolation to station
37827b29 421 AliMUONTrackExtrap::ExtrapToStation(trackParam1, station, trackParam);
79e1e601 422 extrapSegment = new AliMUONSegment(); // empty segment
a9e2aefa 423 // multiple scattering factor corresponding to one chamber
424 // and momentum in bending plane (not total)
425 mcsFactor = 0.0136 * trackParam1->GetInverseBendingMomentum();
426 mcsFactor = fChamberThicknessInX0 * mcsFactor * mcsFactor;
427 // Z difference from previous station
f93b9bd8 428 dZ1 = AliMUONConstants::DefaultChamberZ(2 * station) -
429 AliMUONConstants::DefaultChamberZ(2 * station + 2);
a9e2aefa 430 // Z difference between the two previous stations
f93b9bd8 431 dZ2 = AliMUONConstants::DefaultChamberZ(2 * station + 2) -
432 AliMUONConstants::DefaultChamberZ(2 * station + 4);
04b5ea16 433 // Z difference between the two chambers in the previous station
f93b9bd8 434 dZ3 = AliMUONConstants::DefaultChamberZ(2 * station) -
435 AliMUONConstants::DefaultChamberZ(2 * station + 1);
04b5ea16 436 extrapSegment->SetBendingCoorReso2(fBendingResolution * fBendingResolution);
de2cd600 437 extrapSegment->SetNonBendingCoorReso2(fNonBendingResolution * fNonBendingResolution);
438 extrapSegment->UpdateFromStationTrackParam(trackParam, mcsFactor, dZ1, dZ2, dZ3, station,
439 trackParam1->GetInverseBendingMomentum());
a9e2aefa 440 bestChi2 = 5.0;
441 bestSegment = NULL;
b840996b 442 if (AliLog::GetGlobalDebugLevel() > 2) {
a9e2aefa 443 cout << "FollowTracks: track candidate(0..): " << trackIndex
444 << " Look for segment in station(0..): " << station << endl;
445 }
1a38e749 446
a9e2aefa 447 // Loop over segments in station
448 for (iSegment = 0; iSegment < fNSegments[station]; iSegment++) {
449 // Look for best compatible Segment in station
450 // should consider all possibilities ????
451 // multiple scattering ????
452 // separation in 2 functions: Segment and HitForRec ????
453 segment = (AliMUONSegment*) ((*fSegmentsPtr[station])[iSegment]);
d0bfce8d 454 // correction of corrected segment (fBendingCoor and fNonBendingCoor)
455 // according to real Z value of "segment" and slopes of "extrapSegment"
37827b29 456 AliMUONTrackExtrap::ExtrapToZ(&(trackParam[0]), segment->GetZ());
457 AliMUONTrackExtrap::ExtrapToZ(&(trackParam[1]), segment->GetZ()); // now same as trackParam[0] !?!?!?!?!?!
e516b01d 458 extrapSegment->SetBendingCoor((&(trackParam[0]))->GetBendingCoor());
459 extrapSegment->SetNonBendingCoor((&(trackParam[0]))->GetNonBendingCoor());
460 extrapSegment->SetBendingSlope((&(trackParam[0]))->GetBendingSlope());
461 extrapSegment->SetNonBendingSlope((&(trackParam[0]))->GetNonBendingSlope());
de2cd600 462 chi2 = segment->NormalizedChi2WithSegment(extrapSegment, maxSigma2Distance);
a9e2aefa 463 if (chi2 < bestChi2) {
464 // update best Chi2 and Segment if better found
465 bestSegment = segment;
466 bestChi2 = chi2;
467 }
468 }
469 if (bestSegment) {
470 // best segment found: add it to track candidate
37827b29 471 AliMUONTrackExtrap::ExtrapToZ(&(trackParam[0]), bestSegment->GetZ());
de2cd600 472 track->AddTrackParamAtHit(&(trackParam[0]),bestSegment->GetHitForRec1());
37827b29 473 AliMUONTrackExtrap::ExtrapToZ(&(trackParam[1]), bestSegment->GetZ()); // now same as trackParam[0] !?!?!?!?!?!
de2cd600 474 track->AddTrackParamAtHit(&(trackParam[1]),bestSegment->GetHitForRec2());
b840996b 475 AliDebug(3, Form("FollowTracks: track candidate(0..): %d Added segment in station(0..): %d", trackIndex, station));
476 if (AliLog::GetGlobalDebugLevel()>2) track->RecursiveDump();
de2cd600 477 } else {
a9e2aefa 478 // No best segment found:
479 // Look for best compatible HitForRec in station:
480 // should consider all possibilities ????
481 // multiple scattering ???? do about like for extrapSegment !!!!
79e1e601 482 extrapHit = new AliMUONHitForRec(); // empty hit
a9e2aefa 483 bestChi2 = 3.0;
484 bestHit = NULL;
b840996b 485 AliDebug(3, Form("FollowTracks: track candidate(0..): %d Segment not found, look for hit in station(0..): %d ",
486 trackIndex, station));
487
a9e2aefa 488 // Loop over chambers of the station
d82671a0 489 for (chInStation = 0; chInStation < 2; chInStation++) {
956019b6 490 ch = 2 * station + chInStation;
de2cd600 491 for (iHit = fIndexOfFirstHitForRecPerChamber[ch]; iHit < fIndexOfFirstHitForRecPerChamber[ch]+fNHitsForRecPerChamber[ch]; iHit++) {
a9e2aefa 492 hit = (AliMUONHitForRec*) ((*fHitsForRecPtr)[iHit]);
e516b01d 493 // coordinates of extrapolated hit
37827b29 494 AliMUONTrackExtrap::ExtrapToZ(&(trackParam[chInStation]), hit->GetZ());
de2cd600 495 extrapHit->SetBendingCoor((&(trackParam[chInStation]))->GetBendingCoor());
496 extrapHit->SetNonBendingCoor((&(trackParam[chInStation]))->GetNonBendingCoor());
e516b01d 497 // resolutions from "extrapSegment"
498 extrapHit->SetBendingReso2(extrapSegment->GetBendingCoorReso2());
499 extrapHit->SetNonBendingReso2(extrapSegment->GetNonBendingCoorReso2());
500 // Loop over hits in the chamber
a9e2aefa 501 // condition for hit not already in segment ????
502 chi2 = hit->NormalizedChi2WithHitForRec(extrapHit, maxSigma2Distance);
503 if (chi2 < bestChi2) {
504 // update best Chi2 and HitForRec if better found
505 bestHit = hit;
506 bestChi2 = chi2;
d82671a0 507 chBestHit = chInStation;
a9e2aefa 508 }
509 }
510 }
511 if (bestHit) {
512 // best hit found: add it to track candidate
37827b29 513 AliMUONTrackExtrap::ExtrapToZ(&(trackParam[chBestHit]), bestHit->GetZ());
de2cd600 514 track->AddTrackParamAtHit(&(trackParam[chBestHit]),bestHit);
b840996b 515 if (AliLog::GetGlobalDebugLevel() > 2) {
a9e2aefa 516 cout << "FollowTracks: track candidate(0..): " << trackIndex
517 << " Added hit in station(0..): " << station << endl;
518 track->RecursiveDump();
519 }
de2cd600 520 } else {
8429a5e4 521 // Remove current track candidate
522 // and corresponding TrackHit's, ...
de2cd600 523 fRecTracksPtr->Remove(track);
524 fNRecTracks--;
04b5ea16 525 delete extrapSegment;
d0bfce8d 526 delete extrapHit;
a9e2aefa 527 break; // stop the search for this candidate:
528 // exit from the loop over station
529 }
d0bfce8d 530 delete extrapHit;
a9e2aefa 531 }
04b5ea16 532 delete extrapSegment;
de2cd600 533 // Sort TrackParamAtHit according to increasing Z
534 track->GetTrackParamAtHit()->Sort();
a9e2aefa 535 // Update track parameters at first track hit (smallest Z)
de2cd600 536 trackParam1 = (AliMUONTrackParam*) (track->GetTrackParamAtHit()->First());
d0bfce8d 537 bendingMomentum = 0.;
538 if (TMath::Abs(trackParam1->GetInverseBendingMomentum()) > 0.)
539 bendingMomentum = TMath::Abs(1/(trackParam1->GetInverseBendingMomentum()));
540 // Track removed if bendingMomentum not in window [min, max]
541 if ((bendingMomentum < fMinBendingMomentum) || (bendingMomentum > fMaxBendingMomentum)) {
de2cd600 542 fRecTracksPtr->Remove(track);
543 fNRecTracks--;
d0bfce8d 544 break; // stop the search for this candidate:
545 // exit from the loop over station
546 }
de2cd600 547 // Track fit from parameters at first hit
548 // -> with 5 parameters (momentum and position)
04b5ea16 549 // with multiple Coulomb scattering if all stations
de2cd600 550 if (station == 0) Fit(track,1,1);
956019b6 551 // without multiple Coulomb scattering if not all stations
de2cd600 552 else Fit(track,1,0);
d0bfce8d 553 Double_t numberOfDegFree = (2.0 * track->GetNTrackHits() - 5);
554 if (numberOfDegFree > 0) {
555 chi2Norm = track->GetFitFMin() / numberOfDegFree;
556 } else {
557 chi2Norm = 1.e10;
558 }
559 // Track removed if normalized chi2 too high
560 if (chi2Norm > fMaxChi2) {
de2cd600 561 fRecTracksPtr->Remove(track);
562 fNRecTracks--;
d0bfce8d 563 break; // stop the search for this candidate:
564 // exit from the loop over station
565 }
b840996b 566 if (AliLog::GetGlobalDebugLevel() > 2) {
a9e2aefa 567 cout << "FollowTracks: track candidate(0..): " << trackIndex
568 << " after fit from station(0..): " << station << " to 4" << endl;
569 track->RecursiveDump();
570 }
04b5ea16 571 // Track extrapolation to the vertex through the absorber (Branson)
956019b6 572 // after going through the first station
04b5ea16 573 if (station == 0) {
de2cd600 574 trackParamVertex = *((AliMUONTrackParam*) (track->GetTrackParamAtHit()->First()));
37827b29 575 AliMUONTrackExtrap::ExtrapToVertex(&trackParamVertex, 0., 0., 0.);
04b5ea16 576 track->SetTrackParamAtVertex(&trackParamVertex);
b840996b 577 if (AliLog::GetGlobalDebugLevel() > 0) {
6ae236ee 578 cout << "FollowTracks: track candidate(0..): " << trackIndex
579 << " after extrapolation to vertex" << endl;
580 track->RecursiveDump();
581 }
04b5ea16 582 }
a9e2aefa 583 } // for (station = 2;...
04b5ea16 584 // go really to next track
585 track = nextTrack;
586 } // while (track)
de2cd600 587 // Compression of track array (necessary after Remove)
04b5ea16 588 fRecTracksPtr->Compress();
589 return;
590}
591
de2cd600 592 //__________________________________________________________________________
593void AliMUONTrackReconstructor::Fit(AliMUONTrack *Track, Int_t FitStart, Int_t FitMCS)
594{
2457f726 595 /// Fit the track "Track",
596 /// with or without multiple Coulomb scattering according to "FitMCS",
597 /// starting, according to "FitStart",
598 /// with track parameters at vertex or at the first TrackHit.
de2cd600 599
600 if ((FitStart != 0) && (FitStart != 1)) {
601 cout << "ERROR in AliMUONTrackReconstructor::Fit(...)" << endl;
602 cout << "FitStart = " << FitStart << " is neither 0 nor 1" << endl;
603 exit(0);
604 }
605 if ((FitMCS != 0) && (FitMCS != 1)) {
606 cout << "ERROR in AliMUONTrackReconstructor::Fit(...)" << endl;
607 cout << "FitMCS = " << FitMCS << " is neither 0 nor 1" << endl;
608 exit(0);
609 }
610
611 Double_t arg[1], benC, errorParam, invBenP, lower, nonBenC, upper, x, y;
612 char parName[50];
613 AliMUONTrackParam *trackParam;
614 // Check if Minuit is initialized...
615 fgFitter = TVirtualFitter::Fitter(Track,5);
616 fgFitter->Clear(); // necessary ???? probably yes
617 // how to reset the printout number at every fit ????
618 // is there any risk to leave it like that ????
619 // how to go faster ???? choice of Minuit parameters like EDM ????
620 // choice of function to be minimized according to fFitMCS
621 if (FitMCS == 0) fgFitter->SetFCN(TrackChi2);
622 else fgFitter->SetFCN(TrackChi2MCS);
623 // Switch off printout
624 arg[0] = -1;
625 fgFitter->ExecuteCommand("SET PRINT", arg, 1); // More printing !!!!
626 // No warnings
627 fgFitter->ExecuteCommand("SET NOW", arg, 0);
628 // Parameters according to "fFitStart"
629 // (should be a function to be used at every place where needed ????)
630 if (FitStart == 0) trackParam = Track->GetTrackParamAtVertex();
631 else trackParam = (AliMUONTrackParam*) (Track->GetTrackParamAtHit()->First());
632 // set first 3 Minuit parameters
633 // could be tried with no limits for the search (min=max=0) ????
634 fgFitter->SetParameter(0, "InvBenP",
635 trackParam->GetInverseBendingMomentum(),
636 0.003, -0.4, 0.4);
637 fgFitter->SetParameter(1, "BenS",
638 trackParam->GetBendingSlope(),
639 0.001, -0.5, 0.5);
640 fgFitter->SetParameter(2, "NonBenS",
641 trackParam->GetNonBendingSlope(),
642 0.001, -0.5, 0.5);
643 if (FitStart == 1) {
644 // set last 2 Minuit parameters when we start from first track hit
645 // mandatory limits in Bending to avoid NaN values of parameters
646 fgFitter->SetParameter(3, "X",
647 trackParam->GetNonBendingCoor(),
648 0.03, -500.0, 500.0);
649 // mandatory limits in non Bending to avoid NaN values of parameters
650 fgFitter->SetParameter(4, "Y",
651 trackParam->GetBendingCoor(),
652 0.10, -500.0, 500.0);
653 }
654 // search without gradient calculation in the function
655 fgFitter->ExecuteCommand("SET NOGRADIENT", arg, 0);
656 // minimization
657 fgFitter->ExecuteCommand("MINIMIZE", arg, 0);
658 // exit from Minuit
659 // fgFitter->ExecuteCommand("EXIT", arg, 0); // necessary ????
660 // get results into "invBenP", "benC", "nonBenC" ("x", "y")
661 fgFitter->GetParameter(0, parName, invBenP, errorParam, lower, upper);
662 trackParam->SetInverseBendingMomentum(invBenP);
663 fgFitter->GetParameter(1, parName, benC, errorParam, lower, upper);
664 trackParam->SetBendingSlope(benC);
665 fgFitter->GetParameter(2, parName, nonBenC, errorParam, lower, upper);
666 trackParam->SetNonBendingSlope(nonBenC);
667 if (FitStart == 1) {
668 fgFitter->GetParameter(3, parName, x, errorParam, lower, upper);
669 trackParam->SetNonBendingCoor(x);
670 fgFitter->GetParameter(4, parName, y, errorParam, lower, upper);
671 trackParam->SetBendingCoor(y);
672 }
673 // global result of the fit
2457f726 674 Double_t fedm, errdef, fitFMin;
de2cd600 675 Int_t npari, nparx;
2457f726 676 fgFitter->GetStats(fitFMin, fedm, errdef, npari, nparx);
677 Track->SetFitFMin(fitFMin);
de2cd600 678}
679
680 //__________________________________________________________________________
681void TrackChi2(Int_t &NParam, Double_t * /*Gradient*/, Double_t &Chi2, Double_t *Param, Int_t /*Flag*/)
682{
2457f726 683 /// Return the "Chi2" to be minimized with Minuit for track fitting,
684 /// with "NParam" parameters
685 /// and their current values in array pointed to by "Param".
686 /// Assumes that the track hits are sorted according to increasing Z.
687 /// Track parameters at each TrackHit are updated accordingly.
688 /// Multiple Coulomb scattering is not taken into account
de2cd600 689 AliMUONTrack *trackBeingFitted;
690 AliMUONTrackParam param1;
2457f726 691 AliMUONTrackParam* trackParamAtHit;
692 AliMUONHitForRec* hitForRec;
de2cd600 693 Chi2 = 0.0; // initialize Chi2
694 // copy of track parameters to be fitted
695 trackBeingFitted = (AliMUONTrack*) AliMUONTrackReconstructor::Fitter()->GetObjectFit();
696 // 3 parameters means fit track candidate from parameters at vertex (X_vertex and Y_vertex are fixed)
697 if (NParam == 3) param1 = *(trackBeingFitted->GetTrackParamAtVertex());
698 else param1 = *((AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->First()));
699 // Minuit parameters to be fitted into this copy
700 param1.SetInverseBendingMomentum(Param[0]);
701 param1.SetBendingSlope(Param[1]);
702 param1.SetNonBendingSlope(Param[2]);
703 if (NParam == 5) {
704 param1.SetNonBendingCoor(Param[3]);
705 param1.SetBendingCoor(Param[4]);
706 }
707 // Follow track through all planes of track hits
2457f726 708 trackParamAtHit = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->First());
709 while (trackParamAtHit) {
710 hitForRec = trackParamAtHit->GetHitForRecPtr();
711 // extrapolation to the plane of the hitForRec attached to the current trackParamAtHit
37827b29 712 AliMUONTrackExtrap::ExtrapToZ(&param1, hitForRec->GetZ());
de2cd600 713 // update track parameters of the current hit
2457f726 714 trackParamAtHit->SetTrackParam(param1);
de2cd600 715 // Increment Chi2
716 // done hit per hit, with hit resolution,
717 // and not with point and angle like in "reco_muon.F" !!!!
718 // Needs to add multiple scattering contribution ????
2457f726 719 Double_t dX = hitForRec->GetNonBendingCoor() - param1.GetNonBendingCoor();
720 Double_t dY = hitForRec->GetBendingCoor() - param1.GetBendingCoor();
721 Chi2 = Chi2 + dX * dX / hitForRec->GetNonBendingReso2() + dY * dY / hitForRec->GetBendingReso2();
722 trackParamAtHit = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->After(trackParamAtHit));
de2cd600 723 }
724}
725
726 //__________________________________________________________________________
727void TrackChi2MCS(Int_t &NParam, Double_t * /*Gradient*/, Double_t &Chi2, Double_t *Param, Int_t /*Flag*/)
728{
2457f726 729 /// Return the "Chi2" to be minimized with Minuit for track fitting,
730 /// with "NParam" parameters
731 /// and their current values in array pointed to by "Param".
732 /// Assumes that the track hits are sorted according to increasing Z.
733 /// Track parameters at each TrackHit are updated accordingly.
734 /// Multiple Coulomb scattering is taken into account with covariance matrix.
de2cd600 735 AliMUONTrack *trackBeingFitted;
736 AliMUONTrackParam param1;
2457f726 737 AliMUONTrackParam* trackParamAtHit;
738 AliMUONHitForRec* hitForRec;
de2cd600 739 Chi2 = 0.0; // initialize Chi2
740 // copy of track parameters to be fitted
741 trackBeingFitted = (AliMUONTrack*) AliMUONTrackReconstructor::Fitter()->GetObjectFit();
742 // 3 parameters means fit track candidate from parameters at vertex (X_vertex and Y_vertex are fixed)
743 if (NParam == 3) param1 = *(trackBeingFitted->GetTrackParamAtVertex());
744 else param1 = *((AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->First()));
745 // Minuit parameters to be fitted into this copy
746 param1.SetInverseBendingMomentum(Param[0]);
747 param1.SetBendingSlope(Param[1]);
748 param1.SetNonBendingSlope(Param[2]);
749 if (NParam == 5) {
750 param1.SetNonBendingCoor(Param[3]);
751 param1.SetBendingCoor(Param[4]);
752 }
753
754 Int_t chCurrent, chPrev = 0, hitNumber, hitNumber1, hitNumber2, hitNumber3;
755 Double_t z1, z2, z3;
2457f726 756 AliMUONTrackParam *trackParamAtHit1, *trackParamAtHit2, *trackParamAtHit3;
757 AliMUONHitForRec *hitForRec1, *hitForRec2;
de2cd600 758 Double_t hbc1, hbc2, pbc1, pbc2;
759 Double_t hnbc1, hnbc2, pnbc1, pnbc2;
760 Int_t numberOfHit = trackBeingFitted->GetNTrackHits();
761 TMatrixD *covBending = new TMatrixD(numberOfHit, numberOfHit);
762 TMatrixD *covNonBending = new TMatrixD(numberOfHit, numberOfHit);
763 Double_t *msa2 = new Double_t[numberOfHit];
764
765 // Predicted coordinates and multiple scattering angles are first calculated
766 for (hitNumber = 0; hitNumber < numberOfHit; hitNumber++) {
2457f726 767 trackParamAtHit = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->UncheckedAt(hitNumber));
768 hitForRec = trackParamAtHit->GetHitForRecPtr();
769 // extrapolation to the plane of the hitForRec attached to the current trackParamAtHit
37827b29 770 AliMUONTrackExtrap::ExtrapToZ(&param1, hitForRec->GetZ());
de2cd600 771 // update track parameters of the current hit
2457f726 772 trackParamAtHit->SetTrackParam(param1);
de2cd600 773 // square of multiple scattering angle at current hit, with one chamber
774 msa2[hitNumber] = MultipleScatteringAngle2(&param1);
775 // correction for eventual missing hits or multiple hits in a chamber,
776 // according to the number of chambers
777 // between the current hit and the previous one
2457f726 778 chCurrent = hitForRec->GetChamberNumber();
de2cd600 779 if (hitNumber > 0) msa2[hitNumber] = msa2[hitNumber] * (chCurrent - chPrev);
780 chPrev = chCurrent;
781 }
782
783 // Calculates the covariance matrix
784 for (hitNumber1 = 0; hitNumber1 < numberOfHit; hitNumber1++) {
2457f726 785 trackParamAtHit1 = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->UncheckedAt(hitNumber1));
786 hitForRec1 = trackParamAtHit1->GetHitForRecPtr();
787 z1 = hitForRec1->GetZ();
de2cd600 788 for (hitNumber2 = hitNumber1; hitNumber2 < numberOfHit; hitNumber2++) {
2457f726 789 trackParamAtHit2 = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->UncheckedAt(hitNumber2));
790 z2 = trackParamAtHit2->GetHitForRecPtr()->GetZ();
de2cd600 791 // initialization to 0 (diagonal plus upper triangular part)
792 (*covBending)(hitNumber2, hitNumber1) = 0.0;
793 // contribution from multiple scattering in bending plane:
794 // loop over upstream hits
795 for (hitNumber3 = 0; hitNumber3 < hitNumber1; hitNumber3++) {
2457f726 796 trackParamAtHit3 = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->UncheckedAt(hitNumber3));
797 z3 = trackParamAtHit3->GetHitForRecPtr()->GetZ();
de2cd600 798 (*covBending)(hitNumber2, hitNumber1) = (*covBending)(hitNumber2, hitNumber1) + ((z1 - z3) * (z2 - z3) * msa2[hitNumber3]);
799 }
800 // equal contribution from multiple scattering in non bending plane
801 (*covNonBending)(hitNumber2, hitNumber1) = (*covBending)(hitNumber2, hitNumber1);
802 if (hitNumber1 == hitNumber2) {
803 // Diagonal elements: add contribution from position measurements
804 // in bending plane
2457f726 805 (*covBending)(hitNumber2, hitNumber1) = (*covBending)(hitNumber2, hitNumber1) + hitForRec1->GetBendingReso2();
de2cd600 806 // and in non bending plane
2457f726 807 (*covNonBending)(hitNumber2, hitNumber1) = (*covNonBending)(hitNumber2, hitNumber1) + hitForRec1->GetNonBendingReso2();
de2cd600 808 } else {
809 // Non diagonal elements: symmetrization
810 // for bending plane
811 (*covBending)(hitNumber1, hitNumber2) = (*covBending)(hitNumber2, hitNumber1);
812 // and non bending plane
813 (*covNonBending)(hitNumber1, hitNumber2) = (*covNonBending)(hitNumber2, hitNumber1);
814 }
815 } // for (hitNumber2 = hitNumber1;...
816 } // for (hitNumber1 = 0;...
817
818 // Inversion of covariance matrices
819 // with "mnvertLocal", local "mnvert" function of Minuit.
820 // One cannot use directly "mnvert" since "TVirtualFitter" does not know it.
821 // One will have to replace this local function by the right inversion function
822 // from a specialized Root package for symmetric positive definite matrices,
823 // when available!!!!
824 Int_t ifailBending;
825 mnvertLocal(&((*covBending)(0,0)), numberOfHit, numberOfHit, numberOfHit, ifailBending);
826 Int_t ifailNonBending;
827 mnvertLocal(&((*covNonBending)(0,0)), numberOfHit, numberOfHit, numberOfHit, ifailNonBending);
828
829 // It would be worth trying to calculate the inverse of the covariance matrix
830 // only once per fit, since it cannot change much in principle,
831 // and it would save a lot of computing time !!!!
832
833 // Calculates Chi2
834 if ((ifailBending == 0) && (ifailNonBending == 0)) {
835 // with Multiple Scattering if inversion correct
836 for (hitNumber1 = 0; hitNumber1 < numberOfHit ; hitNumber1++) {
2457f726 837 trackParamAtHit1 = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->UncheckedAt(hitNumber1));
838 hitForRec1 = trackParamAtHit1->GetHitForRecPtr();
839 hbc1 = hitForRec1->GetBendingCoor();
840 pbc1 = trackParamAtHit1->GetBendingCoor();
841 hnbc1 = hitForRec1->GetNonBendingCoor();
842 pnbc1 = trackParamAtHit1->GetNonBendingCoor();
de2cd600 843 for (hitNumber2 = 0; hitNumber2 < numberOfHit; hitNumber2++) {
2457f726 844 trackParamAtHit2 = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->UncheckedAt(hitNumber2));
845 hitForRec2 = trackParamAtHit2->GetHitForRecPtr();
846 hbc2 = hitForRec2->GetBendingCoor();
847 pbc2 = trackParamAtHit2->GetBendingCoor();
848 hnbc2 = hitForRec2->GetNonBendingCoor();
849 pnbc2 = trackParamAtHit2->GetNonBendingCoor();
de2cd600 850 Chi2 += ((*covBending)(hitNumber2, hitNumber1) * (hbc1 - pbc1) * (hbc2 - pbc2)) +
851 ((*covNonBending)(hitNumber2, hitNumber1) * (hnbc1 - pnbc1) * (hnbc2 - pnbc2));
852 }
853 }
854 } else {
855 // without Multiple Scattering if inversion impossible
856 for (hitNumber1 = 0; hitNumber1 < numberOfHit ; hitNumber1++) {
2457f726 857 trackParamAtHit1 = (AliMUONTrackParam*) (trackBeingFitted->GetTrackParamAtHit()->UncheckedAt(hitNumber1));
858 hitForRec1 = trackParamAtHit1->GetHitForRecPtr();
859 hbc1 = hitForRec1->GetBendingCoor();
860 pbc1 = trackParamAtHit1->GetBendingCoor();
861 hnbc1 = hitForRec1->GetNonBendingCoor();
862 pnbc1 = trackParamAtHit1->GetNonBendingCoor();
863 Chi2 += ((hbc1 - pbc1) * (hbc1 - pbc1) / hitForRec1->GetBendingReso2()) +
864 ((hnbc1 - pnbc1) * (hnbc1 - pnbc1) / hitForRec1->GetNonBendingReso2());
de2cd600 865 }
866 }
867
868 delete covBending;
869 delete covNonBending;
870 delete [] msa2;
871}
872
873Double_t MultipleScatteringAngle2(AliMUONTrackParam *param)
874{
2457f726 875 /// Returns square of multiple Coulomb scattering angle
876 /// from TrackParamAtHit pointed to by "param"
de2cd600 877 Double_t slopeBending, slopeNonBending, radiationLength, inverseBendingMomentum2, inverseTotalMomentum2;
878 Double_t varMultipleScatteringAngle;
879 // Better implementation in AliMUONTrack class ????
880 slopeBending = param->GetBendingSlope();
881 slopeNonBending = param->GetNonBendingSlope();
882 // thickness in radiation length for the current track,
883 // taking local angle into account
884 radiationLength = AliMUONConstants::DefaultChamberThicknessInX0() *
885 TMath::Sqrt(1.0 + slopeBending*slopeBending + slopeNonBending*slopeNonBending);
886 inverseBendingMomentum2 = param->GetInverseBendingMomentum() * param->GetInverseBendingMomentum();
887 inverseTotalMomentum2 = inverseBendingMomentum2 * (1.0 + slopeBending * slopeBending) /
888 (1.0 + slopeBending *slopeBending + slopeNonBending * slopeNonBending);
889 varMultipleScatteringAngle = 0.0136 * (1.0 + 0.038 * TMath::Log(radiationLength));
890 // The velocity is assumed to be 1 !!!!
891 varMultipleScatteringAngle = inverseTotalMomentum2 * radiationLength * varMultipleScatteringAngle * varMultipleScatteringAngle;
892 return varMultipleScatteringAngle;
893}
894
895//______________________________________________________________________________
896 void mnvertLocal(Double_t *a, Int_t l, Int_t, Int_t n, Int_t &ifail)
897{
2457f726 898///*-*-*-*-*-*-*-*-*-*-*-*Inverts a symmetric matrix*-*-*-*-*-*-*-*-*-*-*-*-*
899///*-* ==========================
900///*-* inverts a symmetric matrix. matrix is first scaled to
901///*-* have all ones on the diagonal (equivalent to change of units)
902///*-* but no pivoting is done since matrix is positive-definite.
903///*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
de2cd600 904
905 // taken from TMinuit package of Root (l>=n)
906 // fVERTs, fVERTq and fVERTpp changed to localVERTs, localVERTq and localVERTpp
907 // Double_t localVERTs[n], localVERTq[n], localVERTpp[n];
908 Double_t * localVERTs = new Double_t[n];
909 Double_t * localVERTq = new Double_t[n];
910 Double_t * localVERTpp = new Double_t[n];
911 // fMaxint changed to localMaxint
912 Int_t localMaxint = n;
913
914 /* System generated locals */
915 Int_t aOffset;
916
917 /* Local variables */
918 Double_t si;
919 Int_t i, j, k, kp1, km1;
920
921 /* Parameter adjustments */
922 aOffset = l + 1;
923 a -= aOffset;
924
925 /* Function Body */
926 ifail = 0;
927 if (n < 1) goto L100;
928 if (n > localMaxint) goto L100;
929//*-*- scale matrix by sqrt of diag elements
930 for (i = 1; i <= n; ++i) {
931 si = a[i + i*l];
932 if (si <= 0) goto L100;
933 localVERTs[i-1] = 1 / TMath::Sqrt(si);
934 }
935 for (i = 1; i <= n; ++i) {
936 for (j = 1; j <= n; ++j) {
937 a[i + j*l] = a[i + j*l]*localVERTs[i-1]*localVERTs[j-1];
938 }
939 }
940//*-*- . . . start main loop . . . .
941 for (i = 1; i <= n; ++i) {
942 k = i;
943//*-*- preparation for elimination step1
944 if (a[k + k*l] != 0) localVERTq[k-1] = 1 / a[k + k*l];
945 else goto L100;
946 localVERTpp[k-1] = 1;
947 a[k + k*l] = 0;
948 kp1 = k + 1;
949 km1 = k - 1;
950 if (km1 < 0) goto L100;
951 else if (km1 == 0) goto L50;
952 else goto L40;
953L40:
954 for (j = 1; j <= km1; ++j) {
955 localVERTpp[j-1] = a[j + k*l];
956 localVERTq[j-1] = a[j + k*l]*localVERTq[k-1];
957 a[j + k*l] = 0;
958 }
959L50:
960 if (k - n < 0) goto L51;
961 else if (k - n == 0) goto L60;
962 else goto L100;
963L51:
964 for (j = kp1; j <= n; ++j) {
965 localVERTpp[j-1] = a[k + j*l];
966 localVERTq[j-1] = -a[k + j*l]*localVERTq[k-1];
967 a[k + j*l] = 0;
968 }
969//*-*- elimination proper
970L60:
971 for (j = 1; j <= n; ++j) {
972 for (k = j; k <= n; ++k) { a[j + k*l] += localVERTpp[j-1]*localVERTq[k-1]; }
973 }
974 }
975//*-*- elements of left diagonal and unscaling
976 for (j = 1; j <= n; ++j) {
977 for (k = 1; k <= j; ++k) {
978 a[k + j*l] = a[k + j*l]*localVERTs[k-1]*localVERTs[j-1];
979 a[j + k*l] = a[k + j*l];
980 }
981 }
982 delete [] localVERTs;
983 delete [] localVERTq;
984 delete [] localVERTpp;
985 return;
986//*-*- failure return
987L100:
988 delete [] localVERTs;
989 delete [] localVERTq;
990 delete [] localVERTpp;
991 ifail = 1;
992} /* mnvertLocal */
993
04b5ea16 994 //__________________________________________________________________________
29f1b13a 995void AliMUONTrackReconstructor::RemoveDoubleTracks(void)
8429a5e4 996{
2457f726 997 /// To remove double tracks.
998 /// Tracks are considered identical
999 /// if they have at least half of their hits in common.
1000 /// Among two identical tracks, one keeps the track with the larger number of hits
1001 /// or, if these numbers are equal, the track with the minimum Chi2.
8429a5e4 1002 AliMUONTrack *track1, *track2, *trackToRemove;
8429a5e4 1003 Int_t hitsInCommon, nHits1, nHits2;
de2cd600 1004 Bool_t removedTrack1;
1005 // Loop over first track of the pair
1006 track1 = (AliMUONTrack*) fRecTracksPtr->First();
1007 while (track1) {
1008 removedTrack1 = kFALSE;
1009 nHits1 = track1->GetNTrackHits();
1010 // Loop over second track of the pair
1011 track2 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1012 while (track2) {
1013 nHits2 = track2->GetNTrackHits();
1014 // number of hits in common between two tracks
1015 hitsInCommon = track1->HitsInCommon(track2);
1016 // check for identical tracks
1017 if ((4 * hitsInCommon) >= (nHits1 + nHits2)) {
1018 // decide which track to remove
1019 if ((nHits1 > nHits2) || ((nHits1 == nHits2) && (track1->GetFitFMin() < track2->GetFitFMin()))) {
1020 // remove track2 and continue the second loop with the track next to track2
1021 trackToRemove = track2;
1022 track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
1023 fRecTracksPtr->Remove(trackToRemove);
1024 fNRecTracks--;
1025 fRecTracksPtr->Compress(); // this is essential to retrieve the TClonesArray afterwards
1026 } else {
1027 // else remove track1 and continue the first loop with the track next to track1
1028 trackToRemove = track1;
1029 track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1030 fRecTracksPtr->Remove(trackToRemove);
1031 fNRecTracks--;
1032 fRecTracksPtr->Compress(); // this is essential to retrieve the TClonesArray afterwards
1033 removedTrack1 = kTRUE;
1034 break;
1035 }
1036 } else track2 = (AliMUONTrack*) fRecTracksPtr->After(track2);
1037 } // track2
1038 if (removedTrack1) continue;
1039 track1 = (AliMUONTrack*) fRecTracksPtr->After(track1);
1040 } // track1
d837040f 1041 return;
1042}
1043
b8dc484b 1044 //__________________________________________________________________________
29f1b13a 1045void AliMUONTrackReconstructor::UpdateHitForRecAtHit()
b8dc484b 1046{
2457f726 1047 /// Set cluster parameters after track fitting. Fill fHitForRecAtHit of AliMUONTrack's
b8dc484b 1048 AliMUONTrack *track;
de2cd600 1049 AliMUONTrackParam *trackParamAtHit;
b8dc484b 1050 track = (AliMUONTrack*) fRecTracksPtr->First();
1051 while (track) {
de2cd600 1052 trackParamAtHit = (AliMUONTrackParam*) (track->GetTrackParamAtHit()->First());
1053 while (trackParamAtHit) {
1054 track->AddHitForRecAtHit(trackParamAtHit->GetHitForRecPtr());
1055 trackParamAtHit = (AliMUONTrackParam*) (track->GetTrackParamAtHit()->After(trackParamAtHit));
1056 }
b8dc484b 1057 track = (AliMUONTrack*) fRecTracksPtr->After(track);
de2cd600 1058 }
b8dc484b 1059 return;
1060}
1061
8429a5e4 1062 //__________________________________________________________________________
29f1b13a 1063void AliMUONTrackReconstructor::EventDump(void)
04b5ea16 1064{
2457f726 1065 /// Dump reconstructed event (track parameters at vertex and at first hit),
1066 /// and the particle parameters
04b5ea16 1067 AliMUONTrack *track;
1068 AliMUONTrackParam *trackParam, *trackParam1;
04b5ea16 1069 Double_t bendingSlope, nonBendingSlope, pYZ;
1070 Double_t pX, pY, pZ, x, y, z, c;
5e671e06 1071 Int_t trackIndex, nTrackHits;
04b5ea16 1072
b840996b 1073 AliDebug(1,"****** enter EventDump ******");
1074 AliDebug(1, Form("Number of Reconstructed tracks : %d", fNRecTracks));
1075
04b5ea16 1076 fRecTracksPtr->Compress(); // for simple loop without "Next" since no hole
1077 // Loop over reconstructed tracks
1078 for (trackIndex = 0; trackIndex < fNRecTracks; trackIndex++) {
b840996b 1079 AliDebug(1, Form("track number: %d", trackIndex));
04b5ea16 1080 // function for each track for modularity ????
1081 track = (AliMUONTrack*) ((*fRecTracksPtr)[trackIndex]);
1082 nTrackHits = track->GetNTrackHits();
b840996b 1083 AliDebug(1, Form("Number of track hits: %d ", nTrackHits));
04b5ea16 1084 // track parameters at Vertex
1085 trackParam = track->GetTrackParamAtVertex();
1086 x = trackParam->GetNonBendingCoor();
1087 y = trackParam->GetBendingCoor();
1088 z = trackParam->GetZ();
1089 bendingSlope = trackParam->GetBendingSlope();
1090 nonBendingSlope = trackParam->GetNonBendingSlope();
1091 pYZ = 1/TMath::Abs(trackParam->GetInverseBendingMomentum());
1092 pZ = pYZ/TMath::Sqrt(1+bendingSlope*bendingSlope);
1093 pX = pZ * nonBendingSlope;
1094 pY = pZ * bendingSlope;
1095 c = TMath::Sign(1.0, trackParam->GetInverseBendingMomentum());
b840996b 1096 AliDebug(1, Form("Track parameters at Vertex z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1097 z, x, y, pX, pY, pZ, c));
04b5ea16 1098
1099 // track parameters at first hit
de2cd600 1100 trackParam1 = (AliMUONTrackParam*) track->GetTrackParamAtHit()->First();
04b5ea16 1101 x = trackParam1->GetNonBendingCoor();
1102 y = trackParam1->GetBendingCoor();
1103 z = trackParam1->GetZ();
1104 bendingSlope = trackParam1->GetBendingSlope();
1105 nonBendingSlope = trackParam1->GetNonBendingSlope();
1106 pYZ = 1/TMath::Abs(trackParam1->GetInverseBendingMomentum());
1107 pZ = pYZ/TMath::Sqrt(1.0 + bendingSlope * bendingSlope);
1108 pX = pZ * nonBendingSlope;
1109 pY = pZ * bendingSlope;
1110 c = TMath::Sign(1.0, trackParam1->GetInverseBendingMomentum());
b840996b 1111 AliDebug(1, Form("track parameters at z= %f: X= %f Y= %f pX= %f pY= %f pZ= %f c= %f\n",
1112 z, x, y, pX, pY, pZ, c));
04b5ea16 1113 }
8f373020 1114 // informations about generated particles NO !!!!!!!!
04b5ea16 1115
88cb7938 1116// for (Int_t iPart = 0; iPart < np; iPart++) {
1117// p = gAlice->Particle(iPart);
1118// printf(" particle %d: type= %d px= %f py= %f pz= %f pdg= %d\n",
1119// iPart, p->GetPdgCode(), p->Px(), p->Py(), p->Pz(), p->GetPdgCode());
1120// }
a9e2aefa 1121 return;
1122}
1123
276c44b7 1124