LDMX Software
PhotoNuclear.h
1
2#ifndef SIMCORE_BIASOPERATORS_PHOTONUCLEAR_H_
3#define SIMCORE_BIASOPERATORS_PHOTONUCLEAR_H_
4
5#include "SimCore/BiasOperators/XsecBiasingOperator.h"
6#include "SimCore/G4User/PtrRetrieval.h"
7
8namespace simcore {
9namespace biasoperators {
10
15 public:
17 PhotoNuclear(std::string name, const framework::config::Parameters& p);
18
20 ~PhotoNuclear() override { delete em_xsec_operation_; }
21
23 void StartRun() override;
24
29 G4VBiasingOperation* ProposeOccurenceBiasingOperation(
30 const G4Track* track,
31 const G4BiasingProcessInterface* callingProcess) override;
32
34 std::string getProcessToBias() const override { return "photonNuclear"; }
35
37 std::string getParticleToBias() const override { return "gamma"; }
38
40 std::string getVolumeToBias() const override { return volume_; }
41
43 void RecordConfig(ldmx::RunHeader& h) const override {
44 h.setStringParameter("BiasOperators::PhotoNuclear::Volume", volume_);
45 h.setFloatParameter("BiasOperators::PhotoNuclear::Threshold", threshold_);
46 h.setFloatParameter("BiasOperators::PhotoNuclear::Factor", factor_);
47 h.setIntParameter("BiasOperators::PhotoNuclear::Bias Conv Down",
49 h.setIntParameter("BiasOperators::PhotoNuclear::Only Children Of Primary",
51 }
52
53 private:
55 static const std::string CONVERSION_PROCESS;
56
58 G4BOptnChangeCrossSection* em_xsec_operation_{nullptr};
59
62
64 double pn_xsec_biased_{0};
65
67 std::string volume_;
68
70 double threshold_;
71
73 double factor_;
74
77
80
81}; // PhotoNuclear
82} // namespace biasoperators
83} // namespace simcore
84
85#endif // SIMCORE_BIASOPERATORS_PHOTONUCLEAR_H_
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:29
Run-specific configuration and data stored in its own output TTree alongside the event TTree in the o...
Definition RunHeader.h:57
void setFloatParameter(const std::string &name, float value)
Set a float parameter value.
Definition RunHeader.h:197
void setStringParameter(const std::string &name, std::string value)
Set a string parameter value.
Definition RunHeader.h:222
void setIntParameter(const std::string &name, int value)
Set an int parameter value.
Definition RunHeader.h:172
Our specialization of the biasing operator used with Geant4.
Bias the Photon-Nuclear process.
double factor_
factor to bias PN by
bool only_children_of_primary_
Should we restrict biasing to only children of primary?
~PhotoNuclear() override
Destructor.
std::string getVolumeToBias() const override
return the volume we want to bias within
bool down_bias_conv_
Should we down-bias the gamma conversion process?
G4BOptnChangeCrossSection * em_xsec_operation_
Cross-section biasing operation for conversion process.
std::string getProcessToBias() const override
return the process we want to bias
std::string getParticleToBias() const override
return the particle that we want to bias
G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess) override
void StartRun() override
Method called at the beginning of a run.
double pn_xsec_unbiased_
Unbiased photonuclear xsec.
double threshold_
minimum kinetic energy [MeV] for a track to be biased
static const std::string CONVERSION_PROCESS
Geant4 gamma conversion process name.
PhotoNuclear(std::string name, const framework::config::Parameters &p)
Constructor.
double pn_xsec_biased_
Biased photonuclear xsec.
std::string volume_
Volume we are going to bias within.
void RecordConfig(ldmx::RunHeader &h) const override
record the configuration into the run header
Dynamically loadable photonuclear models either from SimCore or external libraries implementing this ...