LDMX Software
simcore::biasoperators::APrimeToFCPPair Class Reference

Bias the A' -> fcp+ fcp- process. More...

#include <APrimeToFCPPair.h>

Public Member Functions

 APrimeToFCPPair (std::string name, const framework::config::Parameters &p)
 Constructor.
 
virtual ~APrimeToFCPPair ()=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 (A' short name in Geant4)
 
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.
 
 DECLARE_FACTORY_WITH_WAREHOUSE (XsecBiasingOperator, std::shared_ptr< XsecBiasingOperator >, std::string, const framework::config::Parameters &)
 The BiasingOperator factory.
 
virtual ~XsecBiasingOperator ()
 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

- 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 * 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.
 

Detailed Description

Bias the A' -> fcp+ fcp- process.

Modeled after the GammaToMuPair biasing operator.

Definition at line 21 of file APrimeToFCPPair.h.

Constructor & Destructor Documentation

◆ APrimeToFCPPair()

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

Constructor.

Calls parent constructor and allows access to configuration parameters.

Definition at line 12 of file APrimeToFCPPair.cxx.

14 : XsecBiasingOperator(name, p) {
15 volume_ = p.get<std::string>("volume");
16 factor_ = p.get<double>("factor");
17 threshold_ = p.get<double>("threshold");
18}
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78
XsecBiasingOperator(std::string name, const framework::config::Parameters &parameters)
Constructor.
double threshold_
Minimum kinetic energy [MeV] to allow a track to be biased.
std::string volume_
The volume to bias in.

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

Member Function Documentation

◆ getParticleToBias()

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

Return the particle to bias (A' short name in Geant4)

Implements simcore::XsecBiasingOperator.

Definition at line 48 of file APrimeToFCPPair.h.

48{ return "A^1"; }

◆ getProcessToBias()

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

Return the process to bias.

Implements simcore::XsecBiasingOperator.

Definition at line 43 of file APrimeToFCPPair.h.

43 {
45 }
static const std::string PROCESS_NAME
The name of this process.

References simcore::APrimeConversionToFCPs::PROCESS_NAME.

◆ getVolumeToBias()

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

Return the volume to bias in.

Implements simcore::XsecBiasingOperator.

Definition at line 51 of file APrimeToFCPPair.h.

51{ return volume_; }

References volume_.

◆ ProposeOccurenceBiasingOperation()

G4VBiasingOperation * simcore::biasoperators::APrimeToFCPPair::ProposeOccurenceBiasingOperation ( const G4Track * track,
const G4BiasingProcessInterface * callingProcess )
overridevirtual
Returns
Method that returns the biasing operation that will be used to bias the conversion of A' to fcp pairs.

Implements simcore::XsecBiasingOperator.

Definition at line 20 of file APrimeToFCPPair.cxx.

21 {
22 if (track->GetKineticEnergy() < threshold_) {
23 return nullptr;
24 }
25
26 std::string current_process =
27 callingProcess->GetWrappedProcess()->GetProcessName();
28 if (current_process.compare(this->getProcessToBias()) == 0) {
29 G4double interaction_length =
30 callingProcess->GetWrappedProcess()->GetCurrentInteractionLength();
31
32 double xsec_unbiased = 1. / interaction_length;
33 double xsec_biased = xsec_unbiased * factor_;
34
35 return BiasedXsec(xsec_biased);
36 }
37 return nullptr;
38}
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::APrimeToFCPPair::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 58 of file APrimeToFCPPair.h.

58 {
59 header.setStringParameter("BiasOperator::APrimeToFCPPair::Volume", volume_);
60 header.setFloatParameter("BiasOperator::APrimeToFCPPair::Factor", factor_);
61 header.setFloatParameter("BiasOperator::APrimeToFCPPair::Threshold",
63 }
void setFloatParameter(const std::string &name, float value)
Set a float parameter value.
Definition RunHeader.h:197
void setStringParameter(const std::string &name, std::string value)
Set a string parameter value.
Definition RunHeader.h:222

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

Member Data Documentation

◆ factor_

double simcore::biasoperators::APrimeToFCPPair::factor_
private

The biasing factor.

Definition at line 70 of file APrimeToFCPPair.h.

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

◆ threshold_

double simcore::biasoperators::APrimeToFCPPair::threshold_
private

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

Definition at line 73 of file APrimeToFCPPair.h.

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

◆ volume_

std::string simcore::biasoperators::APrimeToFCPPair::volume_
private

The volume to bias in.

Definition at line 67 of file APrimeToFCPPair.h.

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


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