8#include <tracking/trackFindingCDC/eventdata/utils/EnergyLossEstimator.h>
10#include <tracking/trackingUtilities/eventdata/trajectories/CDCBFieldUtil.h>
12#include <tracking/trackingUtilities/numerics/ESign.h>
14#include <geometry/GeometryManager.h>
16#include <framework/gearbox/Unit.h>
17#include <framework/gearbox/Const.h>
18#include <framework/logging/Logger.h>
20#include "G4ThreeVector.hh"
21#include "G4Navigator.hh"
22#include "G4VPhysicalVolume.hh"
23#include "G4Material.hh"
26using namespace TrackFindingCDC;
27using namespace TrackingUtilities;
35 g4Nav.SetWorldVolume(g4World);
42 G4ThreeVector g4Pos(posX * CLHEP::cm, posY * CLHEP::cm, posZ * CLHEP::cm);
43 const G4VPhysicalVolume* g4Volume = g4Nav.LocateGlobalPointAndSetup(g4Pos);
44 const G4Material* g4Mat = g4Volume->GetLogicalVolume()->GetMaterial();
49 if (g4Mat->GetNumberOfElements() == 1) {
50 A = g4Mat->GetA() * CLHEP::mole / CLHEP::g;
54 for (
unsigned i = 0; i < g4Mat->GetNumberOfElements(); ++i) {
55 const G4Element* element = (*g4Mat->GetElementVector())[i];
56 const double elementA = element->GetA() * CLHEP::mole / CLHEP::g;
57 const double elementZ = element->GetZ();
58 const double frac = g4Mat->GetFractionVector()[i];
59 B2RESULT(
"Part " << i <<
" Z=" << elementZ <<
" A=" << elementA <<
" (" << frac <<
")");
67 const double density = g4Mat->GetDensity() / CLHEP::g * CLHEP::cm3;
68 const double radiationLength = g4Mat->GetRadlen() / CLHEP::cm;
69 const double mEE = g4Mat->GetIonisation()->GetMeanExcitationEnergy() / CLHEP::GeV;
71 B2RESULT(
"Received Z " << Z);
72 B2RESULT(
"Received A " << A);
73 B2RESULT(
"Received density " << density);
74 B2RESULT(
"Received radiation length " << radiationLength);
75 B2RESULT(
"Received mean excitation energy " << mEE);
77 const double eDensity = Z * density / A;
78 B2RESULT(
"Received electron density " << eDensity);
88 const double electronMass = 0.511 *
Unit::MeV;
90 }
else if (std::abs(pdgCode) ==
Const::muon.getPDGCode()) {
91 const double muonMass = 105.658 *
Unit::MeV;
93 }
else if (std::abs(pdgCode) ==
Const::kaon.getPDGCode()) {
94 const double kaonMass = 0.493677;
96 }
else if (std::abs(pdgCode) ==
Const::pion.getPDGCode()) {
97 const double pionMass = 0.13957;
100 const double protonMass = 0.938272;
109 return -sign(pdgCode);
110 }
else if (std::abs(pdgCode) ==
Const::muon.getPDGCode()) {
111 return -sign(pdgCode);
112 }
else if (std::abs(pdgCode) ==
Const::kaon.getPDGCode()) {
113 return sign(pdgCode);
114 }
else if (std::abs(pdgCode) ==
Const::pion.getPDGCode()) {
115 return sign(pdgCode);
116 }
else if (std::abs(pdgCode) ==
Const::proton.getPDGCode()) {
117 return sign(pdgCode);
131 static const double eMass =
getMass(11);
133 const double M =
getMass(pdgCode);
135 const double E = std::sqrt(p * p + M * M);
136 const double gamma =
E / M;
137 const double beta = p /
E;
139 const double beta2 = beta * beta;
140 const double gamma2 = gamma * gamma;
142 const double Wmax = 2 * eMass * beta2 * gamma2 / (1 + 2 * gamma * eMass / M);
143 const double I2 =
m_I *
m_I;
147 (1.0 / 2.0 * std::log(2 * eMass * beta2 * gamma2 * Wmax / I2) - beta2);
154 const double eLoss = arcLength * dEdx;
161 return (p - eLoss) / p;
163 const double mass =
getMass(pdgCode);
165 const double eBefore = std::sqrt(p * p + mass * mass);
166 const double eAfter = eBefore - eLoss;
167 const double pAfter = std::sqrt(eAfter * eAfter - mass * mass);
173 const double eLoss =
getEnergyLoss(pt, pdgCode, arcLength2D);
177 return radius * eLoss / (pt - eLoss);
static const ChargedStable muon
muon particle
static const ChargedStable pion
charged pion particle
static const ChargedStable proton
proton particle
static const ChargedStable kaon
charged kaon particle
static const ChargedStable electron
electron particle
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...
double getMomentumLossFactor(double p, int pdgCode, double arcLength) const
Calculates a factor applicable scaling the current momentum to the momentum after traveling given arc...
static EnergyLossEstimator forCDC()
Create an energy loss estimator with the material properties of the CDC.
double m_eDensity
Electron density in mol / cm^3.
double m_bZ
B field to be used for the distance translation.
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.
static double getMass(int pdgCode)
Lookup the mass for the given pdg code.
double m_I
Mean excitation energy in GeV.
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.
static double getBFieldZ()
Getter for the signed magnetic field strength in z direction at the origin ( in Tesla )
static double absMom2DToBendRadius(double absMom2D, double bZ)
Conversion helper for momenta to two dimensional (absolute) bend radius.
static const double MeV
[megaelectronvolt]
static const double cm2
[square centimeters]
G4VPhysicalVolume * getTopVolume()
Return a pointer to the top volume.
static GeometryManager & getInstance()
Return a reference to the instance.
Abstract base class for different kinds of events.