Geant4-11
Public Member Functions | Static Public Member Functions | Protected Member Functions | Private Attributes | Static Private Attributes | Friends
G4VBiasingOperator Class Referenceabstract

#include <G4VBiasingOperator.hh>

Inheritance diagram for G4VBiasingOperator:
G4BOptrForceCollision G4ChannelingOptrChangeCrossSection G4ChannelingOptrMultiParticleChangeCrossSection

Public Member Functions

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

Static Public Member Functions

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

Protected Member Functions

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

Private Attributes

std::map< const G4LogicalVolume *, G4intfDepthInTree
 
G4VBiasingOperationfFinalStateBiasingOperation
 
const G4String fName
 
G4VBiasingOperationfNonPhysicsBiasingOperation
 
G4VBiasingOperationfOccurenceBiasingOperation
 
const G4VBiasingOperationfPreviousAppliedFinalStateBiasingOperation
 
const G4VBiasingOperationfPreviousAppliedNonPhysicsBiasingOperation
 
const G4VBiasingOperationfPreviousAppliedOccurenceBiasingOperation
 
G4BiasingAppliedCase fPreviousBiasingAppliedCase
 
const G4VBiasingOperationfPreviousProposedFinalStateBiasingOperation
 
const G4VBiasingOperationfPreviousProposedNonPhysicsBiasingOperation
 
const G4VBiasingOperationfPreviousProposedOccurenceBiasingOperation
 
std::vector< const G4LogicalVolume * > fRootVolumes
 

Static Private Attributes

static G4MapCache< const G4LogicalVolume *, G4VBiasingOperator * > fLogicalToSetupMap
 
static G4VectorCache< G4VBiasingOperator * > fOperators
 
static G4Cache< G4BiasingOperatorStateNotifier * > fStateNotifier
 

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 ),
48{
50
52
53}
value_type & Get() const
Definition: G4Cache.hh:315
void Put(const value_type &val) const
Definition: G4Cache.hh:321
friend class G4BiasingOperatorStateNotifier
const G4VBiasingOperation * fPreviousProposedNonPhysicsBiasingOperation
static G4VectorCache< G4VBiasingOperator * > fOperators
G4VBiasingOperation * fNonPhysicsBiasingOperation
G4BiasingAppliedCase fPreviousBiasingAppliedCase
G4VBiasingOperation * fOccurenceBiasingOperation
G4VBiasingOperation * fFinalStateBiasingOperation
const G4VBiasingOperation * fPreviousProposedOccurenceBiasingOperation
const G4VBiasingOperation * fPreviousAppliedOccurenceBiasingOperation
const G4VBiasingOperation * fPreviousAppliedNonPhysicsBiasingOperation
const G4VBiasingOperation * fPreviousProposedFinalStateBiasingOperation
const G4VBiasingOperation * fPreviousAppliedFinalStateBiasingOperation
static G4Cache< G4BiasingOperatorStateNotifier * > fStateNotifier
void Push_back(const value_type &val)
Definition: G4Cache.hh:375
const char * name(G4int ptype)

References fOperators, fStateNotifier, G4BiasingOperatorStateNotifier, G4Cache< VALTYPE >::Get(), G4VectorCache< VALTYPE >::Push_back(), and G4Cache< VALTYPE >::Put().

◆ ~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:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4endl
Definition: G4ios.hh:57
const G4String & GetName() const
iterator Find(const key_type &k)
Definition: G4Cache.hh:459
iterator End()
Definition: G4Cache.hh:453
static G4MapCache< const G4LogicalVolume *, G4VBiasingOperator * > fLogicalToSetupMap
const G4String GetName() const

References G4MapCache< KEYTYPE, VALTYPE >::End(), G4MapCache< KEYTYPE, VALTYPE >::Find(), fLogicalToSetupMap, G4endl, G4Exception(), G4LogicalVolume::GetName(), GetName(), and JustWarning.

◆ 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 
)

◆ GetBiasingOperator()

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

◆ GetBiasingOperators()

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

◆ GetName()

const G4String G4VBiasingOperator::GetName ( ) const
inline

◆ GetPreviousBiasingAppliedCase()

G4BiasingAppliedCase G4VBiasingOperator::GetPreviousBiasingAppliedCase ( ) const
inline

Definition at line 291 of file G4VBiasingOperator.hh.

References fPreviousBiasingAppliedCase.

◆ GetPreviousNonPhysicsAppliedOperation()

const G4VBiasingOperation * G4VBiasingOperator::GetPreviousNonPhysicsAppliedOperation ( )
inline

◆ GetProposedFinalStateBiasingOperation()

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

Definition at line 92 of file G4VBiasingOperator.cc.

93{
96}
virtual G4VBiasingOperation * ProposeFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

References fFinalStateBiasingOperation, and ProposeFinalStateBiasingOperation().

◆ GetProposedNonPhysicsBiasingOperation()

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

