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

#include <GB05BOptnSplitAndKillByCrossSection.hh>

Inheritance diagram for GB05BOptnSplitAndKillByCrossSection:
Collaboration diagram for GB05BOptnSplitAndKillByCrossSection:

Public Member Functions

 GB05BOptnSplitAndKillByCrossSection (G4String name)
 
virtual ~GB05BOptnSplitAndKillByCrossSection ()
 
virtual const
G4VBiasingInteractionLaw
ProvideOccurenceBiasingInteractionLaw (const G4BiasingProcessInterface *, G4ForceCondition &) final
 
virtual G4VParticleChangeApplyFinalStateBiasing (const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &) final
 
virtual G4double DistanceToApplyOperation (const G4Track *, G4double, G4ForceCondition *condition) final
 
virtual G4VParticleChangeGenerateBiasingFinalState (const G4Track *, const G4Step *) final
 
void SetInteractionLength (G4double interactionLength)
 
- Public Member Functions inherited from G4VBiasingOperation
 G4VBiasingOperation (G4String name)
 
virtual ~G4VBiasingOperation ()
 
virtual G4double ProposeAlongStepLimit (const G4BiasingProcessInterface *)
 
virtual G4GPILSelection ProposeGPILSelection (const G4GPILSelection wrappedProcessSelection)
 
virtual void AlongMoveBy (const G4BiasingProcessInterface *, const G4Step *, G4double)
 
const G4StringGetName () const
 
std::size_t GetUniqueID () const
 

Detailed Description

Definition at line 39 of file GB05BOptnSplitAndKillByCrossSection.hh.

Constructor & Destructor Documentation

GB05BOptnSplitAndKillByCrossSection::GB05BOptnSplitAndKillByCrossSection ( G4String  name)

Definition at line 36 of file GB05BOptnSplitAndKillByCrossSection.cc.

37 : G4VBiasingOperation(name),
38  fParticleChange(),
39  fInteractionLength(-1.0)
40 {}
G4VBiasingOperation(G4String name)
GB05BOptnSplitAndKillByCrossSection::~GB05BOptnSplitAndKillByCrossSection ( )
virtual

Definition at line 44 of file GB05BOptnSplitAndKillByCrossSection.cc.

45 {}

Member Function Documentation

virtual G4VParticleChange* GB05BOptnSplitAndKillByCrossSection::ApplyFinalStateBiasing ( const G4BiasingProcessInterface ,
const G4Track ,
const G4Step ,
G4bool  
)
inlinefinalvirtual

Implements G4VBiasingOperation.

Definition at line 56 of file GB05BOptnSplitAndKillByCrossSection.hh.

60  {return 0;}
G4double GB05BOptnSplitAndKillByCrossSection::DistanceToApplyOperation ( const G4Track ,
G4double  ,
G4ForceCondition condition 
)
finalvirtual

Implements G4VBiasingOperation.

Definition at line 50 of file GB05BOptnSplitAndKillByCrossSection.cc.

53 {
55 
56  // -- Sample the exponential law using the total interaction length of processes
57  // -- to couterbalance for:
58  G4double proposedStepLength = -std::log( G4UniformRand() ) * fInteractionLength;
59 
60  return proposedStepLength;
61 }
G4double condition(const G4ErrorSymMatrix &m)
#define G4UniformRand()
Definition: Randomize.hh:97
double G4double
Definition: G4Types.hh:76
G4VParticleChange * GB05BOptnSplitAndKillByCrossSection::GenerateBiasingFinalState ( const G4Track track,
const G4Step  
)
finalvirtual

Implements G4VBiasingOperation.

Definition at line 66 of file GB05BOptnSplitAndKillByCrossSection.cc.

68 {
69 
70  // -- This method is called if we have limited the step.
71  // -- We hence make the splitting or killing.
72 
73  // Get track weight:
74  G4double initialWeight = track->GetWeight();
75 
76  // The "particle change" is the object to be used to communicate to
77  // the tracking the update of the primary state and/or creation
78  // secondary tracks.
79  fParticleChange.Initialize(*track);
80 
81  // -- Splitting and killing factors.
82  // -- They are taken the same, but the killing factor can be make bigger.
83  G4double splittingFactor = 2.0;
84  G4double killingFactor = 2.0;
85 
86 
87  if ( track->GetMomentumDirection().z() > 0 )
88  {
89  // -- We split if the track is moving forward:
90 
91  // Define the tracks weight:
92  G4double weightOfTrack = initialWeight/splittingFactor;
93 
94  // Ask currect track weight to be changed to new value:
95  fParticleChange.ProposeParentWeight( weightOfTrack );
96  // Now make clones of this track (this is the actual splitting):
97  // we will then have the primary and clone of it, hence the
98  // splitting by a factor 2:
99  G4Track* clone = new G4Track( *track );
100  clone->SetWeight( weightOfTrack );
101  fParticleChange.AddSecondary( clone );
102  // -- Below's call added for safety & illustration : inform particle change to not
103  // -- modify the clone (ie : daughter) weight to male it that of the
104  // -- primary. Here call not mandatory and both tracks have same weights.
105  fParticleChange.SetSecondaryWeightByProcess(true);
106  }
107  else
108  {
109  // -- We apply Russian roulette if the track is moving backward:
110 
111  // Shoot a random number (in ]0,1[ segment):
112  G4double random = G4UniformRand();
113  G4double killingProbability = 1.0 - 1.0/killingFactor;
114  if ( random < killingProbability )
115  {
116  // We ask for the the track to be killed:
117  fParticleChange.ProposeTrackStatus(fStopAndKill);
118  }
119  else
120  {
121  // In this case, the track survives. We change its weight
122  // to conserve weight among killed and survival tracks:
123  fParticleChange.ProposeParentWeight( initialWeight*killingFactor );
124  }
125  }
126 
127  return &fParticleChange;
128 
129 }
void ProposeParentWeight(G4double finalWeight)
double z() const
void SetWeight(G4double aValue)
void SetSecondaryWeightByProcess(G4bool)
#define G4UniformRand()
Definition: Randomize.hh:97
virtual void Initialize(const G4Track &)
const G4ThreeVector & GetMomentumDirection() const
void AddSecondary(G4Track *aSecondary)
G4double GetWeight() const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)

Here is the call graph for this function:

virtual const G4VBiasingInteractionLaw* GB05BOptnSplitAndKillByCrossSection::ProvideOccurenceBiasingInteractionLaw ( const G4BiasingProcessInterface ,
G4ForceCondition  
)
inlinefinalvirtual

Implements G4VBiasingOperation.

Definition at line 52 of file GB05BOptnSplitAndKillByCrossSection.hh.

54  {return 0;}
void GB05BOptnSplitAndKillByCrossSection::SetInteractionLength ( G4double  interactionLength)
inline

Definition at line 81 of file GB05BOptnSplitAndKillByCrossSection.hh.

82  {
83  fInteractionLength = interactionLength;
84  }

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