G4DarkBreM  v2.1.0
Geant4 Dark Bremmstrahlung from MadGraph
g4db::G4DarkBreMModel Class Reference

Geant4 implementation of the model for a particle undergoing a dark brem where we use an imported event library to decide the outgoing kinematics. More...

#include <G4DarkBreMModel.h>

Inheritance diagram for g4db::G4DarkBreMModel:
[legend]
Collaboration diagram for g4db::G4DarkBreMModel:
[legend]

Public Types

enum class  ScalingMethod { ForwardOnly = 1 , CMScaling = 2 , Undefined = 3 }
 Possible methods to use the dark brem events from the imported library inside of this model. More...
 
enum class  XsecMethod { Full = 1 , Improved = 2 , HyperImproved = 3 , Auto = 4 }
 Possible methods for calculating the xsec using the WW approximation technique. More...
 

Public Member Functions

 G4DarkBreMModel (const std::string &library_path, bool muons, double threshold=0.0, double epsilon=1.0, ScalingMethod sm=ScalingMethod::ForwardOnly, XsecMethod xm=XsecMethod::Auto, double max_R_for_full=50.0, int aprime_lhe_id=622, bool load_library=true, bool scale_APrime=false, double dist_decay_min=0.0, double dist_decay_max=1.0)
 Set the parameters for this model. More...
 
virtual ~G4DarkBreMModel ()=default
 Destructor.
 
virtual void PrintInfo () const
 Print the configuration of this model.
 
virtual G4double ComputeCrossSectionPerAtom (G4double lepton_ke, G4double atomicA, G4double atomicZ)
 Calculates the cross section per atom in GEANT4 internal units. More...
 
std::pair< G4ThreeVector, G4ThreeVector > scale (double target_Z, double incident_energy, double lepton_mass)
 Scale one of the MG events in our library to the input incident lepton energy. More...
 
virtual void GenerateChange (G4ParticleChange &particleChange, const G4Track &track, const G4Step &step, const G4Element &element)
 Simulates the emission of a dark photon + lepton. More...
 
- Public Member Functions inherited from g4db::PrototypeModel
 PrototypeModel (bool muons)
 Constructor. More...
 
virtual ~PrototypeModel ()=default
 Destructor, nothing on purpose.
 
void SetVerboseLevel (int l)
 Define the verbose level. More...
 
int GetVerboseLevel () const
 Get the verbose level of the model. More...
 
bool DarkBremOffMuons () const
 Check if we are bremming off muons. More...
 

Private Member Functions

void SetMadGraphDataLibrary (const std::string &path)
 Set the library of dark brem events to be scaled. More...
 
void MakePlaceholders ()
 Fill vector of currentDataPoints_ with the same number of items as the madgraph data. More...
 
OutgoingKinematics sample (double target_Z, double incident_energy)
 Returns MadGraph data given an energy [GeV]. More...
 

Private Attributes

unsigned int maxIterations_ {10000}
 maximum number of iterations to check before giving up on an event More...
 
double threshold_
 Threshold for non-zero xsec [GeV]. More...
 
double epsilon_
 Epsilon value to plug into xsec calculation. More...
 
int aprime_lhe_id_
 PDG ID number for the A' (dark photon) as written in the LHE files being loaded as the dark brem event library. More...
 
ScalingMethod scaling_method_
 method to scale events for this model
 
XsecMethod xsec_method_
 method for calculating total cross section
 
std::string method_name_
 Name of method for persisting into the RunHeader.
 
std::string library_path_
 Full path to the vertex library used for persisting into the RunHeader.
 
bool alwaysCreateNewLepton_ {true}
 should we always create a totally new lepton when we dark brem? More...
 
bool scale_APrime_ {false}
 whether to scale the outgoing A' momentum based on the MadGraph library. More...
 
double dist_decay_min_ {0.0}
 Minimum flight distance [mm] at which to decay the A' if G4APrime::DecayMode is set to FlatDecay.
 
