Geant4
10.03
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
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
// $Id$
27
//
30
31
#include "
Randomize.hh
"
32
#include "
GB05BOptnSplitAndKillByCrossSection.hh
"
33
34
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
35
36
GB05BOptnSplitAndKillByCrossSection::GB05BOptnSplitAndKillByCrossSection
(
G4String
name
)
37
:
G4VBiasingOperation
(name),
38
fParticleChange(),
39
fInteractionLength(-1.0)
40
{}
41
42
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43
44
GB05BOptnSplitAndKillByCrossSection::~GB05BOptnSplitAndKillByCrossSection
()
45
{}
46
47
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48
49
G4double
GB05BOptnSplitAndKillByCrossSection::
50
DistanceToApplyOperation
(
const
G4Track
*,
51
G4double
,
52
G4ForceCondition
*
condition
)
53
{
54
*condition =
NotForced
;
55
56
// -- Sample the exponential law using the total interaction length of processes
57
// -- to couterbalance for:
58
G4double
proposedStepLength = -std::log(
G4UniformRand
() ) *
fInteractionLength
;
59
60
return
proposedStepLength;
61
}
62
63
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64
65
G4VParticleChange
*
66
GB05BOptnSplitAndKillByCrossSection::GenerateBiasingFinalState
(
const
G4Track
* track,
67
const
G4Step
* )
68
{
69
70
// -- This method is called if we have limited the step.
71
// -- We hence make the splitting or killing.
72
73
// Get track weight:
74
G4double
initialWeight = track->
GetWeight
();
75
76
// The "particle change" is the object to be used to communicate to
77
// the tracking the update of the primary state and/or creation
78
// secondary tracks.
79
fParticleChange
.
Initialize
(*track);
80
81
// -- Splitting and killing factors.
82
// -- They are taken the same, but the killing factor can be make bigger.
83
G4double
splittingFactor = 2.0;
84
G4double
killingFactor = 2.0;
85
86
87
if
( track->
GetMomentumDirection
().z() > 0 )
88
{
89
// -- We split if the track is moving forward:
90
91
// Define the tracks weight:
92
G4double
weightOfTrack = initialWeight/splittingFactor;
93
94
// Ask currect track weight to be changed to new value:
95
fParticleChange
.
ProposeParentWeight
( weightOfTrack );
96
// Now make clones of this track (this is the actual splitting):
97
// we will then have the primary and clone of it, hence the
98
// splitting by a factor 2:
99
G4Track
* clone =
new
G4Track
( *track );
100
clone->
SetWeight
( weightOfTrack );
101
fParticleChange
.
AddSecondary
( clone );
102
// -- Below's call added for safety & illustration : inform particle change to not
103
// -- modify the clone (ie : daughter) weight to male it that of the
104
// -- primary. Here call not mandatory and both tracks have same weights.
105
fParticleChange
.
SetSecondaryWeightByProcess
(
true
);
106
}
107
else
108
{
109
// -- We apply Russian roulette if the track is moving backward:
110
111
// Shoot a random number (in ]0,1[ segment):
112
G4double
random =
G4UniformRand
();
113
G4double
killingProbability = 1.0 - 1.0/killingFactor;
114
if
( random < killingProbability )
115
{
116
// We ask for the the track to be killed:
117
fParticleChange
.
ProposeTrackStatus
(
fStopAndKill
);
118
}
119
else
120
{
121
// In this case, the track survives. We change its weight
122
// to conserve weight among killed and survival tracks:
123
fParticleChange
.
ProposeParentWeight
( initialWeight*killingFactor );
124
}
125
}
126
127
return
&
fParticleChange
;
128
129
}
130
131
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
condition
G4double condition(const G4ErrorSymMatrix &m)
GB05BOptnSplitAndKillByCrossSection::fParticleChange
G4ParticleChange fParticleChange
Definition:
GB05BOptnSplitAndKillByCrossSection.hh:88
G4VParticleChange::ProposeParentWeight
void ProposeParentWeight(G4double finalWeight)
GB05BOptnSplitAndKillByCrossSection.hh
Definition of the GB05BOptnSplitAndKillByCrossSection class.
G4InuclParticleNames::name
const char * name(G4int ptype)
Definition:
G4InuclParticleNames.hh:77
G4Track::SetWeight
void SetWeight(G4double aValue)
G4VParticleChange::SetSecondaryWeightByProcess
void SetSecondaryWeightByProcess(G4bool)
GB05BOptnSplitAndKillByCrossSection::GenerateBiasingFinalState
virtual G4VParticleChange * GenerateBiasingFinalState(const G4Track *, const G4Step *) final
Definition:
GB05BOptnSplitAndKillByCrossSection.cc:66
GB05BOptnSplitAndKillByCrossSection::fInteractionLength
G4double fInteractionLength
Definition:
GB05BOptnSplitAndKillByCrossSection.hh:89
G4UniformRand
#define G4UniformRand()
Definition:
Randomize.hh:97
G4VBiasingOperation
Definition:
G4VBiasingOperation.hh:77
G4Step
Definition:
G4Step.hh:76
Randomize.hh
GB05BOptnSplitAndKillByCrossSection::DistanceToApplyOperation
virtual G4double DistanceToApplyOperation(const G4Track *, G4double, G4ForceCondition *condition) final
Definition:
GB05BOptnSplitAndKillByCrossSection.cc:50
G4ParticleChange::Initialize
virtual void Initialize(const G4Track &)
Definition:
G4ParticleChange.cc:228
G4Track::GetMomentumDirection
const G4ThreeVector & GetMomentumDirection() const
NotForced
Definition:
G4ForceCondition.hh:56
GB05BOptnSplitAndKillByCrossSection::GB05BOptnSplitAndKillByCrossSection
GB05BOptnSplitAndKillByCrossSection(G4String name)
Definition:
GB05BOptnSplitAndKillByCrossSection.cc:36
fStopAndKill
Definition:
G4TrackStatus.hh:56
G4Track
Definition:
G4Track.hh:76
G4ParticleChange::AddSecondary
void AddSecondary(G4Track *aSecondary)
Definition:
G4ParticleChange.cc:218
G4Track::GetWeight
G4double GetWeight() const
G4double
double G4double
Definition:
G4Types.hh:76
G4VParticleChange::ProposeTrackStatus
void ProposeTrackStatus(G4TrackStatus status)
G4ForceCondition
G4ForceCondition
Definition:
G4ForceCondition.hh:49
GB05BOptnSplitAndKillByCrossSection::~GB05BOptnSplitAndKillByCrossSection
virtual ~GB05BOptnSplitAndKillByCrossSection()
Definition:
GB05BOptnSplitAndKillByCrossSection.cc:44
G4VParticleChange
Definition:
G4VParticleChange.hh:94
G4String
Definition:
G4String.hh:45
geant4.10.03
examples
extended
biasing
GB05
src
GB05BOptnSplitAndKillByCrossSection.cc
Generated on Thu Feb 14 2002 02:27:43 for Geant4 by
1.8.8