Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Private Attributes | List of all members
GB03BOptrGeometryBasedBiasing Class Reference

Biasing operator class. More...

#include <Doxymodules_biasing.h>

Inheritance diagram for GB03BOptrGeometryBasedBiasing:
G4VBiasingOperator

Public Member Functions

 GB03BOptrGeometryBasedBiasing ()
 
virtual ~GB03BOptrGeometryBasedBiasing ()
 
GB03BOptnSplitOrKillOnBoundaryGetSplitAndKillOperation () const
 
void StartRun ()
 

Private Member Functions

virtual G4VBiasingOperationProposeNonPhysicsBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
virtual G4VBiasingOperationProposeOccurenceBiasingOperation (const G4Track *, const G4BiasingProcessInterface *)
 
virtual G4VBiasingOperationProposeFinalStateBiasingOperation (const G4Track *, const G4BiasingProcessInterface *)
 

Private Attributes

GB03BOptnSplitOrKillOnBoundaryfSplitAndKillOperation
 
G4int fSplittingFactor
 
G4double fApplyProbability
 
G4GenericMessengerfSplittingFactorMessenger
 
G4GenericMessengerfApplyProbabilityMessenger
 

Detailed Description

Biasing operator class.

Definition at line 95 of file Doxymodules_biasing.h.

Constructor & Destructor Documentation

◆ GB03BOptrGeometryBasedBiasing()

GB03BOptrGeometryBasedBiasing::GB03BOptrGeometryBasedBiasing ( )

Definition at line 37 of file GB03BOptrGeometryBasedBiasing.cc.

38: G4VBiasingOperator("GB03BOptrGeometryBasedBiasing"),
41{
43
44 // -- Define messengers:
46 new G4GenericMessenger(this, "/GB03/biasing/","Biasing control" );
47
48 G4GenericMessenger::Command& splittingFactorCmd =
49 fSplittingFactorMessenger->DeclareProperty("setSplittingFactor", fSplittingFactor,
50 "Define the splitting factor." );
51 splittingFactorCmd.SetStates(G4State_Idle);
52
54 new G4GenericMessenger(this, "/GB03/biasing/","Biasing control" );
55
56 G4GenericMessenger::Command& applyProbCmd =
57 fApplyProbabilityMessenger->DeclareProperty("setApplyProbability", fApplyProbability,
58 "Define the probability to apply the splitting/killing." );
59 applyProbCmd.SetStates(G4State_Idle);
60
61}
GB03BOptnSplitOrKillOnBoundary * fSplitAndKillOperation

◆ ~GB03BOptrGeometryBasedBiasing()

GB03BOptrGeometryBasedBiasing::~GB03BOptrGeometryBasedBiasing ( )
virtual

Definition at line 65 of file GB03BOptrGeometryBasedBiasing.cc.

66{
68}

Member Function Documentation

◆ GetSplitAndKillOperation()

GB03BOptnSplitOrKillOnBoundary * GB03BOptrGeometryBasedBiasing::GetSplitAndKillOperation ( ) const
inline

Definition at line 49 of file GB03BOptrGeometryBasedBiasing.hh.

50 { return fSplitAndKillOperation; }

◆ StartRun()

void GB03BOptrGeometryBasedBiasing::StartRun ( )

Definition at line 72 of file GB03BOptrGeometryBasedBiasing.cc.

73{
76 G4cout << GetName() << " : starting run with splitting factor = " << fSplittingFactor
77 << ", and probability for applying the technique " << fApplyProbability
78 << " . " << G4endl;
79}

◆ ProposeNonPhysicsBiasingOperation()

G4VBiasingOperation * GB03BOptrGeometryBasedBiasing::ProposeNonPhysicsBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
privatevirtual

Definition at line 84 of file GB03BOptrGeometryBasedBiasing.cc.

87{
88 // Here, we always return the split and kill biasing operation:
90}

◆ ProposeOccurenceBiasingOperation()

virtual G4VBiasingOperation * GB03BOptrGeometryBasedBiasing::ProposeOccurenceBiasingOperation ( const G4Track ,
const G4BiasingProcessInterface  
)
inlineprivatevirtual

Definition at line 68 of file GB03BOptrGeometryBasedBiasing.hh.

69 { return 0; }

◆ ProposeFinalStateBiasingOperation()

virtual G4VBiasingOperation * GB03BOptrGeometryBasedBiasing::ProposeFinalStateBiasingOperation ( const G4Track ,
const G4BiasingProcessInterface  
)
inlineprivatevirtual

Definition at line 72 of file GB03BOptrGeometryBasedBiasing.hh.

73 { return 0; }

Member Data Documentation

◆ fSplitAndKillOperation

GB03BOptnSplitOrKillOnBoundary* GB03BOptrGeometryBasedBiasing::fSplitAndKillOperation
private

Definition at line 76 of file GB03BOptrGeometryBasedBiasing.hh.

◆ fSplittingFactor

G4int GB03BOptrGeometryBasedBiasing::fSplittingFactor
private

Definition at line 77 of file GB03BOptrGeometryBasedBiasing.hh.

◆ fApplyProbability

G4double GB03BOptrGeometryBasedBiasing::fApplyProbability
private

Definition at line 78 of file GB03BOptrGeometryBasedBiasing.hh.

◆ fSplittingFactorMessenger

G4GenericMessenger* GB03BOptrGeometryBasedBiasing::fSplittingFactorMessenger
private

Definition at line 80 of file GB03BOptrGeometryBasedBiasing.hh.

◆ fApplyProbabilityMessenger

G4GenericMessenger* GB03BOptrGeometryBasedBiasing::fApplyProbabilityMessenger
private

Definition at line 81 of file GB03BOptrGeometryBasedBiasing.hh.


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

Applications | User Support | Publications | Collaboration