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

#include <Doxymodules_biasing.h>

Inheritance diagram for GB03BOptnSplitOrKillOnBoundary:
G4VBiasingOperation

Public Member Functions

 GB03BOptnSplitOrKillOnBoundary (G4String name)
 
virtual ~GB03BOptnSplitOrKillOnBoundary ()
 
virtual const G4VBiasingInteractionLawProvideOccurenceBiasingInteractionLaw (const G4BiasingProcessInterface *, G4ForceCondition &)
 
virtual G4VParticleChangeApplyFinalStateBiasing (const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)
 
virtual G4double DistanceToApplyOperation (const G4Track *, G4double, G4ForceCondition *condition)
 
virtual G4VParticleChangeGenerateBiasingFinalState (const G4Track *, const G4Step *)
 
void SetSplittingFactor (G4int splittingFactor)
 
void SetApplyProbability (G4double proba)
 
G4int GetSplittingFactor () const
 
G4double GetApplyProbability () const
 

Private Attributes

G4ParticleChange fParticleChange
 
G4ParticleChangeForNothing fParticleChangeForNothing
 
G4int fSplittingFactor
 
G4double fApplyProbability
 

Detailed Description

Definition at line 94 of file Doxymodules_biasing.h.

Constructor & Destructor Documentation

◆ GB03BOptnSplitOrKillOnBoundary()

GB03BOptnSplitOrKillOnBoundary::GB03BOptnSplitOrKillOnBoundary ( G4String  name)

◆ ~GB03BOptnSplitOrKillOnBoundary()

GB03BOptnSplitOrKillOnBoundary::~GB03BOptnSplitOrKillOnBoundary ( )
virtual

Definition at line 43 of file GB03BOptnSplitOrKillOnBoundary.cc.

44{}

Member Function Documentation

◆ ProvideOccurenceBiasingInteractionLaw()

virtual const G4VBiasingInteractionLaw * GB03BOptnSplitOrKillOnBoundary::ProvideOccurenceBiasingInteractionLaw ( const G4BiasingProcessInterface ,
G4ForceCondition &   
)
inlinevirtual

Definition at line 52 of file GB03BOptnSplitOrKillOnBoundary.hh.

53 {return 0;}

◆ ApplyFinalStateBiasing()

virtual G4VParticleChange * GB03BOptnSplitOrKillOnBoundary::ApplyFinalStateBiasing ( const G4BiasingProcessInterface ,
const G4Track ,
const G4Step ,
G4bool &   
)
inlinevirtual

Definition at line 55 of file GB03BOptnSplitOrKillOnBoundary.hh.

58 {return 0;}

◆ DistanceToApplyOperation()

G4double GB03BOptnSplitOrKillOnBoundary::DistanceToApplyOperation ( const G4Track ,
G4double  ,
G4ForceCondition *  condition 
)
virtual

Definition at line 48 of file GB03BOptnSplitOrKillOnBoundary.cc.

52{
53 // -- return "infinite" distance for interaction, but asks for GenerateBiasingFinalState
54 // -- being called anyway at the end of the step, by returning the "Forced" condition
55 // -- flag.
56 *condition = Forced;
57 return DBL_MAX;
58}

◆ GenerateBiasingFinalState()

G4VParticleChange * GB03BOptnSplitOrKillOnBoundary::GenerateBiasingFinalState ( const G4Track track,
const G4Step step 
)
virtual

Definition at line 63 of file GB03BOptnSplitOrKillOnBoundary.cc.

