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

#include <Doxymodules_biasing.h>

Inheritance diagram for GB06BOptnSplitAndKillByImportance:
G4VBiasingOperation

Public Member Functions

 GB06BOptnSplitAndKillByImportance (G4String name)
 
virtual ~GB06BOptnSplitAndKillByImportance ()
 
virtual const G4VBiasingInteractionLawProvideOccurenceBiasingInteractionLaw (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 SetParallelWorldIndex (G4int parallelWorldIndex)
 
G4int GetParallelWorldIndex () const
 
void SetBiasingSharedData (const G4BiasingProcessSharedData *sharedData)
 
void SetImportanceMap (std::map< G4int, G4int > *importanceMap)
 

Private Attributes

G4int fParallelWorldIndex
 
const G4BiasingProcessSharedDatafBiasingSharedData
 
G4TouchableHandle fPreStepTouchableHistory
 
G4TouchableHandle fPostStepTouchableHistory
 
G4ParticleChange fParticleChange
 
G4ParticleChangeForNothing fDummyParticleChange
 
std::map< G4int, G4int > * fImportanceMap
 

Detailed Description

Definition at line 138 of file Doxymodules_biasing.h.

Constructor & Destructor Documentation

◆ GB06BOptnSplitAndKillByImportance()

GB06BOptnSplitAndKillByImportance::GB06BOptnSplitAndKillByImportance ( G4String  name)

◆ ~GB06BOptnSplitAndKillByImportance()

GB06BOptnSplitAndKillByImportance::~GB06BOptnSplitAndKillByImportance ( )
virtual

Definition at line 52 of file GB06BOptnSplitAndKillByImportance.cc.

53{}

Member Function Documentation

◆ ProvideOccurenceBiasingInteractionLaw()

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

Definition at line 54 of file GB06BOptnSplitAndKillByImportance.hh.

56 {return 0;}

◆ ApplyFinalStateBiasing()

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

Definition at line 58 of file GB06BOptnSplitAndKillByImportance.hh.

62 {return 0;}

◆ DistanceToApplyOperation()

G4double GB06BOptnSplitAndKillByImportance::DistanceToApplyOperation ( const G4Track ,
G4double  ,
G4ForceCondition *  condition 
)
finalvirtual

Definition at line 57 of file GB06BOptnSplitAndKillByImportance.cc.

61{
62
63 // -- Remember the touchable history (ie geometry state) at the beginning of the step:
64 // -- Start by getting the process handling the step limitation in parallel
65 // -- geometries for the generic biasing:
66 auto biasingLimiterProcess = fBiasingSharedData->GetParallelGeometriesLimiterProcess();
68 biasingLimiterProcess
69 ->GetNavigator( fParallelWorldIndex ) // -- get the navigator of the geometry...
70 ->CreateTouchableHistoryHandle(); // -- ...and collect the geometry state.
71
72 // -- We request to be "forced" : meaning GenerateBiasingFinalState(...) will be called
73 // -- anyway at the PostStepDoIt(...) stage (ie, when the track will have been moved to
74 // -- its end point position) and, there, we will have to handle the decision to
75 // -- split/kill if we are on a volume boundary, or do nothing, if we're not:
76 *condition = Forced;
77
78 // -- As we're forced, we make this returned length void:
79 return DBL_MAX;
80}

◆ GenerateBiasingFinalState()

G4VParticleChange * GB06BOptnSplitAndKillByImportance::GenerateBiasingFinalState ( const G4Track track,
const G4Step  
)
finalvirtual

Definition at line 85 of file GB06BOptnSplitAndKillByImportance.cc.

87{
88 // -- Given we used the "Forced" condition, this method is always called.
89 // -- We check the status of the step in the parallel geometry, and apply
90 // -- splitting/killing if the step has been limited by the geometry.
91
92 // -- We first check if we limit the step in the parallel geometry:
93 G4bool isLimiting = fBiasingSharedData
94 ->GetParallelGeometriesLimiterProcess()
95 ->GetIsLimiting( fParallelWorldIndex );
96
97 // -- if not limiting, we leave the track unchanged:
98 if ( !isLimiting )
99 {
100 fDummyParticleChange.Initialize( *track );
101 return &fDummyParticleChange;
102 }
103
104 // -- We collect the geometry state at the end point step:
105 // -- Note this is the same call than in the DistanceToApplyOperation(...) for the
106 // -- fPreStepTouchableHistory, but the navigator has pushed the track in the next
107 // -- volume since then (even if the track is still on the boundary), and hence the
108 // -- geometry state has changed.
109 auto biasingLimiterProcess = fBiasingSharedData->GetParallelGeometriesLimiterProcess();
111 biasingLimiterProcess
112 ->GetNavigator( fParallelWorldIndex )
113 ->CreateTouchableHistoryHandle();
114
115
116 // -- We verify we are still in the "same" physical volume:
117 // -- remember : using a replica, we have "one" physical volume
118 // -- but representing many different placements, distinguished
119 // -- by replica number. By checking we are in the same physical
120 // -- volume, we verify the end step point has not reached the
121 // -- world volume of the parallel geometry:
122 if ( fPreStepTouchableHistory ->GetVolume() !=
123 fPostStepTouchableHistory->GetVolume() )
124 {
125 // -- the track is leaving the volumes in which biasing is applied,
126 // -- we leave this track unchanged:
127 fDummyParticleChange.Initialize( *track );
128 return &fDummyParticleChange;
129 }
130
131 // -------------------------------------------------------------------------------------
132 // -- At this stage, we know we have a particle crossing a boundary between two slices,
133 // -- each of this slice has a well defined importance : we apply the biasing.
134 // -- We will split if the track is entering a volume with higher importance, and
135 // -- kill (applying Russian roulette) in the other case.
136 // -------------------------------------------------------------------------------------
137
138 // -- We start by getting the replica numbers:
139 G4int preReplicaNumber = fPreStepTouchableHistory->GetReplicaNumber();
140 G4int postReplicaNumber = fPostStepTouchableHistory->GetReplicaNumber();
141
142 // -- and get the related importance we defined in the importance map:
143 G4int preImportance = (*fImportanceMap)[ preReplicaNumber];
144 G4int postImportance = (*fImportanceMap)[postReplicaNumber];
145
146
147 // -- Get track weight:
148 G4double initialWeight = track->GetWeight();
149
150 // -- Initialize the "particle change" (it will communicate the new track state to
151 // -- the tracking):
152 fParticleChange.Initialize(*track);
153
154 if ( postImportance > preImportance )
155 {
156 // -- We split :
157 G4int splittingFactor = postImportance/preImportance;
158
159 // Define the tracks weight:
160 G4double weightOfTrack = initialWeight/splittingFactor;
161
162 // Ask currect track weight to be changed to the new value:
163 fParticleChange.ProposeParentWeight( weightOfTrack );
164 // Now we clone this track (this is the actual splitting):
165 // we will then have the primary and clone of it, hence the
166 // splitting by a factor 2:
167 G4Track* clone = new G4Track( *track );
168 clone->SetWeight( weightOfTrack );
169 fParticleChange.AddSecondary( clone );
170 // -- Below's call added for safety & illustration : inform particle change to not
171 // -- modify the clone (ie : daughter) weight to make it that of the primary.
172 // -- Here this call is not mandatory as both tracks have same weights.
173 fParticleChange.SetSecondaryWeightByProcess(true);
174 }
175 else
176 {
177 // -- We apply Russian roulette:
178 G4double survivingProbability = G4double(postImportance) / G4double(preImportance);
179
180 // Shoot a random number (in ]0,1[ segment):
181 G4double random = G4UniformRand();
182 if ( random > survivingProbability )
183 {
184 // We ask for the the track to be killed:
185 fParticleChange.ProposeTrackStatus(fStopAndKill);
186 }
187 else
188 {
189 // In this case, the track survives. We change its weight
190 // to conserve weight among killed and survival tracks:
191 fParticleChange.ProposeParentWeight( initialWeight/survivingProbability );
192 }
193 }
194
195 return &fParticleChange;
196
197}

◆ SetParallelWorldIndex()

void GB06BOptnSplitAndKillByImportance::SetParallelWorldIndex ( G4int  parallelWorldIndex)
inline

Definition at line 86 of file GB06BOptnSplitAndKillByImportance.hh.

87 { fParallelWorldIndex = parallelWorldIndex; }

◆ GetParallelWorldIndex()

G4int GB06BOptnSplitAndKillByImportance::GetParallelWorldIndex ( ) const
inline

Definition at line 88 of file GB06BOptnSplitAndKillByImportance.hh.

89 { return fParallelWorldIndex; }

◆ SetBiasingSharedData()

void GB06BOptnSplitAndKillByImportance::SetBiasingSharedData ( const G4BiasingProcessSharedData sharedData)
inline

Definition at line 92 of file GB06BOptnSplitAndKillByImportance.hh.

93 { fBiasingSharedData = sharedData; }

◆ SetImportanceMap()

void GB06BOptnSplitAndKillByImportance::SetImportanceMap ( std::map< G4int, G4int > *  importanceMap)
inline

Definition at line 95 of file GB06BOptnSplitAndKillByImportance.hh.

96 { fImportanceMap = importanceMap; }

Member Data Documentation

◆ fParallelWorldIndex

G4int GB06BOptnSplitAndKillByImportance::fParallelWorldIndex
private

Definition at line 100 of file GB06BOptnSplitAndKillByImportance.hh.

◆ fBiasingSharedData

const G4BiasingProcessSharedData* GB06BOptnSplitAndKillByImportance::fBiasingSharedData
private

Definition at line 101 of file GB06BOptnSplitAndKillByImportance.hh.

◆ fPreStepTouchableHistory

G4TouchableHandle GB06BOptnSplitAndKillByImportance::fPreStepTouchableHistory
private

Definition at line 102 of file GB06BOptnSplitAndKillByImportance.hh.

◆ fPostStepTouchableHistory

G4TouchableHandle GB06BOptnSplitAndKillByImportance::fPostStepTouchableHistory
private

Definition at line 103 of file GB06BOptnSplitAndKillByImportance.hh.

◆ fParticleChange

G4ParticleChange GB06BOptnSplitAndKillByImportance::fParticleChange
private

Definition at line 104 of file GB06BOptnSplitAndKillByImportance.hh.

◆ fDummyParticleChange

G4ParticleChangeForNothing GB06BOptnSplitAndKillByImportance::fDummyParticleChange
private

Definition at line 105 of file GB06BOptnSplitAndKillByImportance.hh.

◆ fImportanceMap

std::map< G4int, G4int >* GB06BOptnSplitAndKillByImportance::fImportanceMap
private

Definition at line 106 of file GB06BOptnSplitAndKillByImportance.hh.


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

Applications | User Support | Publications | Collaboration