LDMX Software
XsecBiasingOperator.h
1#ifndef SIMCORE_XSECBIASINGOPERATOR_H_
2#define SIMCORE_XSECBIASINGOPERATOR_H_
3
4#include "Framework/Configure/Parameters.h"
6#include "Framework/Factory.h"
7#include "Framework/Logger.h"
8#include "Framework/RunHeader.h"
9
10//------------//
11// Geant4 //
12//------------//
13#include "G4BOptnChangeCrossSection.hh"
14#include "G4BiasingProcessInterface.hh"
15#include "G4BiasingProcessSharedData.hh"
16#include "G4Electron.hh"
17#include "G4Gamma.hh"
18#include "G4KaonZeroLong.hh"
19#include "G4Neutron.hh"
20#include "G4ParticleDefinition.hh"
21#include "G4ParticleTable.hh"
22#include "G4ProcessManager.hh"
23#include "G4RunManager.hh"
24#include "G4Track.hh"
25#include "G4VBiasingOperator.hh"
26
27namespace simcore {
28
41class XsecBiasingOperator : public G4VBiasingOperator {
42 public:
53 XsecBiasingOperator(std::string name,
54 const framework::config::Parameters& parameters);
55
60 std::shared_ptr<XsecBiasingOperator>,
61 std::string,
63
66
80 virtual G4VBiasingOperation* ProposeOccurenceBiasingOperation(
81 const G4Track* track,
82 const G4BiasingProcessInterface* callingProcess) = 0;
83
94 void StartRun();
95
102 virtual std::string getProcessToBias() const = 0;
103
111 virtual std::string getParticleToBias() const = 0;
112
120 virtual std::string getVolumeToBias() const = 0;
121
128 // NOLINTNEXTLINE
129 virtual void RecordConfig(ldmx::RunHeader& header) const = 0;
130
131 protected:
146 // NOLINTNEXTLINE
147 G4VBiasingOperation* BiasedXsec(double biased_xsec) {
148 xsec_operation_->SetBiasedCrossSection(biased_xsec);
149 xsec_operation_->Sample();
150 return xsec_operation_;
151 }
152
159 bool processIsBiased(std::string process);
160
162 G4BOptnChangeCrossSection* xsec_operation_{nullptr};
163
165 G4ProcessManager* process_manager_{nullptr};
166
171 const G4Track*, const G4BiasingProcessInterface*) {
172 return nullptr;
173 }
174
179 const G4Track*, const G4BiasingProcessInterface*) {
180 return nullptr;
181 }
183 framework::logging::logger the_log_;
184
185}; // XsecBiasingOperator
186} // namespace simcore
187
194#define DECLARE_XSECBIASINGOPERATOR(CLASS) \
195 FACTORY_REGISTRATION(simcore::XsecBiasingOperator, CLASS)
196
197#endif // SIMCORE_XSECBIASINGOPERATOR_H_
Base classes for all user event processing components to extend.
Header holding Factory class and supporting macros.
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
Our specialization of the biasing operator used with Geant4.
G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *, const G4BiasingProcessInterface *)
Do not propose any non-physics biasing.
DECLARE_FACTORY_WITH_WAREHOUSE(XsecBiasingOperator, std::shared_ptr< XsecBiasingOperator >, std::string, const framework::config::Parameters &)
The BiasingOperator factory.
virtual ~XsecBiasingOperator()
Destructor.
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.
XsecBiasingOperator(std::string name, const framework::config::Parameters &parameters)
Constructor.
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.
G4ProcessManager * process_manager_
Process manager associated with the particle of interest.
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.
G4BOptnChangeCrossSection * xsec_operation_
Cross-section biasing operation.
virtual std::string getVolumeToBias() const =0
Return the volume which should be biased.
framework::logging::logger the_log_
Enable logging.
Dynamically loadable photonuclear models either from SimCore or external libraries implementing this ...