]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MFT/AliMFTAnalysisTools.cxx
Merge branch 'master' into TPCdev
[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
35 #include "AliMFTAnalysisTools.h"
36
37 ClassImp(AliMFTAnalysisTools)
38
39 //====================================================================================================================================================
40
41 Bool_t AliMFTAnalysisTools::ExtrapAODMuonToZ(AliAODTrack *muon, Double_t z, Double_t xy[2]) {
42
43   if (!(muon->Pz()!=0)) return kFALSE;
44
45   AliMUONTrackParam *param = new AliMUONTrackParam();
46
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))) );
53
54   AliMUONTrackExtrap::ExtrapToZ(param, z);
55   xy[0] = param->GetNonBendingCoor();
56   xy[1] = param->GetBendingCoor();
57
58   delete param;
59
60   return kTRUE;
61
62 }
63
64 //====================================================================================================================================================
65
66 Bool_t AliMFTAnalysisTools::ExtrapAODMuonToZ(AliAODTrack *muon, Double_t z, Double_t xy[2], TLorentzVector &kinem) {
67
68   if (!(muon->Pz()!=0)) return kFALSE;
69
70   AliMUONTrackParam *param = new AliMUONTrackParam();
71
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))) );
78
79   AliMUONTrackExtrap::ExtrapToZ(param, z);
80   xy[0] = param->GetNonBendingCoor();
81   xy[1] = param->GetBendingCoor();
82
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());
85   
86   kinem.SetPxPyPzE(param->Px(), param->Py(), param->Pz(), energy);
87
88   delete param;
89
90   return kTRUE;
91
92 }
93
94 //====================================================================================================================================================
95
96 Bool_t AliMFTAnalysisTools::ExtrapAODMuonToZ(AliAODTrack *muon, Double_t z, Double_t xy[2], TLorentzVector &kinem, TMatrixD &cov) {
97
98   if (!(muon->Pz()!=0)) return kFALSE;
99
100   AliMUONTrackParam *param = new AliMUONTrackParam();
101
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))) );
108
109   param -> SetCovariances(ConvertCovMatrixAOD2MUON(muon));
110
111   AliMUONTrackExtrap::ExtrapToZCov(param, z);
112   xy[0] = param->GetNonBendingCoor();
113   xy[1] = param->GetBendingCoor();
114
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());
117   
118   kinem.SetPxPyPzE(param->Px(), param->Py(), param->Pz(), energy);
119
120   cov = param->GetCovariances();
121
122   delete param;
123
124   return kTRUE;
125
126 }
127
128 //====================================================================================================================================================
129
130 Double_t AliMFTAnalysisTools::GetAODMuonOffset(AliAODTrack *muon, Double_t xv, Double_t yv, Double_t zv) {
131
132   Double_t xy[2] = {0};
133   ExtrapAODMuonToZ(muon, zv, xy);
134   
135   return TMath::Sqrt((xv-xy[0])*(xv-xy[0]) + (yv-xy[1])*(yv-xy[1]));
136
137 }
138
139 //====================================================================================================================================================
140
141 Double_t AliMFTAnalysisTools::GetAODMuonWeightedOffset(AliAODTrack *muon, Double_t xv, Double_t yv, Double_t zv) {
142
143   Double_t xy[2] = {0};
144   TLorentzVector kinem(0,0,0,0);
145   TMatrixD cov(5,5);
146
147   ExtrapAODMuonToZ(muon, zv, xy, kinem, cov);
148
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);
154   
155   TMatrixD covCoordinatesInverse = covCoordinates.Invert();
156
157   Double_t dX = xy[0] - xv;
158   Double_t dY = xy[1] - yv;
159   
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)));
163
164   return weightedOffset;
165
166 }
167
168 //====================================================================================================================================================
169
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) {
173   
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 
176
177   if (ptDimu != 0) {
178     Double_t decayLengthXY = TMath::Sqrt((xVtx-xDimu)*(xVtx-xDimu)+(yVtx-yDimu)*(yVtx-yDimu));
179     return (decayLengthXY * mDimu/ptDimu)/TMath::Ccgs()*1E12;   // in ps
180   }
181   
182   return -99999999;
183   
184 }
185
186 //====================================================================================================================================================
187
188 Double_t AliMFTAnalysisTools::GetPseudoProperDecayTimeZ(Double_t zVtx, 
189                                                         Double_t zDimu, 
190                                                         Double_t mDimu, Double_t pzDimu) {
191
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 
194
195   if (pzDimu != 0) {
196     Double_t decayLengthZ = zDimu - zVtx;
197     return (decayLengthZ * mDimu/pzDimu)/TMath::Ccgs()*1E12;  // in ps
198   }
199
200   return -99999999;
201
202 }
203
204 //====================================================================================================================================================
205
206 Bool_t AliMFTAnalysisTools::CalculatePCA(AliAODDimuon *dimuon, Double_t *pca, Double_t &pcaQuality, TLorentzVector &kinem) {
207
208   TObjArray *muons = new TObjArray();
209   muons -> Add(dimuon->GetMu(0));
210   muons -> Add(dimuon->GetMu(1));
211
212   return CalculatePCA(muons, pca, pcaQuality, kinem);
213
214 }
215
216 //====================================================================================================================================================
217
218 Bool_t AliMFTAnalysisTools::CalculatePCA(TObjArray *muons, Double_t *pca, Double_t &pcaQuality, TLorentzVector &kinem) {
219   
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");
223     return kFALSE;
224   }
225   
226   Double_t fXPointOfClosestApproach=0, fYPointOfClosestApproach=0, fZPointOfClosestApproach=0;
227   
228   AliAODTrack *muon[nMuons];
229   AliMUONTrackParam *param[nMuons];
230
231   // Finding AliMUONTrackParam objects for each muon
232   
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))) );
243   }
244
245   // here we want to understand in which direction we have to search the minimum...
246   
247   Double_t step = 1.;  // initial step, in cm
248   Double_t startPoint = 0.;
249
250   Double_t r[3]={0}, z[3]={startPoint, startPoint+step, startPoint+2*step};
251   
252   TVector3 **points = new TVector3*[nMuons];
253   
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]);
258     }
259     r[i] = GetDistanceBetweenPoints(points,nMuons);
260   }
261
262   Int_t researchDirection = 0;
263   
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");
268     return kFALSE;
269   }
270
271   while (TMath::Abs(researchDirection)>0.5) {
272
273     if (researchDirection>0) {
274       z[0] = z[1];
275       z[1] = z[2];
276       z[2] = z[1]+researchDirection*step;
277     }
278     else {
279       z[2] = z[1];
280       z[1] = z[0];
281       z[0] = z[1]+researchDirection*step;
282     }
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");
285       return kFALSE;
286     }
287     
288     for (Int_t iMu=0; iMu<nMuons; iMu++) delete points[iMu];
289     
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]);
294       }      
295       r[i] = GetDistanceBetweenPoints(points,nMuons);
296     }
297     researchDirection=0;
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
300
301   }
302
303   // now we now that the minimum is between z[0] and z[2] and we search for it
304   
305   step *= 0.5;
306   while (step>AliMFTConstants::fPrecisionPointOfClosestApproach) {
307     z[0] = z[1]-step;
308     z[2] = z[1]+step;
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]);
313       }      
314       r[i] = GetDistanceBetweenPoints(points,nMuons);
315     }
316     if      (r[0]<r[1]) z[1] = z[0];
317     else if (r[2]<r[1]) z[1] = z[2];
318     else step *= 0.5;
319   }
320   
321   delete [] points;
322
323   // Once z of minimum is found, we evaluate the x and y coordinates by averaging over the contributing tracks
324   
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();
332   }
333   fXPointOfClosestApproach /= Double_t(nMuons);
334   fYPointOfClosestApproach /= Double_t(nMuons);
335   
336   pca[0] = fXPointOfClosestApproach;
337   pca[1] = fYPointOfClosestApproach;
338   pca[2] = fZPointOfClosestApproach;
339   
340   // Evaluating the kinematics of the N-muon
341
342   Double_t pTot[3] = {0};
343   Double_t ene = 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());
350   }
351   
352   kinem.SetPxPyPzE(pTot[0], pTot[1], pTot[2], ene);
353
354   // Evaluating the PCA quality of the N-muon
355
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);
360     sum += f;
361     squareSum += f*f;
362   }
363   if (sum > 0.) pcaQuality =  (sum-squareSum/sum) / (nMuons-1);
364   else pcaQuality = 0.;
365
366   return kTRUE;
367   
368 }
369
370 //=========================================================================================================================
371
372 Double_t AliMFTAnalysisTools::GetDistanceBetweenPoints(TVector3 **points, Int_t nPoints) {
373   
374   if (nPoints>AliMFTConstants::fNMaxMuonsForPCA) {
375     printf("W-AliMFTAnalysisTools::GetDistanceBetweenPoints: number of points not valid\n");
376     return 1.e9;
377   }
378   
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()) );
383
384   const Int_t nEdgesMax = ((AliMFTConstants::fNMaxMuonsForPCA) * (AliMFTConstants::fNMaxMuonsForPCA - 1)) / 2;
385   
386   Int_t startID[nEdgesMax]       = {0};
387   Int_t stopID[nEdgesMax]        = {0};
388   Double_t edgeLength[nEdgesMax] = {0};
389
390   Bool_t pointStatus[AliMFTConstants::fNMaxMuonsForPCA] = {0};
391   
392   Int_t nEdges=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()) );
398       stopID[nEdges]  = i;
399       startID[nEdges] = j;
400       nEdges++;
401     }
402   }
403  
404   // Order Edges
405
406   Double_t min = 0;
407   Int_t   iMin = 0;
408
409   for (Int_t iEdge=0; iEdge<nEdges-1; iEdge++) {
410     min  = edgeLength[iEdge];
411     iMin = iEdge;
412     for (Int_t j=iEdge+1; j<nEdges; j++) {
413       if (edgeLength[j]<min) {
414         min  = edgeLength[j];
415         iMin = j;
416       }
417     }
418     
419     if (iMin != iEdge) {
420
421       Double_t edgeLengthMin = edgeLength[iMin];
422       Int_t startIDmin = startID[iMin];
423       Int_t stopIDmin  = stopID[iMin];
424       
425       edgeLength[iMin] = edgeLength[iEdge];
426       startID[iMin]    = startID[iEdge];
427       stopID[iMin]     = stopID[iEdge];
428
429       edgeLength[iEdge] = edgeLengthMin;
430       startID[iEdge]    = startIDmin;
431       stopID[iEdge]     = stopIDmin;
432
433     }
434     
435   }
436   
437   // Connect
438
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];
445     }
446   }
447   
448   return length;
449   
450 }
451
452 //====================================================================================================================================================
453
454 void AliMFTAnalysisTools::ConvertCovMatrixMUON2AOD(const TMatrixD& covMUON, Double_t covAOD[21]) {
455
456   // Converts the cov matrix from the MUON format (TMatrixD) to the AOD one (Double_t[21])
457   // 
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]
464
465   covAOD[0]  = covMUON(0,0);
466
467   covAOD[1]  = covMUON(1,0);
468   covAOD[2]  = covMUON(1,1);
469
470   covAOD[3]  = covMUON(2,0);  
471   covAOD[4]  = covMUON(2,1);  
472   covAOD[5]  = covMUON(2,2);  
473
474   covAOD[6]  = covMUON(3,0);  
475   covAOD[7]  = covMUON(3,1);  
476   covAOD[8]  = covMUON(3,2);  
477   covAOD[9]  = covMUON(3,3);  
478
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);  
484
485   covAOD[15] = 0;  
486   covAOD[16] = 0;  
487   covAOD[17] = 0;  
488   covAOD[18] = 0;  
489   covAOD[19] = 0;  
490   covAOD[20] = 0;  
491
492 }
493
494 //====================================================================================================================================================
495
496 const TMatrixD AliMFTAnalysisTools::ConvertCovMatrixAOD2MUON(AliAODTrack *muon) {
497
498   Double_t covAOD[21] = {0};
499   muon -> GetCovarianceXYZPxPyPz(covAOD);
500
501   TMatrixD covMUON(5,5);
502
503   covMUON(0,0) = covAOD[0];
504                                 
505   covMUON(1,0) = covAOD[1];
506   covMUON(1,1) = covAOD[2];
507                                 
508   covMUON(2,0) = covAOD[3];  
509   covMUON(2,1) = covAOD[4];  
510   covMUON(2,2) = covAOD[5];  
511                                 
512   covMUON(3,0) = covAOD[6];  
513   covMUON(3,1) = covAOD[7];  
514   covMUON(3,2) = covAOD[8];  
515   covMUON(3,3) = covAOD[9];  
516                                 
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]; 
522
523   return covMUON;
524
525 }
526
527 //====================================================================================================================================================
528