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

Bias operator for the dark brem process. More...

#include <DarkBrem.h>

Public Member Functions

 DarkBrem (std::string name, const framework::config::Parameters &p)
 Constructor.
 
virtual ~DarkBrem ()=default
 Destructor.
 
G4VBiasingOperation * ProposeOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess) override
 Calculate the biased cross section given the input process and track.
 
std::string getProcessToBias () const override
 Return the name of the process this operator biases.
 
std::string getParticleToBias () const override
 Return the name of the particle this operator biases.
 
std::string getVolumeToBias () const override
 Return the volume this operator biases.
 
void RecordConfig (ldmx::RunHeader &header) const override
 Record the configuration of this biasing operator 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

std::string volume_
 DEBUG FUNCTION This function is called by the biasing interface class during PostStepDoIt.
 
double factor_
 factor we want to bias by
 
bool bias_all_
 should we bias all electrons? (or only the primary)
 

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 operator for the dark brem process.

Definition at line 16 of file DarkBrem.h.

Constructor & Destructor Documentation

◆ DarkBrem()

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

Constructor.

Calls base class constructor and allows access to configuration parameters.

Definition at line 12 of file DarkBrem.cxx.

13 : XsecBiasingOperator(name, p) {
14 volume_ = p.getParameter<std::string>("volume");
15 factor_ = p.getParameter<double>("factor");
16 bias_all_ = p.getParameter<bool>("bias_all");
17}
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 we want to bias by
Definition DarkBrem.h:100
bool bias_all_
should we bias all electrons? (or only the primary)
Definition DarkBrem.h:103
std::string volume_
DEBUG FUNCTION This function is called by the biasing interface class during PostStepDoIt.
Definition DarkBrem.h:97

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

◆ ~DarkBrem()

virtual simcore::biasoperators::DarkBrem::~DarkBrem ( )
virtualdefault

Destructor.

Blank right now

Member Function Documentation

◆ getParticleToBias()

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

Return the name of the particle this operator biases.

Implements simcore::XsecBiasingOperator.

Definition at line 56 of file DarkBrem.h.

56{ return "e-"; }

◆ getProcessToBias()

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

Return the name of the process this operator biases.

Implements simcore::XsecBiasingOperator.

Definition at line 51 of file DarkBrem.h.

51 {
52 return G4DarkBremsstrahlung::PROCESS_NAME;
53 }

◆ getVolumeToBias()

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

Return the volume this operator biases.

Implements simcore::XsecBiasingOperator.

Definition at line 59 of file DarkBrem.h.

59{ return volume_; }

References volume_.

◆ ProposeOccurenceBiasingOperation()

G4VBiasingOperation * simcore::biasoperators::DarkBrem::ProposeOccurenceBiasingOperation ( const G4Track *  track,
const G4BiasingProcessInterface *  callingProcess 
)
overridevirtual

Calculate the biased cross section given the input process and track.

This allows for us to have access to the current information about the track while calculating the biased cross section.

See also
XsecBiasingOperator::BiasedXsec
Parameters
[in]trackconst pointer to track to Bias
[in]callingProcessprocess that might be biased by this operator
Returns
Method that returns the biasing operation that will be used to bias the occurence of events.

Implements simcore::XsecBiasingOperator.

Definition at line 19 of file DarkBrem.cxx.

20 {
21 std::string currentProcess =
22 callingProcess->GetWrappedProcess()->GetProcessName();
23 if (currentProcess.compare(this->getProcessToBias()) == 0) {
24 // bias only the primary particle if we don't want to bias all particles
25 if (not bias_all_ and track->GetParentID() != 0) return nullptr;
26
27 G4double interactionLength =
28 callingProcess->GetWrappedProcess()->GetCurrentInteractionLength();
29
30 double dbXsecUnbiased = 1. / interactionLength;
31 double dbXsecBiased = dbXsecUnbiased * factor_;
32
33 if (G4RunManager::GetRunManager()->GetVerboseLevel() > 1) {
34 std::cout << "[ DarkBremXsecBiasingOperator ]: "
35 << " Unbiased DBrem xsec: " << dbXsecUnbiased
36 << " -> Biased xsec: " << dbXsecBiased << std::endl;
37 }
38
39 return BiasedXsec(dbXsecBiased);
40 }
41 return nullptr;
42}
G4VBiasingOperation * BiasedXsec(double biased_xsec)
Helper method for passing a biased interaction length to the Geant4 biasing framework.

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

◆ RecordConfig()

void simcore::biasoperators::DarkBrem::RecordConfig ( ldmx::RunHeader header) const
overridevirtual

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

Parameters
[in,out]headerRunHeader to record configuration to

Implements simcore::XsecBiasingOperator.

Definition at line 44 of file DarkBrem.cxx.

44 {
45 h.setIntParameter("BiasOperator::DarkBrem::Bias All Electrons", bias_all_);
46 h.setFloatParameter("BiasOperator::DarkBrem::Factor", factor_);
47 h.setStringParameter("BiasOperator::DarkBrem::Volume", volume_);
48}

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

Member Data Documentation

◆ bias_all_

bool simcore::biasoperators::DarkBrem::bias_all_
private

should we bias all electrons? (or only the primary)

Definition at line 103 of file DarkBrem.h.

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

◆ factor_

double simcore::biasoperators::DarkBrem::factor_
private

factor we want to bias by

Definition at line 100 of file DarkBrem.h.

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

◆ volume_

std::string simcore::biasoperators::DarkBrem::volume_
private

DEBUG FUNCTION This function is called by the biasing interface class during PostStepDoIt.

You can observe the particle change that was produced by the process and the weight that will be multiplied into this particle change.

This is called inside G4VBiasingOperator::ReportOperationApplied which is called inside G4BiasingProcessInterface::PostStepDoIt void OperationApplied(const G4BiasingProcessInterface* callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation* operationApplied, G4double weight, G4VBiasingOperation* finalStateOpApplied, const G4VParticleChange* particleChangeProduced ) { std::string currentProcess = callingProcess->GetWrappedProcess()->GetProcessName(); if (currentProcess.compare(this->getProcessToBias()) == 0) { std::cout << "DB Final State Biasing Operator Applied: " << callingProcess->GetProcessName() << " -> " << weight*particleChangeProduced->GetWeight() << std::endl; } } volume we want to bias in

Definition at line 97 of file DarkBrem.h.

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


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