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#include "LmvmCand.h"
10#include "TLorentzVector.h"
11#include "TMath.h"
12
13#define M2E 2.6112004954086e-7
14
16 public:
17 Double_t fMomentumMag = 0.; // Absolute value of momentum
18 Double_t fPt = 0.; // Transverse momentum
19 Double_t fRapidity = 0.; // Rapidity
20 Double_t fMinv = 0.; // Invariant mass
21 Double_t fAngle = 0.; // Opening angle
22
23 /*
24 * Calculate kinematic parameters for MC tracks.
25 */
26 static LmvmKinePar Create(const CbmMCTrack* mctrackP, const CbmMCTrack* mctrackM)
27 {
28 LmvmKinePar params;
29 if (mctrackP == nullptr || mctrackM == nullptr) return params;
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 LMVM candidates.
60 */
61 static LmvmKinePar Create(const LmvmCand* candP, const LmvmCand* candM)
62 {
63 LmvmKinePar params;
64 if (candP == nullptr || candM == nullptr) return params;
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
Definition LmvmKinePar.h:13
void GetMomentum(TVector3 &momentum) const
Definition CbmMCTrack.h:169
double fEnergy
Definition LmvmCand.h:61
TVector3 fMomentum
Definition LmvmCand.h:57
static LmvmKinePar Create(const LmvmCand *candP, const LmvmCand *candM)
Definition LmvmKinePar.h:61
Double_t fAngle
Definition LmvmKinePar.h:21
Double_t fMinv
Definition LmvmKinePar.h:20
Double_t fPt
Definition LmvmKinePar.h:18
Double_t fMomentumMag
Definition LmvmKinePar.h:17
Double_t fRapidity
Definition LmvmKinePar.h:19
static LmvmKinePar Create(const CbmMCTrack *mctrackP, const CbmMCTrack *mctrackM)
Definition LmvmKinePar.h:26