LDMX Software
XsecBiasingOperator.h
1#ifndef SIMCORE_XSECBIASINGOPERATOR_H_
2#define SIMCORE_XSECBIASINGOPERATOR_H_
3
4#include "Framework/Configure/Parameters.h"
5#include "Framework/RunHeader.h"
6#include "SimCore/Factory.h"
7
8//------------//
9// Geant4 //
10//------------//
11#include "G4BOptnChangeCrossSection.hh"
12#include "G4BiasingProcessInterface.hh"
13#include "G4BiasingProcessSharedData.hh"
14#include "G4Electron.hh"
15#include "G4Gamma.hh"
16#include "G4KaonZeroLong.hh"
17#include "G4Neutron.hh"
18#include "G4ParticleDefinition.hh"
19#include "G4ParticleTable.hh"
20#include "G4ProcessManager.hh"
21#include "G4RunManager.hh"
22#include "G4Track.hh"
23#include "G4VBiasingOperator.hh"
24
25namespace simcore {
26
39class XsecBiasingOperator : public G4VBiasingOperator {
40 public:
51 XsecBiasingOperator(std::string name,
52 const framework::config::Parameters& parameters);
53
57 using Factory =
59 std::shared_ptr<XsecBiasingOperator>, std::string,
61
63 virtual ~XsecBiasingOperator() = default;
64
78 virtual G4VBiasingOperation* ProposeOccurenceBiasingOperation(
79 const G4Track* track,
80 const G4BiasingProcessInterface* callingProcess) = 0;
81
92 void StartRun();
93
100 virtual std::string getProcessToBias() const = 0;
101
109 virtual std::string getParticleToBias() const = 0;
110
118 virtual std::string getVolumeToBias() const = 0;
119
126 virtual void RecordConfig(ldmx::RunHeader& header) const = 0;
127
128 protected:
143 G4VBiasingOperation* BiasedXsec(double biased_xsec) {
144 xsecOperation_->SetBiasedCrossSection(biased_xsec);
145 xsecOperation_->Sample();
146 return xsecOperation_;
147 }
148
155 bool processIsBiased(std::string process);
156
158 G4BOptnChangeCrossSection* xsecOperation_{nullptr};
159
161 G4ProcessManager* processManager_{nullptr};
162
167 const G4Track*, const G4BiasingProcessInterface*) {
168 return nullptr;
169 }
170
175 const G4Track*, const G4BiasingProcessInterface*) {
176 return nullptr;
177 }
178
179}; // XsecBiasingOperator
180} // namespace simcore
181
188#define DECLARE_XSECBIASINGOPERATOR(CLASS) \
189 namespace { \
190 auto v = ::simcore::XsecBiasingOperator::Factory::get().declare<CLASS>(); \
191 }
192
193#endif // SIMCORE_XSECBIASINGOPERATOR_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
Factory to dynamically create objects derived from a specific prototype class.
Definition Factory.h:197
Our specialization of the biasing operator used with Geant4.
G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *, const G4BiasingProcessInterface *)
Do not propose any non-physics biasing.
virtual std::string getProcessToBias() const =0
Return the process whose cross-section will be biased.
G4VBiasingOperation * ProposeFinalStateBiasingOperation(const G4Track *, const G4BiasingProcessInterface *)
Do not propose any biasing on final states.
void StartRun()
Method called at the beginning of a run.
virtual void RecordConfig(ldmx::RunHeader &header) const =0
Record the configuration of this biasing operator into the run header.
bool processIsBiased(std::string process)
Check if the given processed is being biased.
virtual G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
Propose a biasing operation for the current track and calling process.
G4VBiasingOperation * BiasedXsec(double biased_xsec)
Helper method for passing a biased interaction length to the Geant4 biasing framework.
virtual std::string getParticleToBias() const =0
Return the particle which should be biased.
G4ProcessManager * processManager_
Process manager associated with the particle of interest.
virtual ~XsecBiasingOperator()=default
Destructor.
virtual std::string getVolumeToBias() const =0
Return the volume which should be biased.
G4BOptnChangeCrossSection * xsecOperation_
Cross-section biasing operation.