CbmRoot
Loading...
Searching...
No Matches
CbmAnaJpsiKinematicParams.h
Go to the documentation of this file.
1/* Copyright (C) 2015 UGiessen, JINR-LIT
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Semen Lebedev [committer] */
4
5#ifndef CBM_ANAJPSI_KINEMATIC_PARAMS_H
6#define CBM_ANAJPSI_KINEMATIC_PARAMS_H
7
9#include "CbmMCTrack.h"
10
11#include "TLorentzVector.h"
12#include "TMath.h"
13
14#define M2E 2.6112004954086e-7
15
17public:
18 Double_t fMomentumMag; // Absolute value of momentum
19 Double_t fPt; // Transverse momentum
20 Double_t fRapidity; // Rapidity
21 Double_t fMinv; // Invariant mass
22 Double_t fAngle; // Opening angle
23
24 /*
25 * Calculate kinematic parameters for MC tracks.
26 */
28 {
30
31 TVector3 momP; //momentum e+
32 mctrackP->GetMomentum(momP);
33 Double_t energyP = TMath::Sqrt(momP.Mag2() + M2E);
34 TLorentzVector lorVecP(momP, energyP);
35
36 TVector3 momM; //momentum e-
37 mctrackM->GetMomentum(momM);
38 Double_t energyM = TMath::Sqrt(momM.Mag2() + M2E);
39 TLorentzVector lorVecM(momM, energyM);
40
41 TVector3 momPair = momP + momM;
42 Double_t energyPair = energyP + energyM;
43 Double_t ptPair = momPair.Perp();
44 Double_t pzPair = momPair.Pz();
45 Double_t yPair = 0.5 * TMath::Log((energyPair + pzPair) / (energyPair - pzPair));
46 Double_t anglePair = lorVecM.Angle(lorVecP.Vect());
47 Double_t theta = 180. * anglePair / TMath::Pi();
48 Double_t minv = 2. * TMath::Sin(anglePair / 2.) * TMath::Sqrt(momM.Mag() * momP.Mag());
49
50 params.fMomentumMag = momPair.Mag();
51 params.fPt = ptPair;
52 params.fRapidity = yPair;
53 params.fMinv = minv;
54 params.fAngle = theta;
55 return params;
56 }
57
58 /*
59 * Calculate kinematic parameters for JPSI candidates.
60 */
62 const CbmAnaJpsiCandidate* candM)
63 {
65
66 TLorentzVector lorVecP(candP->fMomentum, candP->fEnergy);
67 TLorentzVector lorVecM(candM->fMomentum, candM->fEnergy);
68
69 TVector3 momPair = candP->fMomentum + candM->fMomentum;
70 Double_t energyPair = candP->fEnergy + candM->fEnergy;
71 Double_t ptPair = momPair.Perp();
72 Double_t pzPair = momPair.Pz();
73 Double_t yPair = 0.5 * TMath::Log((energyPair + pzPair) / (energyPair - pzPair));
74 Double_t anglePair = lorVecM.Angle(lorVecP.Vect());
75 Double_t theta = 180. * anglePair / TMath::Pi();
76 Double_t minv = 2. * TMath::Sin(anglePair / 2.) * TMath::Sqrt(candM->fMomentum.Mag() * candP->fMomentum.Mag());
77
78 params.fMomentumMag = momPair.Mag();
79 params.fPt = ptPair;
80 params.fRapidity = yPair;
81 params.fMinv = minv;
82 params.fAngle = theta;
83 return params;
84 }
85};
86
87#endif
#define M2E
static CbmAnaJpsiKinematicParams KinematicParamsWithCandidates(const CbmAnaJpsiCandidate *candP, const CbmAnaJpsiCandidate *candM)
static CbmAnaJpsiKinematicParams KinematicParamsWithMcTracks(const CbmMCTrack *mctrackP, const CbmMCTrack *mctrackM)
void GetMomentum(TVector3 &momentum) const
Definition CbmMCTrack.h:170