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

Bias the Photon-Nuclear process. More...

#include <PhotoNuclear.h>

Public Member Functions

 PhotoNuclear (std::string name, const framework::config::Parameters &p)
 Constructor.
 
void StartRun () override
 Method called at the beginning of a run.
 
G4VBiasingOperation * ProposeOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess) override
 
std::string getProcessToBias () const override
 return the process we want to bias
 
std::string getParticleToBias () const override
 return the particle that we want to bias
 
std::string getVolumeToBias () const override
 return the volume we want to bias within
 
void RecordConfig (ldmx::RunHeader &h) const override
 record the configuration into 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

G4BOptnChangeCrossSection * emXsecOperation {nullptr}
 Cross-section biasing operation for conversion process.
 
double pnXsecUnbiased_ {0}
 Unbiased photonuclear xsec.
 
double pnXsecBiased_ {0}
 Biased photonuclear xsec.
 
std::string volume_
 Volume we are going to bias within.
 
double threshold_
 minimum kinetic energy [MeV] for a track to be biased
 
double factor_
 factor to bias PN by
 
bool down_bias_conv_
 Should we down-bias the gamma conversion process?
 
bool only_children_of_primary_
 Should we restrict biasing to only children of primary?
 

Static Private Attributes

static const std::string CONVERSION_PROCESS = "conv"
 Geant4 gamma conversion process name.
 

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 Photon-Nuclear process.

Definition at line 13 of file PhotoNuclear.h.

Constructor & Destructor Documentation

◆ PhotoNuclear()

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

Constructor.

Definition at line 8 of file PhotoNuclear.cxx.

10 : XsecBiasingOperator(name, p) {
11 volume_ = p.getParameter<std::string>("volume");
12 threshold_ = p.getParameter<double>("threshold");
13 factor_ = p.getParameter<double>("factor");
14 down_bias_conv_ = p.getParameter<bool>("down_bias_conv");
15 only_children_of_primary_ = p.getParameter<bool>("only_children_of_primary");
16}
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.
double factor_
factor to bias PN by
bool only_children_of_primary_
Should we restrict biasing to only children of primary?
bool down_bias_conv_
Should we down-bias the gamma conversion process?
double threshold_
minimum kinetic energy [MeV] for a track to be biased
std::string volume_
Volume we are going to bias within.

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

Member Function Documentation

◆ getParticleToBias()

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

return the particle that we want to bias

Implements simcore::XsecBiasingOperator.

Definition at line 33 of file PhotoNuclear.h.

33{ return "gamma"; }

◆ getProcessToBias()

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

return the process we want to bias

Implements simcore::XsecBiasingOperator.

Definition at line 30 of file PhotoNuclear.h.

30{ return "photonNuclear"; }

◆ getVolumeToBias()

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

return the volume we want to bias within

Implements simcore::XsecBiasingOperator.

Definition at line 36 of file PhotoNuclear.h.

36{ return volume_; }

References volume_.

◆ ProposeOccurenceBiasingOperation()

G4VBiasingOperation * simcore::biasoperators::PhotoNuclear::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 30 of file PhotoNuclear.cxx.