double dist_decay_max_ {1.0}
 Maximum flight distance [mm] at which to decay the A' if G4APrime::DecayMode is set to FlatDecay.
 
std::map< int, std::map< double, std::vector< OutgoingKinematics > > > madGraphData_
 Storage of data from mad graph. More...
 
std::map< int, std::map< double, unsigned int > > currentDataPoints_
 Stores a map of current access points to mad graph data. More...
 

Additional Inherited Members

- Protected Attributes inherited from g4db::PrototypeModel
bool muons_
 whether muons (true) or electrons (false) are dark bremming
 
int verbose_level_ {0}
 verbose level for this model
 

Detailed Description

Geant4 implementation of the model for a particle undergoing a dark brem where we use an imported event library to decide the outgoing kinematics.

This is where all the heavy lifting in terms of calculating cross sections and actually having a lepton do a dark brem occurs. This model depends on several configurable parameters which should be studied for your specific use case in order to better tune this model to your lepton/incident energy/ target material situation.

Member Enumeration Documentation

◆ ScalingMethod

Possible methods to use the dark brem events from the imported library inside of this model.

Enumerator
ForwardOnly 

Use actual lepton energy and get pT from LHE (such that \(p_T^2+m_l^2 < E_{acc}^2\))

Scales the energy so that the fraction of kinetic energy is constant, keeping the \(p_T\) constant.

If the \(p_T\) is larger than the new energy, that event is skipped, and a new one is taken from the file. If the loaded library does not fully represent the range of incident energies being seen by the simulation, this will occur frequently.

With only the kinetic energy fraction and \(p_T\), the sign of the longitudinal momentum \(p_z\) is undetermined. This method simply chooses the \(p_z\) of the recoil lepton to always be positive.

CMScaling 

Boost LHE vertex momenta to the actual lepton energy.

The scaling is done via two boosts.

  1. Boost out of the center-of-momentum (CoM) frame read in along with the MadGraph event library.
  2. Boost into approximately) the incident lepton energy frame by constructing a "new" CoM frame using the actual CoM frame's transverse momentum and lowering the \(p_z\) and energy of the CoM by the difference between the input incident energy and the sampled incident energy.

After these boosts, the energy of the recoil and its \(p_T\) are extracted.

Undefined 

Use LHE vertex as is.

We simply copy the read-in recoil energy's energy, momentum, and \(p_T\).

◆ XsecMethod

Possible methods for calculating the xsec using the WW approximation technique.

See also
flux_factor_chi_numerical for how the flux factor is calculated Here, we refer to the flux factor as \(\chi(x,\theta)\).
integrate for how the numerical integrations are performed

As the methods are related, they share many variable definitions.

\begin{equation} \tilde{u} = -xE_0^2\theta^2 - m_A^2\frac{1-x}{x} - m_\ell^2x \end{equation}

  • \(x_{max} = 1 - \max(m_A,m_\ell)/E_0\) is the upper limit on the x integration
  • \(E_0\) is the incoming lepton's energy in GeV
  • \(m_\ell\) is the lepton's mass in GeV
  • \(m_A\) is the dark photon's mass in GeV
  • \(\alpha_{EW} = 1/137\) is the fine-structure constant
  • \(\epsilon\) is the dark photon mixsing strength with the standard photon
  • \(pb/GeV = 3.894\times10^8\) is the conversion from GeV to pico-barns

The integrand limits for \(\chi\) are calculated as

\begin{equation} t_{min} = \left(\frac{\tilde{u}}{2E_0(1-x)}\right)^2 \qquad t_{max} = E_0^2 \end{equation}

Enumerator
Full 

use the full WW approximation

This is the 2D integration of a differential cross section including evaluating the flux factor chi numerically at each sampled x and theta. This effectively means we perform a 3D numerical integration when using this method.

It is slow but very faithful to the sampled MadGraph total cross section for high (> 200 GeV) incident energies.