65{
66
67 // Check if step is limited by the geometry: as we attached the biasing operator
68 // to the absorber layer, this volume boundary is the one of the absorber.
69 // (check of current step # of track is inelegant, but is to fix a "feature"
70 // that a cloned track can wrongly be seen in the wrong volume, because of numerical
71 // precision issue. In this case it makes a tiny step, which should be disregarded).
72 if ( ( step->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ) &&
73 ( track->GetCurrentStepNumber() != 1 ) )
74 {
75
76 // -- Before deciding for killing or splitting, we make decision on applying
77 // -- the technique or not:
78 G4double trial = G4UniformRand(); // -- Note: G4UniformRand() is thread-safe
79 // -- engine for random numbers
80 if ( trial <= fApplyProbability )
81 {
82 // -- Get z component of track, to see if it moves forward or backward:
83 G4double pz = track->GetMomentum().z();
84
85 if ( pz > 0.0 )
86 {
87 // -------------------------------------------------
88 // Here, we are moving "forward". We do "splitting":
89 // -------------------------------------------------
90
91 // Get track weight:
92 G4double initialWeight = track->GetWeight();
93 // Define the tracks weight:
94 G4double weightOfTrack = initialWeight/fSplittingFactor;
95
96 // The "particle change" is the object to be used to communicate to
97 // the tracking the update of the primary state and/or creation
98 // secondary tracks.
99 fParticleChange.Initialize(*track);
100
101 // ask currect track weight to be changed to new value:
102 fParticleChange.ProposeParentWeight( weightOfTrack );
103
104 // Now make clones of this track (this is the actual splitting):
105 // we will then have the primary and N-1 clones of it, hence the
106 // splitting by a factor N:
107 fParticleChange.SetNumberOfSecondaries( fSplittingFactor-1 );
108 for ( G4int iSplit = 1 ; iSplit < fSplittingFactor ; iSplit++ )
109 {
110 G4Track* clone = new G4Track( *track );
111 clone->SetWeight( weightOfTrack );
112 fParticleChange.AddSecondary( clone );
113 }
114 fParticleChange.SetSecondaryWeightByProcess(true); // -- tricky
115 // -- take it as is ;) [though not necessary here, put for safety]
116
117 // this new final state is returned to the tracking;
118 return &fParticleChange;
119
120 }
121
122 else
123
124 {
125 // --------------------------------------------------------------
126 // Here, we are moving backward. We do killing, playing a russian
127 // roulette, killing 1/fSplittingFactor of the tracks in average:
128 // --------------------------------------------------------------
129
130 // Get track weight:
131 G4double initialWeight = track->GetWeight();
132
133 // The "particle change" is the object to be used to communicate to
134 // the tracking the update of the primary state and/or creation
135 // secondary tracks.
136 fParticleChange.Initialize(*track);
137
138 // Shoot a random number (in ]0,1[ segment):
139 G4double random = G4UniformRand();
140
141 // Decide to kill, keeping 1/fSplittingFactor of tracks:
142 G4double survivingProbability = 1.0/fSplittingFactor;
143 if ( random > survivingProbability )
144 {
145 // We ask for the the track to be killed:
146 fParticleChange.ProposeTrackStatus(fStopAndKill);
147 }
148 else
149 {
150 // In this case, the track survives. We change its weight
151 // to conserve weight among killed and survival tracks:
152 fParticleChange.ProposeParentWeight( initialWeight*fSplittingFactor );
153 }
154
155 // this new final state is returned to the tracking;
156 return &fParticleChange;
157 }
158 } // -- end of : if ( trial > probaForApplying )
159 } // -- end of : if ( ( step->GetPostStepPoint()->GetStepStatus() ==
160 // fGeomBoundary ) ...
161
162
163 // Here, the step was not limited by the geometry (but certainly by a physics
164 // process). We do nothing: ie we make no change to the current track.
165 fParticleChangeForNothing.Initialize(*track);
167
168}

◆ SetSplittingFactor()

void GB03BOptnSplitOrKillOnBoundary::SetSplittingFactor ( G4int  splittingFactor)
inline

Definition at line 79 of file GB03BOptnSplitOrKillOnBoundary.hh.

80 { fSplittingFactor = splittingFactor; }

◆ SetApplyProbability()

void GB03BOptnSplitOrKillOnBoundary::SetApplyProbability ( G4double  proba)
inline

Definition at line 93 of file GB03BOptnSplitOrKillOnBoundary.hh.

94 { fApplyProbability = proba; }

◆ GetSplittingFactor()

G4int GB03BOptnSplitOrKillOnBoundary::GetSplittingFactor ( ) const
inline

Definition at line 96 of file GB03BOptnSplitOrKillOnBoundary.hh.

96{ return fSplittingFactor; }

◆ GetApplyProbability()

G4double GB03BOptnSplitOrKillOnBoundary::GetApplyProbability ( ) const
inline

Definition at line 97 of file GB03BOptnSplitOrKillOnBoundary.hh.

97{ return fApplyProbability; }

Member Data Documentation

◆ fParticleChange

G4ParticleChange GB03BOptnSplitOrKillOnBoundary::fParticleChange
private

Definition at line 100 of file GB03BOptnSplitOrKillOnBoundary.hh.

◆ fParticleChangeForNothing

G4ParticleChangeForNothing GB03BOptnSplitOrKillOnBoundary::fParticleChangeForNothing
private

Definition at line 101 of file GB03BOptnSplitOrKillOnBoundary.hh.

◆ fSplittingFactor

G4int GB03BOptnSplitOrKillOnBoundary::fSplittingFactor
private

Definition at line 102 of file GB03BOptnSplitOrKillOnBoundary.hh.

◆ fApplyProbability

G4double GB03BOptnSplitOrKillOnBoundary::fApplyProbability
private

Definition at line 103 of file GB03BOptnSplitOrKillOnBoundary.hh.


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

Applications | User Support | Publications | Collaboration