Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VBiasingOperator Class Referenceabstract

#include <G4VBiasingOperator.hh>

+ Inheritance diagram for G4VBiasingOperator:

Public Member Functions

 G4VBiasingOperator (const G4String &name)
 
virtual ~G4VBiasingOperator ()=default
 
virtual void Configure ()
 
virtual void ConfigureForWorker ()
 
virtual void StartRun ()
 
virtual void StartTracking (const G4Track *)
 
virtual void EndTracking ()
 
const G4StringGetName () const
 
void AttachTo (const G4LogicalVolume *)
 
G4BiasingAppliedCase GetPreviousBiasingAppliedCase () const
 
G4VBiasingOperationGetProposedOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
G4VBiasingOperationGetProposedFinalStateBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
G4VBiasingOperationGetProposedNonPhysicsBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
void ExitingBiasing (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
void ReportOperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
 
void ReportOperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
 
const G4VBiasingOperationGetPreviousNonPhysicsAppliedOperation ()
 

Static Public Member Functions

static const std::vector< G4VBiasingOperator * > & GetBiasingOperators ()
 
static G4VBiasingOperatorGetBiasingOperator (const G4LogicalVolume *)
 

Protected Member Functions

virtual G4VBiasingOperationProposeNonPhysicsBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
 
virtual G4VBiasingOperationProposeOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
 
virtual G4VBiasingOperationProposeFinalStateBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
 
virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual void ExitBiasing (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 

Friends

class G4BiasingOperatorStateNotifier
 

Detailed Description

Definition at line 173 of file G4VBiasingOperator.hh.

Constructor & Destructor Documentation

◆ G4VBiasingOperator()

G4VBiasingOperator::G4VBiasingOperator ( const G4String & name)

◆ ~G4VBiasingOperator()

virtual G4VBiasingOperator::~G4VBiasingOperator ( )
virtualdefault

Member Function Documentation

◆ AttachTo()

void G4VBiasingOperator::AttachTo ( const G4LogicalVolume * logical)

Definition at line 46 of file G4VBiasingOperator.cc.

47{
48 auto it = fLogicalToSetupMap.Find(logical);
49 if ( it == fLogicalToSetupMap.End() )
50 {
51 fLogicalToSetupMap[logical] = this;
52 }
53 else if ( (*it).second != this )
54 {
56 ed << "Biasing operator `" << GetName()
57 << "' can not be attached to Logical volume `"
58 << logical->GetName() << "' which is already used by another operator !"
59 << G4endl;
60 G4Exception("G4VBiasingOperator::AttachTo(...)",
61 "BIAS.MNG.01", JustWarning, ed);
62 }
63}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
#define G4endl
Definition G4ios.hh:67
const G4String & GetName() const
const G4String & GetName() const

◆ Configure()

virtual void G4VBiasingOperator::Configure ( )
inlinevirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 190 of file G4VBiasingOperator.hh.

190{}

◆ ConfigureForWorker()

virtual void G4VBiasingOperator::ConfigureForWorker ( )
inlinevirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 193 of file G4VBiasingOperator.hh.

193{}

◆ EndTracking()

virtual void G4VBiasingOperator::EndTracking ( )
inlinevirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 198 of file G4VBiasingOperator.hh.

198{}

◆ ExitBiasing()

void G4VBiasingOperator::ExitBiasing ( const G4Track * track,
const G4BiasingProcessInterface * callingProcess )
protectedvirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 164 of file G4VBiasingOperator.cc.

166{}

Referenced by ExitingBiasing().

◆ ExitingBiasing()

void G4VBiasingOperator::ExitingBiasing ( const G4Track * track,
const G4BiasingProcessInterface * callingProcess )

Definition at line 143 of file G4VBiasingOperator.cc.

145{
146 ExitBiasing( track, callingProcess );
147
148 // -- reset all data members:
149 fOccurenceBiasingOperation = nullptr ;
150 fFinalStateBiasingOperation = nullptr ;
151 fNonPhysicsBiasingOperation = nullptr ;
152 fPreviousProposedOccurenceBiasingOperation = nullptr ;
153 fPreviousProposedFinalStateBiasingOperation = nullptr ;
154 fPreviousProposedNonPhysicsBiasingOperation = nullptr ;
155 fPreviousAppliedOccurenceBiasingOperation = nullptr ;
156 fPreviousAppliedFinalStateBiasingOperation = nullptr ;
157 fPreviousAppliedNonPhysicsBiasingOperation = nullptr ;
158 fPreviousBiasingAppliedCase = BAC_None ;
159}
virtual void ExitBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)

◆ GetBiasingOperator()

G4VBiasingOperator * G4VBiasingOperator::GetBiasingOperator ( const G4LogicalVolume * logical)
static

Definition at line 71 of file G4VBiasingOperator.cc.

72{
73 auto it = fLogicalToSetupMap.Find(logical);
74 if ( it == fLogicalToSetupMap.End() ) { return nullptr; }
75 else { return (*it).second; }
76}

Referenced by G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength().

◆ GetBiasingOperators()

const std::vector< G4VBiasingOperator * > & G4VBiasingOperator::GetBiasingOperators ( )
static

◆ GetName()

const G4String & G4VBiasingOperator::GetName ( ) const
inline

◆ GetPreviousBiasingAppliedCase()

G4BiasingAppliedCase G4VBiasingOperator::GetPreviousBiasingAppliedCase ( ) const
inline

Definition at line 207 of file G4VBiasingOperator.hh.

208 { return fPreviousBiasingAppliedCase; }

◆ GetPreviousNonPhysicsAppliedOperation()

const G4VBiasingOperation * G4VBiasingOperator::GetPreviousNonPhysicsAppliedOperation ( )
inline

Definition at line 236 of file G4VBiasingOperator.hh.

237 { return fPreviousAppliedNonPhysicsBiasingOperation; }

◆ GetProposedFinalStateBiasingOperation()

G4VBiasingOperation * G4VBiasingOperator::GetProposedFinalStateBiasingOperation ( const G4Track * track,
const G4BiasingProcessInterface * callingProcess )

Definition at line 84 of file G4VBiasingOperator.cc.

85{
86 fFinalStateBiasingOperation = ProposeFinalStateBiasingOperation(track, callingProcess);
87 return fFinalStateBiasingOperation;
88}
virtual G4VBiasingOperation * ProposeFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

◆ GetProposedNonPhysicsBiasingOperation()

G4VBiasingOperation * G4VBiasingOperator::GetProposedNonPhysicsBiasingOperation ( const G4Track * track,
const G4BiasingProcessInterface * callingProcess )

Definition at line 90 of file G4VBiasingOperator.cc.

91{
92 fNonPhysicsBiasingOperation = ProposeNonPhysicsBiasingOperation(track, callingProcess);
93 return fNonPhysicsBiasingOperation;
94}
virtual G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

◆ GetProposedOccurenceBiasingOperation()

G4VBiasingOperation * G4VBiasingOperator::GetProposedOccurenceBiasingOperation ( const G4Track * track,
const G4BiasingProcessInterface * callingProcess )

Definition at line 78 of file G4VBiasingOperator.cc.

79{
80 fOccurenceBiasingOperation = ProposeOccurenceBiasingOperation(track, callingProcess);
81 return fOccurenceBiasingOperation;
82}
virtual G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

◆ OperationApplied() [1/2]

void G4VBiasingOperator::OperationApplied ( const G4BiasingProcessInterface * callingProcess,
G4BiasingAppliedCase biasingCase,
G4VBiasingOperation * occurenceOperationApplied,
G4double weightForOccurenceInteraction,
G4VBiasingOperation * finalStateOperationApplied,
const G4VParticleChange * particleChangeProduced )
protectedvirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 172 of file G4VBiasingOperator.cc.

176{
177}

◆ OperationApplied() [2/2]

void G4VBiasingOperator::OperationApplied ( const G4BiasingProcessInterface * callingProcess,
G4BiasingAppliedCase biasingCase,
G4VBiasingOperation * operationApplied,
const G4VParticleChange * particleChangeProduced )
protectedvirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 167 of file G4VBiasingOperator.cc.

170{
171}

Referenced by ReportOperationApplied(), and ReportOperationApplied().

◆ ProposeFinalStateBiasingOperation()

virtual G4VBiasingOperation * G4VBiasingOperator::ProposeFinalStateBiasingOperation ( const G4Track * track,
const G4BiasingProcessInterface * callingProcess )
protectedpure virtual

◆ ProposeNonPhysicsBiasingOperation()

virtual G4VBiasingOperation * G4VBiasingOperator::ProposeNonPhysicsBiasingOperation ( const G4Track * track,
const G4BiasingProcessInterface * callingProcess )
protectedpure virtual

◆ ProposeOccurenceBiasingOperation()

virtual G4VBiasingOperation * G4VBiasingOperator::ProposeOccurenceBiasingOperation ( const G4Track * track,
const G4BiasingProcessInterface * callingProcess )
protectedpure virtual

◆ ReportOperationApplied() [1/2]

void G4VBiasingOperator::ReportOperationApplied ( const G4BiasingProcessInterface * callingProcess,
G4BiasingAppliedCase biasingCase,
G4VBiasingOperation * occurenceOperationApplied,
G4double weightForOccurenceInteraction,
G4VBiasingOperation * finalStateOperationApplied,
const G4VParticleChange * particleChangeProduced )

Definition at line 129 of file G4VBiasingOperator.cc.

136{
137 fPreviousBiasingAppliedCase = biasingCase;
138 fPreviousAppliedOccurenceBiasingOperation = occurenceOperationApplied;
139 fPreviousAppliedFinalStateBiasingOperation = finalStateOperationApplied;
140 OperationApplied( callingProcess, biasingCase, occurenceOperationApplied, weightForOccurenceInteraction, finalStateOperationApplied, particleChangeProduced );
141}
virtual void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)

