Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
GB04BOptrBremSplitting Class Reference

#include <GB04BOptrBremSplitting.hh>

Inheritance diagram for GB04BOptrBremSplitting:
Collaboration diagram for GB04BOptrBremSplitting:

Public Member Functions

 GB04BOptrBremSplitting ()
 
virtual ~GB04BOptrBremSplitting ()
 
virtual void StartRun ()
 
virtual void StartTracking (const G4Track *track)
 
- Public Member Functions inherited from G4VBiasingOperator
 G4VBiasingOperator (G4String name)
 
virtual ~G4VBiasingOperator ()
 
virtual void Configure ()
 
virtual void ConfigureForWorker ()
 
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 ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VBiasingOperator
static const std::vector
< G4VBiasingOperator * > & 
GetBiasingOperators ()
 
static G4VBiasingOperatorGetBiasingOperator (const G4LogicalVolume *)
 
- Protected Member Functions inherited from G4VBiasingOperator
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)
 

Detailed Description

Definition at line 37 of file GB04BOptrBremSplitting.hh.

Constructor & Destructor Documentation

GB04BOptrBremSplitting::GB04BOptrBremSplitting ( )

Definition at line 39 of file GB04BOptrBremSplitting.cc.

40 : G4VBiasingOperator("BremSplittingOperator"),
41  fSplittingFactor(1),
42  fBiasPrimaryOnly(true),
43  fBiasOnlyOnce(true)
44 {
45  fBremSplittingOperation = new GB04BOptnBremSplitting("BremSplittingOperation");
46 
47  // -- Define messengers:
48  // -- Splitting factor:
49  fSplittingFactorMessenger =
50  new G4GenericMessenger(this, "/GB04/biasing/","Biasing control" );
51  G4GenericMessenger::Command& splittingFactorCmd =
52  fSplittingFactorMessenger->DeclareProperty("setSplittingFactor", fSplittingFactor,
53  "Define the brem. splitting factor." );
54  splittingFactorCmd.SetStates(G4State_Idle);
55  // -- Bias ony primary particle:
56  fBiasPrimaryOnlyMessenger =
57  new G4GenericMessenger(this, "/GB04/biasing/","Biasing control" );
58  G4GenericMessenger::Command& biasPrimaryCmd =
59  fBiasPrimaryOnlyMessenger->DeclareProperty("biasPrimaryOnly", fBiasPrimaryOnly,
60  "Chose if brem. splitting applies to primary particles only." );
61  biasPrimaryCmd.SetStates(G4State_Idle);
62  // -- Bias ony primary particle:
63  fBiasOnlyOnceMessenger =
64  new G4GenericMessenger(this, "/GB04/biasing/","Biasing control" );
65  G4GenericMessenger::Command& biasOnlyOnceCmd =
66  fBiasPrimaryOnlyMessenger->DeclareProperty("biasOnlyOnce", fBiasOnlyOnce,
67  "Chose if apply the brem. splitting only once for the track." );
68  biasOnlyOnceCmd.SetStates(G4State_Idle);
69 
70 }
This class is generic messenger.
Command & DeclareProperty(const G4String &name, const G4AnyType &variable, const G4String &doc="")
Declare Methods.
Command & SetStates(G4ApplicationState s0)
G4VBiasingOperator(G4String name)

Here is the call graph for this function:

virtual GB04BOptrBremSplitting::~GB04BOptrBremSplitting ( )
inlinevirtual

Definition at line 40 of file GB04BOptrBremSplitting.hh.

40 {}

Member Function Documentation

void GB04BOptrBremSplitting::StartRun ( )
virtual

Reimplemented from G4VBiasingOperator.

Definition at line 74 of file GB04BOptrBremSplitting.cc.

75 {
76  fBremSplittingOperation->SetSplittingFactor ( fSplittingFactor );
77  G4cout << GetName() << " : starting run with brem. splitting factor = "
78  << fSplittingFactor;
79  if ( fBiasPrimaryOnly ) G4cout << ", biasing only primaries ";
80  else G4cout << ", biasing primary and secondary tracks ";
81  if ( fBiasOnlyOnce ) G4cout << ", biasing only once per track ";
82  else G4cout << ", biasing several times per track ";
83  G4cout << " . " << G4endl;
84 }
void SetSplittingFactor(G4int splittingFactor)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
const G4String GetName() const

Here is the call graph for this function:

void GB04BOptrBremSplitting::StartTracking ( const G4Track track)
virtual

Reimplemented from G4VBiasingOperator.

Definition at line 88 of file GB04BOptrBremSplitting.cc.

89 {
90  // -- reset the number of times the brem. splitting was applied:
91  fNInteractions = 0;
92 }

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