Geant4  10.02.p03
G4VBiasingOperator Class Referenceabstract

#include <G4VBiasingOperator.hh>

Inheritance diagram for G4VBiasingOperator:
Collaboration 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)
 

Private Attributes

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

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 193 of file G4VBiasingOperator.hh.

Constructor & Destructor Documentation

◆ G4VBiasingOperator()

G4VBiasingOperator::G4VBiasingOperator ( G4String  name)

Definition at line 36 of file G4VBiasingOperator.cc.

37  : fName (name),
48 {
49  fOperators.Push_back(this);
50 
52 
53 }
void Put(const value_type &val) const
Definition: G4Cache.hh:286
const G4VBiasingOperation * fPreviousProposedOccurenceBiasingOperation
static G4VectorCache< G4VBiasingOperator *> fOperators
G4BiasingAppliedCase fPreviousBiasingAppliedCase
value_type & Get() const
Definition: G4Cache.hh:282
const G4VBiasingOperation * fPreviousAppliedOccurenceBiasingOperation
const G4VBiasingOperation * fPreviousAppliedNonPhysicsBiasingOperation
const G4VBiasingOperation * fPreviousProposedFinalStateBiasingOperation
const G4VBiasingOperation * fPreviousAppliedFinalStateBiasingOperation
void Push_back(const value_type &val)
Definition: G4Cache.hh:333
G4VBiasingOperation * fFinalStateBiasingOperation
G4VBiasingOperation * fOccurenceBiasingOperation
static G4Cache< G4BiasingOperatorStateNotifier *> fStateNotifier
const G4VBiasingOperation * fPreviousProposedNonPhysicsBiasingOperation
G4VBiasingOperation * fNonPhysicsBiasingOperation
Here is the call graph for this function:

◆ ~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 an other operator !" << G4endl;
70  G4Exception("G4VBiasingOperator::AttachTo(...)",
71  "BIAS.MNG.01",
73  ed);
74  }
75 }
static G4MapCache< const G4LogicalVolume *, G4VBiasingOperator *> fLogicalToSetupMap
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
iterator End()
Definition: G4Cache.hh:411
const G4String & GetName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
iterator Find(const key_type &k)
Definition: G4Cache.hh:417
map_type::iterator iterator
Definition: G4Cache.hh:175
#define G4endl
Definition: G4ios.hh:61
const G4String GetName() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Configure()

virtual void G4VBiasingOperator::Configure ( )
inlinevirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 263 of file G4VBiasingOperator.hh.

263 {}

◆ ConfigureForWorker()

virtual void G4VBiasingOperator::ConfigureForWorker ( )
inlinevirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 266 of file G4VBiasingOperator.hh.

266 {}

◆ EndTracking()

virtual void G4VBiasingOperator::EndTracking ( )
inlinevirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 271 of file G4VBiasingOperator.hh.

271 {}

◆ ExitBiasing()

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

Reimplemented in GB02BOptrMultiParticleForceCollision, and G4BOptrForceCollision.

Definition at line 176 of file G4VBiasingOperator.cc.

177 {}
Here is the caller graph for this function:

◆ ExitingBiasing()

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

Definition at line 155 of file G4VBiasingOperator.cc.

156 {
157  ExitBiasing( track, callingProcess );
158 
159  // -- reset all data members:
170 }
const G4VBiasingOperation * fPreviousProposedOccurenceBiasingOperation
G4BiasingAppliedCase fPreviousBiasingAppliedCase
virtual void ExitBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
const G4VBiasingOperation * fPreviousAppliedOccurenceBiasingOperation
const G4VBiasingOperation * fPreviousAppliedNonPhysicsBiasingOperation
const G4VBiasingOperation * fPreviousProposedFinalStateBiasingOperation
const G4VBiasingOperation * fPreviousAppliedFinalStateBiasingOperation
G4VBiasingOperation * fFinalStateBiasingOperation
G4VBiasingOperation * fOccurenceBiasingOperation
const G4VBiasingOperation * fPreviousProposedNonPhysicsBiasingOperation
G4VBiasingOperation * fNonPhysicsBiasingOperation
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 0;
83  else return (*it).second;
84 }
static G4MapCache< const G4LogicalVolume *, G4VBiasingOperator *> fLogicalToSetupMap
iterator End()
Definition: G4Cache.hh:411
iterator Find(const key_type &k)
Definition: G4Cache.hh:417
map_type::iterator iterator
Definition: G4Cache.hh:175
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetBiasingOperators()

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

