LDMX Software
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
simcore::XsecBiasingOperator Class Referenceabstract

Our specialization of the biasing operator used with Geant4. More...

#include <XsecBiasingOperator.h>

Public Types

using Factory = ::simcore::Factory< XsecBiasingOperator, std::shared_ptr< XsecBiasingOperator >, std::string, const framework::config::Parameters & >
 The BiasingOperator factory.
 

Public Member Functions

 XsecBiasingOperator (std::string name, const framework::config::Parameters &parameters)
 Constructor.
 
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 * xsecOperation_ {nullptr}
 Cross-section biasing operation.
 
G4ProcessManager * processManager_ {nullptr}
 Process manager associated with the particle of interest.
 

Detailed Description

Our specialization of the biasing operator used with Geant4.

This specialization accomplishes three main tasks.

  1. Allows any derived class to be dynamically loaded after using the declaration macro given below.
  2. Interfaces with the derived class using our parameters class.
  3. Pre-defines the necessary biasing operation so the derived class only needs to worry about calculating the biased xsec.

Definition at line 39 of file XsecBiasingOperator.h.

Member Typedef Documentation

◆ Factory

The BiasingOperator factory.

Definition at line 57 of file XsecBiasingOperator.h.

Constructor & Destructor Documentation

◆ XsecBiasingOperator()

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.

Parameters
[in]nameunique instance name for this biasing operator
[in]parameterspython configuration parameters

Definition at line 7 of file XsecBiasingOperator.cxx.

9 : G4VBiasingOperator(name) {}

Member Function Documentation

◆ BiasedXsec()

G4VBiasingOperation * simcore::XsecBiasingOperator::BiasedXsec ( double  biased_xsec)
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.

Parameters
[in]biased_xsecthe biased cross section
Returns
the biasing operation with the input biased cross section

Definition at line 143 of file XsecBiasingOperator.h.

143 {
144 xsecOperation_->SetBiasedCrossSection(biased_xsec);
145 xsecOperation_->Sample();
146 return xsecOperation_;
147 }
G4BOptnChangeCrossSection * xsecOperation_
Cross-section biasing operation.

References xsecOperation_.

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().

◆ getParticleToBias()

virtual std::string simcore::XsecBiasingOperator::getParticleToBias ( ) const
pure virtual

Return the particle which should be biased.

We need this to be able to tell the physics list which particle to bias.

See also
RunManager::setupPhysics

Implemented in simcore::biasoperators::DarkBrem, simcore::biasoperators::ElectroNuclear, simcore::biasoperators::GammaToMuPair, simcore::biasoperators::K0LongInelastic, simcore::biasoperators::NeutronInelastic, and simcore::biasoperators::PhotoNuclear.

Referenced by StartRun().

◆ getProcessToBias()

virtual std::string simcore::XsecBiasingOperator::getProcessToBias ( ) const
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().

◆ getVolumeToBias()

virtual std::string simcore::XsecBiasingOperator::getVolumeToBias ( ) const
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.

◆ processIsBiased()

bool simcore::XsecBiasingOperator::processIsBiased ( std::string  process)
protected

Check if the given processed is being biased.

Parameters
processProcess of interest
Returns
true if the process is being biased, false otherwise

Definition at line 39 of file XsecBiasingOperator.cxx.

39 {
40 // Loop over all processes and check if the given process is being
41 // biased.
42 const G4BiasingProcessSharedData* sharedData =
43 G4BiasingProcessInterface::GetSharedData(processManager_);
44 if (sharedData) {
45 for (size_t iprocess = 0;
46 iprocess < (sharedData->GetPhysicsBiasingProcessInterfaces()).size();
47 ++iprocess) {
48 const G4BiasingProcessInterface* wrapperProcess =
49 (sharedData->GetPhysicsBiasingProcessInterfaces())[iprocess];
50
51 if (wrapperProcess->GetWrappedProcess()->GetProcessName().compareTo(
52 process) == 0) {
53 return true;
54 }
55 }
56 }
57 return false;
58}
G4ProcessManager * processManager_
Process manager associated with the particle of interest.

References processManager_.

Referenced by StartRun(), and simcore::biasoperators::PhotoNuclear::StartRun().

◆ ProposeFinalStateBiasingOperation()

G4VBiasingOperation * simcore::XsecBiasingOperator::ProposeFinalStateBiasingOperation ( const G4Track *  ,
const G4BiasingProcessInterface *   
)
inlineprotected