\begin{equation} \sigma = \frac{pb}{GeV} \int_0^{x_{max}} \int_0^{0.3} \frac{d\sigma}{dx d\theta} d\theta~dx \end{equation}

\begin{equation} \frac{d\sigma}{dx~d\cos\theta} = 2 \alpha_{EW}^3\epsilon^2 \sqrt{x^2E_0^2 - m_A^2}E_0(1-x) \frac{\chi(x,\theta)}{\tilde{u}^2} \mathcal{A}^2 \end{equation}

\begin{equation} \mathcal{A}^2 = 2\frac{2-2x+x^2}{1-x}+\frac{4(m_A^2+2m_\ell^2)}{\tilde{u}^2}(\tilde{u}x + m_A^2(1-x) + m_\ell^2x^2) \end{equation}

Improved 

assume \(\theta=0\) in the integrand

This reduces one of the dimensions of the numerical integration, allowing the theta integration to be done analytically and so we only have a numerical integration over x and t

\begin{equation} \sigma = \frac{pb}{GeV} \int_0^{x_{max}} \chi(x,\theta=0) \frac{d\sigma}{dx} dx \end{equation}

\begin{equation} \frac{d\sigma}{dx}(x) = 4 \alpha_{EW}^3\epsilon^2 \sqrt{1-\frac{m_A^2}{E_0^2}}\frac{1-x+x^2/3}{m_A^2(1-x)/x+m_\ell^2x} \end{equation}

HyperImproved 

only calculate the flux factor once

This simplifies the numerical integration even further by calculating the flux factor chi only at \((x=1, \theta=0)\) and then using that value everywhere in the integration of the differential cross section over x.

\begin{equation} \sigma = \frac{pb}{GeV} \chi(x=1,\theta=0) \int_0^{x_{max}} \frac{d\sigma}{dx} dx \end{equation}

where \(d\sigma/dx\) is the same as Improved. Since this uses the maximum value of \(\chi\) for all values of \(\chi\), we modify the upper limit of integration to avoid making this overestimate too much of an overestimate. Instead of \(t_{max}=E_0^2\) we use \(t_{max}=m_A^2+m_\ell^2\).

Auto 

the default cross section calculation

This requires its own enum variant since it is determined by the \(R = m_A/m_\ell\) ratio as configured by other parameters.

If \(R < R_{max}\), Full is used while Improved is used otherwise. These boarders were determine qualitatively by looking at a series of cross section plots comparing all of these methods.

The parameter \(R_{max}\) is configurable in G4DarkBreMModel::G4DarkBreMModel as max_R_for_full.

Constructor & Destructor Documentation

◆ G4DarkBreMModel()

g4db::G4DarkBreMModel::G4DarkBreMModel ( const std::string &  library_path,
bool  muons,
double  threshold = 0.0,
double  epsilon = 1.0,
ScalingMethod  sm = ScalingMethod::ForwardOnly,
XsecMethod  xm = XsecMethod::Auto,
double  max_R_for_full = 50.0,
int  aprime_lhe_id = 622,
bool  load_library = true,
bool  scale_APrime = false,
double  dist_decay_min = 0.0,
double  dist_decay_max = 1.0 
)

Set the parameters for this model.