31 {
32 // if we want to only bias children of primary, leave if this track is NOT a
33 // child of the primary
34 if (only_children_of_primary_ and track->GetParentID() != 1) {
35 return nullptr;
36 }
37
38 // is this track too low energy to be biased?
39 if (track->GetKineticEnergy() < threshold_) {
40 return nullptr;
41 }
42
43 std::string currentProcess =
44 callingProcess->GetWrappedProcess()->GetProcessName();
45 if (currentProcess.compare(this->getProcessToBias()) == 0) {
46 G4double interactionLength =
47 callingProcess->GetWrappedProcess()->GetCurrentInteractionLength();
48
49 pnXsecUnbiased_ = 1. / interactionLength;
50
52
54 }
55 if ((currentProcess.compare(CONVERSION_PROCESS) == 0) and down_bias_conv_) {
56 G4double interactionLength =
57 callingProcess->GetWrappedProcess()->GetCurrentInteractionLength();
58
59 double emXsecUnbiased = 1. / interactionLength;
60
61 double emXsecBiased = std::max(
63 if (emXsecBiased == pnXsecUnbiased_) {
64 G4cout << "[ PhotoNuclearXsecBiasingOperator ]: [ WARNING ]: "
65 << "Biasing factor is too large." << std::endl;
66 }
67
68 emXsecOperation->SetBiasedCrossSection(emXsecBiased);
69 emXsecOperation->Sample();
70
71 return emXsecOperation;
72 }
73 return nullptr;
74}
G4VBiasingOperation * BiasedXsec(double biased_xsec)
Helper method for passing a biased interaction length to the Geant4 biasing framework.
G4BOptnChangeCrossSection * emXsecOperation
Cross-section biasing operation for conversion process.
static const std::string CONVERSION_PROCESS
Geant4 gamma conversion process name.
double pnXsecUnbiased_
Unbiased photonuclear xsec.
double pnXsecBiased_
Biased photonuclear xsec.

References simcore::XsecBiasingOperator::BiasedXsec(), CONVERSION_PROCESS, down_bias_conv_, emXsecOperation, factor_, only_children_of_primary_, pnXsecBiased_, pnXsecUnbiased_, and threshold_.

◆ RecordConfig()

void simcore::biasoperators::PhotoNuclear::RecordConfig ( ldmx::RunHeader h) const
inlineoverridevirtual

record the configuration into the run header

Implements simcore::XsecBiasingOperator.

Definition at line 39 of file PhotoNuclear.h.

39 {
40 h.setStringParameter("BiasOperators::PhotoNuclear::Volume", volume_);
41 h.setFloatParameter("BiasOperators::PhotoNuclear::Threshold", threshold_);
42 h.setFloatParameter("BiasOperators::PhotoNuclear::Factor", factor_);
43 h.setIntParameter("BiasOperators::PhotoNuclear::Bias Conv Down",
45 h.setIntParameter("BiasOperators::PhotoNuclear::Only Children Of Primary",
47 }
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
void setIntParameter(const std::string &name, int value)
Set an int parameter value.
Definition RunHeader.h:164

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

◆ StartRun()

void simcore::biasoperators::PhotoNuclear::StartRun ( )
override

Method called at the beginning of a run.

Definition at line 18 of file PhotoNuclear.cxx.

18 {
20
22 emXsecOperation = new G4BOptnChangeCrossSection("changeXsec-conv");
23 } else if (down_bias_conv_) {
24 EXCEPTION_RAISE(
25 "PhotoNuclearBiasing",
26 "Gamma Conversion process '" + CONVERSION_PROCESS + "' is not biased!");
27 }
28}
void StartRun()
Method called at the beginning of a run.
bool processIsBiased(std::string process)
Check if the given processed is being biased.

References CONVERSION_PROCESS, down_bias_conv_, emXsecOperation, simcore::XsecBiasingOperator::processIsBiased(), and simcore::XsecBiasingOperator::StartRun().

Member Data Documentation

◆ CONVERSION_PROCESS

const std::string simcore::biasoperators::PhotoNuclear::CONVERSION_PROCESS = "conv"
staticprivate

Geant4 gamma conversion process name.

Definition at line 51 of file PhotoNuclear.h.

Referenced by ProposeOccurenceBiasingOperation(), and StartRun().

◆ down_bias_conv_

bool simcore::biasoperators::PhotoNuclear::down_bias_conv_
private

Should we down-bias the gamma conversion process?

Definition at line 72 of file PhotoNuclear.h.

Referenced by PhotoNuclear(), ProposeOccurenceBiasingOperation(), RecordConfig(), and StartRun().

◆ emXsecOperation

G4BOptnChangeCrossSection* simcore::biasoperators::PhotoNuclear::emXsecOperation {nullptr}
private

Cross-section biasing operation for conversion process.

Definition at line 54 of file PhotoNuclear.h.

54{nullptr};

Referenced by ProposeOccurenceBiasingOperation(), and StartRun().

◆ factor_

double simcore::biasoperators::PhotoNuclear::factor_
private

factor to bias PN by

Definition at line 69 of file PhotoNuclear.h.

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

◆ only_children_of_primary_

bool simcore::biasoperators::PhotoNuclear::only_children_of_primary_
private

Should we restrict biasing to only children of primary?

Definition at line 75 of file PhotoNuclear.h.

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

◆ pnXsecBiased_

double simcore::biasoperators::PhotoNuclear::pnXsecBiased_ {0}
private

Biased photonuclear xsec.

Definition at line 60 of file PhotoNuclear.h.

60{0};

Referenced by ProposeOccurenceBiasingOperation().

◆ pnXsecUnbiased_

double simcore::biasoperators::PhotoNuclear::pnXsecUnbiased_ {0}
private

Unbiased photonuclear xsec.

Definition at line 57 of file PhotoNuclear.h.

57{0};

Referenced by ProposeOccurenceBiasingOperation().

◆ threshold_

double simcore::biasoperators::PhotoNuclear::threshold_
private

minimum kinetic energy [MeV] for a track to be biased

Definition at line 66 of file PhotoNuclear.h.

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

◆ volume_

std::string simcore::biasoperators::PhotoNuclear::volume_
private

Volume we are going to bias within.

Definition at line 63 of file PhotoNuclear.h.

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


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