]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MFT/AliMFTAnalysisTools.cxx
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / MFT / AliMFTAnalysisTools.cxx
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
16 //====================================================================================================================================================
17 //
18 //      Support class for various common operation on MFT objects
19 //
20 //      Contact author: antonio.uras@cern.ch
21 //
22 //====================================================================================================================================================
23
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"
31 #include "TMath.h"
32 #include "AliLog.h"
33 #include "TObjArray.h"
34 #include "TDecompLU.h"
35
36 #include "AliMFTAnalysisTools.h"
37
38 ClassImp(AliMFTAnalysisTools)
39
40 //====================================================================================================================================================
41
42 Bool_t AliMFTAnalysisTools::ExtrapAODMuonToZ(AliAODTrack *muon, Double_t z, Double_t xy[2]) {
43
44   if (!(muon->Pz()!=0)) return kFALSE;
45
46   AliMUONTrackParam *param = new AliMUONTrackParam();
47
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))) );
54
55   AliMUONTrackExtrap::ExtrapToZ(param, z);
56   xy[0] = param->GetNonBendingCoor();
57   xy[1] = param->GetBendingCoor();
58
59   delete param;
60
61   return kTRUE;
62
63 }
64
65 //====================================================================================================================================================
66
67 Bool_t AliMFTAnalysisTools::ExtrapAODMuonToZ(AliAODTrack *muon, Double_t z, Double_t xy[2], TLorentzVector &kinem) {
68
69   if (!(muon->Pz()!=0)) return kFALSE;
70
71   AliMUONTrackParam *param = new AliMUONTrackParam();
72
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))) );
79
80   AliMUONTrackExtrap::ExtrapToZ(param, z);
81   xy[0] = param->GetNonBendingCoor();
82   xy[1] = param->GetBendingCoor();
83
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());
86   
87   kinem.SetPxPyPzE(param->Px(), param->Py(), param->Pz(), energy);
88
89   delete param;
90
91   return kTRUE;
92
93 }
94
95 //====================================================================================================================================================
96
97 Bool_t AliMFTAnalysisTools::ExtrapAODMuonToZ(AliAODTrack *muon, Double_t z, Double_t xy[2], TLorentzVector &kinem, TMatrixD &cov) {
98
99   if (!(muon->Pz()!=0)) return kFALSE;
100
101   AliMUONTrackParam *param = new AliMUONTrackParam();
102
103   param -> SetNonBendingCoor(muon->XAtDCA());
104   param -> SetBendingCoor(muon->YAtDCA());
105   param -> SetZ(AliMFTConstants::fZEvalKinem);
106   param -> SetNonBendingSlope(muon->Px()/muon->Pz());
107   param -> SetBendingSlope(muon->Py()/muon->Pz());
108   param -> SetInverseBendingMomentum( muon->Charge() * (1./muon->Pz()) / (TMath::Sqrt(1+TMath::Power(muon->Py()/muon->Pz(),2))) );
109
110   param -> SetCovariances(ConvertCovMatrixAOD2MUON(muon));
111
112   AliMUONTrackExtrap::ExtrapToZCov(param, z);
113   xy[0] = param->GetNonBendingCoor();
114   xy[1] = param->GetBendingCoor();
115
116   Double_t massMu = TDatabasePDG::Instance()->GetParticle("mu-")->Mass();
117   Double_t energy = TMath::Sqrt(massMu*massMu + param->Px()*param->Px() + param->Py()*param->Py() + param->Pz()*param->Pz());
118   
119   kinem.SetPxPyPzE(param->Px(), param->Py(), param->Pz(), energy);
120
121   cov = param->GetCovariances();
122
123   delete param;
124
125   return kTRUE;
126
127 }
128
129 //====================================================================================================================================================
130
131 Double_t AliMFTAnalysisTools::GetAODMuonOffset(AliAODTrack *muon, Double_t xv, Double_t yv, Double_t zv) {
132
133   Double_t xy[2] = {0};
134   ExtrapAODMuonToZ(muon, zv, xy);
135   
136   return TMath::Sqrt((xv-xy[0])*(xv-xy[0]) + (yv-xy[1])*(yv-xy[1]));
137
138 }
139
140 //====================================================================================================================================================
141
142 Double_t AliMFTAnalysisTools::GetAODMuonWeightedOffset(AliAODTrack *muon, Double_t xv, Double_t yv, Double_t zv) {
143
144   Double_t xy[2] = {0};
145   TLorentzVector kinem(0,0,0,0);
146   TMatrixD cov(5,5);
147
148   ExtrapAODMuonToZ(muon, zv, xy, kinem, cov);
149
150   TMatrixD covCoordinates(2,2);
151   covCoordinates(0,0) = cov(0,0);
152   covCoordinates(0,1) = cov(0,2);
153   covCoordinates(1,0) = cov(2,0);
154   covCoordinates(1,1) = cov(2,2);
155   
156   if (TDecompLU::InvertLU(covCoordinates,covCoordinates.GetTol(),0)) {
157
158     TMatrixD covCoordinatesInverse = covCoordinates;
159     Double_t dX = xy[0] - xv;
160     Double_t dY = xy[1] - yv;
161     
162     Double_t weightedOffset = TMath::Sqrt(0.5*(dX*dX*covCoordinatesInverse(0,0) +
163                                                dY*dY*covCoordinatesInverse(1,1) +
164                                                2.*dX*dY*covCoordinatesInverse(0,1)));
165     
166     return weightedOffset;
167
168   }
169   
170   return -999;
171
172 }
173
174 //====================================================================================================================================================
175
176 Double_t AliMFTAnalysisTools::GetPseudoProperDecayTimeXY(Double_t xVtx,  Double_t yVtx, 
177                                                          Double_t xDimu, Double_t yDimu, 
178                                                          Double_t mDimu, Double_t ptDimu) {
179   
180   // pseudo-proper decay time of a particle produced in the primary vertex and decaying into a dimuon (+ X)
181   // evaluated using the transverse degree of freedom of the decay topology 
182
183   if (ptDimu != 0) {
184     Double_t decayLengthXY = TMath::Sqrt((xVtx-xDimu)*(xVtx-xDimu)+(yVtx-yDimu)*(yVtx-yDimu));
185     return (decayLengthXY * mDimu/ptDimu)/TMath::Ccgs()*1E12;   // in ps
186   }
187   
188   return -99999999;
189   
190 }
191
192 //====================================================================================================================================================
193
194 Double_t AliMFTAnalysisTools::GetPseudoProperDecayTimeZ(Double_t zVtx, 
195                                                         Double_t zDimu, 
196                                                         Double_t mDimu, Double_t pzDimu) {
197
198   // pseudo-proper decay time of a particle produced in the primary vertex and decaying into a dimuon (+ X)
199   // evaluated using the longitudinal degree of freedom of the decay topology 
200
201   if (pzDimu != 0) {
202     Double_t decayLengthZ = zDimu - zVtx;
203     return (decayLengthZ * mDimu/pzDimu)/TMath::Ccgs()*1E12;  // in ps
204   }
205
206   return -99999999;
207
208 }
209
210 //====================================================================================================================================================
211
212 Bool_t AliMFTAnalysisTools::CalculatePCA(AliAODDimuon *dimuon, Double_t *pca, Double_t &pcaQuality, TLorentzVector &kinem) {
213
214   TObjArray *muons = new TObjArray();
215   muons -> Add(dimuon->GetMu(0));
216   muons -> Add(dimuon->GetMu(1));
217   
218   Bool_t result = CalculatePCA(muons, pca, pcaQuality, kinem);
219   delete muons;
220   return result;
221
222 }
223
224 //====================================================================================================================================================
225
226 Bool_t AliMFTAnalysisTools::CalculatePCA(TObjArray *muons, Double_t *pca, Double_t &pcaQuality, TLorentzVector &kinem) {
227   
228   const Int_t nMuons = muons->GetEntriesFast();
229   if (nMuons<2 || nMuons>AliMFTConstants::fNMaxMuonsForPCA) {
230     printf("W-AliMFTAnalysisTools::CalculatePCA: number of muons not valid\n");
231     return kFALSE;
232   }
233   
234   Double_t fXPointOfClosestApproach=0, fYPointOfClosestApproach=0, fZPointOfClosestApproach=0;
235
236   AliAODTrack *muon[AliMFTConstants::fNMaxMuonsForPCA]        = {0};
237   AliMUONTrackParam *param[AliMFTConstants::fNMaxMuonsForPCA] = {0};
238
239   // Finding AliMUONTrackParam objects for each muon
240   
241   for (Int_t iMu=0; iMu<nMuons; iMu++) {
242     muon[iMu] = (AliAODTrack*) muons->At(iMu);
243     if (TMath::Abs(muon[iMu]->Pz())<1.e-6) {
244       for(Int_t i=0;i<iMu;i++) delete param[i];
245       return kFALSE;
246     }
247     param[iMu] = new AliMUONTrackParam();
248     param[iMu] -> SetNonBendingCoor(muon[iMu]->XAtDCA());
249     param[iMu] -> SetBendingCoor(muon[iMu]->YAtDCA());
250     param[iMu] -> SetZ(0.);
251     param[iMu] -> SetNonBendingSlope(muon[iMu]->Px()/muon[iMu]->Pz());
252     param[iMu] -> SetBendingSlope(muon[iMu]->Py()/muon[iMu]->Pz());
253     param[iMu] -> SetInverseBendingMomentum( muon[iMu]->Charge() * (1./muon[iMu]->Pz()) / (TMath::Sqrt(1+TMath::Power(muon[iMu]->Py()/muon[iMu]->Pz(),2))) );
254   }
255   
256   // here we want to understand in which direction we have to search the minimum...
257   
258   Double_t step = 1.;  // initial step, in cm
259   Double_t startPoint = 0.;
260   
261   Double_t r[3]={0}, z[3]={startPoint, startPoint+step, startPoint+2*step};
262   
263   TVector3 **points = new TVector3*[AliMFTConstants::fNMaxMuonsForPCA];
264   
265   for (Int_t i=0; i<3; i++) {
266     for (Int_t iMu=0; iMu<nMuons; iMu++) {
267       AliMUONTrackExtrap::ExtrapToZ(param[iMu], z[i]);
268       points[iMu] = new TVector3(param[iMu]->GetNonBendingCoor(),param[iMu]->GetBendingCoor(),z[i]);
269                 }
270     r[i] = GetDistanceBetweenPoints(points,nMuons);
271     for (Int_t iMu=0; iMu<nMuons; iMu++) delete points[iMu]; 
272   }
273   
274   Int_t researchDirection = 0;
275   
276   if      (r[0]>r[1] && r[1]>r[2]) researchDirection = +1;   // towards z positive
277   else if (r[0]<r[1] && r[1]<r[2]) researchDirection = -1;   // towards z negative
278   else if (r[0]<r[1] && r[1]>r[2]) {
279     printf("E-AliMFTAnalysisTools::CalculatePCA: Point of closest approach cannot be found for dimuon (no minima)\n");
280     for (Int_t iMu=0;iMu<nMuons;iMu++) delete param[iMu];
281     delete points;
282     return kFALSE;
283   }     
284   
285   while (TMath::Abs(researchDirection)>0.5) {
286       
287     if (researchDirection>0) {
288       z[0] = z[1];
289       z[1] = z[2];
290       z[2] = z[1]+researchDirection*step;
291     }
292     else {
293       z[2] = z[1];
294       z[1] = z[0];
295       z[0] = z[1]+researchDirection*step;
296     }
297     if (TMath::Abs(z[0])>900.) {
298       printf("E-AliMFTAnalysisTools::CalculatePCA: Point of closest approach cannot be found for dimuon (no minima in the fiducial region)\n");
299       for (Int_t iMu=0;iMu<nMuons;iMu++) delete param[iMu];
300       delete points;
301       return kFALSE;
302     }
303     
304     for (Int_t i=0; i<3; i++) {
305       for (Int_t iMu=0; iMu<nMuons; iMu++) {
306         AliMUONTrackExtrap::ExtrapToZ(param[iMu], z[i]);
307         points[iMu] = new TVector3(param[iMu]->GetNonBendingCoor(),param[iMu]->GetBendingCoor(),z[i]);
308       }      
309       r[i] = GetDistanceBetweenPoints(points,nMuons);
310       for (Int_t iMu=0;iMu<nMuons;iMu++) delete points[iMu];
311     }
312     researchDirection=0;
313     if      (r[0]>r[1] && r[1]>r[2]) researchDirection = +1;   // towards z positive
314     else if (r[0]<r[1] && r[1]<r[2]) researchDirection = -1;   // towards z negative
315     
316   }
317   
318   // now we now that the minimum is between z[0] and z[2] and we search for it
319   
320   step *= 0.5;
321   while (step>AliMFTConstants::fPrecisionPointOfClosestApproach) {
322     z[0] = z[1]-step;
323     z[2] = z[1]+step;
324     for (Int_t i=0; i<3; i++) {
325       for (Int_t iMu=0; iMu<nMuons; iMu++) {
326         AliMUONTrackExtrap::ExtrapToZ(param[iMu], z[i]);
327         points[iMu] = new TVector3(param[iMu]->GetNonBendingCoor(),param[iMu]->GetBendingCoor(),z[i]);
328       }      
329       r[i] = GetDistanceBetweenPoints(points,nMuons);
330       for (Int_t iMu=0;iMu<nMuons;iMu++)        delete points[iMu];
331     }
332     if      (r[0]<r[1]) z[1] = z[0];
333     else if (r[2]<r[1]) z[1] = z[2];
334     else step *= 0.5;
335   }
336   
337   // Once z of minimum is found, we evaluate the x and y coordinates by averaging over the contributing tracks
338   
339   fZPointOfClosestApproach = z[1];
340   fXPointOfClosestApproach = 0.;
341   fYPointOfClosestApproach = 0.;
342   for (Int_t iMu=0; iMu<nMuons; iMu++) {
343     AliMUONTrackExtrap::ExtrapToZ(param[iMu], fZPointOfClosestApproach);
344     fXPointOfClosestApproach += param[iMu]->GetNonBendingCoor();
345     fYPointOfClosestApproach += param[iMu]->GetBendingCoor();
346   }
347   fXPointOfClosestApproach /= Double_t(nMuons);
348   fYPointOfClosestApproach /= Double_t(nMuons);
349   
350   pca[0] = fXPointOfClosestApproach;
351   pca[1] = fYPointOfClosestApproach;
352   pca[2] = fZPointOfClosestApproach;
353   
354   // Evaluating the kinematics of the N-muon
355   
356   Double_t pTot[3] = {0};
357   Double_t ene = 0.;
358   Double_t massMu = TDatabasePDG::Instance()->GetParticle("mu-")->Mass();
359   for (Int_t iMu=0; iMu<nMuons; iMu++) {
360     pTot[0] += param[iMu]->Px();
361     pTot[1] += param[iMu]->Py();
362     pTot[2] += param[iMu]->Pz();
363     ene += TMath::Sqrt(massMu*massMu + param[iMu]->Px()*param[iMu]->Px() + param[iMu]->Py()*param[iMu]->Py() + param[iMu]->Pz()*param[iMu]->Pz());
364   }
365   
366   kinem.SetPxPyPzE(pTot[0], pTot[1], pTot[2], ene);
367   
368   // Evaluating the PCA quality of the N-muon
369   
370   Double_t sum=0.,squareSum=0.;
371   for (Int_t iMu=0; iMu<nMuons; iMu++) {
372     Double_t wOffset = AliMFTAnalysisTools::GetAODMuonWeightedOffset(muon[iMu],fXPointOfClosestApproach, fYPointOfClosestApproach, fZPointOfClosestApproach);
373     Double_t f = TMath::Exp(-0.5 * wOffset);
374     sum += f;
375     squareSum += f*f;
376   }
377   if (sum > 0.) pcaQuality =  (sum-squareSum/sum) / (nMuons-1);
378   else pcaQuality = 0.;
379   
380   for(Int_t iMu=0;iMu<nMuons;iMu++) delete param[iMu];
381   delete points;
382   return kTRUE;
383   
384 }
385
386 //=========================================================================================================================
387
388 Double_t AliMFTAnalysisTools::GetDistanceBetweenPoints(TVector3 **points, Int_t nPoints) {
389   
390   if (nPoints>AliMFTConstants::fNMaxMuonsForPCA) {
391     printf("W-AliMFTAnalysisTools::GetDistanceBetweenPoints: number of points not valid\n");
392     return 1.e9;
393   }
394   
395   if (nPoints<2) return 0.;
396   if (nPoints<3) return TMath::Sqrt( (points[0]->X()-points[1]->X()) * (points[0]->X()-points[1]->X()) +
397                                      (points[0]->Y()-points[1]->Y()) * (points[0]->Y()-points[1]->Y()) +
398                                      (points[0]->Z()-points[1]->Z()) * (points[0]->Z()-points[1]->Z()) );
399
400   const Int_t nEdgesMax = ((AliMFTConstants::fNMaxMuonsForPCA) * (AliMFTConstants::fNMaxMuonsForPCA - 1)) / 2;
401   
402   Int_t startID[nEdgesMax]       = {0};
403   Int_t stopID[nEdgesMax]        = {0};
404   Double_t edgeLength[nEdgesMax] = {0};
405
406   Bool_t pointStatus[AliMFTConstants::fNMaxMuonsForPCA] = {0};
407   
408   Int_t nEdges=0;
409   for (Int_t i=0; i<nPoints-1; i++) {
410     for (Int_t j=i+1; j<nPoints; j++) {
411       edgeLength[nEdges] = TMath::Sqrt( (points[i]->X()-points[j]->X()) * (points[i]->X()-points[j]->X()) +
412                                         (points[i]->Y()-points[j]->Y()) * (points[i]->Y()-points[j]->Y()) +
413                                         (points[i]->Z()-points[j]->Z()) * (points[i]->Z()-points[j]->Z()) );
414       stopID[nEdges]  = i;
415       startID[nEdges] = j;
416       nEdges++;
417     }
418   }
419  
420   // Order Edges
421
422   Double_t min = 0;
423   Int_t   iMin = 0;
424
425   for (Int_t iEdge=0; iEdge<nEdges-1; iEdge++) {
426     min  = edgeLength[iEdge];
427     iMin = iEdge;
428     for (Int_t j=iEdge+1; j<nEdges; j++) {
429       if (edgeLength[j]<min) {
430         min  = edgeLength[j];
431         iMin = j;
432       }
433     }
434     
435     if (iMin != iEdge) {
436
437       Double_t edgeLengthMin = edgeLength[iMin];
438       Int_t startIDmin = startID[iMin];
439       Int_t stopIDmin  = stopID[iMin];
440       
441       edgeLength[iMin] = edgeLength[iEdge];
442       startID[iMin]    = startID[iEdge];
443       stopID[iMin]     = stopID[iEdge];
444
445       edgeLength[iEdge] = edgeLengthMin;
446       startID[iEdge]    = startIDmin;
447       stopID[iEdge]     = stopIDmin;
448
449     }
450     
451   }
452   
453   // Connect
454
455   Double_t length = 0.;
456   for (Int_t i=0; i<nEdges; i++) {
457     if (!(pointStatus[startID[i]] && pointStatus[stopID[i]])) {
458       pointStatus[startID[i]] = kTRUE;
459       pointStatus[stopID[i]]  = kTRUE;
460       length += edgeLength[i];
461     }
462   }
463   
464   return length;
465   
466 }
467
468 //====================================================================================================================================================
469
470 void AliMFTAnalysisTools::ConvertCovMatrixMUON2AOD(const TMatrixD& covMUON, Double_t covAOD[21]) {
471
472   // Converts the cov matrix from the MUON format (TMatrixD) to the AOD one (Double_t[21])
473   // 
474   // Cov(x,x)       ... :   cv[0]
475   // Cov(x,slopeX)  ... :   cv[1]  cv[2]
476   // Cov(x,y)       ... :   cv[3]  cv[4]  cv[5]
477   // Cov(x,slopeY)  ... :   cv[6]  cv[7]  cv[8]  cv[9]
478   // Cov(x,invP_yz) ... :   cv[10] cv[11] cv[12] cv[13] cv[14]
479   // not-used       ... :   cv[15] cv[16] cv[17] cv[18] cv[19] cv[20]
480
481   covAOD[0]  = covMUON(0,0);
482
483   covAOD[1]  = covMUON(1,0);
484   covAOD[2]  = covMUON(1,1);
485
486   covAOD[3]  = covMUON(2,0);  
487   covAOD[4]  = covMUON(2,1);  
488   covAOD[5]  = covMUON(2,2);  
489
490   covAOD[6]  = covMUON(3,0);  
491   covAOD[7]  = covMUON(3,1);  
492   covAOD[8]  = covMUON(3,2);  
493   covAOD[9]  = covMUON(3,3);  
494
495   covAOD[10] = covMUON(4,0);  
496   covAOD[11] = covMUON(4,1);  
497   covAOD[12] = covMUON(4,2);  
498   covAOD[13] = covMUON(4,3);  
499   covAOD[14] = covMUON(4,4);  
500
501   covAOD[15] = 0;  
502   covAOD[16] = 0;  
503   covAOD[17] = 0;  
504   covAOD[18] = 0;  
505   covAOD[19] = 0;  
506   covAOD[20] = 0;  
507
508 }
509
510 //====================================================================================================================================================
511
512 const TMatrixD AliMFTAnalysisTools::ConvertCovMatrixAOD2MUON(AliAODTrack *muon) {
513
514   Double_t covAOD[21] = {0};
515   muon -> GetCovarianceXYZPxPyPz(covAOD);
516
517   TMatrixD covMUON(5,5);
518
519   covMUON(0,0) = covAOD[0];
520                                 
521   covMUON(1,0) = covAOD[1];
522   covMUON(1,1) = covAOD[2];
523                                 
524   covMUON(2,0) = covAOD[3];  
525   covMUON(2,1) = covAOD[4];  
526   covMUON(2,2) = covAOD[5];  
527                                 
528   covMUON(3,0) = covAOD[6];  
529   covMUON(3,1) = covAOD[7];  
530   covMUON(3,2) = covAOD[8];  
531   covMUON(3,3) = covAOD[9];  
532                                 
533   covMUON(4,0) = covAOD[10];  
534   covMUON(4,1) = covAOD[11];  
535   covMUON(4,2) = covAOD[12];  
536   covMUON(4,3) = covAOD[13];  
537   covMUON(4,4) = covAOD[14]; 
538
539   return covMUON;
540
541 }
542
543 //====================================================================================================================================================
544