G4DarkBreM  v2.1.0
Geant4 Dark Bremmstrahlung from MadGraph
g4db Namespace Reference

G4DarkBreM internal namespace. More...

Namespaces

 example
 example simulation application required classes
 
 parse
 namespace holding implementation of library parsing
 

Classes

class  ElementXsecCache
 The cache of already computed cross sections. More...
 
class  ElementXsecInterpolation
 Interpolation of cross sections by element. More...
 
class  G4DarkBreMModel
 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...
 
struct  OutgoingKinematics
 Data frame to store necessary information from LHE files. More...
 
class  PrototypeModel
 Abstract class representing a model for dark brem. More...
 

Functions

template<typename IntegrandFunction >
static double integrate (IntegrandFunction F, double low, double high)
 The integration method we will be using for our numerical integrals. More...
 
static double flux_factor_chi_numerical (G4double A, G4double Z, double tmin, double tmax)
 numerically integrate the value of the flux factor chi More...
 
bool hasEnding (const std::string &str, const std::string &end)
 Check if the input string has the input ending. More...
 
void parseLibrary (const std::string &path, int aprime_lhe_id, std::map< int, std::map< double, std::vector< OutgoingKinematics >>> &lib)
 parse the input library and return the in-memory kinematics library More...
 
void dumpLibrary (std::ostream &o, const std::map< int, std::map< double, std::vector< OutgoingKinematics >>> &lib)
 Dump the input library to the input output stream. More...
 

Detailed Description

G4DarkBreM internal namespace.

Function Documentation

◆ dumpLibrary()

void g4db::dumpLibrary ( std::ostream &  o,
const std::map< int, std::map< double, std::vector< OutgoingKinematics >>> &  lib 
)

Dump the input library to the input output stream.

Note
This is only helpful for our extraction executable and for testing that the parsing is performing correctly.
Parameters
[in,out]ooutput stream to write CSV to
[in]liblibrary to write out

This function writes out the input library as CSV to the input output stream in the same format as expected by parse::csv.

◆ flux_factor_chi_numerical()

static double g4db::flux_factor_chi_numerical ( G4double  A,
G4double  Z,
double  tmin,
double  tmax 
)
static

numerically integrate the value of the flux factor chi

The integration of the form factor into the flux factor can be done analytically with a tool like mathematica, but when including the inelastic term, it produces such a complicated result that the numerical integration is actually faster than the analytical one.

The form factors are copied from Appendix A (Eq A18 and A19) of https://journals.aps.org/prd/pdf/10.1103/PhysRevD.80.075018

Here, the equations are listed for reference.

\begin{equation} \chi(x,\theta) = \int^{t_{max}}_{t_{min}} dt \left( \frac{Z^2a^4t^2}{(1+a^2t)^2(1+t/d)^2}+\frac{Za_p^4t^2}{(1+a_p^2t)^2(1+t/0.71)^8}\left(1+\frac{t(\mu_p^2-1)}{4m_p^2}\right)^2\right)\frac{t-t_{min}}{t^2} \end{equation}

where

\begin{equation} a = \frac{111.0}{m_e Z^{1/3}} \quad a_p = \frac{773.0}{m_e Z^{2/3}} \quad d = \frac{0.164}{A^{2/3}} \end{equation}

  • \(m_e\) is the mass of the electron in GeV
  • \(m_p = 0.938\) is the mass of the proton in GeV
  • \(\mu_p = 2.79\) is the proton \(\mu\)
  • \(A\) is the atomic mass of the target nucleus in amu
  • \(Z\) is the atomic number of the target nucleus
Parameters
[in]Aatomic mass of the target nucleus in amu
[in]Zatomic number of target nucleus
[in]tminlower limit of integration over t
[in]tmaxupper limit of integration over t

We've manually expanded the integrand to cancel out the 1/t^2 factor from the differential, this helps the numerical integration converge because we aren't teetering on the edge of division by zero

The auto used in the integrand definition represents a function whose return value is a double and which has a single input lnt. This lambda expression saves us the time of having to re-calculate the form factor constants that do not depend on t because it can inherit their values from the environment. The return value is a double since it is calculated by simple arithmetic operations on doubles.

The integrand is so sharply peaked at t close to tmin, it is very helpful to do the integration in the variable u = ln(t) rather than t itself.

◆ hasEnding()

bool g4db::hasEnding ( const std::string &  str,
const std::string &  end 
)