Definition at line 98 of file G4VBiasingOperator.cc.

99{
102}
virtual G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

References fNonPhysicsBiasingOperation, and ProposeNonPhysicsBiasingOperation().

◆ GetProposedOccurenceBiasingOperation()

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

Definition at line 86 of file G4VBiasingOperator.cc.

87{
90}
virtual G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

References fOccurenceBiasingOperation, and ProposeOccurenceBiasingOperation().

Referenced by G4ChannelingOptrMultiParticleChangeCrossSection::ProposeOccurenceBiasingOperation().

◆ OperationApplied() [1/2]

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

◆ OperationApplied() [2/2]

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

◆ 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)

References fPreviousAppliedFinalStateBiasingOperation, fPreviousAppliedOccurenceBiasingOperation, fPreviousBiasingAppliedCase, and OperationApplied().

◆ 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;
113 switch ( biasingCase )
114 {
115 case BAC_None:
116 break;
117 case BAC_NonPhysics:
119 break;
120 case BAC_FinalState:
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

References BAC_FinalState, BAC_None, BAC_NonPhysics, BAC_Occurence, fPreviousAppliedFinalStateBiasingOperation, fPreviousAppliedNonPhysicsBiasingOperation, fPreviousAppliedOccurenceBiasingOperation, fPreviousBiasingAppliedCase, G4Exception(), JustWarning, and OperationApplied().

Referenced by G4ChannelingOptrMultiParticleChangeCrossSection::OperationApplied().

◆ 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().

Field Documentation

◆ fDepthInTree

std::map< const G4LogicalVolume*, G4int > G4VBiasingOperator::fDepthInTree
private

Definition at line 333 of file G4VBiasingOperator.hh.

◆ fFinalStateBiasingOperation

G4VBiasingOperation* G4VBiasingOperator::fFinalStateBiasingOperation
private

Definition at line 337 of file G4VBiasingOperator.hh.

Referenced by ExitingBiasing(), and GetProposedFinalStateBiasingOperation().

◆ fLogicalToSetupMap

G4MapCache< const G4LogicalVolume *, G4VBiasingOperator * > G4VBiasingOperator::fLogicalToSetupMap
staticprivate

Definition at line 322 of file G4VBiasingOperator.hh.

Referenced by AttachTo(), and GetBiasingOperator().

◆ fName

const G4String G4VBiasingOperator::fName
private

Definition at line 319 of file G4VBiasingOperator.hh.

Referenced by GetName().

◆ fNonPhysicsBiasingOperation

G4VBiasingOperation* G4VBiasingOperator::fNonPhysicsBiasingOperation
private

Definition at line 338 of file G4VBiasingOperator.hh.

Referenced by ExitingBiasing(), and GetProposedNonPhysicsBiasingOperation().

◆ fOccurenceBiasingOperation

G4VBiasingOperation* G4VBiasingOperator::fOccurenceBiasingOperation
private

Definition at line 336 of file G4VBiasingOperator.hh.

Referenced by ExitingBiasing(), and GetProposedOccurenceBiasingOperation().

◆ fOperators

G4VectorCache< G4VBiasingOperator * > G4VBiasingOperator::fOperators
staticprivate

◆ fPreviousAppliedFinalStateBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousAppliedFinalStateBiasingOperation
private

Definition at line 345 of file G4VBiasingOperator.hh.

Referenced by ExitingBiasing(), and ReportOperationApplied().

◆ fPreviousAppliedNonPhysicsBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousAppliedNonPhysicsBiasingOperation
private

◆ fPreviousAppliedOccurenceBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousAppliedOccurenceBiasingOperation
private

Definition at line 344 of file G4VBiasingOperator.hh.

Referenced by ExitingBiasing(), and ReportOperationApplied().

◆ fPreviousBiasingAppliedCase

G4BiasingAppliedCase G4VBiasingOperator::fPreviousBiasingAppliedCase
private

◆ fPreviousProposedFinalStateBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousProposedFinalStateBiasingOperation
private

Definition at line 342 of file G4VBiasingOperator.hh.

Referenced by ExitingBiasing().

◆ fPreviousProposedNonPhysicsBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousProposedNonPhysicsBiasingOperation
private

Definition at line 343 of file G4VBiasingOperator.hh.

Referenced by ExitingBiasing().

◆ fPreviousProposedOccurenceBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousProposedOccurenceBiasingOperation
private

Definition at line 341 of file G4VBiasingOperator.hh.

Referenced by ExitingBiasing().

◆ fRootVolumes

std::vector< const G4LogicalVolume* > G4VBiasingOperator::fRootVolumes
private

Definition at line 332 of file G4VBiasingOperator.hh.

◆ fStateNotifier

G4Cache< G4BiasingOperatorStateNotifier * > G4VBiasingOperator::fStateNotifier
staticprivate

Definition at line 328 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator().


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