Geant4  10.02
GB03BOptnSplitOrKillOnBoundary.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id$
27 //
30 
31 #include "Randomize.hh"
33 
34 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
35 
37 : G4VBiasingOperation(name),
38  fParticleChange(),
39  fParticleChangeForNothing()
40 {}
41 
42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43 
45 {}
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48 
51  G4double,
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 }
60 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62 
65  const G4Step* step )
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:
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:
148  }
149  else
150  {
151  // In this case, the track survives. We change its weight
152  // to conserve weight among killed and survival tracks:
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.
168 
169 }
170 
171 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
virtual void Initialize(const G4Track &track)
G4double condition(const G4ErrorSymMatrix &m)
G4String name
Definition: TRTMaterials.hh:40
G4StepStatus GetStepStatus() const
void ProposeParentWeight(G4double finalWeight)
int G4int
Definition: G4Types.hh:78
void SetWeight(G4double aValue)
void SetSecondaryWeightByProcess(G4bool)
#define G4UniformRand()
Definition: Randomize.hh:97
G4int GetCurrentStepNumber() const
Definition: G4Step.hh:76
G4ParticleChangeForNothing fParticleChangeForNothing
G4ThreeVector GetMomentum() const
virtual void Initialize(const G4Track &)
virtual G4double DistanceToApplyOperation(const G4Track *, G4double, G4ForceCondition *condition)
void SetNumberOfSecondaries(G4int totSecondaries)
G4StepPoint * GetPostStepPoint() const
virtual G4VParticleChange * GenerateBiasingFinalState(const G4Track *, const G4Step *)
void AddSecondary(G4Track *aSecondary)
G4double GetWeight() const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ForceCondition
Definition of the GB03BOptnSplitOrKillOnBoundary class.
#define DBL_MAX
Definition: templates.hh:83