]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - PWGLF/RESONANCES/AliRsnMother.cxx
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / AliRsnMother.cxx
index 7be7ab05acaa1378c8413d8abefd1945d962dce9..09e45bc6ee186993cb4f1add5d4621763680e080 100644 (file)
@@ -44,7 +44,8 @@ AliRsnMother::AliRsnMother(const AliRsnMother &obj) :
    fSum(obj.fSum),
    fSumMC(obj.fSumMC),
    fRef(obj.fRef),
-   fRefMC(obj.fRefMC)
+   fRefMC(obj.fRefMC),
+   fDCAproduct(obj.fDCAproduct)
 {
 //
 // Copy constructor.
@@ -72,6 +73,7 @@ AliRsnMother &AliRsnMother::operator=(const AliRsnMother &obj)
    fRefEvent = obj.fRefEvent;
    fDaughter[0] = obj.fDaughter[0];
    fDaughter[1] = obj.fDaughter[1];
+   fDCAproduct = obj.fDCAproduct;
 
    return (*this);
 }
@@ -196,19 +198,19 @@ Double_t AliRsnMother::CosThetaStar(Bool_t first, Bool_t useMC)
 
    TLorentzVector &mother    = Sum(useMC);
    TLorentzVector &daughter0 = (first ? fDaughter[0]->P(useMC) : fDaughter[1]->P(useMC));
-   TLorentzVector &daughter1 = (first ? fDaughter[1]->P(useMC) : fDaughter[0]->P(useMC));
+//    TLorentzVector &daughter1 = (first ? fDaughter[1]->P(useMC) : fDaughter[0]->P(useMC));
    TVector3 momentumM(mother.Vect());
    TVector3 normal(mother.Y() / momentumM.Mag(), -mother.X() / momentumM.Mag(), 0.0);
 
-   // Computes first the invariant mass of the mother
-   Double_t mass0      = fDaughter[0]->P(useMC).M();
-   Double_t mass1      = fDaughter[1]->P(useMC).M();
-   Double_t p0         = daughter0.Vect().Mag();
-   Double_t p1         = daughter1.Vect().Mag();
-   Double_t E0         = TMath::Sqrt(mass0 * mass0 + p0 * p0);
-   Double_t E1         = TMath::Sqrt(mass1 * mass1 + p1 * p1);
-   Double_t MotherMass = TMath::Sqrt((E0 + E1) * (E0 + E1) - (p0 * p0 + 2.0 * daughter0.Vect().Dot(daughter1.Vect()) + p1 * p1));
-   MotherMass = fSum.M();
+//    // Computes first the invariant mass of the mother
+//    Double_t mass0      = fDaughter[0]->P(useMC).M();
+//    Double_t mass1      = fDaughter[1]->P(useMC).M();
+//    Double_t p0         = daughter0.Vect().Mag();
+//    Double_t p1         = daughter1.Vect().Mag();
+//    Double_t E0         = TMath::Sqrt(mass0 * mass0 + p0 * p0);
+//    Double_t E1         = TMath::Sqrt(mass1 * mass1 + p1 * p1);
+//    Double_t MotherMass = TMath::Sqrt((E0 + E1) * (E0 + E1) - (p0 * p0 + 2.0 * daughter0.Vect().Dot(daughter1.Vect()) + p1 * p1));
+//    MotherMass = fSum.M();
 
    // Computes components of beta
    Double_t betaX = -mother.X() / mother.E();
@@ -225,6 +227,50 @@ Double_t AliRsnMother::CosThetaStar(Bool_t first, Bool_t useMC)
    return cosThetaStar;
 }
 
+//__________________________________________________________________________________________________
+Double_t AliRsnMother::DCAproduct()
+{
+  //
+  // returns product of DCA of the two daughters
+  //
+  AliRsnEvent *event1 = fDaughter[0]->GetOwnerEvent();
+  AliRsnEvent *event2 = fDaughter[1]->GetOwnerEvent();
+  
+  if(event1 != event2){  
+    AliError("Attempting to build pair with tracks coming from different events");
+    return 0.0;
+  }
+   
+  if (event1->IsAOD()) {
+    AliAODEvent *aodEvent = (AliAODEvent*)event1->GetRefAOD();  
+    if (!aodEvent) return 0.0;
+    AliAODTrack *track1 = (AliAODTrack*)fDaughter[0]->Ref2AODtrack();
+    AliAODTrack *track2 = (AliAODTrack*)fDaughter[1]->Ref2AODtrack();
+    AliVVertex *vertex = aodEvent->GetPrimaryVertex();   
+    if (!vertex || !track1 || !track2) return 0.0;
+     
+    Double_t b1[2], cov1[3], b2[2], cov2[3];
+    if ( track1->PropagateToDCA(vertex, aodEvent->GetMagneticField(), kVeryBig, b1, cov1) &&   
+        track2->PropagateToDCA(vertex, aodEvent->GetMagneticField(), kVeryBig, b2, cov2) )
+      fDCAproduct = b1[0]*b2[0]; 
+    else fDCAproduct = 0.0;
+  } else {
+    AliESDEvent *esdEvent = (AliESDEvent*)event1->GetRefESD();  
+    if (!esdEvent) return 0.0;
+    AliESDtrack *track1 = (AliESDtrack*)fDaughter[0]->Ref2ESDtrack();
+    AliESDtrack *track2 = (AliESDtrack*)fDaughter[1]->Ref2ESDtrack();
+    const AliVVertex *vertex = esdEvent->GetPrimaryVertex();
+    if (!vertex || !track1 || !track2) return 0.0;
+     
+    Double_t b1[2], cov1[3], b2[2], cov2[3];
+    if ( track1->PropagateToDCA(vertex, esdEvent->GetMagneticField(), kVeryBig, b1, cov1) &&    
+        track2->PropagateToDCA(vertex, esdEvent->GetMagneticField(), kVeryBig, b2, cov2) )
+      fDCAproduct = b1[0]*b2[0];  
+    else fDCAproduct = 0.0;  
+  }
+  return fDCAproduct;
+}
+
 //__________________________________________________________________________________________________
 void AliRsnMother::PrintInfo(const Option_t * /*option*/) const
 {