Definition at line 285 of file G4VBiasingOperator.hh.

285 {return fOperators.Get();}
static G4VectorCache< G4VBiasingOperator *> fOperators
value_type & Get() const
Definition: G4Cache.hh:282
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetName()

const G4String G4VBiasingOperator::GetName ( void  ) const
inline

Definition at line 280 of file G4VBiasingOperator.hh.

280 {return fName;}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetPreviousBiasingAppliedCase()

G4BiasingAppliedCase G4VBiasingOperator::GetPreviousBiasingAppliedCase ( ) const
inline

Definition at line 283 of file G4VBiasingOperator.hh.

G4BiasingAppliedCase fPreviousBiasingAppliedCase

◆ GetPreviousNonPhysicsAppliedOperation()

const G4VBiasingOperation* G4VBiasingOperator::GetPreviousNonPhysicsAppliedOperation ( )
inline

Definition at line 307 of file G4VBiasingOperator.hh.

const G4VBiasingOperation * fPreviousAppliedNonPhysicsBiasingOperation

◆ 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
G4VBiasingOperation * fFinalStateBiasingOperation
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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
G4VBiasingOperation * fNonPhysicsBiasingOperation
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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
G4VBiasingOperation * fOccurenceBiasingOperation
Here is the call graph for this function:
Here is the caller graph for this function:

◆ OperationApplied() [1/2]

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

Reimplemented in G4BOptrForceCollision, and GB02BOptrMultiParticleForceCollision.

Definition at line 178 of file G4VBiasingOperator.cc.

180 {
181 }
Here is the caller graph for this function:

◆ OperationApplied() [2/2]

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

◆ ProposeFinalStateBiasingOperation()

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

Implemented in GB01BOptrMultiParticleChangeCrossSection, GB01BOptrChangeCrossSection, GB03BOptrGeometryBasedBiasing, G4BOptrForceCollision, GB04BOptrBremSplitting, and GB02BOptrMultiParticleForceCollision.

Here is the caller graph for this function:

◆ ProposeNonPhysicsBiasingOperation()

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

Implemented in GB01BOptrMultiParticleChangeCrossSection, GB01BOptrChangeCrossSection, G4BOptrForceCollision, GB03BOptrGeometryBasedBiasing, GB04BOptrBremSplitting, and GB02BOptrMultiParticleForceCollision.

Here is the caller graph for this function:

◆ ProposeOccurenceBiasingOperation()

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

Implemented in GB01BOptrMultiParticleChangeCrossSection, GB01BOptrChangeCrossSection, GB03BOptrGeometryBasedBiasing, G4BOptrForceCollision, GB04BOptrBremSplitting, and GB02BOptrMultiParticleForceCollision.

Here is the caller graph for this function:

◆ ReportOperationApplied() [1/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:
118  fPreviousAppliedNonPhysicsBiasingOperation = operationApplied ;
119  break;
120  case BAC_DenyInteraction: // -- ยงยง will b deprecated
122  break;
123  case BAC_FinalState:
125  break;
126  case BAC_Occurence:
127  G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
128  "BIAS.MNG.02",
129  JustWarning,
130  "Internal logic error, please report !");
131  break;
132  default:
133  G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
134  "BIAS.MNG.03",
135  JustWarning,
136  "Internal logic error, please report !");
137  }
138  OperationApplied( callingProcess, biasingCase, operationApplied, particleChangeProduced );
139 }
virtual void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
G4BiasingAppliedCase fPreviousBiasingAppliedCase
const G4VBiasingOperation * fPreviousAppliedOccurenceBiasingOperation
const G4VBiasingOperation * fPreviousAppliedNonPhysicsBiasingOperation
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4VBiasingOperation * fPreviousAppliedFinalStateBiasingOperation
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ReportOperationApplied() [2/2]

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

