Geant4 11.1.1
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 (G4String name)
 
virtual ~G4VBiasingOperator ()
 
virtual void Configure ()
 
virtual void ConfigureForWorker ()
 
virtual void StartRun ()
 
virtual void StartTracking (const G4Track *)
 
virtual void EndTracking ()
 
const G4String GetName () 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 180 of file G4VBiasingOperator.hh.

Constructor & Destructor Documentation

◆ G4VBiasingOperator()

G4VBiasingOperator::G4VBiasingOperator ( G4String  name)

Definition at line 36 of file G4VBiasingOperator.cc.

37 : fName ( name ),
38 fOccurenceBiasingOperation ( nullptr ),
39 fFinalStateBiasingOperation ( nullptr ),
40 fNonPhysicsBiasingOperation ( nullptr ),
41 fPreviousProposedOccurenceBiasingOperation ( nullptr ),
42 fPreviousProposedFinalStateBiasingOperation( nullptr ),
43 fPreviousProposedNonPhysicsBiasingOperation( nullptr ),
44 fPreviousAppliedOccurenceBiasingOperation ( nullptr ),
45 fPreviousAppliedFinalStateBiasingOperation ( nullptr ),
46 fPreviousAppliedNonPhysicsBiasingOperation ( nullptr ),
47 fPreviousBiasingAppliedCase ( BAC_None )
48{
49 fOperators.Push_back(this);
50
51 if ( fStateNotifier.Get() == 0 ) fStateNotifier.Put( new G4BiasingOperatorStateNotifier() );
52
53}
value_type & Get() const
Definition: G4Cache.hh:315
void Put(const value_type &val) const
Definition: G4Cache.hh:321
friend class G4BiasingOperatorStateNotifier
void Push_back(const value_type &val)
Definition: G4Cache.hh:374

◆ ~G4VBiasingOperator()

G4VBiasingOperator::~G4VBiasingOperator ( )
virtual

Definition at line 55 of file G4VBiasingOperator.cc.

56{
57}

Member Function Documentation

◆ AttachTo()

void G4VBiasingOperator::AttachTo ( const G4LogicalVolume logical)

Definition at line 59 of file G4VBiasingOperator.cc.

60{
62 it = fLogicalToSetupMap.Find(logical);
63 if ( it == fLogicalToSetupMap.End() ) fLogicalToSetupMap[logical] = this;
64 else if ( (*it).second != this )
65 {
67 ed << "Biasing operator `" << GetName()
68 << "' can not be attached to Logical volume `"
69 << logical->GetName() << "' which is already used by another operator !" << G4endl;
70 G4Exception("G4VBiasingOperator::AttachTo(...)",
71 "BIAS.MNG.01",
73 ed);
74 }
75}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4endl
Definition: G4ios.hh:57
const G4String & GetName() const
typename map_type::iterator iterator
Definition: G4Cache.hh:184
iterator Find(const key_type &k)
Definition: G4Cache.hh:458
iterator End()
Definition: G4Cache.hh:452
const G4String GetName() const

◆ Configure()

virtual void G4VBiasingOperator::Configure ( )
inlinevirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 271 of file G4VBiasingOperator.hh.

271{}

◆ ConfigureForWorker()

virtual void G4VBiasingOperator::ConfigureForWorker ( )
inlinevirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 274 of file G4VBiasingOperator.hh.

274{}

◆ EndTracking()

virtual void G4VBiasingOperator::EndTracking ( )
inlinevirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 279 of file G4VBiasingOperator.hh.

279{}

◆ ExitBiasing()

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

Reimplemented in G4BOptrForceCollision.

Definition at line 173 of file G4VBiasingOperator.cc.

174{}

Referenced by ExitingBiasing().

◆ ExitingBiasing()

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

Definition at line 152 of file G4VBiasingOperator.cc.

153{
154 ExitBiasing( track, callingProcess );
155
156 // -- reset all data members:
157 fOccurenceBiasingOperation = nullptr ;
158 fFinalStateBiasingOperation = nullptr ;
159 fNonPhysicsBiasingOperation = nullptr ;
160 fPreviousProposedOccurenceBiasingOperation = nullptr ;
161 fPreviousProposedFinalStateBiasingOperation = nullptr ;
162 fPreviousProposedNonPhysicsBiasingOperation = nullptr ;
163 fPreviousAppliedOccurenceBiasingOperation = nullptr ;
164 fPreviousAppliedFinalStateBiasingOperation = nullptr ;
165 fPreviousAppliedNonPhysicsBiasingOperation = nullptr ;
166 fPreviousBiasingAppliedCase = BAC_None ;
167}
virtual void ExitBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)