Do not propose any biasing on final states.

Definition at line 166 of file XsecBiasingOperator.h.

167 {
168 return nullptr;
169 }

◆ ProposeNonPhysicsBiasingOperation()

G4VBiasingOperation * simcore::XsecBiasingOperator::ProposeNonPhysicsBiasingOperation ( const G4Track *  ,
const G4BiasingProcessInterface *   
)
inlineprotected

Do not propose any non-physics biasing.

Definition at line 174 of file XsecBiasingOperator.h.

175 {
176 return nullptr;
177 }

◆ ProposeOccurenceBiasingOperation()

virtual G4VBiasingOperation * simcore::XsecBiasingOperator::ProposeOccurenceBiasingOperation ( const G4Track *  track,
const G4BiasingProcessInterface *  callingProcess 
)
pure virtual

Propose a biasing operation for the current track and calling process.

Note
Returning 0 from this function will mean that the current track and process will not be biased.
See also
BiasedXsec for a method that allows the derived class to not interact with the biasing operation itself.
Parameters
[in]trackhandle to current track that could be biased
[in]callingProcesshandle to process asking if it should be biased
Returns
the biasing operation with the biased xsec

Implemented in simcore::biasoperators::DarkBrem, simcore::biasoperators::ElectroNuclear, simcore::biasoperators::GammaToMuPair, simcore::biasoperators::K0LongInelastic, simcore::biasoperators::NeutronInelastic, and simcore::biasoperators::PhotoNuclear.

◆ RecordConfig()

virtual void simcore::XsecBiasingOperator::RecordConfig ( ldmx::RunHeader header) const
pure virtual

Record the configuration of this biasing operator into the run header.

Parameters
[in,out]headerRunHeader to write configuration to

Implemented in simcore::biasoperators::PhotoNuclear, simcore::biasoperators::DarkBrem, simcore::biasoperators::ElectroNuclear, simcore::biasoperators::GammaToMuPair, simcore::biasoperators::K0LongInelastic, and simcore::biasoperators::NeutronInelastic.

◆ StartRun()

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 11 of file XsecBiasingOperator.cxx.

11 {
12 if (this->getParticleToBias().compare("gamma") == 0) {
13 processManager_ = G4Gamma::GammaDefinition()->GetProcessManager();
14 } else if (this->getParticleToBias().compare("e-") == 0) {
15 processManager_ = G4Electron::ElectronDefinition()->GetProcessManager();
16 } else if (this->getParticleToBias().compare("neutron") == 0) {
17 processManager_ = G4Neutron::NeutronDefinition()->GetProcessManager();
18 } else if (this->getParticleToBias().compare("kaon0L") == 0) {
20 G4KaonZeroLong::KaonZeroLongDefinition()->GetProcessManager();
21 } else {
22 EXCEPTION_RAISE("BiasSetup", "Invalid particle type '" +
23 this->getParticleToBias() + "'.");
24 }
25
26 std::cout << "[ XsecBiasingOperator ]: Biasing particles of type "
27 << this->getParticleToBias() << std::endl;
28
29 if (processIsBiased(this->getProcessToBias())) {
31 new G4BOptnChangeCrossSection("changeXsec-" + this->getProcessToBias());
32 } else {
33 EXCEPTION_RAISE("BiasSetup",
34 this->getProcessToBias() +
35 " is not found in list of biased processes!");
36 }
37}
virtual std::string getProcessToBias() const =0
Return the process whose cross-section will be biased.
bool processIsBiased(std::string process)
Check if the given processed is being biased.
virtual std::string getParticleToBias() const =0
Return the particle which should be biased.

References getParticleToBias(), getProcessToBias(), processIsBiased(), processManager_, and xsecOperation_.

Referenced by simcore::biasoperators::PhotoNuclear::StartRun().

Member Data Documentation

◆ processManager_

G4ProcessManager* simcore::XsecBiasingOperator::processManager_ {nullptr}
protected

Process manager associated with the particle of interest.

Definition at line 161 of file XsecBiasingOperator.h.

161{nullptr};

Referenced by processIsBiased(), and StartRun().

◆ xsecOperation_

G4BOptnChangeCrossSection* simcore::XsecBiasingOperator::xsecOperation_ {nullptr}
protected

Cross-section biasing operation.

Definition at line 158 of file XsecBiasingOperator.h.

158{nullptr};

Referenced by BiasedXsec(), and StartRun().


The documentation for this class was generated from the following files: