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

#include <GB06BOptnSplitAndKillByImportance.hh>

Inheritance diagram for GB06BOptnSplitAndKillByImportance:
Collaboration diagram for GB06BOptnSplitAndKillByImportance:

Public Member Functions

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

Constructor & Destructor Documentation

GB06BOptnSplitAndKillByImportance::GB06BOptnSplitAndKillByImportance ( G4String  name)

Definition at line 43 of file GB06BOptnSplitAndKillByImportance.cc.

44 : G4VBiasingOperation (name),
45  fParallelWorldIndex ( -1 ),
46  fBiasingSharedData ( nullptr ),
47  fParticleChange(),
48  fDummyParticleChange()
49 {}
G4VBiasingOperation(G4String name)
GB06BOptnSplitAndKillByImportance::~GB06BOptnSplitAndKillByImportance ( )
virtual

Definition at line 53 of file GB06BOptnSplitAndKillByImportance.cc.

54 {}

Member Function Documentation

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

Implements G4VBiasingOperation.

Definition at line 59 of file GB06BOptnSplitAndKillByImportance.hh.

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

Implements G4VBiasingOperation.

Definition at line 59 of file GB06BOptnSplitAndKillByImportance.cc.

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

Here is the call graph for this function:

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

Implements G4VBiasingOperation.

Definition at line 86 of file GB06BOptnSplitAndKillByImportance.cc.

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

G4int GB06BOptnSplitAndKillByImportance::GetParallelWorldIndex ( ) const
inline

Definition at line 89 of file GB06BOptnSplitAndKillByImportance.hh.

90  { return fParallelWorldIndex; }
virtual const G4VBiasingInteractionLaw* GB06BOptnSplitAndKillByImportance::ProvideOccurenceBiasingInteractionLaw ( const G4BiasingProcessInterface ,
G4ForceCondition  
)
inlinefinalvirtual

Implements G4VBiasingOperation.

Definition at line 55 of file GB06BOptnSplitAndKillByImportance.hh.

57  {return 0;}
void GB06BOptnSplitAndKillByImportance::SetBiasingSharedData ( const G4BiasingProcessSharedData sharedData)
inline

Definition at line 93 of file GB06BOptnSplitAndKillByImportance.hh.

94  { fBiasingSharedData = sharedData; }
void GB06BOptnSplitAndKillByImportance::SetImportanceMap ( std::map< G4int, G4int > *  importanceMap)
inline

Definition at line 96 of file GB06BOptnSplitAndKillByImportance.hh.

97  { fImportanceMap = importanceMap; }
void GB06BOptnSplitAndKillByImportance::SetParallelWorldIndex ( G4int  parallelWorldIndex)
inline

Definition at line 87 of file GB06BOptnSplitAndKillByImportance.hh.

88  { fParallelWorldIndex = parallelWorldIndex; }

Here is the caller graph for this function:


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