Belle II Software development
EnergyLossEstimator Class Reference

Helper struct to provide consistent energy loss estimation throughout the CDC track finding. More...

#include <EnergyLossEstimator.h>

Public Member Functions

 EnergyLossEstimator (double eDensity, double I, double bZ=NAN)
 Constructor from the material properties.
 
double getBetheStoppingPower (double p, int pdgCode) const
 Stopping power aka energy loss / arc length.
 
double getEnergyLoss (double p, int pdgCode, double arcLength) const
 Calculates the total energy loss after travelling the given distance.
 
double getMomentumLossFactor (double p, int pdgCode, double arcLength) const
 Calculates a factor applicable scaling the current momentum to the momentum after traveling given arc length.
 
double getLossDist2D (double pt, int pdgCode, double arcLength2D) const
 Calculates a correction term for the two dimensional distance undoing the energy loss after the given arc length.
 

Static Public Member Functions

static EnergyLossEstimator forCDC ()
 Create an energy loss estimator with the material properties of the CDC.
 
static double getMass (int pdgCode)
 Lookup the mass for the given pdg code.
 
static int getCharge (int pdgCode)
 Lookup the charge for the given pdg code.
 

Private Attributes

double m_eDensity
 Electron density in mol / cm^3.
 
double m_I
 Mean excitation energy in GeV.
 
double m_bZ
 B field to be used for the distance translation.
 

Detailed Description

Helper struct to provide consistent energy loss estimation throughout the CDC track finding.

Definition at line 22 of file EnergyLossEstimator.h.

Constructor & Destructor Documentation

◆ EnergyLossEstimator()

EnergyLossEstimator ( double  eDensity,
double  I,
double  bZ = NAN 
)

Constructor from the material properties.

Parameters
eDensityElectron density in mol / cm ^ 3, e.g. from Z * density / A
IExcitation energy (must be given in standard GeV units)
bZMean magentic field in the z direction

Definition at line 121 of file EnergyLossEstimator.cc.

122 : m_eDensity(eDensity)
123 , m_I(I)
124 , m_bZ(bZ)
125{
126}
double m_eDensity
Electron density in mol / cm^3.
double m_bZ
B field to be used for the distance translation.
double m_I
Mean excitation energy in GeV.

Member Function Documentation

◆ forCDC()

EnergyLossEstimator forCDC ( )
static

Create an energy loss estimator with the material properties of the CDC.

Definition at line 28 of file EnergyLossEstimator.cc.

29{
30
31 // Look up material properties from the Belle2 geometry
32 G4Navigator g4Nav;
33 G4VPhysicalVolume* g4World = geometry::GeometryManager::getInstance().getTopVolume();
34 g4Nav.SetWorldVolume(g4World);
35
36 // Choose position well inside the CDC
37 double posX = 0;
38 double posY = 50;
39 double posZ = 0;
40
41 G4ThreeVector g4Pos(posX * CLHEP::cm, posY * CLHEP::cm, posZ * CLHEP::cm);
42 const G4VPhysicalVolume* g4Volume = g4Nav.LocateGlobalPointAndSetup(g4Pos);
43 const G4Material* g4Mat = g4Volume->GetLogicalVolume()->GetMaterial();
44
45 double A = 0;
46 double Z = 0;
47
48 if (g4Mat->GetNumberOfElements() == 1) {
49 A = g4Mat->GetA() * CLHEP::mole / CLHEP::g;
50 Z = g4Mat->GetZ();
51 } else {
52 // Calculate weight-averaged A, Z
53 for (unsigned i = 0; i < g4Mat->GetNumberOfElements(); ++i) {
54 const G4Element* element = (*g4Mat->GetElementVector())[i];
55 const double elementA = element->GetA() * CLHEP::mole / CLHEP::g;
56 const double elementZ = element->GetZ();
57 const double frac = g4Mat->GetFractionVector()[i];
58 B2RESULT("Part " << i << " Z=" << elementZ << " A=" << elementA << " (" << frac << ")");
59 Z += elementZ * frac;
60 A += elementA * frac;
61 }
62 }
63
64 // Make sure to translate to Belle units from the Geant4 / CLHEP units
65 // Z has correct units (is unitless)
66 const double density = g4Mat->GetDensity() / CLHEP::g * CLHEP::cm3;
67 const double radiationLength = g4Mat->GetRadlen() / CLHEP::cm;
68 const double mEE = g4Mat->GetIonisation()->GetMeanExcitationEnergy() / CLHEP::GeV;
69
70 B2RESULT("Received Z " << Z);
71 B2RESULT("Received A " << A);
72 B2RESULT("Received density " << density);
73 B2RESULT("Received radiation length " << radiationLength);
74 B2RESULT("Received mean excitation energy " << mEE);
75
76 const double eDensity = Z * density / A;
77 B2RESULT("Received electron density " << eDensity);
78
79 const double bZ = CDCBFieldUtil::getBFieldZ();
80
81 return EnergyLossEstimator(eDensity, mEE, bZ);
82}
static double getBFieldZ()
Getter for the signed magnetic field strength in z direction at the origin ( in Tesla )
Helper struct to provide consistent energy loss estimation throughout the CDC track finding.
G4VPhysicalVolume * getTopVolume()
Return a pointer to the top volume.
static GeometryManager & getInstance()
Return a reference to the instance.

◆ getBetheStoppingPower()

double getBetheStoppingPower ( double  p,
int  pdgCode 
) const

Stopping power aka energy loss / arc length.

Parameters
pabsolute momentum
pdgCode

Definition at line 128 of file EnergyLossEstimator.cc.