◆ GetBiasingOperator()

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

Definition at line 78 of file G4VBiasingOperator.cc.

79{
81 it = fLogicalToSetupMap.Find(logical);
82 if ( it == fLogicalToSetupMap.End() ) return nullptr;
83 else return (*it).second;
84}

Referenced by G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength().

◆ GetBiasingOperators()

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

◆ GetName()

◆ GetPreviousBiasingAppliedCase()

G4BiasingAppliedCase G4VBiasingOperator::GetPreviousBiasingAppliedCase ( ) const
inline

Definition at line 291 of file G4VBiasingOperator.hh.

291{return fPreviousBiasingAppliedCase;}

◆ GetPreviousNonPhysicsAppliedOperation()

const G4VBiasingOperation * G4VBiasingOperator::GetPreviousNonPhysicsAppliedOperation ( )
inline

Definition at line 315 of file G4VBiasingOperator.hh.

315{return fPreviousAppliedNonPhysicsBiasingOperation;}

◆ GetProposedFinalStateBiasingOperation()

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

Definition at line 92 of file G4VBiasingOperator.cc.

93{
94 fFinalStateBiasingOperation = ProposeFinalStateBiasingOperation(track, callingProcess);
95 return fFinalStateBiasingOperation;
96}
virtual G4VBiasingOperation * ProposeFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

◆ GetProposedNonPhysicsBiasingOperation()

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

Definition at line 98 of file G4VBiasingOperator.cc.

99{
100 fNonPhysicsBiasingOperation = ProposeNonPhysicsBiasingOperation(track, callingProcess);
101 return fNonPhysicsBiasingOperation;
102}
virtual G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

◆ GetProposedOccurenceBiasingOperation()

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

Definition at line 86 of file G4VBiasingOperator.cc.

87{
88 fOccurenceBiasingOperation = ProposeOccurenceBiasingOperation(track, callingProcess);
89 return fOccurenceBiasingOperation;
90}
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 179 of file G4VBiasingOperator.cc.

182{
183}

◆ OperationApplied() [2/2]

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

Reimplemented in G4BOptrForceCollision.

Definition at line 175 of file G4VBiasingOperator.cc.

177{
178}

Referenced by 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 138 of file G4VBiasingOperator.cc.

144{
145 fPreviousBiasingAppliedCase = biasingCase;
146 fPreviousAppliedOccurenceBiasingOperation = occurenceOperationApplied;
147 fPreviousAppliedFinalStateBiasingOperation = finalStateOperationApplied;
148 OperationApplied( callingProcess, biasingCase, occurenceOperationApplied, weightForOccurenceInteraction, finalStateOperationApplied, particleChangeProduced );
149}
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 104 of file G4VBiasingOperator.cc.

108{
109 fPreviousBiasingAppliedCase = biasingCase;
110 fPreviousAppliedOccurenceBiasingOperation = nullptr;
111 fPreviousAppliedFinalStateBiasingOperation = nullptr;
112 fPreviousAppliedNonPhysicsBiasingOperation = nullptr;
113 switch ( biasingCase )
114 {
115 case BAC_None:
116 break;
117 case BAC_NonPhysics:
118 fPreviousAppliedNonPhysicsBiasingOperation = operationApplied ;
119 break;
120 case BAC_FinalState:
121 fPreviousAppliedFinalStateBiasingOperation = operationApplied;
122 break;
123 case BAC_Occurence:
124 G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
125 "BIAS.MNG.02",
127 "Internal logic error, please report !");
128 break;
129 default:
130 G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
131 "BIAS.MNG.03",
133 "Internal logic error, please report !");
134 }
135 OperationApplied( callingProcess, biasingCase, operationApplied, particleChangeProduced );
136}
@ BAC_Occurence
@ BAC_FinalState
@ BAC_NonPhysics

◆ StartRun()

virtual void G4VBiasingOperator::StartRun ( )
inlinevirtual

Reimplemented in G4ChannelingOptrChangeCrossSection, and G4BOptrForceCollision.

Definition at line 276 of file G4VBiasingOperator.hh.

276{}

◆ StartTracking()

virtual void G4VBiasingOperator::StartTracking ( const G4Track )
inlinevirtual

Reimplemented in G4ChannelingOptrMultiParticleChangeCrossSection, and G4BOptrForceCollision.

Definition at line 278 of file G4VBiasingOperator.hh.

278{}

Friends And Related Function Documentation

◆ G4BiasingOperatorStateNotifier

friend class G4BiasingOperatorStateNotifier
friend

Definition at line 185 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator().


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