Check if the input string has the input ending.

Parameters
[in]strstring to test
[in]endending to have
Returns
true if the last characters of str are end

◆ integrate()

template<typename IntegrandFunction >
static double g4db::integrate ( IntegrandFunction  F,
double  low,
double  high 
)
static

The integration method we will be using for our numerical integrals.

The Gauss-Kronrod method was chosen due to its ability to limit the number of calls to the function representing the integrand which should help improve performance for us due to the complexity of our integrand. The order of the GK method was chosen after some experimentation, starting at a high value (61) and then lowering it to achieve better performance while checking the accuracy of the results.

As explained in the Boost GK Docs, generally the error estimation technique for this method is overly pessimistic, so we can confidently set the maximum depth low and the desired relative error high compared to other methods. We have followed the examples in the docs where we use max_depth to 5 and relative error to 1e-9.

Template Parameters
IntegrandFunctionfunction type whose signature is double(*)(double).
Parameters
[in]Ffunction to integrate
[in]lowlower limit of integration
[in]highupper limit of integration

◆ parseLibrary()

void g4db::parseLibrary ( const std::string &  path,
int  aprime_lhe_id,
std::map< int, std::map< double, std::vector< OutgoingKinematics >>> &  lib 
)

parse the input library and return the in-memory kinematics library

Library Specification

A "library" includes a single file or a directory holding multiple files. No subdirectories are inspected, so the directory must be "flat". Each file in a library can either be a CSV file or an LHE file, both of which can be compressed by gzip if desired.

CSV

The CSV file is expected to have a single header line which names the columns. These column names have no requirements (besides the existence of this line).

The CSV is required to have 10 columns on all non-empty lines of the file. The 10 columns of the CSV all are in MeV and in order are

  1. The target Z
  2. The incident lepton energy
  3. The total energy of the recoil
  4. The x-component of the recoil momentum
  5. The y-component of the recoil momentum
  6. The z-component of the recoil momentum
  7. The total energy of the CoM frame
  8. The x-component of the CoM
  9. The y-component of the CoM
  10. The z-component of the CoM

LHE

The LHE files must have dark brem events in it where "dark brem event" in this context is defined below.

lepton_id -1 <skip> <skip> <skip> <skip> px py pz E m
<skip-line>
<skip-line>
lepton_id 1 <skip> <skip> <skip> <skip> px py pz E m
<skip-line>
<skip-line>
aprime_id 1 <skip> <skip> <skip> <skip> px py pz E m

We also check each line if it contains a 'Znuc' string. Looking for the target Z value from this line matching the following format.

<num> <Z> # Znuc <other comments>

This is supposed to match a line of the param_card dumped into the header of the LHE. In other words, this parser should function properly if your MG/ME model has a configurable parameter in param_card.dat for the target Z and that parameter is called 'Znuc'.

This matches a subcomponent of the LHE scheme written by MadGraph/MadEvent (hence the reason this is the "lhe" parser); however, a lot of information is skipped and additional assumptions are made in order to increase the parsing speed.

The lepton_id is allowed to be either 11 or 13 everywhere. No consistency checking is done.

The E from the first line is used as the incident lepton energy. The four-momentum from the middle line is the recoil lepton's four momentum, and the four-momentum from the last line (A') is used in conjunction with the recoil four-momentum to calculate the center of momentum vector.

Parameters
[in]pathpath to library to parse
[in]aprime_lhe_idthe ID number for the A' (aka dark photon) in the library being parsed
[in,out]libmap of target Z and incident energy keys to set of outgoing kinematics

If the input path has one of the four file extensions below, we assume it is a file to be parsed into the library.

  • '.csv'
  • '.csv.gz'
  • '.lhe'
  • '.lhe.gz'

If the extension ends with '.gz', then a decompression step is added to the input stream. Boost.Iostream provides the gzip_decompressor "filter" which does the decompression on the data stream as it is being read in.

See also
parse::csv for files ending with '.csv' or '.csv.gz'
parse::lhe for files ending with '.lhe' or '.lhe.gz'

If the input path does not match one of the four accepted extensions, then we assume it is a directory

Note
We do not recursively enter subdirectories.

If any of the directory entries has one of the acceptable extensions, we recursively call this function on that file path so that it can be parsed into the library.

If we can't open the path that we assumed was a directory as a directory, we end processing.