G4DarkBreM  v2.1.0
Geant4 Dark Bremmstrahlung from MadGraph
G4DarkBreMModel.h
1 #ifndef G4DARKBREM_G4DARKBREMMODEL_H
2 #define G4DARKBREM_G4DARKBREMMODEL_H
3 
4 #include <map>
5 #include <memory>
6 
9 
10 namespace g4db {
11 
25  public:
32  enum class ScalingMethod {
49  ForwardOnly = 1,
50 
66  CMScaling = 2,
73  Undefined = 3
74  };
75 
107  enum class XsecMethod {
130  Full = 1,
131 
145  Improved = 2,
146 
165  HyperImproved = 3,
166 
180  Auto = 4
181  };
182 
216  G4DarkBreMModel(const std::string& library_path, bool muons,
217  double threshold = 0.0, double epsilon = 1.0,
220  double max_R_for_full = 50.0, int aprime_lhe_id = 622,
221  bool load_library = true, bool scale_APrime = false,
222  double dist_decay_min = 0.0, double dist_decay_max = 1.0);
223 
227  virtual ~G4DarkBreMModel() = default;
228 
232  virtual void PrintInfo() const;
233 
249  virtual G4double ComputeCrossSectionPerAtom(G4double lepton_ke,
250  G4double atomicA,
251  G4double atomicZ);
252 
278  std::pair<G4ThreeVector, G4ThreeVector> scale(double target_Z,
279  double incident_energy,
280  double lepton_mass);
281 
298  virtual void GenerateChange(G4ParticleChange& particleChange,
299  const G4Track& track, const G4Step& step,
300  const G4Element& element);
301 
302  private:
314  void SetMadGraphDataLibrary(const std::string& path);
315 
332  void MakePlaceholders();
333 
347  OutgoingKinematics sample(double target_Z, double incident_energy);
348 
349  private:
357  unsigned int maxIterations_{10000};
358 
365  double threshold_;
366 
374  double epsilon_;
375 
384 
389 
394 
398  std::string method_name_;
399 
403  std::string library_path_;
404 
414 
437  bool scale_APrime_{false};
438 
443  double dist_decay_min_{0.0};
444 
449  double dist_decay_max_{1.0};
450 
458  std::map<int, std::map<double, std::vector<OutgoingKinematics>>>
460 
470  std::map<int, std::map<double, unsigned int>> currentDataPoints_;
471 };
472 
473 } // namespace g4db
474 
475 #endif // G4DARKBREM_G4DARKBREMMODEL_H
Declaration of library parsing function.
Class providing a prototype model for dark brem.
Geant4 implementation of the model for a particle undergoing a dark brem where we use an imported eve...
Definition: G4DarkBreMModel.h:24
double threshold_
Threshold for non-zero xsec [GeV].
Definition: G4DarkBreMModel.h:365
ScalingMethod
Possible methods to use the dark brem events from the imported library inside of this model.
Definition: G4DarkBreMModel.h:32
@ CMScaling
Boost LHE vertex momenta to the actual lepton energy.
@ ForwardOnly
Use actual lepton energy and get pT from LHE (such that )
@ Undefined
Use LHE vertex as is.
virtual void PrintInfo() const
Print the configuration of this model.
Definition: G4DarkBreMModel.cxx:179
std::string method_name_
Name of method for persisting into the RunHeader.
Definition: G4DarkBreMModel.h:398
OutgoingKinematics sample(double target_Z, double incident_energy)
Returns MadGraph data given an energy [GeV].
Definition: G4DarkBreMModel.cxx:652
int aprime_lhe_id_
PDG ID number for the A' (dark photon) as written in the LHE files being loaded as the dark brem even...
Definition: G4DarkBreMModel.h:383
std::map< int, std::map< double, std::vector< OutgoingKinematics > > > madGraphData_
Storage of data from mad graph.
Definition: G4DarkBreMModel.h:459
bool alwaysCreateNewLepton_
should we always create a totally new lepton when we dark brem?
Definition: G4DarkBreMModel.h:413
void SetMadGraphDataLibrary(const std::string &path)
Set the library of dark brem events to be scaled.
Definition: G4DarkBreMModel.cxx:597
XsecMethod xsec_method_
method for calculating total cross section
Definition: G4DarkBreMModel.h:393
double epsilon_
Epsilon value to plug into xsec calculation.
Definition: G4DarkBreMModel.h:374
bool scale_APrime_
whether to scale the outgoing A' momentum based on the MadGraph library.
Definition: G4DarkBreMModel.h:437
XsecMethod
Possible methods for calculating the xsec using the WW approximation technique.
Definition: G4DarkBreMModel.h:107
@ Auto
the default cross section calculation
@ Full
use the full WW approximation
@ Improved
assume in the integrand
@ HyperImproved
only calculate the flux factor once
virtual void GenerateChange(G4ParticleChange &particleChange, const G4Track &track, const G4Step &step, const G4Element &element)
Simulates the emission of a dark photon + lepton.
Definition: G4DarkBreMModel.cxx:536
void MakePlaceholders()
Fill vector of currentDataPoints_ with the same number of items as the madgraph data.
Definition: G4DarkBreMModel.cxx:639
std::string library_path_
Full path to the vertex library used for persisting into the RunHeader.
Definition: G4DarkBreMModel.h:403
virtual ~G4DarkBreMModel()=default
Destructor.
unsigned int maxIterations_
maximum number of iterations to check before giving up on an event
Definition: G4DarkBreMModel.h:357
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.
Definition: G4DarkBreMModel.cxx:145
std::map< int, std::map< double, unsigned int > > currentDataPoints_
Stores a map of current access points to mad graph data.
Definition: G4DarkBreMModel.h:470
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.
Definition: G4DarkBreMModel.cxx:404
double dist_decay_max_
Maximum flight distance [mm] at which to decay the A' if G4APrime::DecayMode is set to FlatDecay.
Definition: G4DarkBreMModel.h:449
double dist_decay_min_
Minimum flight distance [mm] at which to decay the A' if G4APrime::DecayMode is set to FlatDecay.
Definition: G4DarkBreMModel.h:443
virtual G4double ComputeCrossSectionPerAtom(G4double lepton_ke, G4double atomicA, G4double atomicZ)
Calculates the cross section per atom in GEANT4 internal units.
Definition: G4DarkBreMModel.cxx:204
ScalingMethod scaling_method_
method to scale events for this model
Definition: G4DarkBreMModel.h:388
Abstract class representing a model for dark brem.
Definition: PrototypeModel.h:30
G4DarkBreM internal namespace.
Definition: ElementXsecCache.cxx:3
Data frame to store necessary information from LHE files.
Definition: ParseLibrary.h:21