Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LFission.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//
28// G4 Model: Low Energy Fission
29// F.W. Jones, TRIUMF, 03-DEC-96
30//
31// This is a prototype of a low-energy fission process.
32// Currently it is based on the GHEISHA routine FISSIO,
33// and conforms fairly closely to the original Fortran.
34// Note: energy is in MeV and momentum is in MeV/c.
35//
36// use -scheme for elastic scattering: HPW, 20th June 1997
37// the code comes mostly from the old Low-energy Fission class
38//
39// 25-JUN-98 FWJ: replaced missing Initialize for ParticleChange.
40
41#include <iostream>
42
43#include "G4LFission.hh"
44#include "globals.hh"
45#include "G4Exp.hh"
46#include "G4Log.hh"
47#include "G4Pow.hh"
49#include "G4SystemOfUnits.hh"
50#include "Randomize.hh"
51
54{
55 init();
56 SetMinEnergy(0.0*GeV);
58}
59
60
62{
64}
65
66
67void G4LFission::ModelDescription(std::ostream& outFile) const
68{
69 outFile << "G4LFission is one of the Low Energy Parameterized\n"
70 << "(LEP) models used to implement neutron-induced fission of\n"
71 << "nuclei. It is a re-engineered version of the GHEISHA code\n"
72 << "of H. Fesefeldt which emits neutrons and gammas but no\n"
73 << "nuclear fragments. The model is applicable to all incident\n"
74 << "neutron energies.\n";
75}
76
77void G4LFission::init()
78{
79 G4int i;
80 G4double xx = 1. - 0.5;
81 G4double xxx = std::sqrt(2.29*xx);
82 spneut[0] = G4Exp(-xx/0.965)*(G4Exp(xxx) - G4Exp(-xxx))/2.;
83 for (i = 2; i <= 10; i++) {
84 xx = i*1. - 0.5;
85 xxx = std::sqrt(2.29*xx);
86 spneut[i-1] = spneut[i-2] + G4Exp(-xx/0.965)*(G4Exp(xxx) - G4Exp(-xxx))/2.;
87 }
88 for (i = 1; i <= 10; i++) {
89 spneut[i-1] = spneut[i-1]/spneut[9];
90 if (verboseLevel > 1) G4cout << "G4LFission::init: i=" << i <<
91 " spneut=" << spneut[i-1] << G4endl;
92 }
93}
94
95
97 G4Nucleus& targetNucleus)
98{
100 const G4HadProjectile* aParticle = &aTrack;
101
102 G4double N = targetNucleus.GetA_asInt();
103 G4double Z = targetNucleus.GetZ_asInt();
105
106 G4double P = aParticle->GetTotalMomentum()/MeV;
107 G4double Px = aParticle->Get4Momentum().vect().x();
108 G4double Py = aParticle->Get4Momentum().vect().y();
109 G4double Pz = aParticle->Get4Momentum().vect().z();
110 G4double E = aParticle->GetTotalEnergy()/MeV;
111 G4double E0 = aParticle->GetDefinition()->GetPDGMass()/MeV;
112 G4double Q = aParticle->GetDefinition()->GetPDGCharge();
113 if (verboseLevel > 1) {
114 G4cout << "G4LFission:ApplyYourself: incident particle:" << G4endl;
115 G4cout << "P " << P << " MeV/c" << G4endl;
116 G4cout << "Px " << Px << " MeV/c" << G4endl;
117 G4cout << "Py " << Py << " MeV/c" << G4endl;
118 G4cout << "Pz " << Pz << " MeV/c" << G4endl;
119 G4cout << "E " << E << " MeV" << G4endl;
120 G4cout << "mass " << E0 << " MeV" << G4endl;
121 G4cout << "charge " << Q << G4endl;
122 }
123 // GHEISHA ADD operation to get total energy, mass, charge:
124 if (verboseLevel > 1) {
125 G4cout << "G4LFission:ApplyYourself: material:" << G4endl;
126 G4cout << "A " << N << G4endl;
127 G4cout << "Z " << Z << G4endl;
128 G4cout << "atomic mass " <<
129 Atomas(N, Z) << "MeV" << G4endl;
130 }
131 E = E + Atomas(N, Z);
132 G4double E02 = E*E - P*P;
133 E0 = std::sqrt(std::abs(E02));
134 if (E02 < 0) E0 = -E0;
135 Q = Q + Z;
136 if (verboseLevel > 1) {
137 G4cout << "G4LFission:ApplyYourself: total:" << G4endl;
138 G4cout << "E " << E << " MeV" << G4endl;
139 G4cout << "mass " << E0 << " MeV" << G4endl;
140 G4cout << "charge " << Q << G4endl;
141 }
142 Px = -Px;
143 Py = -Py;
144 Pz = -Pz;
145
146 G4double e1 = aParticle->GetKineticEnergy()/MeV;
147 if (e1 < 1.) e1 = 1.;
148
149// Average number of neutrons
150 G4double avern = 2.569 + 0.559*G4Log(e1);
151 G4bool photofission = 0; // For now
152// Take the following value if photofission is not included
153 if (!photofission) avern = 2.569 + 0.900*G4Log(e1);
154
155// Average number of gammas
156 G4double averg = 9.500 + 0.600*G4Log(e1);
157
158 G4double ran = G4RandGauss::shoot();
159// Number of neutrons
160 G4int nn = static_cast<G4int>(avern + ran*1.23 + 0.5);
161 ran = G4RandGauss::shoot();
162// Number of gammas
163 G4int ng = static_cast<G4int>(averg + ran*3. + 0.5);
164 if (nn < 1) nn = 1;
165 if (ng < 1) ng = 1;
166 G4double exn = 0.;
167 G4double exg = 0.;
168
169// Make secondary neutrons and distribute kinetic energy
170 G4DynamicParticle* aNeutron;
171 G4int i;
172 for (i = 1; i <= nn; i++) {
173 ran = G4UniformRand();
174 G4int j;
175 for (j = 1; j <= 10; j++) {
176 if (ran < spneut[j-1]) goto label12;
177 }
178 j = 10;
179 label12:
180 ran = G4UniformRand();
181 G4double ekin = (j - 1)*1. + ran;
182 exn = exn + ekin;
184 G4ParticleMomentum(1.,0.,0.),
185 ekin*MeV);
187 }
188
189// Make secondary gammas and distribute kinetic energy
190 G4DynamicParticle* aGamma;
191 for (i = 1; i <= ng; i++) {
192 ran = G4UniformRand();
193 G4double ekin = -0.87*G4Log(ran);
194 exg = exg + ekin;
196 G4ParticleMomentum(1.,0.,0.),
197 ekin*MeV);
199 }
200
201// Distribute momentum vectors and do Lorentz transformation
202
203 G4HadSecondary* theSecondary;
204
205 for (i = 1; i <= nn + ng; i++) {
206 G4double ran1 = G4UniformRand();
207 G4double ran2 = G4UniformRand();
208 G4double cost = -1. + 2.*ran1;
209 G4double sint = std::sqrt(std::abs(1. - cost*cost));
210 G4double phi = ran2*twopi;
211 // G4cout << ran1 << " " << ran2 << G4endl;
212 // G4cout << cost << " " << sint << " " << phi << G4endl;
213 theSecondary = theParticleChange.GetSecondary(i - 1);
214 G4double pp = theSecondary->GetParticle()->GetTotalMomentum()/MeV;
215 G4double px = pp*sint*std::sin(phi);
216 G4double py = pp*sint*std::cos(phi);
217 G4double pz = pp*cost;
218 // G4cout << pp << G4endl;
219 // G4cout << px << " " << py << " " << pz << G4endl;
220 G4double e = theSecondary->GetParticle()->GetTotalEnergy()/MeV;
221 G4double e0 = theSecondary->GetParticle()->GetDefinition()->GetPDGMass()/MeV;
222
223 G4double a = px*Px + py*Py + pz*Pz;
224 a = (a/(E + E0) - e)/E0;
225
226 px = px + a*Px;
227 py = py + a*Py;
228 pz = pz + a*Pz;
229 G4double p2 = px*px + py*py + pz*pz;
230 pp = std::sqrt(p2);
231 e = std::sqrt(e0*e0 + p2);
232 G4double ekin = e - theSecondary->GetParticle()->GetDefinition()->GetPDGMass()/MeV;
234 py/pp,
235 pz/pp));
236 theSecondary->GetParticle()->SetKineticEnergy(ekin*MeV);
237 }
238
239 return &theParticleChange;
240}
241
242// Computes atomic mass in MeV (translation of GHEISHA routine ATOMAS)
243// Not optimized: conforms closely to original Fortran.
244
246{
252
253 G4int ia = static_cast<G4int>(A + 0.5);
254 if (ia < 1) return 0;
255 G4int iz = static_cast<G4int>(Z + 0.5);
256 if (iz < 0) return 0;
257 if (iz > ia) return 0;
258
259 if (ia == 1) {
260 if (iz == 0) return rmn; //neutron
261 if (iz == 1) return rmp + rmel; //Hydrogen
262 }
263 else if (ia == 2 && iz == 1) {
264 return rmd; //Deuteron
265 }
266 else if (ia == 4 && iz == 2) {
267 return rma; //Alpha
268 }
269
271 G4double mass = (A - Z)*rmn + Z*rmp + Z*rmel - 15.67*A
272 + 17.23*Pow->A23(A)
273 + 93.15*(A/2. - Z)*(A/2. - Z)/A
274 + 0.6984523*Z*Z/Pow->A13(A);
275 G4int ipp = (ia - iz)%2;
276 G4int izz = iz%2;
277 if (ipp == izz) mass = mass + (ipp + izz -1)*12.*Pow->powA(A, -0.5);
278
279 return mass;
280}
281
282const std::pair<G4double, G4double> G4LFission::GetFatalEnergyCheckLevels() const
283{
284 // max energy non-conservation is mass of heavy nucleus
285 return std::pair<G4double, G4double>(5*perCent,250*GeV);
286}
double A(double temperature)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
@ stopAndKill
G4double G4Log(G4double x)
Definition: G4Log.hh:226
G4ThreeVector G4ParticleMomentum
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
double x() const
double y() const
Hep3Vector vect() const
static G4Alpha * AlphaDefinition()
Definition: G4Alpha.cc:83
static G4Deuteron * DeuteronDefinition()
Definition: G4Deuteron.cc:88
void SetMomentumDirection(const G4ThreeVector &aDirection)
G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
G4double GetTotalMomentum() const
void SetKineticEnergy(G4double aEnergy)
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:80
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
G4HadSecondary * GetSecondary(size_t i)
G4double GetTotalMomentum() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetTotalEnergy() const
G4DynamicParticle * GetParticle()
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
G4LFission(const G4String &name="G4LFission")
Definition: G4LFission.cc:52
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
Definition: G4LFission.cc:282
virtual void ModelDescription(std::ostream &outFile) const
Definition: G4LFission.cc:67
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
Definition: G4LFission.cc:96
static G4double Atomas(const G4double A, const G4double Z)
Definition: G4LFission.cc:245
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:98
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4double GetPDGCharge() const
Definition: G4Pow.hh:49
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double A13(G4double A) const
Definition: G4Pow.cc:120
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230
G4double A23(G4double A) const
Definition: G4Pow.hh:131
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:87
#define DBL_MAX
Definition: templates.hh:62