|
LDMX Software
|
Our specialization of the biasing operator used with Geant4. More...
#include <XsecBiasingOperator.h>
Public Member Functions | |
| XsecBiasingOperator (std::string name, const framework::config::Parameters ¶meters) | |
| Constructor. | |
| DECLARE_FACTORY_WITH_WAREHOUSE (XsecBiasingOperator, std::shared_ptr< XsecBiasingOperator >, std::string, const framework::config::Parameters &) | |
| The BiasingOperator factory. | |
| virtual | ~XsecBiasingOperator ()=default |
| Destructor. | |
| virtual G4VBiasingOperation * | ProposeOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0 |
| Propose a biasing operation for the current track and calling process. | |
| void | StartRun () |
| Method called at the beginning of a run. | |
| virtual std::string | getProcessToBias () const =0 |
| Return the process whose cross-section will be biased. | |
| virtual std::string | getParticleToBias () const =0 |
| Return the particle which should be biased. | |
| virtual std::string | getVolumeToBias () const =0 |
| Return the volume which should be biased. | |
| virtual void | RecordConfig (ldmx::RunHeader &header) const =0 |
| Record the configuration of this biasing operator into the run header. | |
Protected Member Functions | |
| G4VBiasingOperation * | BiasedXsec (double biased_xsec) |
| Helper method for passing a biased interaction length to the Geant4 biasing framework. | |
| bool | processIsBiased (std::string process) |
| Check if the given processed is being biased. | |
| G4VBiasingOperation * | ProposeFinalStateBiasingOperation (const G4Track *, const G4BiasingProcessInterface *) |
| Do not propose any biasing on final states. | |
| G4VBiasingOperation * | ProposeNonPhysicsBiasingOperation (const G4Track *, const G4BiasingProcessInterface *) |
| Do not propose any non-physics biasing. | |
Protected Attributes | |
| G4BOptnChangeCrossSection * | xsec_operation_ {nullptr} |
| Cross-section biasing operation. | |
| G4ProcessManager * | process_manager_ {nullptr} |
| Process manager associated with the particle of interest. | |
| framework::logging::logger | the_log_ |
| Enable logging. | |
Our specialization of the biasing operator used with Geant4.
This specialization accomplishes three main tasks.
Definition at line 41 of file XsecBiasingOperator.h.
| simcore::XsecBiasingOperator::XsecBiasingOperator | ( | std::string | name, |
| const framework::config::Parameters & | parameters ) |
Constructor.
Here, we define a unique name for this biasing operator and are given the configuration parameters loaded from the python script.
| [in] | name | unique instance name for this biasing operator |
| [in] | parameters | python configuration parameters |
Definition at line 7 of file XsecBiasingOperator.cxx.
|
inlineprotected |
Helper method for passing a biased interaction length to the Geant4 biasing framework.
Use like:
return BiasedXsec(biased_xsec);
inside of ProposeOccurenceBiasingOperation when you want to update the biased cross section.
| [in] | biased_xsec | the biased cross section |
Definition at line 145 of file XsecBiasingOperator.h.
References xsec_operation_.
Referenced by simcore::biasoperators::DarkBrem::ProposeOccurenceBiasingOperation(), simcore::biasoperators::ElectroNuclear::ProposeOccurenceBiasingOperation(), simcore::biasoperators::GammaToMuPair::ProposeOccurenceBiasingOperation(), simcore::biasoperators::K0LongInelastic::ProposeOccurenceBiasingOperation(), simcore::biasoperators::NeutronInelastic::ProposeOccurenceBiasingOperation(), and simcore::biasoperators::PhotoNuclear::ProposeOccurenceBiasingOperation().
|
pure virtual |
Return the particle which should be biased.
We need this to be able to tell the physics list which particle to bias.
Implemented in simcore::biasoperators::DarkBrem, simcore::biasoperators::ElectroNuclear, simcore::biasoperators::GammaToMuPair, simcore::biasoperators::K0LongInelastic, simcore::biasoperators::NeutronInelastic, and simcore::biasoperators::PhotoNuclear.
Referenced by StartRun().
|
pure virtual |
Return the process whose cross-section will be biased.
We need this to be able to check that the process was biased before creating the biasing operator.
Implemented in simcore::biasoperators::DarkBrem, simcore::biasoperators::ElectroNuclear, simcore::biasoperators::GammaToMuPair, simcore::biasoperators::K0LongInelastic, simcore::biasoperators::NeutronInelastic, and simcore::biasoperators::PhotoNuclear.
Referenced by StartRun().
|
pure virtual |
Return the volume which should be biased.
We need this to be able to tell the detector construction which volumes to attach this operator to.
Implemented in simcore::biasoperators::DarkBrem, simcore::biasoperators::ElectroNuclear, simcore::biasoperators::GammaToMuPair, simcore::biasoperators::K0LongInelastic, simcore::biasoperators::NeutronInelastic, and simcore::biasoperators::PhotoNuclear.
|
protected |
Check if the given processed is being biased.
| process | Process of interest |
Definition at line 39 of file XsecBiasingOperator.cxx.
References process_manager_.
Referenced by simcore::biasoperators::PhotoNuclear::StartRun(), and StartRun().
|
inlineprotected |
Do not propose any biasing on final states.
Definition at line 168 of file XsecBiasingOperator.h.
|
inlineprotected |
Do not propose any non-physics biasing.
Definition at line 176 of file XsecBiasingOperator.h.
|
pure virtual |
Propose a biasing operation for the current track and calling process.
0 from this function will mean that the current track and process will not be biased.| [in] | track | handle to current track that could be biased |
| [in] | callingProcess | handle to process asking if it should be biased |
Implemented in simcore::biasoperators::DarkBrem, simcore::biasoperators::ElectroNuclear, simcore::biasoperators::GammaToMuPair, simcore::biasoperators::K0LongInelastic, simcore::biasoperators::NeutronInelastic, and simcore::biasoperators::PhotoNuclear.
|
pure virtual |
Record the configuration of this biasing operator into the run header.
| [in,out] | header | RunHeader to write configuration to |
Implemented in simcore::biasoperators::DarkBrem, simcore::biasoperators::ElectroNuclear, simcore::biasoperators::GammaToMuPair, simcore::biasoperators::K0LongInelastic, simcore::biasoperators::NeutronInelastic, and simcore::biasoperators::PhotoNuclear.
| void simcore::XsecBiasingOperator::StartRun | ( | ) |
Method called at the beginning of a run.
This makes sure that the process we want to bias can be biased and constructs a corresponding biasing operation.
It can be over-written, but then the derived class should call XsecBiasingOperator::StartRun() at the beginning of their own StartRun.
Definition at line 12 of file XsecBiasingOperator.cxx.
References getParticleToBias(), getProcessToBias(), process_manager_, processIsBiased(), and xsec_operation_.
Referenced by simcore::biasoperators::PhotoNuclear::StartRun().
|
protected |
Process manager associated with the particle of interest.
Definition at line 163 of file XsecBiasingOperator.h.
Referenced by processIsBiased(), and StartRun().
|
protected |
Enable logging.
Definition at line 181 of file XsecBiasingOperator.h.
|
protected |
Cross-section biasing operation.
Definition at line 160 of file XsecBiasingOperator.h.
Referenced by BiasedXsec(), and StartRun().