◆ ReportOperationApplied() [2/2]

void G4VBiasingOperator::ReportOperationApplied ( const G4BiasingProcessInterface * callingProcess,
G4BiasingAppliedCase biasingCase,
G4VBiasingOperation * operationApplied,
const G4VParticleChange * particleChangeProduced )

Definition at line 96 of file G4VBiasingOperator.cc.

101{
102 fPreviousBiasingAppliedCase = biasingCase;
103 fPreviousAppliedOccurenceBiasingOperation = nullptr;
104 fPreviousAppliedFinalStateBiasingOperation = nullptr;
105 fPreviousAppliedNonPhysicsBiasingOperation = nullptr;
106 switch ( biasingCase )
107 {
108 case BAC_None:
109 break;
110 case BAC_NonPhysics:
111 fPreviousAppliedNonPhysicsBiasingOperation = operationApplied ;
112 break;
113 case BAC_FinalState:
114 fPreviousAppliedFinalStateBiasingOperation = operationApplied;
115 break;
116 case BAC_Occurence:
117 G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
118 "BIAS.MNG.02", JustWarning,
119 "Internal logic error, please report !");
120 break;
121 default:
122 G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
123 "BIAS.MNG.03", JustWarning,
124 "Internal logic error, please report !");
125 }
126 OperationApplied( callingProcess, biasingCase, operationApplied, particleChangeProduced );
127}
@ BAC_Occurence
@ BAC_FinalState
@ BAC_NonPhysics

◆ StartRun()

virtual void G4VBiasingOperator::StartRun ( )
inlinevirtual

Reimplemented in G4BOptrForceCollision, and G4ChannelingOptrChangeCrossSection.

Definition at line 195 of file G4VBiasingOperator.hh.

195{}

◆ StartTracking()

virtual void G4VBiasingOperator::StartTracking ( const G4Track * )
inlinevirtual

Reimplemented in G4BOptrForceCollision, and G4ChannelingOptrMultiParticleChangeCrossSection.

Definition at line 197 of file G4VBiasingOperator.hh.

197{}

Friends And Related Symbol Documentation

◆ G4BiasingOperatorStateNotifier

friend class G4BiasingOperatorStateNotifier
friend

Definition at line 178 of file G4VBiasingOperator.hh.

Referenced by G4BiasingOperatorStateNotifier, and G4VBiasingOperator().


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