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"
35 #include "AliMFTAnalysisTools.h"
37 ClassImp(AliMFTAnalysisTools)
39 //====================================================================================================================================================
41 Bool_t AliMFTAnalysisTools::ExtrapAODMuonToZ(AliAODTrack *muon, Double_t z, Double_t xy[2]) {
43 if (!(muon->Pz()!=0)) return kFALSE;
45 AliMUONTrackParam *param = new AliMUONTrackParam();
47 param -> SetNonBendingCoor(muon->XAtDCA());
48 param -> SetBendingCoor(muon->YAtDCA());
49 param -> SetZ(AliMFTConstants::fZEvalKinem);
50 param -> SetNonBendingSlope(muon->Px()/muon->Pz());
51 param -> SetBendingSlope(muon->Py()/muon->Pz());
52 param -> SetInverseBendingMomentum( muon->Charge() * (1./muon->Pz()) / (TMath::Sqrt(1+TMath::Power(muon->Py()/muon->Pz(),2))) );
54 AliMUONTrackExtrap::ExtrapToZ(param, z);
55 xy[0] = param->GetNonBendingCoor();
56 xy[1] = param->GetBendingCoor();
64 //====================================================================================================================================================
66 Bool_t AliMFTAnalysisTools::ExtrapAODMuonToZ(AliAODTrack *muon, Double_t z, Double_t xy[2], TLorentzVector &kinem) {
68 if (!(muon->Pz()!=0)) return kFALSE;
70 AliMUONTrackParam *param = new AliMUONTrackParam();
72 param -> SetNonBendingCoor(muon->XAtDCA());
73 param -> SetBendingCoor(muon->YAtDCA());
74 param -> SetZ(AliMFTConstants::fZEvalKinem);
75 param -> SetNonBendingSlope(muon->Px()/muon->Pz());
76 param -> SetBendingSlope(muon->Py()/muon->Pz());
77 param -> SetInverseBendingMomentum( muon->Charge() * (1./muon->Pz()) / (TMath::Sqrt(1+TMath::Power(muon->Py()/muon->Pz(),2))) );
79 AliMUONTrackExtrap::ExtrapToZ(param, z);
80 xy[0] = param->GetNonBendingCoor();
81 xy[1] = param->GetBendingCoor();
83 Double_t massMu = TDatabasePDG::Instance()->GetParticle("mu-")->Mass();
84 Double_t energy = TMath::Sqrt(massMu*massMu + param->Px()*param->Px() + param->Py()*param->Py() + param->Pz()*param->Pz());
86 kinem.SetPxPyPzE(param->Px(), param->Py(), param->Pz(), energy);
94 //====================================================================================================================================================
96 Bool_t AliMFTAnalysisTools::ExtrapAODMuonToZ(AliAODTrack *muon, Double_t z, Double_t xy[2], TLorentzVector &kinem, TMatrixD &cov) {
98 if (!(muon->Pz()!=0)) return kFALSE;
100 AliMUONTrackParam *param = new AliMUONTrackParam();
102 param -> SetNonBendingCoor(muon->XAtDCA());
103 param -> SetBendingCoor(muon->YAtDCA());
104 param -> SetZ(AliMFTConstants::fZEvalKinem);
105 param -> SetNonBendingSlope(muon->Px()/muon->Pz());
106 param -> SetBendingSlope(muon->Py()/muon->Pz());
107 param -> SetInverseBendingMomentum( muon->Charge() * (1./muon->Pz()) / (TMath::Sqrt(1+TMath::Power(muon->Py()/muon->Pz(),2))) );
109 param -> SetCovariances(ConvertCovMatrixAOD2MUON(muon));
111 AliMUONTrackExtrap::ExtrapToZCov(param, z);
112 xy[0] = param->GetNonBendingCoor();
113 xy[1] = param->GetBendingCoor();
115 Double_t massMu = TDatabasePDG::Instance()->GetParticle("mu-")->Mass();
116 Double_t energy = TMath::Sqrt(massMu*massMu + param->Px()*param->Px() + param->Py()*param->Py() + param->Pz()*param->Pz());
118 kinem.SetPxPyPzE(param->Px(), param->Py(), param->Pz(), energy);
120 cov = param->GetCovariances();
128 //====================================================================================================================================================
130 Double_t AliMFTAnalysisTools::GetAODMuonOffset(AliAODTrack *muon, Double_t xv, Double_t yv, Double_t zv) {
132 Double_t xy[2] = {0};
133 ExtrapAODMuonToZ(muon, zv, xy);
135 return TMath::Sqrt((xv-xy[0])*(xv-xy[0]) + (yv-xy[1])*(yv-xy[1]));
139 //====================================================================================================================================================
141 Double_t AliMFTAnalysisTools::GetAODMuonWeightedOffset(AliAODTrack *muon, Double_t xv, Double_t yv, Double_t zv) {
143 Double_t xy[2] = {0};
144 TLorentzVector kinem(0,0,0,0);
147 ExtrapAODMuonToZ(muon, zv, xy, kinem, cov);
149 TMatrixD covCoordinates(2,2);
150 covCoordinates(0,0) = cov(0,0);
151 covCoordinates(0,1) = cov(0,2);
152 covCoordinates(1,0) = cov(2,0);
153 covCoordinates(1,1) = cov(2,2);
155 TMatrixD covCoordinatesInverse = covCoordinates.Invert();
157 Double_t dX = xy[0] - xv;
158 Double_t dY = xy[1] - yv;
160 Double_t weightedOffset = TMath::Sqrt(0.5*(dX*dX*covCoordinatesInverse(0,0) +
161 dY*dY*covCoordinatesInverse(1,1) +
162 2.*dX*dY*covCoordinatesInverse(0,1)));
164 return weightedOffset;
168 //====================================================================================================================================================
170 Double_t AliMFTAnalysisTools::GetPseudoProperDecayTimeXY(Double_t xVtx, Double_t yVtx,
171 Double_t xDimu, Double_t yDimu,
172 Double_t mDimu, Double_t ptDimu) {
174 // pseudo-proper decay time of a particle produced in the primary vertex and decaying into a dimuon (+ X)
175 // evaluated using the transverse degree of freedom of the decay topology
178 Double_t decayLengthXY = TMath::Sqrt((xVtx-xDimu)*(xVtx-xDimu)+(yVtx-yDimu)*(yVtx-yDimu));
179 return (decayLengthXY * mDimu/ptDimu)/TMath::Ccgs()*1E12; // in ps
186 //====================================================================================================================================================
188 Double_t AliMFTAnalysisTools::GetPseudoProperDecayTimeZ(Double_t zVtx,
190 Double_t mDimu, Double_t pzDimu) {
192 // pseudo-proper decay time of a particle produced in the primary vertex and decaying into a dimuon (+ X)
193 // evaluated using the longitudinal degree of freedom of the decay topology
196 Double_t decayLengthZ = zDimu - zVtx;
197 return (decayLengthZ * mDimu/pzDimu)/TMath::Ccgs()*1E12; // in ps
204 //====================================================================================================================================================
206 Bool_t AliMFTAnalysisTools::CalculatePCA(AliAODDimuon *dimuon, Double_t *pca, Double_t &pcaQuality, TLorentzVector &kinem) {
208 TObjArray *muons = new TObjArray();
209 muons -> Add(dimuon->GetMu(0));
210 muons -> Add(dimuon->GetMu(1));
212 return CalculatePCA(muons, pca, pcaQuality, kinem);
216 //====================================================================================================================================================
218 Bool_t AliMFTAnalysisTools::CalculatePCA(TObjArray *muons, Double_t *pca, Double_t &pcaQuality, TLorentzVector &kinem) {
220 const Int_t nMuons = muons->GetEntriesFast();
221 if (nMuons<2 || nMuons>AliMFTConstants::fNMaxMuonsForPCA) {
222 printf("W-AliMFTAnalysisTools::CalculatePCA: number of muons not valid\n");
226 Double_t fXPointOfClosestApproach=0, fYPointOfClosestApproach=0, fZPointOfClosestApproach=0;
228 AliAODTrack *muon[nMuons];
229 AliMUONTrackParam *param[nMuons];
231 // Finding AliMUONTrackParam objects for each muon
233 for (Int_t iMu=0; iMu<nMuons; iMu++) {
234 muon[iMu] = (AliAODTrack*) muons->At(iMu);
235 if (TMath::Abs(muon[iMu]->Pz())<1.e-6) return kFALSE;
236 param[iMu] = new AliMUONTrackParam();
237 param[iMu] -> SetNonBendingCoor(muon[iMu]->XAtDCA());
238 param[iMu] -> SetBendingCoor(muon[iMu]->YAtDCA());
239 param[iMu] -> SetZ(0.);
240 param[iMu] -> SetNonBendingSlope(muon[iMu]->Px()/muon[iMu]->Pz());
241 param[iMu] -> SetBendingSlope(muon[iMu]->Py()/muon[iMu]->Pz());
242 param[iMu] -> SetInverseBendingMomentum( muon[iMu]->Charge() * (1./muon[iMu]->Pz()) / (TMath::Sqrt(1+TMath::Power(muon[iMu]->Py()/muon[iMu]->Pz(),2))) );
245 // here we want to understand in which direction we have to search the minimum...
247 Double_t step = 1.; // initial step, in cm
248 Double_t startPoint = 0.;
250 Double_t r[3]={0}, z[3]={startPoint, startPoint+step, startPoint+2*step};
252 TVector3 **points = new TVector3*[nMuons];
254 for (Int_t i=0; i<3; i++) {
255 for (Int_t iMu=0; iMu<nMuons; iMu++) {
256 AliMUONTrackExtrap::ExtrapToZ(param[iMu], z[i]);
257 points[iMu] = new TVector3(param[iMu]->GetNonBendingCoor(),param[iMu]->GetBendingCoor(),z[i]);
259 r[i] = GetDistanceBetweenPoints(points,nMuons);
262 Int_t researchDirection = 0;
264 if (r[0]>r[1] && r[1]>r[2]) researchDirection = +1; // towards z positive
265 else if (r[0]<r[1] && r[1]<r[2]) researchDirection = -1; // towards z negative
266 else if (r[0]<r[1] && r[1]>r[2]) {
267 printf("E-AliMFTAnalysisTools::CalculatePCA: Point of closest approach cannot be found for dimuon (no minima)\n");
271 while (TMath::Abs(researchDirection)>0.5) {
273 if (researchDirection>0) {
276 z[2] = z[1]+researchDirection*step;
281 z[0] = z[1]+researchDirection*step;
283 if (TMath::Abs(z[0])>900.) {
284 printf("E-AliMFTAnalysisTools::CalculatePCA: Point of closest approach cannot be found for dimuon (no minima in the fiducial region)\n");
288 for (Int_t iMu=0; iMu<nMuons; iMu++) delete points[iMu];
290 for (Int_t i=0; i<3; i++) {
291 for (Int_t iMu=0; iMu<nMuons; iMu++) {
292 AliMUONTrackExtrap::ExtrapToZ(param[iMu], z[i]);
293 points[iMu] = new TVector3(param[iMu]->GetNonBendingCoor(),param[iMu]->GetBendingCoor(),z[i]);
295 r[i] = GetDistanceBetweenPoints(points,nMuons);
298 if (r[0]>r[1] && r[1]>r[2]) researchDirection = +1; // towards z positive
299 else if (r[0]<r[1] && r[1]<r[2]) researchDirection = -1; // towards z negative
303 // now we now that the minimum is between z[0] and z[2] and we search for it
306 while (step>AliMFTConstants::fPrecisionPointOfClosestApproach) {
309 for (Int_t i=0; i<3; i++) {
310 for (Int_t iMu=0; iMu<nMuons; iMu++) {
311 AliMUONTrackExtrap::ExtrapToZ(param[iMu], z[i]);
312 points[iMu] = new TVector3(param[iMu]->GetNonBendingCoor(),param[iMu]->GetBendingCoor(),z[i]);
314 r[i] = GetDistanceBetweenPoints(points,nMuons);
316 if (r[0]<r[1]) z[1] = z[0];
317 else if (r[2]<r[1]) z[1] = z[2];
323 // Once z of minimum is found, we evaluate the x and y coordinates by averaging over the contributing tracks
325 fZPointOfClosestApproach = z[1];
326 fXPointOfClosestApproach = 0.;
327 fYPointOfClosestApproach = 0.;
328 for (Int_t iMu=0; iMu<nMuons; iMu++) {
329 AliMUONTrackExtrap::ExtrapToZ(param[iMu], fZPointOfClosestApproach);
330 fXPointOfClosestApproach += param[iMu]->GetNonBendingCoor();
331 fYPointOfClosestApproach += param[iMu]->GetBendingCoor();
333 fXPointOfClosestApproach /= Double_t(nMuons);
334 fYPointOfClosestApproach /= Double_t(nMuons);
336 pca[0] = fXPointOfClosestApproach;
337 pca[1] = fYPointOfClosestApproach;
338 pca[2] = fZPointOfClosestApproach;
340 // Evaluating the kinematics of the N-muon
342 Double_t pTot[3] = {0};
344 Double_t massMu = TDatabasePDG::Instance()->GetParticle("mu-")->Mass();
345 for (Int_t iMu=0; iMu<nMuons; iMu++) {
346 pTot[0] += param[iMu]->Px();
347 pTot[1] += param[iMu]->Py();
348 pTot[2] += param[iMu]->Pz();
349 ene += TMath::Sqrt(massMu*massMu + param[iMu]->Px()*param[iMu]->Px() + param[iMu]->Py()*param[iMu]->Py() + param[iMu]->Pz()*param[iMu]->Pz());
352 kinem.SetPxPyPzE(pTot[0], pTot[1], pTot[2], ene);
354 // Evaluating the PCA quality of the N-muon
356 Double_t sum=0.,squareSum=0.;
357 for (Int_t iMu=0; iMu<nMuons; iMu++) {
358 Double_t wOffset = AliMFTAnalysisTools::GetAODMuonWeightedOffset(muon[iMu],fXPointOfClosestApproach, fYPointOfClosestApproach, fZPointOfClosestApproach);
359 Double_t f = TMath::Exp(-0.5 * wOffset);
363 if (sum > 0.) pcaQuality = (sum-squareSum/sum) / (nMuons-1);
364 else pcaQuality = 0.;
370 //=========================================================================================================================
372 Double_t AliMFTAnalysisTools::GetDistanceBetweenPoints(TVector3 **points, Int_t nPoints) {
374 if (nPoints>AliMFTConstants::fNMaxMuonsForPCA) {
375 printf("W-AliMFTAnalysisTools::GetDistanceBetweenPoints: number of points not valid\n");
379 if (nPoints<2) return 0.;
380 if (nPoints<3) return TMath::Sqrt( (points[0]->X()-points[1]->X()) * (points[0]->X()-points[1]->X()) +
381 (points[0]->Y()-points[1]->Y()) * (points[0]->Y()-points[1]->Y()) +
382 (points[0]->Z()-points[1]->Z()) * (points[0]->Z()-points[1]->Z()) );
384 const Int_t nEdgesMax = ((AliMFTConstants::fNMaxMuonsForPCA) * (AliMFTConstants::fNMaxMuonsForPCA - 1)) / 2;
386 Int_t startID[nEdgesMax] = {0};
387 Int_t stopID[nEdgesMax] = {0};
388 Double_t edgeLength[nEdgesMax] = {0};
390 Bool_t pointStatus[AliMFTConstants::fNMaxMuonsForPCA] = {0};
393 for (Int_t i=0; i<nPoints-1; i++) {
394 for (Int_t j=i+1; j<nPoints; j++) {
395 edgeLength[nEdges] = TMath::Sqrt( (points[i]->X()-points[j]->X()) * (points[i]->X()-points[j]->X()) +
396 (points[i]->Y()-points[j]->Y()) * (points[i]->Y()-points[j]->Y()) +
397 (points[i]->Z()-points[j]->Z()) * (points[i]->Z()-points[j]->Z()) );
409 for (Int_t iEdge=0; iEdge<nEdges-1; iEdge++) {
410 min = edgeLength[iEdge];
412 for (Int_t j=iEdge+1; j<nEdges; j++) {
413 if (edgeLength[j]<min) {
421 Double_t edgeLengthMin = edgeLength[iMin];
422 Int_t startIDmin = startID[iMin];
423 Int_t stopIDmin = stopID[iMin];
425 edgeLength[iMin] = edgeLength[iEdge];
426 startID[iMin] = startID[iEdge];
427 stopID[iMin] = stopID[iEdge];
429 edgeLength[iEdge] = edgeLengthMin;
430 startID[iEdge] = startIDmin;
431 stopID[iEdge] = stopIDmin;
439 Double_t length = 0.;
440 for (Int_t i=0; i<nEdges; i++) {
441 if (!(pointStatus[startID[i]] && pointStatus[stopID[i]])) {
442 pointStatus[startID[i]] = kTRUE;
443 pointStatus[stopID[i]] = kTRUE;
444 length += edgeLength[i];
452 //====================================================================================================================================================
454 void AliMFTAnalysisTools::ConvertCovMatrixMUON2AOD(const TMatrixD& covMUON, Double_t covAOD[21]) {
456 // Converts the cov matrix from the MUON format (TMatrixD) to the AOD one (Double_t[21])
458 // Cov(x,x) ... : cv[0]
459 // Cov(x,slopeX) ... : cv[1] cv[2]
460 // Cov(x,y) ... : cv[3] cv[4] cv[5]
461 // Cov(x,slopeY) ... : cv[6] cv[7] cv[8] cv[9]
462 // Cov(x,invP_yz) ... : cv[10] cv[11] cv[12] cv[13] cv[14]
463 // not-used ... : cv[15] cv[16] cv[17] cv[18] cv[19] cv[20]
465 covAOD[0] = covMUON(0,0);
467 covAOD[1] = covMUON(1,0);
468 covAOD[2] = covMUON(1,1);
470 covAOD[3] = covMUON(2,0);
471 covAOD[4] = covMUON(2,1);
472 covAOD[5] = covMUON(2,2);
474 covAOD[6] = covMUON(3,0);
475 covAOD[7] = covMUON(3,1);
476 covAOD[8] = covMUON(3,2);
477 covAOD[9] = covMUON(3,3);
479 covAOD[10] = covMUON(4,0);
480 covAOD[11] = covMUON(4,1);
481 covAOD[12] = covMUON(4,2);
482 covAOD[13] = covMUON(4,3);
483 covAOD[14] = covMUON(4,4);
494 //====================================================================================================================================================
496 const TMatrixD AliMFTAnalysisTools::ConvertCovMatrixAOD2MUON(AliAODTrack *muon) {
498 Double_t covAOD[21] = {0};
499 muon -> GetCovarianceXYZPxPyPz(covAOD);
501 TMatrixD covMUON(5,5);
503 covMUON(0,0) = covAOD[0];
505 covMUON(1,0) = covAOD[1];
506 covMUON(1,1) = covAOD[2];
508 covMUON(2,0) = covAOD[3];
509 covMUON(2,1) = covAOD[4];
510 covMUON(2,2) = covAOD[5];
512 covMUON(3,0) = covAOD[6];
513 covMUON(3,1) = covAOD[7];
514 covMUON(3,2) = covAOD[8];
515 covMUON(3,3) = covAOD[9];
517 covMUON(4,0) = covAOD[10];
518 covMUON(4,1) = covAOD[11];
519 covMUON(4,2) = covAOD[12];
520 covMUON(4,3) = covAOD[13];
521 covMUON(4,4) = covAOD[14];
527 //====================================================================================================================================================