Parameters
[in]library_pathfull path to MG library The library path is immediately passed to G4DarkBreMModel::SetMadGraphDataLibrary which simply calls g4db::parseLibrary (where you should look for the specification options of a dark brem library).
[in]muonstrue if using muons, false for electrons
[in]thresholdminimum energy lepton needs to have to dark brem [GeV] (overwritten by twice the A' mass if it is less so that it kinematically makes sense).
[in]epsilondark photon mixing strength
[in]smG4DarkBreMModel::ScalingMethod on how to scale the events in the library
[in]xmG4DarkBreMModel::XsecMethod on calculating the total cross section
[in]max_R_for_fullmaximum value of R where Full will still be used Only used with G4DarkBreMModel::XsecMethod::Auto (go there for context)
[in]aprime_lhe_idPDG ID number for the dark photon in the LHE files being loaded for the library. Only used if parsing LHE files.
[in]load_libraryonly used in cross section executable where it is known that the library will not be used during program run
[in]scale_APrimewhether to scale the A' kinematics based on the library. If false, define A' momentum to conserve momentum without taking nuclear recoil into account.
[in]dist_decay_minminimum flight distance [mm] at which to decay the A' if G4APrime::DecayMode is set to FlatDecay
[in]dist_decay_maxmaximum flight distance [mm] at which to decay the A' if G4APrime::DecayMode is set to FlatDecay

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double g4db::G4DarkBreMModel::ComputeCrossSectionPerAtom ( G4double  lepton_ke,
G4double  atomicA,
G4double  atomicZ 
)
virtual

Calculates the cross section per atom in GEANT4 internal units.

The estimate for the total cross section given the material and the lepton's energy is done using an implementation of the WW approximation using Boost's Math Quadrature library to numerically calculate the integrals.

See also
XsecMethod for the different methods of calculating the cross section.
Parameters
lepton_kekinetic energy of incoming particle
atomicZatomic number of atom
atomicAatomic mass of atom
Returns
cross section (0. if outside energy cuts)

Implements g4db::PrototypeModel.

◆ GenerateChange()

void g4db::G4DarkBreMModel::GenerateChange ( G4ParticleChange &  particleChange,
const G4Track &  track,
const G4Step &  step,
const G4Element &  element 
)
virtual

Simulates the emission of a dark photon + lepton.

See also
scale for how the event library is sampled and scaled to the incident lepton's actual energy

After calling scale, we rotate the outgoing lepton's momentum to the frame of the incident particle and then calculate the dark photon's momentum such that three-momentum is conserved.

Parameters
[in,out]particleChangestructure holding changes to make to particle track
[in]trackcurrent track being processesed
[in]stepcurrent step of the track
[in]elementelement we are going to dark brem off of

Implements g4db::PrototypeModel.

◆ MakePlaceholders()

void g4db::G4DarkBreMModel::MakePlaceholders ( )
private

Fill vector of currentDataPoints_ with the same number of items as the madgraph data.

Randomly choose a starting point so that the simulation run isn't dependent on the order of the events as written in the LHE library. The random starting position is uniformly chosen using G4Uniform() so. The sample function will loop from the last event parsed back to the first event parsed so that the starting position does not matter.

In this function, we also update maxIterations_ so that it is equal to the smallest entry in the library (with a maximum of 10k). This saves time in the situation where an incorrect library was accidentally used and the simulation is looping through events attempting to find one that can fit its criteria.

◆ sample()

OutgoingKinematics g4db::G4DarkBreMModel::sample ( double  target_Z,
double  incident_energy 
)
private

Returns MadGraph data given an energy [GeV].

Gets the energy fraction and \(p_T\) from the imported LHE data. incident_energy should be in GeV, returns the sample outgoing kinematics.

Samples from the closest imported incident energy above the given value (this helps avoid biasing issues).

Parameters
target_Zatomic Z of target nucleus
incident_energyenergy of particle undergoing dark brem [GeV]
Returns
sample outgoing kinematics

◆ scale()

std::pair< G4ThreeVector, G4ThreeVector > g4db::G4DarkBreMModel::scale ( double  target_Z,
double  incident_energy,
double  lepton_mass 
)

Scale one of the MG events in our library to the input incident lepton energy.

This is also helpful for testing the scaling procedure in its own executable separate from the Geant4 infrastructure.

Note
The vector returned is relative to the incident lepton as if it came in along the z-axis.

Gets an energy fraction and transverse momentum ( \(p_T\)) from the loaded library of MadGraph events using the entry in the library with the nearest incident energy above the actual input incident energy.

The scaling of this energy fraction and \(p_T\) to the actual lepton energy depends on the input G4DarkBreMModel::ScalingMethod. In all cases, the azimuthal angle is chosen uniformly between 0 and \(2\pi\).

Parameters
[in]target_Zatomic Z of target nucleus
[in]incident_energyincident total energy of the lepton [GeV]
[in]lepton_massmass of incident lepton [GeV]
Returns
G4ThreeVectors representing the outgoing momenta of the recoil lepton and the A', respectively

◆ SetMadGraphDataLibrary()

void g4db::G4DarkBreMModel::SetMadGraphDataLibrary ( const std::string &  path)
private

Set the library of dark brem events to be scaled.

This function loads the library from the on-disk file (or files).

See also
parseLibrary for how libraries are read and what the structure of a dark brem library is
Parameters
pathpath to the on-disk library

Member Data Documentation

◆ alwaysCreateNewLepton_

bool g4db::G4DarkBreMModel::alwaysCreateNewLepton_ {true}
private

should we always create a totally new lepton when we dark brem?

Note
make this configurable? I (Tom E) can't think of a reason NOT to have it... The alternative is to allow Geant4 to decide when to make a new particle by checking if the resulting kinetic energy is below some threshold.

◆ aprime_lhe_id_

int g4db::G4DarkBreMModel::aprime_lhe_id_
private

PDG ID number for the A' (dark photon) as written in the LHE files being loaded as the dark brem event library.

This is only used if parsing LHE files and can be ignored if the library is being loaded from an already-constructed CSV.

◆ currentDataPoints_

std::map<int, std::map<double, unsigned int> > g4db::G4DarkBreMModel::currentDataPoints_
private

Stores a map of current access points to mad graph data.

Maps target Z and incoming lepton energy to the index of the data vector that we will get the data from.

Also sorts the incoming lepton energy so that we can find the sampling energy that is closest above the actual incoming energy.

◆ epsilon_

double g4db::G4DarkBreMModel::epsilon_
private

Epsilon value to plug into xsec calculation.

See also
ComputeCrossSectionPerAtom for how this is used

Configurable with 'epsilon'

◆ madGraphData_

std::map<int, std::map<double, std::vector<OutgoingKinematics> > > g4db::G4DarkBreMModel::madGraphData_
private

Storage of data from mad graph.

Maps target Z and incoming lepton energy to various options for outgoing kinematics. This is a hefty map and is what stores all of the events imported from the LHE library of dark brem events.

◆ maxIterations_

unsigned int g4db::G4DarkBreMModel::maxIterations_ {10000}
private

maximum number of iterations to check before giving up on an event

This is only used in the ForwardOnly scaling method and is only reached if the event library energies are not appropriately matched with the energy range of particles that are existing in the simulation.

◆ scale_APrime_

bool g4db::G4DarkBreMModel::scale_APrime_ {false}
private

whether to scale the outgoing A' momentum based on the MadGraph library.

The same scaling method is used as for the recoil lepton. The difference between the azimuthal angle of the A' and the recoil lepton is designed to be preserved from MadGraph. The preservation of this angle and the overall scaling of the A' has only been validated for electrons with incident energies ranging from 1GeV to 10GeV and dark photon masses ranging from 1MeV to 100MeV. While the scaling is expected to be well behaved for muons and at other energy scales, users are encouraged to double-check this behavior in their situation and report issues to the repository. https://github.com/LDMX-Software/G4DarkBreM

If false, will define the 3-momentum of the A' to conserve 3-momentum with primary and recoil lepton, not taking into account the nuclear recoil. This is the default because then the user is able to "reconstruct" the true incident lepton's momentum using the recoil lepton and produced dark photon's three-momenta. Scaling the A' momentum is only advisable if the A' momentum itself is important to the search (for example, when the A' decays visibly and the decay is observed in which case the A' momentum deciding where the decay occurs is very important).

◆ threshold_

double g4db::G4DarkBreMModel::threshold_
private

Threshold for non-zero xsec [GeV].

Configurable with 'threshold'. At minimum, it is always at least twice the dark photon mass.


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