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"
35 #include "TLorentzVector.h"
36 #include "TDatabasePDG.h"
38 ClassImp(AliMuonForwardTrack)
40 //====================================================================================================================================================
42 AliMuonForwardTrack::AliMuonForwardTrack():
47 fNWrongClustersMC(-1),
53 // default constructor
55 for (Int_t iPlane=0; iPlane<AliMFTConstants::fNMaxPlanes; iPlane++) fPlaneExists[iPlane] = kFALSE;
56 for (Int_t iParent=0; iParent<fgkNParentsMax; iParent++) {
57 fParentMCLabel[iParent] = -1;
58 fParentPDGCode[iParent] = 0;
60 fMFTClusters = new TClonesArray("AliMFTCluster");
61 fMFTClusters -> SetOwner(kTRUE);
65 //====================================================================================================================================================
67 AliMuonForwardTrack::AliMuonForwardTrack(AliMUONTrack *MUONTrack):
72 fNWrongClustersMC(-1),
78 SetMUONTrack(MUONTrack);
79 for (Int_t iPlane=0; iPlane<AliMFTConstants::fNMaxPlanes; iPlane++) fPlaneExists[iPlane] = kFALSE;
80 for (Int_t iParent=0; iParent<fgkNParentsMax; iParent++) {
81 fParentMCLabel[iParent] = -1;
82 fParentPDGCode[iParent] = 0;
84 fMFTClusters = new TClonesArray("AliMFTCluster");
85 fMFTClusters -> SetOwner(kTRUE);
89 //====================================================================================================================================================
91 AliMuonForwardTrack::AliMuonForwardTrack(const AliMuonForwardTrack& track):
96 fNWrongClustersMC(track.fNWrongClustersMC),
97 fTrackMCId(track.fTrackMCId),
99 fParamCovMatrix(track.fParamCovMatrix)
103 fMUONTrack = new AliMUONTrack(*(track.fMUONTrack));
104 if (track.fMCTrackRef) fMCTrackRef = new TParticle(*(track.fMCTrackRef));
105 fMFTClusters = new TClonesArray(*(track.fMFTClusters));
106 fMFTClusters->SetOwner(kTRUE);
107 for (Int_t iPlane=0; iPlane<AliMFTConstants::fNMaxPlanes; iPlane++) fPlaneExists[iPlane] = (track.fPlaneExists)[iPlane];
108 for (Int_t iParent=0; iParent<fgkNParentsMax; iParent++) {
109 fParentMCLabel[iParent] = (track.fParentMCLabel)[iParent];
110 fParentPDGCode[iParent] = (track.fParentPDGCode)[iParent];
115 //====================================================================================================================================================
117 AliMuonForwardTrack& AliMuonForwardTrack::operator=(const AliMuonForwardTrack& track) {
119 // Asignment operator
121 // check assignement to self
122 if (this == &track) return *this;
124 // base class assignement
125 AliMUONTrack::operator=(track);
130 fMUONTrack = new AliMUONTrack(*(track.fMUONTrack));
131 if (track.fMCTrackRef) fMCTrackRef = new TParticle(*(track.fMCTrackRef));
132 fMFTClusters = new TClonesArray(*(track.fMFTClusters));
133 fMFTClusters->SetOwner(kTRUE);
134 fNWrongClustersMC = track.fNWrongClustersMC;
135 fTrackMCId = track.fTrackMCId;
136 fKinem = track.fKinem;
137 fParamCovMatrix = track.fParamCovMatrix;
139 for (Int_t iPlane=0; iPlane<AliMFTConstants::fNMaxPlanes; iPlane++) fPlaneExists[iPlane] = (track.fPlaneExists)[iPlane];
140 for (Int_t iParent=0; iParent<fgkNParentsMax; iParent++) {
141 fParentMCLabel[iParent] = (track.fParentMCLabel)[iParent];
142 fParentPDGCode[iParent] = (track.fParentPDGCode)[iParent];
149 //====================================================================================================================================================
151 void AliMuonForwardTrack::Clear(const Option_t* /*opt*/) {
154 fMFTClusters -> Delete();
155 delete fMFTClusters; fMFTClusters = 0x0;
156 delete fMUONTrack; fMUONTrack = 0x0;
157 delete fMCTrackRef; fMCTrackRef = 0x0;
161 //====================================================================================================================================================
163 AliMuonForwardTrack::~AliMuonForwardTrack() {
167 fMFTClusters -> Delete();
172 //====================================================================================================================================================
174 void AliMuonForwardTrack::SetMUONTrack(AliMUONTrack *MUONTrack) {
177 AliInfo("fMUONTrack already exists, nothing will be done");
181 fMUONTrack = MUONTrack;
182 SetGlobalChi2(fMUONTrack->GetGlobalChi2());
186 //====================================================================================================================================================
188 void AliMuonForwardTrack::SetMCTrackRef(TParticle *MCTrackRef) {
191 AliInfo("fMCTrackRef already exists, nothing will be done");
195 fMCTrackRef = MCTrackRef;
199 //====================================================================================================================================================
201 void AliMuonForwardTrack::AddTrackParamAtMFTCluster(AliMUONTrackParam &trackParam, AliMFTCluster &mftCluster) {
203 AliDebug(1, Form("Before adding: this->fMFTClusters=%p has %d entries", this->fMFTClusters, this->fMFTClusters->GetEntries()));
204 Int_t iMFTCluster = this->fMFTClusters->GetEntries();
205 AliDebug(1, Form("mftCluster->GetX() = %f mftCluster->GetY() = %f mftCluster->GetErrX() = %f cmftCluster->GetErrY() = %f",
206 mftCluster.GetX(), mftCluster.GetY(), mftCluster.GetErrX(), mftCluster.GetErrY()));
207 AliMUONVCluster *muonCluster = (AliMUONVCluster*) mftCluster.CreateMUONCluster();
208 AliDebug(1, Form("Created MUON cluster %p", muonCluster));
209 trackParam.SetUniqueID(iMFTCluster); // we profit of this slot to store the reference to the corresponding MFTCluster
210 AliDebug(1, Form("Now adding muonCluster %p and trackParam %p",muonCluster, &trackParam));
211 AddTrackParamAtCluster(trackParam, *muonCluster, kTRUE);
212 // we pass the parameters this->GetTrackParamAtCluster()->First() to the Kalman Filter algorithm: they will be updated!!
213 Double_t chi2Kalman = RunKalmanFilter(*(AliMUONTrackParam*)(GetTrackParamAtCluster()->First()));
214 AliDebug(1, Form("Adding Kalman chi2 = %f to global chi2 = %f", chi2Kalman, GetGlobalChi2()));
215 Double_t newGlobalChi2 = GetGlobalChi2() + chi2Kalman;
216 mftCluster.SetLocalChi2(chi2Kalman);
217 mftCluster.SetTrackChi2(newGlobalChi2);
218 new ((*(this->fMFTClusters))[iMFTCluster]) AliMFTCluster(mftCluster);
219 AliDebug(1, Form("GetTrackParamAtCluster() = %p has %d entries while this->fMFTClusters=%p has %d entries",
220 GetTrackParamAtCluster(), GetTrackParamAtCluster()->GetEntries(), this->fMFTClusters, this->fMFTClusters->GetEntries()));
221 AliDebug(1, Form("muonCluster->GetZ() = %f, trackParam->GetZ() = %f",muonCluster->GetZ(), trackParam.GetZ()));
222 SetGlobalChi2(newGlobalChi2);
223 AliDebug(1, Form("New global chi2 = %f", GetGlobalChi2()));
224 ((AliMUONTrackParam*) GetTrackParamAtCluster()->First())->SetTrackChi2(newGlobalChi2);
225 for (Int_t iPar=0; iPar<GetTrackParamAtCluster()->GetEntries(); iPar++) {
226 AliDebug(1, Form("GetTrackParamAtCluster()->At(%d)->GetClusterPtr() = %p",
227 iPar, ((AliMUONTrackParam*)(GetTrackParamAtCluster()->At(iPar)))->GetClusterPtr()));
232 //====================================================================================================================================================
234 AliMUONTrackParam* AliMuonForwardTrack::GetTrackParamAtMUONCluster(Int_t iMUONCluster) {
236 if (iMUONCluster<0 || iMUONCluster>=GetNMUONClusters()) {
237 AliError("Invalid MUON cluster index. NULL pointer will be returned");
241 AliMUONTrackParam *trackParam = (AliMUONTrackParam*) fMUONTrack->GetTrackParamAtCluster()->At(iMUONCluster);
247 //====================================================================================================================================================
249 AliMUONTrackParam* AliMuonForwardTrack::GetTrackParamAtMFTCluster(Int_t iMFTCluster) {
251 if (iMFTCluster<0 || iMFTCluster>=GetNMFTClusters()) {
252 AliError("Invalid MFT cluster index. NULL pointer will be returned");
256 AliMUONTrackParam *trackParam = (AliMUONTrackParam*) GetTrackParamAtCluster()->At(iMFTCluster);
262 //====================================================================================================================================================
264 AliMUONVCluster* AliMuonForwardTrack::GetMUONCluster(Int_t iMUONCluster) {
266 if (iMUONCluster<0 || iMUONCluster>=GetNMUONClusters()) {
267 AliError("Invalid MUON cluster index. NULL pointer will be returned");
271 AliMUONTrackParam *trackParam = GetTrackParamAtMUONCluster(iMUONCluster);
272 AliMUONVCluster *muonCluster = trackParam->GetClusterPtr();
278 //====================================================================================================================================================
280 AliMFTCluster* AliMuonForwardTrack::GetMFTCluster(Int_t iMFTCluster) {
282 if (iMFTCluster<0 || iMFTCluster>=GetNMFTClusters()) {
283 AliError(Form("Invalid MFT cluster index (%d). GetNMFTClusters()=%d. NULL pointer will be returned", iMFTCluster, GetNMFTClusters()));
287 AliMUONTrackParam *trackParam = GetTrackParamAtMFTCluster(iMFTCluster);
288 AliMFTCluster *mftCluster = (AliMFTCluster*) this->fMFTClusters->At(trackParam->GetUniqueID());
294 //====================================================================================================================================================
296 Double_t AliMuonForwardTrack::RunKalmanFilter(AliMUONTrackParam &trackParamAtCluster) {
298 AliDebug(1, Form("Running Kalman filter for parameters %p (z = %f) and cluster %p (z = %f)",
299 &trackParamAtCluster, trackParamAtCluster.GetZ(), trackParamAtCluster.GetClusterPtr(), trackParamAtCluster.GetClusterPtr()->GetZ()));
301 // Compute new track parameters and their covariances including new cluster using kalman filter
302 // return the *additional* track chi2
304 // Get actual track parameters (p)
305 TMatrixD param(trackParamAtCluster.GetParameters());
307 // Get new cluster parameters (m)
308 AliMUONVCluster *cluster = trackParamAtCluster.GetClusterPtr();
309 AliDebug(1, Form("cluster->GetX() = %f cluster->GetY() = %f cluster->GetErrX() = %f cluster->GetErrY() = %f",
310 cluster->GetX(), cluster->GetY(), cluster->GetErrX(), cluster->GetErrY()));
311 TMatrixD clusterParam(5,1);
313 clusterParam(0,0) = cluster->GetX();
314 clusterParam(2,0) = cluster->GetY();
316 // Compute the current parameter weight (W)
317 TMatrixD paramWeight(trackParamAtCluster.GetCovariances());
318 if (paramWeight.Determinant() != 0) {
319 paramWeight.Invert();
321 Warning("RunKalmanFilter"," Determinant = 0");
325 // Compute the new cluster weight (U)
326 TMatrixD clusterWeight(5,5);
327 clusterWeight.Zero();
328 clusterWeight(0,0) = 1. / cluster->GetErrX2();
329 clusterWeight(2,2) = 1. / cluster->GetErrY2();
331 // Compute the new parameters covariance matrix ( (W+U)^-1 )
332 TMatrixD newParamCov(paramWeight,TMatrixD::kPlus,clusterWeight);
333 if (newParamCov.Determinant() != 0) {
334 newParamCov.Invert();
337 Warning("RunKalmanFilter"," Determinant = 0");
341 // Save the new parameters covariance matrix
342 trackParamAtCluster.SetCovariances(newParamCov);
344 // Compute the new parameters (p' = ((W+U)^-1)U(m-p) + p)
345 TMatrixD tmp(clusterParam,TMatrixD::kMinus,param);
346 TMatrixD tmp2(clusterWeight,TMatrixD::kMult,tmp); // U(m-p)
347 TMatrixD newParam(newParamCov,TMatrixD::kMult,tmp2); // ((W+U)^-1)U(m-p)
348 newParam += param; // ((W+U)^-1)U(m-p) + p
350 // Save the new parameters
351 trackParamAtCluster.SetParameters(newParam);
353 // Compute the additional chi2 (= ((p'-p)^-1)W(p'-p) + ((p'-m)^-1)U(p'-m))
354 tmp = newParam; // p'
355 tmp -= param; // (p'-p)
356 TMatrixD tmp3(paramWeight,TMatrixD::kMult,tmp); // W(p'-p)
357 TMatrixD addChi2Track(tmp,TMatrixD::kTransposeMult,tmp3); // ((p'-p)^-1)W(p'-p)
358 tmp = newParam; // p'
359 tmp -= clusterParam; // (p'-m)
360 TMatrixD tmp4(clusterWeight,TMatrixD::kMult,tmp); // U(p'-m)
361 addChi2Track += TMatrixD(tmp,TMatrixD::kTransposeMult,tmp4); // ((p'-p)^-1)W(p'-p) + ((p'-m)^-1)U(p'-m)
363 AliDebug(1,Form("Adding Kalman chi2 = %f",addChi2Track(0,0)));
365 return addChi2Track(0,0);
369 //====================================================================================================================================================
371 Double_t AliMuonForwardTrack::GetWeightedOffset(Double_t x, Double_t y, Double_t z) {
373 AliMUONTrackParam *param = GetTrackParamAtMFTCluster(0);
374 AliMUONTrackExtrap::ExtrapToZCov(param, z);
377 cov = param->GetCovariances();
379 TMatrixD covCoordinates(2,2);
380 covCoordinates(0,0) = cov(0,0);
381 covCoordinates(0,1) = cov(0,2);
382 covCoordinates(1,0) = cov(2,0);
383 covCoordinates(1,1) = cov(2,2);
385 TMatrixD covCoordinatesInverse = covCoordinates.Invert();
387 Double_t dX = param->GetNonBendingCoor() - x;
388 Double_t dY = param->GetBendingCoor() - y;
390 Double_t weightedOffset = TMath::Sqrt(0.5*(dX*dX*covCoordinatesInverse(0,0) +
391 dY*dY*covCoordinatesInverse(1,1) +
392 2.*dX*dY*covCoordinatesInverse(0,1)));
394 return weightedOffset;
398 //====================================================================================================================================================
400 Double_t AliMuonForwardTrack::GetOffsetX(Double_t x, Double_t z) {
402 AliMUONTrackParam *param = GetTrackParamAtMFTCluster(0);
403 AliMUONTrackExtrap::ExtrapToZCov(param, z);
404 Double_t dX = param->GetNonBendingCoor() - x;
409 //====================================================================================================================================================
411 Double_t AliMuonForwardTrack::GetOffsetY(Double_t y, Double_t z) {
413 AliMUONTrackParam *param = GetTrackParamAtMFTCluster(0);
414 AliMUONTrackExtrap::ExtrapToZCov(param, z);
415 Double_t dY = param->GetBendingCoor() - y;
420 //====================================================================================================================================================
422 Double_t AliMuonForwardTrack::GetOffset(Double_t x, Double_t y, Double_t z) {
424 AliMUONTrackParam *param = GetTrackParamAtMFTCluster(0);
425 AliMUONTrackExtrap::ExtrapToZCov(param, z);
426 Double_t dX = param->GetNonBendingCoor() - x;
427 Double_t dY = param->GetBendingCoor() - y;
428 return TMath::Sqrt(dX*dX + dY*dY);
432 //====================================================================================================================================================
434 Bool_t AliMuonForwardTrack::IsFromResonance() {
436 Bool_t result = kFALSE;
438 if ( GetParentPDGCode(0) == 113 ||
439 GetParentPDGCode(0) == 221 ||
440 GetParentPDGCode(0) == 223 ||
441 GetParentPDGCode(0) == 331 ||
442 GetParentPDGCode(0) == 333 ||
443 GetParentPDGCode(0) == 443 ||
444 GetParentPDGCode(0) == 100443 ||
445 GetParentPDGCode(0) == 553 ||
446 GetParentPDGCode(0) == 100553 ) result = kTRUE;
448 if (result) AliDebug(1, Form("Muon comes from a resonance %d", GetParentPDGCode(0)));
454 //====================================================================================================================================================
456 Bool_t AliMuonForwardTrack::IsFromCharm() {
458 Bool_t result = kFALSE;
460 if (IsPDGCharm(GetParentPDGCode(0)) && !IsPDGBeauty(GetParentPDGCode(1))) result = kTRUE;
462 if (result) AliDebug(1, Form("Muon comes from a charmed hadron %d", GetParentPDGCode(0)));
468 //====================================================================================================================================================
470 Bool_t AliMuonForwardTrack::IsFromBeauty() {
472 Bool_t result = kFALSE;
474 if (IsPDGBeauty(GetParentPDGCode(0))) result = kTRUE;
476 if (result) AliDebug(1, Form("Muon comes from a beauty hadron %d", GetParentPDGCode(0)));
482 //====================================================================================================================================================
484 Bool_t AliMuonForwardTrack::IsPDGCharm(Int_t pdg) {
486 Bool_t result = kFALSE;
488 if ( TMath::Abs(pdg) == 411 ||
489 TMath::Abs(pdg) == 421 ||
490 TMath::Abs(pdg) == 10411 ||
491 TMath::Abs(pdg) == 10421 ||
492 TMath::Abs(pdg) == 413 ||
493 TMath::Abs(pdg) == 423 ||
494 TMath::Abs(pdg) == 10413 ||
495 TMath::Abs(pdg) == 10423 ||
496 TMath::Abs(pdg) == 20413 ||
497 TMath::Abs(pdg) == 20423 ||
498 TMath::Abs(pdg) == 415 ||
499 TMath::Abs(pdg) == 425 ||
500 TMath::Abs(pdg) == 431 ||
501 TMath::Abs(pdg) == 10431 ||
502 TMath::Abs(pdg) == 433 ||
503 TMath::Abs(pdg) == 10433 ||
504 TMath::Abs(pdg) == 20433 ||
505 TMath::Abs(pdg) == 435 ) result = kTRUE;
511 //====================================================================================================================================================
513 Bool_t AliMuonForwardTrack::IsPDGBeauty(Int_t pdg) {
515 Bool_t result = kFALSE;
517 if ( TMath::Abs(pdg) == 511 ||
518 TMath::Abs(pdg) == 521 ||
519 TMath::Abs(pdg) == 10511 ||
520 TMath::Abs(pdg) == 10521 ||
521 TMath::Abs(pdg) == 513 ||
522 TMath::Abs(pdg) == 523 ||
523 TMath::Abs(pdg) == 10513 ||
524 TMath::Abs(pdg) == 10523 ||
525 TMath::Abs(pdg) == 20513 ||
526 TMath::Abs(pdg) == 20523 ||
527 TMath::Abs(pdg) == 515 ||
528 TMath::Abs(pdg) == 525 ||
529 TMath::Abs(pdg) == 531 ||
530 TMath::Abs(pdg) == 10531 ||
531 TMath::Abs(pdg) == 533 ||
532 TMath::Abs(pdg) == 10533 ||
533 TMath::Abs(pdg) == 20533 ||
534 TMath::Abs(pdg) == 535 ||
535 TMath::Abs(pdg) == 541 ||
536 TMath::Abs(pdg) == 10541 ||
537 TMath::Abs(pdg) == 543 ||
538 TMath::Abs(pdg) == 10543 ||
539 TMath::Abs(pdg) == 20543 ||
540 TMath::Abs(pdg) == 545 ) result = kTRUE;
546 //====================================================================================================================================================
548 Bool_t AliMuonForwardTrack::IsFromBackground() {
550 Bool_t result = kFALSE;
552 if (!IsFromResonance() && !IsFromCharm() && !IsFromBeauty()) result = kTRUE;
554 if (result) AliDebug(1, Form("Muon comes from a background source %d", GetParentPDGCode(0)));
560 //====================================================================================================================================================
562 void AliMuonForwardTrack::EvalKinem(Double_t z) {
564 AliMUONTrackParam *param = GetTrackParamAtMFTCluster(0);
565 AliMUONTrackExtrap::ExtrapToZCov(param, z);
567 Double_t mMu = TDatabasePDG::Instance()->GetParticle("mu-")->Mass();
568 Double_t energy = TMath::Sqrt(param->P()*param->P() + mMu*mMu);
569 fKinem.SetPxPyPzE(param->Px(), param->Py(), param->Pz(), energy);
572 fParamCovMatrix = param->GetCovariances();
576 //====================================================================================================================================================