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