Definition at line 141 of file G4VBiasingOperator.cc.

147 {
148  fPreviousBiasingAppliedCase = biasingCase;
149  fPreviousAppliedOccurenceBiasingOperation = occurenceOperationApplied;
150  fPreviousAppliedFinalStateBiasingOperation = finalStateOperationApplied;
151  OperationApplied( callingProcess, biasingCase, occurenceOperationApplied, weightForOccurenceInteraction, finalStateOperationApplied, particleChangeProduced );
152 }
virtual void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
G4BiasingAppliedCase fPreviousBiasingAppliedCase
const G4VBiasingOperation * fPreviousAppliedOccurenceBiasingOperation
const G4VBiasingOperation * fPreviousAppliedFinalStateBiasingOperation
Here is the call graph for this function:

◆ StartRun()

virtual void G4VBiasingOperator::StartRun ( )
inlinevirtual

◆ StartTracking()

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

Friends And Related Function Documentation

◆ G4BiasingOperatorStateNotifier

friend class G4BiasingOperatorStateNotifier
friend

Definition at line 198 of file G4VBiasingOperator.hh.

Member Data Documentation

◆ fDepthInTree

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

Definition at line 325 of file G4VBiasingOperator.hh.

◆ fFinalStateBiasingOperation

G4VBiasingOperation* G4VBiasingOperator::fFinalStateBiasingOperation
private

Definition at line 329 of file G4VBiasingOperator.hh.

◆ fLogicalToSetupMap

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

Definition at line 314 of file G4VBiasingOperator.hh.

◆ fName

const G4String G4VBiasingOperator::fName
private

Definition at line 311 of file G4VBiasingOperator.hh.

◆ fNonPhysicsBiasingOperation

G4VBiasingOperation* G4VBiasingOperator::fNonPhysicsBiasingOperation
private

Definition at line 330 of file G4VBiasingOperator.hh.

◆ fOccurenceBiasingOperation

G4VBiasingOperation* G4VBiasingOperator::fOccurenceBiasingOperation
private

Definition at line 328 of file G4VBiasingOperator.hh.

◆ fOperators

G4VectorCache< G4VBiasingOperator *> G4VBiasingOperator::fOperators
staticprivate

Definition at line 316 of file G4VBiasingOperator.hh.

◆ fPreviousAppliedFinalStateBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousAppliedFinalStateBiasingOperation
private

Definition at line 337 of file G4VBiasingOperator.hh.

◆ fPreviousAppliedNonPhysicsBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousAppliedNonPhysicsBiasingOperation
private

Definition at line 338 of file G4VBiasingOperator.hh.

◆ fPreviousAppliedOccurenceBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousAppliedOccurenceBiasingOperation
private

Definition at line 336 of file G4VBiasingOperator.hh.

◆ fPreviousBiasingAppliedCase

G4BiasingAppliedCase G4VBiasingOperator::fPreviousBiasingAppliedCase
private

Definition at line 339 of file G4VBiasingOperator.hh.

◆ fPreviousProposedFinalStateBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousProposedFinalStateBiasingOperation
private

Definition at line 334 of file G4VBiasingOperator.hh.

◆ fPreviousProposedNonPhysicsBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousProposedNonPhysicsBiasingOperation
private

Definition at line 335 of file G4VBiasingOperator.hh.

◆ fPreviousProposedOccurenceBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousProposedOccurenceBiasingOperation
private

Definition at line 333 of file G4VBiasingOperator.hh.

◆ fRootVolumes

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

Definition at line 324 of file G4VBiasingOperator.hh.

◆ fStateNotifier

G4Cache< G4BiasingOperatorStateNotifier *> G4VBiasingOperator::fStateNotifier
staticprivate

Definition at line 320 of file G4VBiasingOperator.hh.


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