G4DarkBreM v2.2.2
Geant4 Dark Bremmstrahlung from MadGraph
Loading...
Searching...
No Matches
G4DarkBreMModel.h
1#ifndef G4DARKBREM_G4DARKBREMMODEL_H
2#define G4DARKBREM_G4DARKBREMMODEL_H
3
4#include <map>
5#include <memory>
6
12
13// Geant4
14#include "G4Electron.hh"
15#include "G4EventManager.hh" //for EventID number
16#include "G4MuonMinus.hh"
17#include "G4PhaseSpaceDecayChannel.hh"
18#include "G4PhysicalConstants.hh"
19#include "G4RunManager.hh" //for VerboseLevel
20#include "G4SystemOfUnits.hh"
21#include "Randomize.hh"
22
23// Boost
24#include <boost/math/quadrature/gauss_kronrod.hpp>
25
26// STL
27#include <math.h>
28#include <stdio.h>
29#include <stdlib.h>
30
31namespace g4db {
32
46 public:
53 enum class ScalingMethod {
70 ForwardOnly = 1,
71
87 CMScaling = 2,
94 Undefined = 3
95 };
96
128 enum class XsecMethod {
151 Full = 1,
152
166 Improved = 2,
167
186 HyperImproved = 3,
187
201 Auto = 4
202 };
203
239 G4DarkBreMModel(const std::string& library_path, bool muons,
240 double threshold = 0.0, double epsilon = 1.0,
243 double max_R_for_full = 50.0, int aprime_lhe_id = 622,
244 bool load_library = true, bool scale_APrime = false,
245 double dist_decay_min = 0.0, double dist_decay_max = 1.0,
246 int decay_particle_id = 11);
247
251 virtual ~G4DarkBreMModel() = default;
252
256 virtual void PrintInfo() const;
257
273 virtual G4double ComputeCrossSectionPerAtom(G4double lepton_ke,
274 G4double atomicA,
275 G4double atomicZ);
276
302 std::pair<G4ThreeVector, G4ThreeVector> scale(double target_Z,
303 double incident_energy,
304 double lepton_mass);
305
322 virtual void GenerateChange(G4ParticleChange& particleChange,
323 const G4Track& track, const G4Step& step,
324 const G4Element& element);
325
326 private:
338 void SetMadGraphDataLibrary(const std::string& path);
339
356 void MakePlaceholders();
357
371 OutgoingKinematics sample(double target_Z, double incident_energy);
372
373 private:
381 unsigned int maxIterations_{10000};
382
390
398 double epsilon_;
399
408
413
418
422 std::string method_name_;
423
427 std::string library_path_;
428
438
461 bool scale_APrime_{false};
462
467 double dist_decay_min_{0.0};
468
473 double dist_decay_max_{1.0};
474
482 std::map<int, std::map<double, std::vector<OutgoingKinematics>>>
484
494 std::map<int, std::map<double, unsigned int>> currentDataPoints_;
495
503};
504
505} // namespace g4db
506
507#endif // G4DARKBREM_G4DARKBREMMODEL_H
Class creating the A' particle in Geant.
Class creating the milli-charged particle in Geant.
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:45
double threshold_
Threshold for non-zero xsec [GeV].
Definition G4DarkBreMModel.h:389
ScalingMethod
Possible methods to use the dark brem events from the imported library inside of this model.
Definition G4DarkBreMModel.h:53
@ 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.
int decay_particle_id_
PDG ID of the particle to decay to if G4APrime::DecayMode is set to FlatDecay.
Definition G4DarkBreMModel.h:502
virtual void PrintInfo() const
Print the configuration of this model.
Definition G4DarkBreMModel.cxx:166
std::string method_name_
Name of method for persisting into the RunHeader.
Definition G4DarkBreMModel.h:422
OutgoingKinematics sample(double target_Z, double incident_energy)
Returns MadGraph data given an energy [GeV].
Definition G4DarkBreMModel.cxx:649
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:407
std::map< int, std::map< double, std::vector< OutgoingKinematics > > > madGraphData_
Storage of data from mad graph.
Definition G4DarkBreMModel.h:483
bool alwaysCreateNewLepton_
should we always create a totally new lepton when we dark brem?
Definition G4DarkBreMModel.h:437
void SetMadGraphDataLibrary(const std::string &path)
Set the library of dark brem events to be scaled.
Definition G4DarkBreMModel.cxx:594
XsecMethod xsec_method_
method for calculating total cross section
Definition G4DarkBreMModel.h:417
double epsilon_
Epsilon value to plug into xsec calculation.
Definition G4DarkBreMModel.h:398
bool scale_APrime_
whether to scale the outgoing A' momentum based on the MadGraph library.
Definition G4DarkBreMModel.h:461
XsecMethod
Possible methods for calculating the xsec using the WW approximation technique.
Definition G4DarkBreMModel.h:128
@ 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:523
void MakePlaceholders()
Fill vector of currentDataPoints_ with the same number of items as the madgraph data.
Definition G4DarkBreMModel.cxx:636
std::string library_path_
Full path to the vertex library used for persisting into the RunHeader.
Definition G4DarkBreMModel.h:427
virtual ~G4DarkBreMModel()=default
Destructor.
unsigned int maxIterations_
maximum number of iterations to check before giving up on an event
Definition G4DarkBreMModel.h:381
std::map< int, std::map< double, unsigned int > > currentDataPoints_
Stores a map of current access points to mad graph data.
Definition G4DarkBreMModel.h:494
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:391
double dist_decay_max_
Maximum flight distance [mm] at which to decay the A' if G4APrime::DecayMode is set to FlatDecay.
Definition G4DarkBreMModel.h:473
double dist_decay_min_
Minimum flight distance [mm] at which to decay the A' if G4APrime::DecayMode is set to FlatDecay.
Definition G4DarkBreMModel.h:467
virtual G4double ComputeCrossSectionPerAtom(G4double lepton_ke, G4double atomicA, G4double atomicZ)
Calculates the cross section per atom in GEANT4 internal units.
Definition G4DarkBreMModel.cxx:191
ScalingMethod scaling_method_
method to scale events for this model
Definition G4DarkBreMModel.h:412
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