LDMX Software
DarkBrem.h
1#ifndef SIMCORE_BIASOPERATORS_DARKBREM_H
2#define SIMCORE_BIASOPERATORS_DARKBREM_H
3
4//----------//
5// LDMX //
6//----------//
7#include "G4DarkBreM/G4DarkBremsstrahlung.h"
8#include "SimCore/XsecBiasingOperator.h"
9
10namespace simcore {
11namespace biasoperators {
12
17 public:
24 DarkBrem(std::string name, const framework::config::Parameters& p);
25
31 virtual ~DarkBrem() = default;
32
46 G4VBiasingOperation* ProposeOccurenceBiasingOperation(
47 const G4Track* track,
48 const G4BiasingProcessInterface* callingProcess) override;
49
51 std::string getProcessToBias() const override {
52 return G4DarkBremsstrahlung::PROCESS_NAME;
53 }
54
56 std::string getParticleToBias() const override { return "e-"; }
57
59 std::string getVolumeToBias() const override { return volume_; }
60
66 void RecordConfig(ldmx::RunHeader& header) const override;
67
68 protected:
95 private:
97 std::string volume_;
98
100 double factor_;
101
104
105}; // DarkBrem
106} // namespace biasoperators
107} // namespace simcore
108
109#endif // SIMCORE_DARKBREM_DARKBREMXSECBIASINGOPERATOR_H
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:27
Run-specific configuration and data stored in its own output TTree alongside the event TTree in the o...
Definition RunHeader.h:54
Our specialization of the biasing operator used with Geant4.
Bias operator for the dark brem process.
Definition DarkBrem.h:16
void RecordConfig(ldmx::RunHeader &header) const override
Record the configuration of this biasing operator into the run header.
Definition DarkBrem.cxx:44
std::string getProcessToBias() const override
Return the name of the process this operator biases.
Definition DarkBrem.h:51
std::string getVolumeToBias() const override
Return the volume this operator biases.
Definition DarkBrem.h:59
double factor_
factor we want to bias by
Definition DarkBrem.h:100
bool bias_all_
should we bias all electrons? (or only the primary)
Definition DarkBrem.h:103
std::string getParticleToBias() const override
Return the name of the particle this operator biases.
Definition DarkBrem.h:56
std::string volume_
DEBUG FUNCTION This function is called by the biasing interface class during PostStepDoIt.
Definition DarkBrem.h:97
G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess) override
Calculate the biased cross section given the input process and track.
Definition DarkBrem.cxx:19
virtual ~DarkBrem()=default
Destructor.