1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 //====================================================================================================================================================
18 // Description of an ALICE muon forward track, combining the information of the Muon Spectrometer and the Muon Forward Tracker
20 // Contact author: antonio.uras@cern.ch
22 //====================================================================================================================================================
25 #include "AliMUONTrack.h"
26 #include "AliMFTCluster.h"
27 #include "AliMUONVCluster.h"
28 #include "AliMUONTrackParam.h"
29 #include "AliMUONTrackExtrap.h"
30 #include "TClonesArray.h"
32 #include "TParticle.h"
33 #include "AliMuonForwardTrack.h"
34 #include "AliMFTConstants.h"
36 ClassImp(AliMuonForwardTrack)
38 //====================================================================================================================================================
40 AliMuonForwardTrack::AliMuonForwardTrack():
45 fNWrongClustersMC(-1),
49 // default constructor
51 for (Int_t iPlane=0; iPlane<AliMFTConstants::fNMaxPlanes; iPlane++) fPlaneExists[iPlane] = kFALSE;
52 for (Int_t iParent=0; iParent<fgkNParentsMax; iParent++) {
53 fParentMCLabel[iParent] = -1;
54 fParentPDGCode[iParent] = 0;
56 fMFTClusters = new TClonesArray("AliMFTCluster");
60 //====================================================================================================================================================
62 AliMuonForwardTrack::AliMuonForwardTrack(AliMUONTrack *MUONTrack):
67 fNWrongClustersMC(-1),
71 SetMUONTrack(MUONTrack);
72 for (Int_t iPlane=0; iPlane<AliMFTConstants::fNMaxPlanes; iPlane++) fPlaneExists[iPlane] = kFALSE;
73 for (Int_t iParent=0; iParent<fgkNParentsMax; iParent++) {
74 fParentMCLabel[iParent] = -1;
75 fParentPDGCode[iParent] = 0;
77 fMFTClusters = new TClonesArray("AliMFTCluster");
81 //====================================================================================================================================================
83 AliMuonForwardTrack::AliMuonForwardTrack(const AliMuonForwardTrack& track):
88 fNWrongClustersMC(track.fNWrongClustersMC),
89 fTrackMCId(track.fTrackMCId)
93 fMUONTrack = new AliMUONTrack(*(track.fMUONTrack));
94 if (track.fMCTrackRef) fMCTrackRef = new TParticle(*(track.fMCTrackRef));
95 fMFTClusters = new TClonesArray(*(track.fMFTClusters));
96 for (Int_t iPlane=0; iPlane<AliMFTConstants::fNMaxPlanes; iPlane++) fPlaneExists[iPlane] = (track.fPlaneExists)[iPlane];
97 for (Int_t iParent=0; iParent<fgkNParentsMax; iParent++) {
98 fParentMCLabel[iParent] = (track.fParentMCLabel)[iParent];
99 fParentPDGCode[iParent] = (track.fParentPDGCode)[iParent];
104 //====================================================================================================================================================
106 AliMuonForwardTrack& AliMuonForwardTrack::operator=(const AliMuonForwardTrack& track) {
108 // Asignment operator
110 // check assignement to self
111 if (this == &track) return *this;
113 // base class assignement
114 AliMUONTrack::operator=(track);
119 fMUONTrack = new AliMUONTrack(*(track.fMUONTrack));
120 if (track.fMCTrackRef) fMCTrackRef = new TParticle(*(track.fMCTrackRef));
121 fMFTClusters = new TClonesArray(*(track.fMFTClusters));
122 fNWrongClustersMC = track.fNWrongClustersMC;
123 fTrackMCId = track.fTrackMCId;
125 for (Int_t iPlane=0; iPlane<AliMFTConstants::fNMaxPlanes; iPlane++) fPlaneExists[iPlane] = (track.fPlaneExists)[iPlane];
126 for (Int_t iParent=0; iParent<fgkNParentsMax; iParent++) {
127 fParentMCLabel[iParent] = (track.fParentMCLabel)[iParent];
128 fParentPDGCode[iParent] = (track.fParentPDGCode)[iParent];
135 //====================================================================================================================================================
137 void AliMuonForwardTrack::SetMUONTrack(AliMUONTrack *MUONTrack) {
140 AliInfo("fMUONTrack already exists, nothing will be done");
144 fMUONTrack = MUONTrack;
145 SetGlobalChi2(fMUONTrack->GetGlobalChi2());
149 //====================================================================================================================================================
151 void AliMuonForwardTrack::SetMCTrackRef(TParticle *MCTrackRef) {
154 AliInfo("fMCTrackRef already exists, nothing will be done");
158 fMCTrackRef = MCTrackRef;
162 //====================================================================================================================================================
164 void AliMuonForwardTrack::AddTrackParamAtMFTCluster(AliMUONTrackParam &trackParam, AliMFTCluster &mftCluster) {
166 AliDebug(1, Form("Before adding: this->fMFTClusters=%p has %d entries", this->fMFTClusters, this->fMFTClusters->GetEntries()));
167 Int_t iMFTCluster = this->fMFTClusters->GetEntries();
168 AliDebug(1, Form("mftCluster->GetX() = %f mftCluster->GetY() = %f mftCluster->GetErrX() = %f cmftCluster->GetErrY() = %f",
169 mftCluster.GetX(), mftCluster.GetY(), mftCluster.GetErrX(), mftCluster.GetErrY()));
170 AliMUONVCluster *muonCluster = (AliMUONVCluster*) mftCluster.CreateMUONCluster();
171 AliDebug(1, Form("Created MUON cluster %p", muonCluster));
172 trackParam.SetUniqueID(iMFTCluster); // we profit of this slot to store the reference to the corresponding MFTCluster
173 AliDebug(1, Form("Now adding muonCluster %p and trackParam %p",muonCluster, &trackParam));
174 AddTrackParamAtCluster(trackParam, *muonCluster, kTRUE);
175 // we pass the parameters this->GetTrackParamAtCluster()->First() to the Kalman Filter algorithm: they will be updated!!
176 Double_t chi2Kalman = RunKalmanFilter(*(AliMUONTrackParam*)(GetTrackParamAtCluster()->First()));
177 AliDebug(1, Form("Adding Kalman chi2 = %f to global chi2 = %f", chi2Kalman, GetGlobalChi2()));
178 Double_t newGlobalChi2 = GetGlobalChi2() + chi2Kalman;
179 mftCluster.SetLocalChi2(chi2Kalman);
180 mftCluster.SetTrackChi2(newGlobalChi2);
181 new ((*(this->fMFTClusters))[iMFTCluster]) AliMFTCluster(mftCluster);
182 AliDebug(1, Form("GetTrackParamAtCluster() = %p has %d entries while this->fMFTClusters=%p has %d entries",
183 GetTrackParamAtCluster(), GetTrackParamAtCluster()->GetEntries(), this->fMFTClusters, this->fMFTClusters->GetEntries()));
184 AliDebug(1, Form("muonCluster->GetZ() = %f, trackParam->GetZ() = %f",muonCluster->GetZ(), trackParam.GetZ()));
185 SetGlobalChi2(newGlobalChi2);
186 AliDebug(1, Form("New global chi2 = %f", GetGlobalChi2()));
187 ((AliMUONTrackParam*) GetTrackParamAtCluster()->First())->SetTrackChi2(newGlobalChi2);
188 for (Int_t iPar=0; iPar<GetTrackParamAtCluster()->GetEntries(); iPar++) {
189 AliDebug(1, Form("GetTrackParamAtCluster()->At(%d)->GetClusterPtr() = %p",
190 iPar, ((AliMUONTrackParam*)(GetTrackParamAtCluster()->At(iPar)))->GetClusterPtr()));
195 //====================================================================================================================================================
197 AliMUONTrackParam* AliMuonForwardTrack::GetTrackParamAtMUONCluster(Int_t iMUONCluster) {
199 if (iMUONCluster<0 || iMUONCluster>=GetNMUONClusters()) {
200 AliError("Invalid MUON cluster index. NULL pointer will be returned");
204 AliMUONTrackParam *trackParam = (AliMUONTrackParam*) fMUONTrack->GetTrackParamAtCluster()->At(iMUONCluster);
210 //====================================================================================================================================================
212 AliMUONTrackParam* AliMuonForwardTrack::GetTrackParamAtMFTCluster(Int_t iMFTCluster) {
214 if (iMFTCluster<0 || iMFTCluster>=GetNMFTClusters()) {
215 AliError("Invalid MFT cluster index. NULL pointer will be returned");
219 AliMUONTrackParam *trackParam = (AliMUONTrackParam*) GetTrackParamAtCluster()->At(iMFTCluster);
225 //====================================================================================================================================================
227 AliMUONVCluster* AliMuonForwardTrack::GetMUONCluster(Int_t iMUONCluster) {
229 if (iMUONCluster<0 || iMUONCluster>=GetNMUONClusters()) {
230 AliError("Invalid MUON cluster index. NULL pointer will be returned");
234 AliMUONTrackParam *trackParam = GetTrackParamAtMUONCluster(iMUONCluster);
235 AliMUONVCluster *muonCluster = trackParam->GetClusterPtr();
241 //====================================================================================================================================================
243 AliMFTCluster* AliMuonForwardTrack::GetMFTCluster(Int_t iMFTCluster) {
245 if (iMFTCluster<0 || iMFTCluster>=GetNMFTClusters()) {
246 AliError(Form("Invalid MFT cluster index (%d). GetNMFTClusters()=%d. NULL pointer will be returned", iMFTCluster, GetNMFTClusters()));
250 AliMUONTrackParam *trackParam = GetTrackParamAtMFTCluster(iMFTCluster);
251 AliMFTCluster *mftCluster = (AliMFTCluster*) this->fMFTClusters->At(trackParam->GetUniqueID());
257 //====================================================================================================================================================
259 Double_t AliMuonForwardTrack::RunKalmanFilter(AliMUONTrackParam &trackParamAtCluster) {
261 AliDebug(1, Form("Running Kalman filter for parameters %p (z = %f) and cluster %p (z = %f)",
262 &trackParamAtCluster, trackParamAtCluster.GetZ(), trackParamAtCluster.GetClusterPtr(), trackParamAtCluster.GetClusterPtr()->GetZ()));
264 // Compute new track parameters and their covariances including new cluster using kalman filter
265 // return the *additional* track chi2
267 // Get actual track parameters (p)
268 TMatrixD param(trackParamAtCluster.GetParameters());
270 // Get new cluster parameters (m)
271 AliMUONVCluster *cluster = trackParamAtCluster.GetClusterPtr();
272 AliDebug(1, Form("cluster->GetX() = %f cluster->GetY() = %f cluster->GetErrX() = %f cluster->GetErrY() = %f",
273 cluster->GetX(), cluster->GetY(), cluster->GetErrX(), cluster->GetErrY()));
274 TMatrixD clusterParam(5,1);
276 clusterParam(0,0) = cluster->GetX();
277 clusterParam(2,0) = cluster->GetY();
279 // Compute the current parameter weight (W)
280 TMatrixD paramWeight(trackParamAtCluster.GetCovariances());
281 if (paramWeight.Determinant() != 0) {
282 paramWeight.Invert();
284 Warning("RunKalmanFilter"," Determinant = 0");
288 // Compute the new cluster weight (U)
289 TMatrixD clusterWeight(5,5);
290 clusterWeight.Zero();
291 clusterWeight(0,0) = 1. / cluster->GetErrX2();
292 clusterWeight(2,2) = 1. / cluster->GetErrY2();
294 // Compute the new parameters covariance matrix ( (W+U)^-1 )
295 TMatrixD newParamCov(paramWeight,TMatrixD::kPlus,clusterWeight);
296 if (newParamCov.Determinant() != 0) {
297 newParamCov.Invert();
300 Warning("RunKalmanFilter"," Determinant = 0");
304 // Save the new parameters covariance matrix
305 trackParamAtCluster.SetCovariances(newParamCov);
307 // Compute the new parameters (p' = ((W+U)^-1)U(m-p) + p)
308 TMatrixD tmp(clusterParam,TMatrixD::kMinus,param);
309 TMatrixD tmp2(clusterWeight,TMatrixD::kMult,tmp); // U(m-p)
310 TMatrixD newParam(newParamCov,TMatrixD::kMult,tmp2); // ((W+U)^-1)U(m-p)
311 newParam += param; // ((W+U)^-1)U(m-p) + p
313 // Save the new parameters
314 trackParamAtCluster.SetParameters(newParam);
316 // Compute the additional chi2 (= ((p'-p)^-1)W(p'-p) + ((p'-m)^-1)U(p'-m))
317 tmp = newParam; // p'
318 tmp -= param; // (p'-p)
319 TMatrixD tmp3(paramWeight,TMatrixD::kMult,tmp); // W(p'-p)
320 TMatrixD addChi2Track(tmp,TMatrixD::kTransposeMult,tmp3); // ((p'-p)^-1)W(p'-p)
321 tmp = newParam; // p'
322 tmp -= clusterParam; // (p'-m)
323 TMatrixD tmp4(clusterWeight,TMatrixD::kMult,tmp); // U(p'-m)
324 addChi2Track += TMatrixD(tmp,TMatrixD::kTransposeMult,tmp4); // ((p'-p)^-1)W(p'-p) + ((p'-m)^-1)U(p'-m)
326 AliDebug(1,Form("Adding Kalman chi2 = %f",addChi2Track(0,0)));
328 return addChi2Track(0,0);
332 //====================================================================================================================================================
334 Double_t AliMuonForwardTrack::GetWeightedOffset(Double_t x, Double_t y, Double_t z) {
336 AliMUONTrackParam *param = GetTrackParamAtMFTCluster(0);
337 AliMUONTrackExtrap::ExtrapToZCov(param, z);
340 cov = param->GetCovariances();
342 TMatrixD covCoordinates(2,2);
343 covCoordinates(0,0) = cov(0,0);
344 covCoordinates(0,1) = cov(0,2);
345 covCoordinates(1,0) = cov(2,0);
346 covCoordinates(1,1) = cov(2,2);
348 TMatrixD covCoordinatesInverse = covCoordinates.Invert();
350 Double_t dX = param->GetNonBendingCoor() - x;
351 Double_t dY = param->GetBendingCoor() - y;
353 Double_t weightedOffset = TMath::Sqrt(0.5*(dX*dX*covCoordinatesInverse(0,0) +
354 dY*dY*covCoordinatesInverse(1,1) +
355 2.*dX*dY*covCoordinatesInverse(0,1)));
357 return weightedOffset;
361 //====================================================================================================================================================
363 Double_t AliMuonForwardTrack::GetOffsetX(Double_t x, Double_t z) {
365 AliMUONTrackParam *param = GetTrackParamAtMFTCluster(0);
366 AliMUONTrackExtrap::ExtrapToZCov(param, z);
367 Double_t dX = param->GetNonBendingCoor() - x;
372 //====================================================================================================================================================
374 Double_t AliMuonForwardTrack::GetOffsetY(Double_t y, Double_t z) {
376 AliMUONTrackParam *param = GetTrackParamAtMFTCluster(0);
377 AliMUONTrackExtrap::ExtrapToZCov(param, z);
378 Double_t dY = param->GetBendingCoor() - y;
383 //====================================================================================================================================================
385 Double_t AliMuonForwardTrack::GetOffset(Double_t x, Double_t y, Double_t z) {
387 AliMUONTrackParam *param = GetTrackParamAtMFTCluster(0);
388 AliMUONTrackExtrap::ExtrapToZCov(param, z);
389 Double_t dX = param->GetNonBendingCoor() - x;
390 Double_t dY = param->GetBendingCoor() - y;
391 return TMath::Sqrt(dX*dX + dY*dY);
395 //====================================================================================================================================================