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