]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/trigger/AliHLTD0toKpi.cxx
flat friends update
[u/mrichter/AliRoot.git] / HLT / trigger / AliHLTD0toKpi.cxx
1 // $Id$
2 #include "AliHLTD0toKpi.h"
3 #include "TDatabasePDG.h"
4 #include "TMath.h"
5 #include "AliESDtrack.h"
6 #include "TVector3.h"
7 #include "AliAODVertex.h"
8 #include "AliESDVertex.h"
9 #include "TObjArray.h"
10 #include "AliVertexerTracks.h"
11 #include "AliHLTGlobalBarrelTrack.h"
12 #include "AliExternalTrackParam.h"
13 #include "AliKFVertex.h"
14
15 ClassImp(AliHLTD0toKpi)
16
17 AliHLTD0toKpi::AliHLTD0toKpi() 
18 {
19 }
20
21 Double_t AliHLTD0toKpi::InvMass(const AliExternalTrackParam* d1, const AliExternalTrackParam* d2)
22 {
23   Double_t mpi=TDatabasePDG::Instance()->GetParticle(211)->Mass();
24   Double_t mK=TDatabasePDG::Instance()->GetParticle(321)->Mass();
25
26   Double_t energy[2]; 
27   energy[1] = TMath::Sqrt(mK*mK+d1->GetP()*d1->GetP());
28   energy[0] = TMath::Sqrt(mpi*mpi+d2->GetP()*d2->GetP());
29
30   Double_t p1[3],p2[3];
31   d1->GetPxPyPz(p1);
32   d2->GetPxPyPz(p2);
33   
34   Double_t momTot2 = (p1[0]+p2[0])*(p1[0]+p2[0])+
35                      (p1[1]+p2[1])*(p1[1]+p2[1])+
36                      (p1[2]+p2[2])*(p1[2]+p2[2]);
37
38   return TMath::Sqrt((energy[0]+energy[1])*(energy[0]+energy[1])-momTot2);
39
40 }
41 void AliHLTD0toKpi::CosThetaStar(const AliExternalTrackParam* d1, const AliExternalTrackParam* d2,Double_t &D0,Double_t &D0bar)
42 {
43   Double_t mD0 = TDatabasePDG::Instance()->GetParticle(421)->Mass();
44   Double_t mpi=TDatabasePDG::Instance()->GetParticle(211)->Mass();
45   Double_t mK=TDatabasePDG::Instance()->GetParticle(321)->Mass();
46
47   Double_t pStar = TMath::Sqrt(TMath::Power(mD0*mD0-mK*mK-mpi*mpi,2.)-4.*mK*mK*mpi*mpi)/(2.*mD0);
48  
49   Double_t px = d1->Px()+d2->Px();
50   Double_t py = d1->Py()+d2->Py();
51   Double_t pz = d1->Pz()+d2->Pz();
52   Double_t p = TMath::Sqrt(px*px+py*py+pz*pz);
53   Double_t energy = TMath::Sqrt(p*p+mD0*mD0);
54
55   Double_t beta = p/energy;
56   Double_t gamma = energy/mD0;
57   
58   Double_t qL;
59   TVector3 mom(d1->Px(),d1->Py(),d1->Pz());
60   TVector3 momD(px,py,pz);
61   qL = mom.Dot(momD)/momD.Mag();
62
63   D0 = (qL/gamma-beta*TMath::Sqrt(pStar*pStar+mK*mK))/pStar;
64   
65   TVector3 mom2(d2->Px(),d2->Py(),d2->Pz());
66   TVector3 momD2(px,py,pz);
67   qL = mom2.Dot(momD2)/momD2.Mag();
68
69   D0bar = (qL/gamma-beta*TMath::Sqrt(pStar*pStar+mK*mK))/pStar;
70
71 }
72 Double_t AliHLTD0toKpi::PointingAngle(const AliExternalTrackParam* n, const AliExternalTrackParam* p, const Double_t *pv, const Double_t *sv)
73 {
74
75   TVector3 mom(n->Px()+p->Px(),n->Py()+p->Py(),n->Pz()+p->Pz());
76   TVector3 flight(sv[0]-pv[0],sv[1]-pv[1],sv[2]-pv[2]);
77   
78   double pta = mom.Angle(flight);
79
80   return TMath::Cos(pta); 
81 }
82
83 AliAODVertex* AliHLTD0toKpi::ReconstructSecondaryVertex(TObjArray *trkArray, Double_t b, const AliESDVertex *v, bool useKF)
84 {
85   
86   AliESDVertex *vertexESD = 0;
87   AliAODVertex *vertexAOD = 0;
88   
89   if(!useKF){
90     AliVertexerTracks *vertexer = new AliVertexerTracks(b);
91     AliESDVertex* vertex = new AliESDVertex(*((AliESDVertex*)v));
92     //AliESDVertex* vertex =  const_cast<AliESDVertex*>(v);
93     vertexer->SetVtxStart(vertex);
94     vertexESD = (AliESDVertex*)vertexer->VertexForSelectedESDTracks(trkArray);
95     delete vertexer; vertexer=NULL;
96     delete vertex;
97
98     if(!vertexESD) return vertexAOD;
99     
100     if(vertexESD->GetNContributors()!=trkArray->GetEntriesFast()) { 
101       //AliDebug(2,"vertexing failed"); 
102       delete vertexESD; vertexESD=NULL;
103       return vertexAOD;
104     }
105   }
106   else{
107     AliKFParticle::SetField(b);
108     
109     AliKFVertex vertexKF;
110     
111     Int_t nTrks = trkArray->GetEntriesFast();
112     for(Int_t i=0; i<nTrks; i++) {
113       AliESDtrack *esdTrack = (AliESDtrack*)trkArray->At(i);
114       AliKFParticle daughterKF(*esdTrack,211);
115       vertexKF.AddDaughter(daughterKF);
116     }
117     vertexESD = new AliESDVertex(vertexKF.Parameters(),
118                                  vertexKF.CovarianceMatrix(),
119                                  vertexKF.GetChi2(),
120                                  vertexKF.GetNContributors());
121   }
122   // convert to AliAODVertex
123   Double_t pos[3],cov[6],chi2perNDF;
124   vertexESD->GetXYZ(pos); // position
125   vertexESD->GetCovMatrix(cov); //covariance matrix
126   chi2perNDF = vertexESD->GetChi2toNDF();
127   //dispersion = vertexESD->GetDispersion();
128   delete vertexESD; vertexESD=NULL;
129
130   Int_t nprongs= trkArray->GetEntriesFast();
131   vertexAOD = new AliAODVertex(pos,cov,chi2perNDF,0x0,-1,AliAODVertex::kUndef,nprongs);
132
133   return vertexAOD;
134
135 }
136 Double_t AliHLTD0toKpi::Pt(const AliExternalTrackParam* d1, const AliExternalTrackParam* d2)
137 {
138   Double_t p1[3],p2[3];
139   d1->GetPxPyPz(p1);
140   d2->GetPxPyPz(p2);
141   
142   Double_t pt2 = (p1[0]+p2[0])*(p1[0]+p2[0]) + (p1[1]+p2[1])*(p1[1]+p2[1]);
143
144   return TMath::Sqrt(pt2);
145 }