LDMX Software
Public Member Functions | Private Attributes | List of all members
simcore::biasoperators::ElectroNuclear Class Reference

Bias the Electron-Nuclear process. More...

#include <ElectroNuclear.h>

Public Member Functions

 ElectroNuclear (std::string name, const framework::config::Parameters &p)
 Constructor.
 
virtual ~ElectroNuclear ()=default
 Destructor.
 
G4VBiasingOperation * ProposeOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess) override
 
std::string getProcessToBias () const override
 Return the process to bias.
 
std::string getParticleToBias () const override
 Return the particle to bias.
 
std::string getVolumeToBias () const override
 Return the volume to bias in.
 
void RecordConfig (ldmx::RunHeader &header) const override
 Record the configuration to the run header.
 
- Public Member Functions inherited from simcore::XsecBiasingOperator
 XsecBiasingOperator (std::string name, const framework::config::Parameters &parameters)
 Constructor.
 
virtual ~XsecBiasingOperator ()=default
 Destructor.
 
void StartRun ()
 Method called at the beginning of a run.
 

Private Attributes

std::string volume_
 The volume to bias in.
 
double factor_
 The biasing factor.
 
double threshold_
 Minimum kinetic energy [MeV] to allow a track to be biased.
 

Additional Inherited Members

- Public Types inherited from simcore::XsecBiasingOperator
using Factory = ::simcore::Factory< XsecBiasingOperator, std::shared_ptr< XsecBiasingOperator >, std::string, const framework::config::Parameters & >
 The BiasingOperator factory.
 
- Protected Member Functions inherited from simcore::XsecBiasingOperator
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 inherited from simcore::XsecBiasingOperator
G4BOptnChangeCrossSection * xsecOperation_ {nullptr}
 Cross-section biasing operation.
 
G4ProcessManager * processManager_ {nullptr}
 Process manager associated with the particle of interest.
 

Detailed Description

Bias the Electron-Nuclear process.

Definition at line 12 of file ElectroNuclear.h.

Constructor & Destructor Documentation

◆ ElectroNuclear()

simcore::biasoperators::ElectroNuclear::ElectroNuclear ( std::string  name,
const framework::config::Parameters p 
)

Constructor.

Calls parent constructor and allows accesss to configuration parameters.

Definition at line 7 of file ElectroNuclear.cxx.

9 : XsecBiasingOperator(name, p) {
10 volume_ = p.getParameter<std::string>("volume");
11 factor_ = p.getParameter<double>("factor");
12 threshold_ = p.getParameter<double>("threshold");
13}
T getParameter(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:89
XsecBiasingOperator(std::string name, const framework::config::Parameters &parameters)
Constructor.
std::string volume_
The volume to bias in.
double threshold_
Minimum kinetic energy [MeV] to allow a track to be biased.
double factor_
The biasing factor.

References factor_, framework::config::Parameters::getParameter(), threshold_, and volume_.

Member Function Documentation

◆ getParticleToBias()

std::string simcore::biasoperators::ElectroNuclear::getParticleToBias ( ) const
inlineoverridevirtual

Return the particle to bias.

Implements simcore::XsecBiasingOperator.

Definition at line 37 of file ElectroNuclear.h.

37{ return "e-"; }

◆ getProcessToBias()

std::string simcore::biasoperators::ElectroNuclear::getProcessToBias ( ) const
inlineoverridevirtual

Return the process to bias.

Implements simcore::XsecBiasingOperator.

Definition at line 34 of file ElectroNuclear.h.

34{ return "electronNuclear"; }

◆ getVolumeToBias()

std::string simcore::biasoperators::ElectroNuclear::getVolumeToBias ( ) const
inlineoverridevirtual

Return the volume to bias in.

Implements simcore::XsecBiasingOperator.

Definition at line 40 of file ElectroNuclear.h.

40{ return volume_; }

References volume_.

◆ ProposeOccurenceBiasingOperation()

G4VBiasingOperation * simcore::biasoperators::ElectroNuclear::ProposeOccurenceBiasingOperation ( const G4Track *  track,
const G4BiasingProcessInterface *  callingProcess 
)
overridevirtual
Returns
Method that returns the biasing operation that will be used to bias the occurence of photonuclear events.

Implements simcore::XsecBiasingOperator.

Definition at line 15 of file ElectroNuclear.cxx.

16 {
17 if (track->GetKineticEnergy() < threshold_) {
18 return nullptr;
19 }
20
21 std::string currentProcess =
22 callingProcess->GetWrappedProcess()->GetProcessName();
23 if (currentProcess.compare(this->getProcessToBias()) == 0) {
24 G4double interactionLength =
25 callingProcess->GetWrappedProcess()->GetCurrentInteractionLength();
26 double enXsecUnbiased = 1. / interactionLength;
27
28 double enXsecBiased = enXsecUnbiased * factor_;
29
30 return BiasedXsec(enXsecBiased);
31 }
32 return nullptr;
33}
G4VBiasingOperation * BiasedXsec(double biased_xsec)
Helper method for passing a biased interaction length to the Geant4 biasing framework.

References simcore::XsecBiasingOperator::BiasedXsec(), factor_, and threshold_.

◆ RecordConfig()

void simcore::biasoperators::ElectroNuclear::RecordConfig ( ldmx::RunHeader header) const
inlineoverridevirtual

Record the configuration to the run header.

Parameters
[in,out]headerRunHeader to record to

Implements simcore::XsecBiasingOperator.

Definition at line 47 of file ElectroNuclear.h.

47 {
48 header.setStringParameter("BiasOperator::ElectroNuclear::Volume", volume_);
49 header.setFloatParameter("BiasOperator::ElectroNuclear::Factor", factor_);
50 header.setFloatParameter("BiasOperator::ElectroNuclear::Threshold",
52 }
void setFloatParameter(const std::string &name, float value)
Set a float parameter value.
Definition RunHeader.h:189
void setStringParameter(const std::string &name, std::string value)
Set a string parameter value.
Definition RunHeader.h:214

References factor_, ldmx::RunHeader::setFloatParameter(), ldmx::RunHeader::setStringParameter(), threshold_, and volume_.

Member Data Documentation

◆ factor_

double simcore::biasoperators::ElectroNuclear::factor_
private

The biasing factor.

Definition at line 59 of file ElectroNuclear.h.

Referenced by ElectroNuclear(), ProposeOccurenceBiasingOperation(), and RecordConfig().

◆ threshold_

double simcore::biasoperators::ElectroNuclear::threshold_
private

Minimum kinetic energy [MeV] to allow a track to be biased.

Definition at line 62 of file ElectroNuclear.h.

Referenced by ElectroNuclear(), ProposeOccurenceBiasingOperation(), and RecordConfig().

◆ volume_

std::string simcore::biasoperators::ElectroNuclear::volume_
private

The volume to bias in.

Definition at line 56 of file ElectroNuclear.h.

Referenced by ElectroNuclear(), getVolumeToBias(), and RecordConfig().


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