Loading...
Searching...
No Matches
GB05BOptnSplitAndKillByCrossSection.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//
27/// \file GB05BOptnSplitAndKillByCrossSection.cc
28/// \brief Implementation of the GB05BOptnSplitAndKillByCrossSection class
29
30#include "Randomize.hh"
32
33//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
34
37 fParticleChange(),
38 fInteractionLength(-1.0)
39{}
40
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
42
45
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
47
50 G4double,
51 G4ForceCondition* condition)
52{
53 *condition = NotForced;
54
55 // -- Sample the exponential law using the total interaction length of processes
56 // -- to couterbalance for:
57 G4double proposedStepLength = -std::log( G4UniformRand() ) * fInteractionLength;
58
59 return proposedStepLength;
60}
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63
66 const G4Step* )
67{
68
69 // -- This method is called if we have limited the step.
70 // -- We hence make the splitting or killing.
71
72 // Get track weight:
73 G4double initialWeight = track->GetWeight();
74
75 // The "particle change" is the object to be used to communicate to
76 // the tracking the update of the primary state and/or creation
77 // secondary tracks.
78 fParticleChange.Initialize(*track);
79
80 // -- Splitting and killing factors.
81 // -- They are taken the same, but the killing factor can be make bigger.
82 G4double splittingFactor = 2.0;
83 G4double killingFactor = 2.0;
84
85
86 if ( track->GetMomentumDirection().z() > 0 )
87 {
88 // -- We split if the track is moving forward:
89
90 // Define the tracks weight:
91 G4double weightOfTrack = initialWeight/splittingFactor;
92
93 // Ask currect track weight to be changed to new value:
94 fParticleChange.ProposeParentWeight( weightOfTrack );
95 // Now make clones of this track (this is the actual splitting):
96 // we will then have the primary and clone of it, hence the
97 // splitting by a factor 2:
98 G4Track* clone = new G4Track( *track );
99 clone->SetWeight( weightOfTrack );
100 fParticleChange.AddSecondary( clone );
101 // -- Below's call added for safety & illustration : inform particle change to not
102 // -- modify the clone (ie : daughter) weight to male it that of the
103 // -- primary. Here call not mandatory and both tracks have same weights.
104 fParticleChange.SetSecondaryWeightByProcess(true);
105 }
106 else
107 {
108 // -- We apply Russian roulette if the track is moving backward:
109
110 // Shoot a random number (in ]0,1[ segment):
111 G4double random = G4UniformRand();
112 G4double killingProbability = 1.0 - 1.0/killingFactor;
113 if ( random < killingProbability )
114 {
115 // We ask for the the track to be killed:
116 fParticleChange.ProposeTrackStatus(fStopAndKill);
117 }
118 else
119 {
120 // In this case, the track survives. We change its weight
121 // to conserve weight among killed and survival tracks:
122 fParticleChange.ProposeParentWeight( initialWeight*killingFactor );
123 }
124 }
125
126 return &fParticleChange;
127
128}
129
130//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Definition of the GB05BOptnSplitAndKillByCrossSection class.
virtual G4double DistanceToApplyOperation(const G4Track *, G4double, G4ForceCondition *condition) final
virtual G4VParticleChange * GenerateBiasingFinalState(const G4Track *, const G4Step *) final

Applications | User Support | Publications | Collaboration