]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - MUON/MUONmassPlot_ESD.C
Calculation of pyz corrected
[u/mrichter/AliRoot.git] / MUON / MUONmassPlot_ESD.C
index 454ce595e3acd71e8547123340a04c8acadba70f..f5cd411396ad25db30eb9d06f76b2621dba02bf8 100644 (file)
@@ -181,7 +181,7 @@ TH1F *hInvMassRes;
       thetaY = muonTrack->GetThetaY();
 
       pYZ     =  1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
-      fPzRec1  = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaX));
+      fPzRec1  = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaY));
       fPxRec1  = fPzRec1 * TMath::Tan(thetaX);
       fPyRec1  = fPzRec1 * TMath::Tan(thetaY);
       fCharge = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
@@ -231,7 +231,7 @@ TH1F *hInvMassRes;
          thetaY = muonTrack->GetThetaY();
 
          pYZ     =  1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
-         fPzRec2  = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaX));
+         fPzRec2  = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaY));
          fPxRec2  = fPzRec2 * TMath::Tan(thetaX);
          fPyRec2  = fPzRec2 * TMath::Tan(thetaY);
          fCharge2 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));