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 // Support class for various common operation on MFT objects
20 // Contact author: antonio.uras@cern.ch
22 //====================================================================================================================================================
24 #include "AliMUONTrackParam.h"
25 #include "AliMUONTrackExtrap.h"
26 #include "AliAODTrack.h"
27 #include "AliAODDimuon.h"
28 #include "TLorentzVector.h"
29 #include "AliMFTConstants.h"
30 #include "TDatabasePDG.h"
33 #include "TObjArray.h"
34 #include "TDecompLU.h"
36 #include "AliMFTAnalysisTools.h"
38 ClassImp(AliMFTAnalysisTools)
40 //====================================================================================================================================================
42 Bool_t AliMFTAnalysisTools::ExtrapAODMuonToZ(AliAODTrack *muon, Double_t z, Double_t xy[2]) {
44 if (!(muon->Pz()!=0)) return kFALSE;
46 AliMUONTrackParam *param = new AliMUONTrackParam();
48 param -> SetNonBendingCoor(muon->XAtDCA());
49 param -> SetBendingCoor(muon->YAtDCA());
50 param -> SetZ(AliMFTConstants::fZEvalKinem);
51 param -> SetNonBendingSlope(muon->Px()/muon->Pz());
52 param -> SetBendingSlope(muon->Py()/muon->Pz());
53 param -> SetInverseBendingMomentum( muon->Charge() * (1./muon->Pz()) / (TMath::Sqrt(1+TMath::Power(muon->Py()/muon->Pz(),2))) );
55 AliMUONTrackExtrap::ExtrapToZ(param, z);
56 xy[0] = param->GetNonBendingCoor();
57 xy[1] = param->GetBendingCoor();
65 //====================================================================================================================================================
67 Bool_t AliMFTAnalysisTools::ExtrapAODMuonToZ(AliAODTrack *muon, Double_t z, Double_t xy[2], TLorentzVector &kinem) {
69 if (!(muon->Pz()!=0)) return kFALSE;
71 AliMUONTrackParam *param = new AliMUONTrackParam();
73 param -> SetNonBendingCoor(muon->XAtDCA());
74 param -> SetBendingCoor(muon->YAtDCA());
75 param -> SetZ(AliMFTConstants::fZEvalKinem);
76 param -> SetNonBendingSlope(muon->Px()/muon->Pz());
77 param -> SetBendingSlope(muon->Py()/muon->Pz());
78 param -> SetInverseBendingMomentum( muon->Charge() * (1./muon->Pz()) / (TMath::Sqrt(1+TMath::Power(muon->Py()/muon->Pz(),2))) );
80 AliMUONTrackExtrap::ExtrapToZ(param, z);
81 xy[0] = param->GetNonBendingCoor();
82 xy[1] = param->GetBendingCoor();
84 Double_t massMu = TDatabasePDG::Instance()->GetParticle("mu-")->Mass();
85 Double_t energy = TMath::Sqrt(massMu*massMu + param->Px()*param->Px() + param->Py()*param->Py() + param->Pz()*param->Pz());
87 kinem.SetPxPyPzE(param->Px(), param->Py(), param->Pz(), energy);
95 //====================================================================================================================================================
97 Bool_t AliMFTAnalysisTools::ExtrapAODMuonToZ(AliAODTrack *muon, Double_t z, Double_t xy[2], TLorentzVector &kinem, TMatrixD &cov) {
99 // Extrapolate muon to a given z providing the corresponding (x,y) position and updating kinematics and covariance matrix
101 if (!(muon->Pz()!=0)) return kFALSE;
103 AliMUONTrackParam *param = new AliMUONTrackParam();
105 param -> SetNonBendingCoor(muon->XAtDCA());
106 param -> SetBendingCoor(muon->YAtDCA());
107 param -> SetZ(AliMFTConstants::fZEvalKinem);
108 param -> SetNonBendingSlope(muon->Px()/muon->Pz());
109 param -> SetBendingSlope(muon->Py()/muon->Pz());
110 param -> SetInverseBendingMomentum( muon->Charge() * (1./muon->Pz()) / (TMath::Sqrt(1+TMath::Power(muon->Py()/muon->Pz(),2))) );
112 param -> SetCovariances(ConvertCovMatrixAOD2MUON(muon));
114 AliMUONTrackExtrap::ExtrapToZCov(param, z);
115 xy[0] = param->GetNonBendingCoor();
116 xy[1] = param->GetBendingCoor();
118 Double_t massMu = TDatabasePDG::Instance()->GetParticle("mu-")->Mass();
119 Double_t energy = TMath::Sqrt(massMu*massMu + param->Px()*param->Px() + param->Py()*param->Py() + param->Pz()*param->Pz());
121 kinem.SetPxPyPzE(param->Px(), param->Py(), param->Pz(), energy);
123 cov = param->GetCovariances();
131 //====================================================================================================================================================
133 Bool_t AliMFTAnalysisTools::ExtrapAODMuonToXY(AliAODTrack *muon, Double_t xy[2], Double_t &zFinal, TLorentzVector &kinem, TMatrixD &cov) {
135 // Find the point of closest approach between the muon and the direction parallel to the z-axis defined by the given (x,y)
136 // Provide the z of the above point as weel as the updated kinematics and covariance matrix
138 // We look for the above-defined PCA
140 Double_t xPCA=0, yPCA=0, zPCA=0;
142 AliMUONTrackParam *param = new AliMUONTrackParam();
143 param -> SetNonBendingCoor(muon->XAtDCA());
144 param -> SetBendingCoor(muon->YAtDCA());
145 param -> SetZ(AliMFTConstants::fZEvalKinem);
146 param -> SetNonBendingSlope(muon->Px()/muon->Pz());
147 param -> SetBendingSlope(muon->Py()/muon->Pz());
148 param -> SetInverseBendingMomentum( muon->Charge() * (1./muon->Pz()) / (TMath::Sqrt(1+TMath::Power(muon->Py()/muon->Pz(),2))) );
150 // here we want to understand in which direction we have to search the minimum...
152 Double_t step = 1.; // initial step, in cm
153 Double_t startPoint = 0.;
155 Double_t r[3]={0}, z[3]={startPoint, startPoint+step, startPoint+2*step};
157 TVector3 **points = new TVector3*[2]; // points[0] for the muon, points[1] for the direction parallel to the z-axis defined by the given (x,y)
159 for (Int_t i=0; i<3; i++) {
160 AliMUONTrackExtrap::ExtrapToZ(param, z[i]);
161 points[0] = new TVector3(param->GetNonBendingCoor(),param->GetBendingCoor(),z[i]);
162 points[1] = new TVector3(xy[0],xy[1],z[i]);
163 r[i] = GetDistanceBetweenPoints(points,2);
164 for (Int_t iMu=0; iMu<2; iMu++) delete points[iMu];
167 Int_t researchDirection = 0;
169 if (r[0]>r[1] && r[1]>r[2]) researchDirection = +1; // towards z positive
170 else if (r[0]<r[1] && r[1]<r[2]) researchDirection = -1; // towards z negative
171 else if (r[0]<r[1] && r[1]>r[2]) {
172 printf("E-AliMFTAnalysisTools::ExtrapAODMuonToXY: Point of closest approach cannot be found (no minima)\n");
178 while (TMath::Abs(researchDirection)>0.5) {
180 if (researchDirection>0) {
183 z[2] = z[1]+researchDirection*step;
188 z[0] = z[1]+researchDirection*step;
190 if (TMath::Abs(z[0])>900.) {
191 printf("E-AliMFTAnalysisTools::ExtrapAODMuonToXY: Point of closest approach cannot be found (no minima in the fiducial region)\n");
197 for (Int_t i=0; i<3; i++) {
198 AliMUONTrackExtrap::ExtrapToZ(param, z[i]);
199 points[0] = new TVector3(param->GetNonBendingCoor(),param->GetBendingCoor(),z[i]);
200 points[1] = new TVector3(xy[0],xy[1],z[i]);
201 r[i] = GetDistanceBetweenPoints(points,2);
202 for (Int_t iMu=0; iMu<2; iMu++) delete points[iMu];
206 if (r[0]>r[1] && r[1]>r[2]) researchDirection = +1; // towards z positive
207 else if (r[0]<r[1] && r[1]<r[2]) researchDirection = -1; // towards z negative
211 // now we now that the minimum is between z[0] and z[2] and we search for it
214 while (step>AliMFTConstants::fPrecisionPointOfClosestApproach) {
217 for (Int_t i=0; i<3; i++) {
218 AliMUONTrackExtrap::ExtrapToZ(param, z[i]);
219 points[0] = new TVector3(param->GetNonBendingCoor(),param->GetBendingCoor(),z[i]);
220 points[1] = new TVector3(xy[0],xy[1],z[i]);
221 r[i] = GetDistanceBetweenPoints(points,2);
222 for (Int_t iMu=0; iMu<2; iMu++) delete points[iMu];
224 if (r[0]<r[1]) z[1] = z[0];
225 else if (r[2]<r[1]) z[1] = z[2];
231 Double_t xyMuon[2] = {0};
232 ExtrapAODMuonToZ(muon, zFinal, xyMuon, kinem, cov);
238 //====================================================================================================================================================
240 Bool_t AliMFTAnalysisTools::GetAODMuonOffset(AliAODTrack *muon, Double_t xv, Double_t yv, Double_t zv, Double_t &offset) {
242 Double_t xy[2] = {0};
243 ExtrapAODMuonToZ(muon, zv, xy);
245 offset = TMath::Sqrt((xv-xy[0])*(xv-xy[0]) + (yv-xy[1])*(yv-xy[1]));
251 //====================================================================================================================================================
253 Bool_t AliMFTAnalysisTools::GetAODMuonOffsetZ(AliAODTrack *muon, Double_t xv, Double_t yv, Double_t zv, Double_t &offset) {
255 Double_t xy[2] = {xv, yv};
257 TLorentzVector kinem(0,0,0,0);
260 ExtrapAODMuonToXY(muon, xy, zFinal, kinem, cov);
262 offset = TMath::Abs(zFinal - zv);
268 //====================================================================================================================================================
270 Bool_t AliMFTAnalysisTools::GetAODMuonWeightedOffset(AliAODTrack *muon, Double_t xv, Double_t yv, Double_t zv, Double_t &offset) {
272 Double_t xy[2] = {0};
273 TLorentzVector kinem(0,0,0,0);
276 ExtrapAODMuonToZ(muon, zv, xy, kinem, cov);
278 TMatrixD covCoordinates(2,2);
279 covCoordinates(0,0) = cov(0,0);
280 covCoordinates(0,1) = cov(0,2);
281 covCoordinates(1,0) = cov(2,0);
282 covCoordinates(1,1) = cov(2,2);
284 if (TDecompLU::InvertLU(covCoordinates,covCoordinates.GetTol(),0)) {
286 TMatrixD covCoordinatesInverse = covCoordinates;
287 Double_t dX = xy[0] - xv;
288 Double_t dY = xy[1] - yv;
290 offset = TMath::Sqrt(0.5*(dX*dX*covCoordinatesInverse(0,0) +
291 dY*dY*covCoordinatesInverse(1,1) +
292 2.*dX*dY*covCoordinatesInverse(0,1)));
302 //====================================================================================================================================================
304 Double_t AliMFTAnalysisTools::GetPseudoProperDecayTimeXY(Double_t xVtx, Double_t yVtx,
305 Double_t xDimu, Double_t yDimu,
306 Double_t mDimu, Double_t ptDimu) {
308 // pseudo-proper decay time of a particle produced in the primary vertex and decaying into a dimuon (+ X)
309 // evaluated using the transverse degree of freedom of the decay topology
312 Double_t decayLengthXY = TMath::Sqrt((xVtx-xDimu)*(xVtx-xDimu)+(yVtx-yDimu)*(yVtx-yDimu));
313 return (decayLengthXY * mDimu/ptDimu)/TMath::Ccgs()*1E12; // in ps
320 //====================================================================================================================================================
322 Double_t AliMFTAnalysisTools::GetPseudoProperDecayTimeZ(Double_t zVtx,
324 Double_t mDimu, Double_t pzDimu) {
326 // pseudo-proper decay time of a particle produced in the primary vertex and decaying into a dimuon (+ X)
327 // evaluated using the longitudinal degree of freedom of the decay topology
330 Double_t decayLengthZ = zDimu - zVtx;
331 return (decayLengthZ * mDimu/pzDimu)/TMath::Ccgs()*1E12; // in ps
338 //====================================================================================================================================================
340 Bool_t AliMFTAnalysisTools::CalculatePCA(AliAODDimuon *dimuon, Double_t *pca, Double_t &pcaQuality, TLorentzVector &kinem) {
342 TObjArray *muons = new TObjArray();
343 muons -> Add(dimuon->GetMu(0));
344 muons -> Add(dimuon->GetMu(1));
346 Bool_t result = CalculatePCA(muons, pca, pcaQuality, kinem);
352 //====================================================================================================================================================
354 Bool_t AliMFTAnalysisTools::CalculatePCA(TObjArray *muons, Double_t *pca, Double_t &pcaQuality, TLorentzVector &kinem) {
356 const Int_t nMuons = muons->GetEntriesFast();
357 if (nMuons<2 || nMuons>AliMFTConstants::fNMaxMuonsForPCA) {
358 printf("W-AliMFTAnalysisTools::CalculatePCA: number of muons not valid\n");
362 Double_t fXPointOfClosestApproach=0, fYPointOfClosestApproach=0, fZPointOfClosestApproach=0;
364 AliAODTrack *muon[AliMFTConstants::fNMaxMuonsForPCA] = {0};
365 AliMUONTrackParam *param[AliMFTConstants::fNMaxMuonsForPCA] = {0};
367 // Finding AliMUONTrackParam objects for each muon
369 for (Int_t iMu=0; iMu<nMuons; iMu++) {
370 muon[iMu] = (AliAODTrack*) muons->At(iMu);
371 if (TMath::Abs(muon[iMu]->Pz())<1.e-6) {
372 for(Int_t i=0;i<iMu;i++) delete param[i];
375 param[iMu] = new AliMUONTrackParam();
376 param[iMu] -> SetNonBendingCoor(muon[iMu]->XAtDCA());
377 param[iMu] -> SetBendingCoor(muon[iMu]->YAtDCA());
378 param[iMu] -> SetZ(AliMFTConstants::fZEvalKinem);
379 param[iMu] -> SetNonBendingSlope(muon[iMu]->Px()/muon[iMu]->Pz());
380 param[iMu] -> SetBendingSlope(muon[iMu]->Py()/muon[iMu]->Pz());
381 param[iMu] -> SetInverseBendingMomentum( muon[iMu]->Charge() * (1./muon[iMu]->Pz()) / (TMath::Sqrt(1+TMath::Power(muon[iMu]->Py()/muon[iMu]->Pz(),2))) );
384 // here we want to understand in which direction we have to search the minimum...
386 Double_t step = 1.; // initial step, in cm
387 Double_t startPoint = 0.;
389 Double_t r[3]={0}, z[3]={startPoint, startPoint+step, startPoint+2*step};
391 TVector3 **points = new TVector3*[AliMFTConstants::fNMaxMuonsForPCA];
393 for (Int_t i=0; i<3; i++) {
394 for (Int_t iMu=0; iMu<nMuons; iMu++) {
395 AliMUONTrackExtrap::ExtrapToZ(param[iMu], z[i]);
396 points[iMu] = new TVector3(param[iMu]->GetNonBendingCoor(),param[iMu]->GetBendingCoor(),z[i]);
398 r[i] = GetDistanceBetweenPoints(points,nMuons);
399 for (Int_t iMu=0; iMu<nMuons; iMu++) delete points[iMu];
402 Int_t researchDirection = 0;
404 if (r[0]>r[1] && r[1]>r[2]) researchDirection = +1; // towards z positive
405 else if (r[0]<r[1] && r[1]<r[2]) researchDirection = -1; // towards z negative
406 else if (r[0]<r[1] && r[1]>r[2]) {
407 printf("E-AliMFTAnalysisTools::CalculatePCA: Point of closest approach cannot be found for dimuon (no minima)\n");
408 for (Int_t iMu=0;iMu<nMuons;iMu++) delete param[iMu];
413 while (TMath::Abs(researchDirection)>0.5) {
415 if (researchDirection>0) {
418 z[2] = z[1]+researchDirection*step;
423 z[0] = z[1]+researchDirection*step;
425 if (TMath::Abs(z[0])>900.) {
426 printf("E-AliMFTAnalysisTools::CalculatePCA: Point of closest approach cannot be found for dimuon (no minima in the fiducial region)\n");
427 for (Int_t iMu=0;iMu<nMuons;iMu++) delete param[iMu];
432 for (Int_t i=0; i<3; i++) {
433 for (Int_t iMu=0; iMu<nMuons; iMu++) {
434 AliMUONTrackExtrap::ExtrapToZ(param[iMu], z[i]);
435 points[iMu] = new TVector3(param[iMu]->GetNonBendingCoor(),param[iMu]->GetBendingCoor(),z[i]);
437 r[i] = GetDistanceBetweenPoints(points,nMuons);
438 for (Int_t iMu=0;iMu<nMuons;iMu++) delete points[iMu];
441 if (r[0]>r[1] && r[1]>r[2]) researchDirection = +1; // towards z positive
442 else if (r[0]<r[1] && r[1]<r[2]) researchDirection = -1; // towards z negative
446 // now we now that the minimum is between z[0] and z[2] and we search for it
449 while (step>AliMFTConstants::fPrecisionPointOfClosestApproach) {
452 for (Int_t i=0; i<3; i++) {
453 for (Int_t iMu=0; iMu<nMuons; iMu++) {
454 AliMUONTrackExtrap::ExtrapToZ(param[iMu], z[i]);
455 points[iMu] = new TVector3(param[iMu]->GetNonBendingCoor(),param[iMu]->GetBendingCoor(),z[i]);
457 r[i] = GetDistanceBetweenPoints(points,nMuons);
458 for (Int_t iMu=0;iMu<nMuons;iMu++) delete points[iMu];
460 if (r[0]<r[1]) z[1] = z[0];
461 else if (r[2]<r[1]) z[1] = z[2];
465 // Once z of minimum is found, we evaluate the x and y coordinates by averaging over the contributing tracks
467 fZPointOfClosestApproach = z[1];
468 fXPointOfClosestApproach = 0.;
469 fYPointOfClosestApproach = 0.;
470 for (Int_t iMu=0; iMu<nMuons; iMu++) {
471 AliMUONTrackExtrap::ExtrapToZ(param[iMu], fZPointOfClosestApproach);
472 fXPointOfClosestApproach += param[iMu]->GetNonBendingCoor();
473 fYPointOfClosestApproach += param[iMu]->GetBendingCoor();
475 fXPointOfClosestApproach /= Double_t(nMuons);
476 fYPointOfClosestApproach /= Double_t(nMuons);
478 pca[0] = fXPointOfClosestApproach;
479 pca[1] = fYPointOfClosestApproach;
480 pca[2] = fZPointOfClosestApproach;
482 // Evaluating the kinematics of the N-muon
484 Double_t pTot[3] = {0};
486 Double_t massMu = TDatabasePDG::Instance()->GetParticle("mu-")->Mass();
487 for (Int_t iMu=0; iMu<nMuons; iMu++) {
488 pTot[0] += param[iMu]->Px();
489 pTot[1] += param[iMu]->Py();
490 pTot[2] += param[iMu]->Pz();
491 ene += TMath::Sqrt(massMu*massMu + param[iMu]->Px()*param[iMu]->Px() + param[iMu]->Py()*param[iMu]->Py() + param[iMu]->Pz()*param[iMu]->Pz());
494 kinem.SetPxPyPzE(pTot[0], pTot[1], pTot[2], ene);
496 // Evaluating the PCA quality of the N-muon
498 Double_t sum=0.,squareSum=0.;
499 for (Int_t iMu=0; iMu<nMuons; iMu++) {
500 Double_t wOffset = 0;
501 if (!GetAODMuonWeightedOffset(muon[iMu],fXPointOfClosestApproach, fYPointOfClosestApproach, fZPointOfClosestApproach, wOffset)) {
502 for(Int_t jMu=0;jMu<nMuons;jMu++) delete param[jMu];
506 Double_t f = TMath::Exp(-0.5 * wOffset);
510 if (sum > 0.) pcaQuality = (sum-squareSum/sum) / (nMuons-1);
511 else pcaQuality = 0.;
513 for(Int_t iMu=0;iMu<nMuons;iMu++) delete param[iMu];
519 //=========================================================================================================================
521 Double_t AliMFTAnalysisTools::GetDistanceBetweenPoints(TVector3 **points, Int_t nPoints) {
523 if (nPoints>AliMFTConstants::fNMaxMuonsForPCA) {
524 printf("W-AliMFTAnalysisTools::GetDistanceBetweenPoints: number of points not valid\n");
528 if (nPoints<2) return 0.;
529 if (nPoints<3) return TMath::Sqrt( (points[0]->X()-points[1]->X()) * (points[0]->X()-points[1]->X()) +
530 (points[0]->Y()-points[1]->Y()) * (points[0]->Y()-points[1]->Y()) +
531 (points[0]->Z()-points[1]->Z()) * (points[0]->Z()-points[1]->Z()) );
533 const Int_t nEdgesMax = ((AliMFTConstants::fNMaxMuonsForPCA) * (AliMFTConstants::fNMaxMuonsForPCA - 1)) / 2;
535 Int_t startID[nEdgesMax] = {0};
536 Int_t stopID[nEdgesMax] = {0};
537 Double_t edgeLength[nEdgesMax] = {0};
539 Bool_t pointStatus[AliMFTConstants::fNMaxMuonsForPCA] = {0};
542 for (Int_t i=0; i<nPoints-1; i++) {
543 for (Int_t j=i+1; j<nPoints; j++) {
544 edgeLength[nEdges] = TMath::Sqrt( (points[i]->X()-points[j]->X()) * (points[i]->X()-points[j]->X()) +
545 (points[i]->Y()-points[j]->Y()) * (points[i]->Y()-points[j]->Y()) +
546 (points[i]->Z()-points[j]->Z()) * (points[i]->Z()-points[j]->Z()) );
558 for (Int_t iEdge=0; iEdge<nEdges-1; iEdge++) {
559 min = edgeLength[iEdge];
561 for (Int_t j=iEdge+1; j<nEdges; j++) {
562 if (edgeLength[j]<min) {
570 Double_t edgeLengthMin = edgeLength[iMin];
571 Int_t startIDmin = startID[iMin];
572 Int_t stopIDmin = stopID[iMin];
574 edgeLength[iMin] = edgeLength[iEdge];
575 startID[iMin] = startID[iEdge];
576 stopID[iMin] = stopID[iEdge];
578 edgeLength[iEdge] = edgeLengthMin;
579 startID[iEdge] = startIDmin;
580 stopID[iEdge] = stopIDmin;
588 Double_t length = 0.;
589 for (Int_t i=0; i<nEdges; i++) {
590 if (!(pointStatus[startID[i]] && pointStatus[stopID[i]])) {
591 pointStatus[startID[i]] = kTRUE;
592 pointStatus[stopID[i]] = kTRUE;
593 length += edgeLength[i];
601 //====================================================================================================================================================
603 void AliMFTAnalysisTools::ConvertCovMatrixMUON2AOD(const TMatrixD& covMUON, Double_t covAOD[21]) {
605 // Converts the cov matrix from the MUON format (TMatrixD) to the AOD one (Double_t[21])
607 // Cov(x,x) ... : cv[0]
608 // Cov(x,slopeX) ... : cv[1] cv[2]
609 // Cov(x,y) ... : cv[3] cv[4] cv[5]
610 // Cov(x,slopeY) ... : cv[6] cv[7] cv[8] cv[9]
611 // Cov(x,invP_yz) ... : cv[10] cv[11] cv[12] cv[13] cv[14]
612 // not-used ... : cv[15] cv[16] cv[17] cv[18] cv[19] cv[20]
614 covAOD[0] = covMUON(0,0);
616 covAOD[1] = covMUON(1,0);
617 covAOD[2] = covMUON(1,1);
619 covAOD[3] = covMUON(2,0);
620 covAOD[4] = covMUON(2,1);
621 covAOD[5] = covMUON(2,2);
623 covAOD[6] = covMUON(3,0);
624 covAOD[7] = covMUON(3,1);
625 covAOD[8] = covMUON(3,2);
626 covAOD[9] = covMUON(3,3);
628 covAOD[10] = covMUON(4,0);
629 covAOD[11] = covMUON(4,1);
630 covAOD[12] = covMUON(4,2);
631 covAOD[13] = covMUON(4,3);
632 covAOD[14] = covMUON(4,4);
643 //====================================================================================================================================================
645 const TMatrixD AliMFTAnalysisTools::ConvertCovMatrixAOD2MUON(AliAODTrack *muon) {
647 Double_t covAOD[21] = {0};
648 muon -> GetCovarianceXYZPxPyPz(covAOD);
650 TMatrixD covMUON(5,5);
652 covMUON(0,0) = covAOD[0];
654 covMUON(1,0) = covAOD[1];
655 covMUON(1,1) = covAOD[2];
657 covMUON(2,0) = covAOD[3];
658 covMUON(2,1) = covAOD[4];
659 covMUON(2,2) = covAOD[5];
661 covMUON(3,0) = covAOD[6];
662 covMUON(3,1) = covAOD[7];
663 covMUON(3,2) = covAOD[8];
664 covMUON(3,3) = covAOD[9];
666 covMUON(4,0) = covAOD[10];
667 covMUON(4,1) = covAOD[11];
668 covMUON(4,2) = covAOD[12];
669 covMUON(4,3) = covAOD[13];
670 covMUON(4,4) = covAOD[14];
676 //====================================================================================================================================================