G4DarkBreM  v2.1.0
Geant4 Dark Bremmstrahlung from MadGraph
PrototypeModel.h
Go to the documentation of this file.
1 
8 #ifndef G4DARKBREM_PROTOTYPEMODEL_H
9 #define G4DARKBREM_PROTOTYPEMODEL_H
10 
11 #include "G4ParticleChange.hh"
12 #include "G4Step.hh"
13 #include "G4Track.hh"
14 
18 namespace g4db {
19 
31  public:
39  PrototypeModel(bool muons) : muons_{muons} {}
40 
42  virtual ~PrototypeModel() = default;
43 
54  void SetVerboseLevel(int l) { verbose_level_ = l; }
55 
64  int GetVerboseLevel() const { return verbose_level_; }
65 
73  bool DarkBremOffMuons() const { return muons_; }
74 
81  virtual void PrintInfo() const = 0;
82 
92  virtual G4double ComputeCrossSectionPerAtom(G4double leptonKE,
93  G4double atomicA,
94  G4double atomicZ) = 0;
95 
109  virtual void GenerateChange(G4ParticleChange& particleChange,
110  const G4Track& track, const G4Step& step,
111  const G4Element& element) = 0;
112 
113  protected:
115  bool muons_;
118 }; // PrototypeModel
119 
120 } // namespace g4db
121 
122 #endif
Abstract class representing a model for dark brem.
Definition: PrototypeModel.h:30
void SetVerboseLevel(int l)
Define the verbose level.
Definition: PrototypeModel.h:54
bool muons_
whether muons (true) or electrons (false) are dark bremming
Definition: PrototypeModel.h:115
virtual void PrintInfo() const =0
Print the configuration of this model.
bool DarkBremOffMuons() const
Check if we are bremming off muons.
Definition: PrototypeModel.h:73
int verbose_level_
verbose level for this model
Definition: PrototypeModel.h:117
virtual void GenerateChange(G4ParticleChange &particleChange, const G4Track &track, const G4Step &step, const G4Element &element)=0
Generate the change in the particle now that we can assume the interaction is occuring.
virtual ~PrototypeModel()=default
Destructor, nothing on purpose.
PrototypeModel(bool muons)
Constructor.
Definition: PrototypeModel.h:39
int GetVerboseLevel() const
Get the verbose level of the model.
Definition: PrototypeModel.h:64
virtual G4double ComputeCrossSectionPerAtom(G4double leptonKE, G4double atomicA, G4double atomicZ)=0
Calculate the cross section given the input parameters.
G4DarkBreM internal namespace.
Definition: ElementXsecCache.cxx:3