129{
130 static const double eMass = getMass(11);
131
132 const double M = getMass(pdgCode);
133
134 const double E = std::sqrt(p * p + M * M);
135 const double gamma = E / M;
136 const double beta = p / E;
137
138 const double beta2 = beta * beta;
139 const double gamma2 = gamma * gamma;
140
141 const double Wmax = 2 * eMass * beta2 * gamma2 / (1 + 2 * gamma * eMass / M);
142 const double I2 = m_I * m_I;
143
144 static const double K = 0.307075 * Unit::MeV * Unit::cm2;
145 const double dEdx = K * m_eDensity / beta2 *
146 (1.0 / 2.0 * std::log(2 * eMass * beta2 * gamma2 * Wmax / I2) - beta2);
147 return dEdx;
148}
R E
internal precision of FFTW codelets
#define K(x)
macro autogenerated by FFTW
static double getMass(int pdgCode)
Lookup the mass for the given pdg code.
static const double MeV
[megaelectronvolt]
Definition: Unit.h:114
static const double cm2
[square centimeters]
Definition: Unit.h:78

◆ getCharge()

int getCharge ( int  pdgCode)
static

Lookup the charge for the given pdg code.

Definition at line 105 of file EnergyLossEstimator.cc.

106{
107 if (std::abs(pdgCode) == Const::electron.getPDGCode()) {
108 return -sign(pdgCode);
109 } else if (std::abs(pdgCode) == Const::muon.getPDGCode()) {
110 return -sign(pdgCode);
111 } else if (std::abs(pdgCode) == Const::kaon.getPDGCode()) {
112 return sign(pdgCode);
113 } else if (std::abs(pdgCode) == Const::pion.getPDGCode()) {
114 return sign(pdgCode);
115 } else if (std::abs(pdgCode) == Const::proton.getPDGCode()) {
116 return sign(pdgCode);
117 }
118 return 0;
119}
static const ChargedStable muon
muon particle
Definition: Const.h:660
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
static const ChargedStable proton
proton particle
Definition: Const.h:663
static const ChargedStable kaon
charged kaon particle
Definition: Const.h:662
static const ChargedStable electron
electron particle
Definition: Const.h:659

◆ getEnergyLoss()

double getEnergyLoss ( double  p,
int  pdgCode,
double  arcLength 
) const

Calculates the total energy loss after travelling the given distance.

Definition at line 150 of file EnergyLossEstimator.cc.

151{
152 const double dEdx = getBetheStoppingPower(p, pdgCode);
153 const double eLoss = arcLength * dEdx;
154 return eLoss;
155}
double getBetheStoppingPower(double p, int pdgCode) const
Stopping power aka energy loss / arc length.

◆ getLossDist2D()

double getLossDist2D ( double  pt,
int  pdgCode,
double  arcLength2D 
) const

Calculates a correction term for the two dimensional distance undoing the energy loss after the given arc length.

Definition at line 170 of file EnergyLossEstimator.cc.

171{
172 const double eLoss = getEnergyLoss(pt, pdgCode, arcLength2D);
173
174 const int q = getCharge(pdgCode);
175 const double radius = q * CDCBFieldUtil::absMom2DToBendRadius(pt, m_bZ);;
176 return radius * eLoss / (pt - eLoss);
177
178 // const double pFactor = getMomentumLossFactor(p, pdgCode, arcLength2D);
179 // return radius * (1 / pFactor - 1);
180}
static double absMom2DToBendRadius(double absMom2D, double bZ)
Conversion helper for momenta to two dimensional (absolute) bend radius.
static int getCharge(int pdgCode)
Lookup the charge for the given pdg code.
double getEnergyLoss(double p, int pdgCode, double arcLength) const
Calculates the total energy loss after travelling the given distance.

◆ getMass()

double getMass ( int  pdgCode)
static

Lookup the mass for the given pdg code.

Definition at line 84 of file EnergyLossEstimator.cc.

85{
86 if (std::abs(pdgCode) == Const::electron.getPDGCode()) {
87 const double electronMass = 0.511 * Unit::MeV;
88 return electronMass;
89 } else if (std::abs(pdgCode) == Const::muon.getPDGCode()) {
90 const double muonMass = 105.658 * Unit::MeV;
91 return muonMass;
92 } else if (std::abs(pdgCode) == Const::kaon.getPDGCode()) {
93 const double kaonMass = 0.493677;
94 return kaonMass;
95 } else if (std::abs(pdgCode) == Const::pion.getPDGCode()) {
96 const double pionMass = 0.13957;
97 return pionMass;
98 } else if (std::abs(pdgCode) == Const::proton.getPDGCode()) {
99 const double protonMass = 0.938272;
100 return protonMass;
101 }
102 return NAN;
103}

◆ getMomentumLossFactor()

double getMomentumLossFactor ( double  p,
int  pdgCode,
double  arcLength 
) const

Calculates a factor applicable scaling the current momentum to the momentum after traveling given arc length.

Definition at line 157 of file EnergyLossEstimator.cc.

158{
159 const double eLoss = getEnergyLoss(p, pdgCode, arcLength);
160 return (p - eLoss) / p;
161
162 const double mass = getMass(pdgCode);
163
164 const double eBefore = std::sqrt(p * p + mass * mass);
165 const double eAfter = eBefore - eLoss;
166 const double pAfter = std::sqrt(eAfter * eAfter - mass * mass);
167 return pAfter / p;
168}

Member Data Documentation

◆ m_bZ

double m_bZ
private

B field to be used for the distance translation.

Definition at line 81 of file EnergyLossEstimator.h.

◆ m_eDensity

double m_eDensity
private

Electron density in mol / cm^3.

Definition at line 75 of file EnergyLossEstimator.h.

◆ m_I

double m_I
private

Mean excitation energy in GeV.

Definition at line 78 of file EnergyLossEstimator.h.


The documentation for this class was generated from the following files: