G4DarkBreM  v2.1.0
Geant4 Dark Bremmstrahlung from MadGraph
G4DarkBremsstrahlung.h
Go to the documentation of this file.
1 
8 #ifndef G4DARKBREM_G4DARKBREMSSTRAHLUNG_H_
9 #define G4DARKBREM_G4DARKBREMSSTRAHLUNG_H_
10 
11 // Geant
12 #include "G4DarkBreM/ElementXsecCache.h"
13 #include "G4DarkBreM/ElementXsecInterpolation.h"
15 #include "G4VDiscreteProcess.hh"
16 
17 class G4String;
18 class G4ParticleDefinition;
19 
26 class G4DarkBremsstrahlung : public G4VDiscreteProcess {
27  public:
31  typedef void (*StorageCallback)(const G4Element&);
32 
36  static const std::string PROCESS_NAME;
37 
68  G4DarkBremsstrahlung(std::shared_ptr<g4db::PrototypeModel> the_model,
69  bool only_one_per_event = false, double global_bias = 1.,
70  bool interp_xsec = true, bool cache_xsec = false,
71  int verbose_level = 0, int subtype = 63);
72 
76  virtual ~G4DarkBremsstrahlung() = default;
77 
85 
91  virtual G4bool IsApplicable(const G4ParticleDefinition& p);
92 
98  virtual void PrintInfo();
99 
117  virtual G4VParticleChange* PostStepDoIt(const G4Track& track,
118  const G4Step& step);
119 
128 
129  protected:
167  G4double GetMeanFreePath(const G4Track& track, G4double prevStepSize,
168  G4ForceCondition* condition);
169 
170  private:
173 
176 
188 
192  double global_bias_;
193 
198 
203 
209  std::shared_ptr<g4db::PrototypeModel> model_;
210 
213 
216 
227  std::vector<double> partial_sum_sigma_;
228 
229  StorageCallback storage_func_ = nullptr;
230 }; // G4DarkBremsstrahlung
231 
232 #endif
Class providing a prototype model for dark brem.
Class that represents the dark brem process.
Definition: G4DarkBremsstrahlung.h:26
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.
Definition: G4DarkBremsstrahlung.cxx:21
static const std::string PROCESS_NAME
The name of this process in Geant4.
Definition: G4DarkBremsstrahlung.h:36
G4DarkBremsstrahlung(const G4DarkBremsstrahlung &)
remove ability to copy construct
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)
This is the function actually called by Geant4 that does the dark brem interaction.
Definition: G4DarkBremsstrahlung.cxx:123
G4double GetMeanFreePath(const G4Track &track, G4double prevStepSize, G4ForceCondition *condition)
Calculate the mean free path given the input conditions.
Definition: G4DarkBremsstrahlung.cxx:210
double global_bias_
Bias the dark brem cross section GLOBALLY.
Definition: G4DarkBremsstrahlung.h:192
bool cache_xsec_
Should we have a cache for the computed cross sections?
Definition: G4DarkBremsstrahlung.h:202
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 da...
Definition: G4DarkBremsstrahlung.h:227
virtual ~G4DarkBremsstrahlung()=default
Destructor.
virtual G4bool IsApplicable(const G4ParticleDefinition &p)
Checks if the passed particle should be able to do this process.
Definition: G4DarkBremsstrahlung.cxx:107
g4db::ElementXsecCache & getCache()
Get a reference to the cross section cache.
Definition: G4DarkBremsstrahlung.h:127
void RegisterStorageMechanism(StorageCallback f)
Register a storage mechanism for saving the material in-which the brem occurred.
Definition: G4DarkBremsstrahlung.cxx:102
bool interpolate_xsec_
Should we interpolated over the computed cross sections?
Definition: G4DarkBremsstrahlung.h:197
G4DarkBremsstrahlung & operator=(const G4DarkBremsstrahlung &right)
remove ability to assign this object
void(* StorageCallback)(const G4Element &)
Function signature able to be used as a storage callback.
Definition: G4DarkBremsstrahlung.h:31
std::shared_ptr< g4db::PrototypeModel > model_
The model that we are using in this run.
Definition: G4DarkBremsstrahlung.h:209
g4db::ElementXsecCache element_xsec_cache_
Our instance of a cross section cache.
Definition: G4DarkBremsstrahlung.h:212
virtual void PrintInfo()
Reports the parameters to G4cout.
Definition: G4DarkBremsstrahlung.cxx:114
g4db::ElementXsecInterpolation element_xsec_interpolation_
Our instance of a cross section interpolation.
Definition: G4DarkBremsstrahlung.h:215
bool only_one_per_event_
Only allow the dark brem to happen once per event.
Definition: G4DarkBremsstrahlung.h:187
The cache of already computed cross sections.
Definition: ElementXsecCache.h:17
Interpolation of cross sections by element.
Definition: ElementXsecInterpolation.h:19