1#ifndef SIMCORE_XSECBIASINGOPERATOR_H_
2#define SIMCORE_XSECBIASINGOPERATOR_H_
4#include "Framework/Configure/Parameters.h"
5#include "Framework/RunHeader.h"
6#include "SimCore/Factory.h"
11#include "G4BOptnChangeCrossSection.hh"
12#include "G4BiasingProcessInterface.hh"
13#include "G4BiasingProcessSharedData.hh"
14#include "G4Electron.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"
23#include "G4VBiasingOperator.hh"
59 std::shared_ptr<XsecBiasingOperator>, std::string,
80 const G4BiasingProcessInterface* callingProcess) = 0;
167 const G4Track*,
const G4BiasingProcessInterface*) {
175 const G4Track*,
const G4BiasingProcessInterface*) {
188#define DECLARE_XSECBIASINGOPERATOR(CLASS) \
190 auto v = ::simcore::XsecBiasingOperator::Factory::get().declare<CLASS>(); \
Class encapsulating parameters for configuring a processor.
Factory to dynamically create objects derived from a specific prototype class.
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.