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

#include <GB03BOptnSplitOrKillOnBoundary.hh>

Inheritance diagram for GB03BOptnSplitOrKillOnBoundary:
Collaboration diagram for GB03BOptnSplitOrKillOnBoundary:

Public Member Functions

 GB03BOptnSplitOrKillOnBoundary (G4String name)
 
virtual ~GB03BOptnSplitOrKillOnBoundary ()
 
virtual const
G4VBiasingInteractionLaw
ProvideOccurenceBiasingInteractionLaw (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
 
- 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 40 of file GB03BOptnSplitOrKillOnBoundary.hh.

Constructor & Destructor Documentation

GB03BOptnSplitOrKillOnBoundary::GB03BOptnSplitOrKillOnBoundary ( G4String  name)

Definition at line 36 of file GB03BOptnSplitOrKillOnBoundary.cc.

37 : G4VBiasingOperation(name),
38  fParticleChange(),
39  fParticleChangeForNothing()
40 {}
G4VBiasingOperation(G4String name)
GB03BOptnSplitOrKillOnBoundary::~GB03BOptnSplitOrKillOnBoundary ( )
virtual

Definition at line 44 of file GB03BOptnSplitOrKillOnBoundary.cc.

45 {}

Member Function Documentation

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

Implements G4VBiasingOperation.

Definition at line 56 of file GB03BOptnSplitOrKillOnBoundary.hh.

59  {return 0;}
G4double GB03BOptnSplitOrKillOnBoundary::DistanceToApplyOperation ( const G4Track ,
G4double  ,
G4ForceCondition condition 
)
virtual

Implements G4VBiasingOperation.

Definition at line 50 of file GB03BOptnSplitOrKillOnBoundary.cc.

53 {
54  // -- return "infinite" distance for interaction, but asks for GenerateBiasingFinalState
55  // -- being called anyway at the end of the step, by returning the "Forced" condition
56  // -- flag.
57  *condition = Forced;
58  return DBL_MAX;
59 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
G4VParticleChange * GB03BOptnSplitOrKillOnBoundary::GenerateBiasingFinalState ( const G4Track track,
const G4Step step 
)
virtual

Implements G4VBiasingOperation.

Definition at line 64 of file GB03BOptnSplitOrKillOnBoundary.cc.

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

G4double GB03BOptnSplitOrKillOnBoundary::GetApplyProbability ( ) const
inline

Definition at line 98 of file GB03BOptnSplitOrKillOnBoundary.hh.

98 { return fApplyProbability; }
G4int GB03BOptnSplitOrKillOnBoundary::GetSplittingFactor ( ) const
inline

Definition at line 97 of file GB03BOptnSplitOrKillOnBoundary.hh.

97 { return fSplittingFactor; }
virtual const G4VBiasingInteractionLaw* GB03BOptnSplitOrKillOnBoundary::ProvideOccurenceBiasingInteractionLaw ( const G4BiasingProcessInterface ,
G4ForceCondition  
)
inlinevirtual

Implements G4VBiasingOperation.

Definition at line 53 of file GB03BOptnSplitOrKillOnBoundary.hh.

54  {return 0;}
void GB03BOptnSplitOrKillOnBoundary::SetApplyProbability ( G4double  proba)
inline

Definition at line 94 of file GB03BOptnSplitOrKillOnBoundary.hh.

95  { fApplyProbability = proba; }

Here is the caller graph for this function:

void GB03BOptnSplitOrKillOnBoundary::SetSplittingFactor ( G4int  splittingFactor)
inline

Definition at line 80 of file GB03BOptnSplitOrKillOnBoundary.hh.

81  { fSplittingFactor = splittingFactor; }

Here is the caller graph for this function:


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