G4DarkBreM
v2.1.0
Geant4 Dark Bremmstrahlung from MadGraph
|
Class that represents the dark brem process. More...
#include <G4DarkBremsstrahlung.h>
Public Types | |
typedef void(* | StorageCallback) (const G4Element &) |
Function signature able to be used as a storage callback. | |
Public Member Functions | |
G4DarkBremsstrahlung (std::shared_ptr< g4db::PrototypeModel > the_model, bool only_one_per_event=false, double global_bias=1., bool interp_xsec=true, bool cache_xsec=false, int verbose_level=0, int subtype=63) | |
Constructor. More... | |
virtual | ~G4DarkBremsstrahlung ()=default |
Destructor. | |
void | RegisterStorageMechanism (StorageCallback f) |
Register a storage mechanism for saving the material in-which the brem occurred. More... | |
virtual G4bool | IsApplicable (const G4ParticleDefinition &p) |
Checks if the passed particle should be able to do this process. More... | |
virtual void | PrintInfo () |
Reports the parameters to G4cout. More... | |
virtual G4VParticleChange * | PostStepDoIt (const G4Track &track, const G4Step &step) |
This is the function actually called by Geant4 that does the dark brem interaction. More... | |
g4db::ElementXsecCache & | getCache () |
Get a reference to the cross section cache. More... | |
Static Public Attributes | |
static const std::string | PROCESS_NAME = "DarkBrem" |
The name of this process in Geant4. | |
Protected Member Functions | |
G4double | GetMeanFreePath (const G4Track &track, G4double prevStepSize, G4ForceCondition *condition) |
Calculate the mean free path given the input conditions. More... | |
Private Member Functions | |
G4DarkBremsstrahlung & | operator= (const G4DarkBremsstrahlung &right) |
remove ability to assign this object | |
G4DarkBremsstrahlung (const G4DarkBremsstrahlung &) | |
remove ability to copy construct | |
Private Attributes | |
bool | only_one_per_event_ |
Only allow the dark brem to happen once per event. More... | |
double | global_bias_ |
Bias the dark brem cross section GLOBALLY. | |
bool | interpolate_xsec_ |
Should we interpolated over the computed cross sections? | |
bool | cache_xsec_ |
Should we have a cache for the computed cross sections? | |
std::shared_ptr< g4db::PrototypeModel > | model_ |
The model that we are using in this run. More... | |
g4db::ElementXsecCache | element_xsec_cache_ |
Our instance of a cross section cache. | |
g4db::ElementXsecInterpolation | element_xsec_interpolation_ |
Our instance of a cross section interpolation. | |
std::vector< double > | partial_sum_sigma_ |
Hold the calculations of the cross section per atom here so we can select a random element for the dark brem to go off of when Geant4 decides it is time. More... | |
StorageCallback | storage_func_ = nullptr |
Class that represents the dark brem process.
A muon or electron is allowed to brem a dark photon
G4DarkBremsstrahlung::G4DarkBremsstrahlung | ( | std::shared_ptr< g4db::PrototypeModel > | the_model, |
bool | only_one_per_event = false , |
||
double | global_bias = 1. , |
||
bool | interp_xsec = true , |
||
bool | cache_xsec = false , |
||
int | verbose_level = 0 , |
||
int | subtype = 63 |
||
) |
Constructor.
Configures this process by doing three main things:
We add ourselves to the process table for muons (or electrons), so the caller does not need to do this. Normally, this is done in a physics constructor, but we have chosen to do this differently to avoid mistakenly adding the dark brem process configured for muons to the electron table (and vis-versa).
[in] | the_model | model to use for dark brem simulation |
[in] | only_one_per_event | true if de-activating process after first dark brem |
[in] | global_bias | bias xsec globally by this factor |
[in] | interpolate_xsec | true if we should interpolate cross sections over 10% energy differences |
[in] | cache_xsec | true if we should cache xsecs at the MeV level of precision |
[in] | verbose_level | level of verbosity to print for this process and model |
[in] | subtype | subtype for this process distinct from other EM processes in use (detail below) |
The Geant4 enum G4EmProcessSubType lists the standard subtypes; however, some physics lists introduce more EM subtypes than those listed there. You can use G4ProcessManager::GetProcessList in order to find out what EM subtypes are available in your specific physics list configuration. Deducing an available sub type from the current process list in code is dangerous since there may be more processes constructed after us.
In G4 speak, a "discrete" process is one that only happens at the end of steps. We want the DB to be discrete because it is not a "slow braking" like ionization, the lepton suddenly has the interaction and loses a lot of its energy.
We have our custom dark brem process go first in any process ordering so that we always have the opportunity to dark brem if the cross section allows it.
|
inline |
Get a reference to the cross section cache.
Again, this method is public only to be available to the executable that generates a cross section table and testing. Do not use this unless you really know what you are doing.
|
protected |
Calculate the mean free path given the input conditions.
If the input track's particle definition does not align with how we are configured (checked via IsApplicable), then we return DBL_MAX to signal that this should never happen.
We calculate the total cross section by calculating the total cross section for each element in the current material and weighting those cross sections by the number of atoms per volume in the material. This allows the process to handle the material dependence on the cross section, leaving the detailed elemental cross section calculation to the model.
The global_bias
parameter from the constructor is also used here after-the-calculation in order to allow rudimentary biasing. The use of global_bias
is discouraged in favor of using Geant4's biasing infrastructure.
We maintain a cache for the cross sections calculated by the model so that later in the run it is less likely that the model will need to be called to calculate the cross section. This is done in order to attempt to improve speed of simulation and avoid repetition of the same, deterministic calculations. If you want to turn off the cache-ing behavior, set cache_xsec
to false in the constructor.
If the total cross section is above DBL_MIN, then it is inverted to obtain the mean free path. Otherwise, DBL_MAX is returned.
[in] | track | G4Track that is being stepped |
[in] | prevStepSize | G4double measuring previous step size, unused |
[in] | condition | G4ForceCondition, always NotForced for G4VDiscreteProcess, unused |
|
virtual |
Checks if the passed particle should be able to do this process.
|
virtual |
This is the function actually called by Geant4 that does the dark brem interaction.
aParticleChange is a protected member variable of G4VDiscreteProcess that we should edit here.
If only one per event is set, then we deactivate the dark brem process, ensuring only one dark brem per step and per event. Reactivated in RunManager::TerminateOneEvent.
[in] | track | current G4Track that is being stepped |
[in] | step | current step that just finished |
If configured to do so, we deactivate the process after one dark brem. If this is in the stepping action instead, more than one brem can occur within each step.
Both biased and unbiased process could be in the run (but not at the same time), so we turn off both while silencing the warnings from the process table. Reactivating the process is essentially the same code as here except using true
intead of false
in SetProcessActivation.
|
virtual |
Reports the parameters to G4cout.
void G4DarkBremsstrahlung::RegisterStorageMechanism | ( | StorageCallback | f | ) |
Register a storage mechanism for saving the material in-which the brem occurred.
[in] | f | function that stores information from G4Element as desired |
|
private |
The model that we are using in this run.
Shared with the chaching class.
|
private |
Only allow the dark brem to happen once per event.
This allows for the dark brem process to be de-activated when SampleSecondaries is called.
|
private |
Hold the calculations of the cross section per atom here so we can select a random element for the dark brem to go off of when Geant4 decides it is time.
The elements in this vector correspond to the element vector from the material. It is called "partial sum sigma" because it actually is the running cumulative sum of the cross section as we loop through the different elements in the material.