Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 181 of file G4VBiasingOperator.hh.

Constructor & Destructor Documentation

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:282
void Push_back(const value_type &val)
Definition: G4Cache.hh:333
void Put(const value_type &val) const
Definition: G4Cache.hh:286

Here is the call graph for this function:

G4VBiasingOperator::~G4VBiasingOperator ( )
virtual

Definition at line 55 of file G4VBiasingOperator.cc.

56 {
57 }

Member Function Documentation

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 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
iterator End()
Definition: G4Cache.hh:411
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
#define G4endl
Definition: G4ios.hh:61
const G4String & GetName() const
const G4String GetName() const

Here is the call graph for this function:

virtual void G4VBiasingOperator::Configure ( )
inlinevirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 272 of file G4VBiasingOperator.hh.

272 {}
virtual void G4VBiasingOperator::ConfigureForWorker ( )
inlinevirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 275 of file G4VBiasingOperator.hh.

275 {}
virtual void G4VBiasingOperator::EndTracking ( )
inlinevirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 280 of file G4VBiasingOperator.hh.

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

Reimplemented in G4BOptrForceCollision.

Definition at line 173 of file G4VBiasingOperator.cc.

174 {}

Here is the caller graph for this function:

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)

Here is the call graph for this function:

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 }
iterator End()
Definition: G4Cache.hh:411
iterator Find(const key_type &k)
Definition: G4Cache.hh:417

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 294 of file G4VBiasingOperator.hh.

294 {return fOperators.Get();}
value_type & Get() const
Definition: G4Cache.hh:282

Here is the call graph for this function:

Here is the caller graph for this function:

const G4String G4VBiasingOperator::GetName ( ) const
inline

Definition at line 289 of file G4VBiasingOperator.hh.

289 {return fName;}

Here is the caller graph for this function:

G4BiasingAppliedCase G4VBiasingOperator::GetPreviousBiasingAppliedCase ( ) const
inline

Definition at line 292 of file G4VBiasingOperator.hh.

292 {return fPreviousBiasingAppliedCase;}
const G4VBiasingOperation* G4VBiasingOperator::GetPreviousNonPhysicsAppliedOperation ( )
inline

Definition at line 316 of file G4VBiasingOperator.hh.

316 {return fPreviousAppliedNonPhysicsBiasingOperation;}
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

Here is the call graph for this function:

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

Here is the call graph for this function:

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

Here is the call graph for this function:

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 }

Here is the caller graph for this function:

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 }
virtual G4VBiasingOperation* G4VBiasingOperator::ProposeFinalStateBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
protectedpure virtual

Here is the caller graph for this function:

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

Here is the caller graph for this function:

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

Here is the caller graph for this function:

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",
126  JustWarning,
127  "Internal logic error, please report !");
128  break;
129  default:
130  G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
131  "BIAS.MNG.03",
132  JustWarning,
133  "Internal logic error, please report !");
134  }
135  OperationApplied( callingProcess, biasingCase, operationApplied, particleChangeProduced );
136 }
virtual void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

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)

Here is the call graph for this function:

virtual void G4VBiasingOperator::StartRun ( )
inlinevirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 277 of file G4VBiasingOperator.hh.

277 {}
virtual void G4VBiasingOperator::StartTracking ( const G4Track )
inlinevirtual

Reimplemented in G4BOptrForceCollision.

Definition at line 279 of file G4VBiasingOperator.hh.

279 {}

Friends And Related Function Documentation

friend class G4BiasingOperatorStateNotifier
friend

Definition at line 186 of file G4VBiasingOperator.hh.


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