LDMX Software
K0LongInelastic.h
1#ifndef SIMCORE_BIASOPERATORS_K0LONGINELASTIC_H_
2#define SIMCORE_BIASOPERATORS_K0LONGINELASTIC_H_
3
4#include "SimCore/XsecBiasingOperator.h"
5
6namespace simcore {
7namespace biasoperators {
8
13 public:
20 K0LongInelastic(std::string name, const framework::config::Parameters& p);
21
23 virtual ~K0LongInelastic() = default;
24
29 G4VBiasingOperation* ProposeOccurenceBiasingOperation(
30 const G4Track* track,
31 const G4BiasingProcessInterface* callingProcess) override;
32
34 std::string getProcessToBias() const override { return "kaon0LInelastic"; }
35
37 std::string getParticleToBias() const override { return "kaon0L"; }
38
40 std::string getVolumeToBias() const override { return volume_; }
41
47 void RecordConfig(ldmx::RunHeader& header) const override {
48 header.setStringParameter("BiasOperator::K0LongInelastic::Volume", volume_);
49 header.setFloatParameter("BiasOperator::K0LongInelastic::Factor", factor_);
50 header.setFloatParameter("BiasOperator::K0LongInelastic::Threshold",
52 }
53
54 private:
56 std::string volume_;
57
59 double factor_;
60
62 double threshold_;
63};
64
65} // namespace biasoperators
66} // namespace simcore
67
68#endif // SIMCORE_BIASOPERATORS_K0LONGINELASTIC_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
void setFloatParameter(const std::string &name, float value)
Set a float parameter value.
Definition RunHeader.h:189
void setStringParameter(const std::string &name, std::string value)
Set a string parameter value.
Definition RunHeader.h:214
Our specialization of the biasing operator used with Geant4.
Bias the k0 long inelastic collisions.
virtual ~K0LongInelastic()=default
Destructor.
void RecordConfig(ldmx::RunHeader &header) const override
Record the configuration to the run header.
std::string volume_
The volume to bias in.
G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess) override
std::string getVolumeToBias() const override
Return the volume to bias in.
std::string getProcessToBias() const override
Return the process to bias.
double threshold_
Minimum kinetic energy [MeV] to allow a track to be biased.
std::string getParticleToBias() const override
Return the particle to bias.