Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4LFission Class Reference

#include <G4LFission.hh>

Inheritance diagram for G4LFission:
Collaboration diagram for G4LFission:

Public Member Functions

 G4LFission (const G4String &name="G4LFission")
 
 ~G4LFission ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual const std::pair
< G4double, G4double
GetFatalEnergyCheckLevels () const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
virtual std::pair< G4double,
G4double
GetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void InitialiseModel ()
 

Static Public Member Functions

static G4double Atomas (const G4double A, const G4double Z)
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 65 of file G4LFission.hh.

Constructor & Destructor Documentation

G4LFission::G4LFission ( const G4String name = "G4LFission")

Definition at line 53 of file G4LFission.cc.

54  : G4HadronicInteraction(name)
55 {
56  init();
57  SetMinEnergy(0.0*GeV);
59 }
void SetMinEnergy(G4double anEnergy)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static constexpr double GeV
Definition: G4SIunits.hh:217
void SetMaxEnergy(const G4double anEnergy)
#define DBL_MAX
Definition: templates.hh:83

Here is the call graph for this function:

G4LFission::~G4LFission ( )

Definition at line 62 of file G4LFission.cc.

63 {
65 }

Here is the call graph for this function:

Member Function Documentation

G4HadFinalState * G4LFission::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)
virtual

Implements G4HadronicInteraction.

Definition at line 97 of file G4LFission.cc.

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 }
const int N
Definition: mixmax.h:43
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
ThreeVector shoot(const G4int Ap, const G4int Af)
G4HadSecondary * GetSecondary(size_t i)
G4double GetTotalEnergy() const
double x() const
static double Q[]
void SetMomentumDirection(const G4ThreeVector &aDirection)
G4ParticleDefinition * GetDefinition() const
static G4double Atomas(const G4double A, const G4double Z)
Definition: G4LFission.cc:246
int G4int
Definition: G4Types.hh:78
double z() const
static double P[]
static constexpr double twopi
Definition: G4SIunits.hh:76
G4double GetTotalMomentum() const
void SetStatusChange(G4HadFinalStateStatus aS)
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ParticleDefinition * GetDefinition() const
bool G4bool
Definition: G4Types.hh:79
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
void SetKineticEnergy(G4double aEnergy)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double GetPDGMass() const
G4DynamicParticle * GetParticle()
double y() const
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
#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
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:99
G4double GetTotalMomentum() const
G4double GetTotalEnergy() const
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81

Here is the call graph for this function:

G4double G4LFission::Atomas ( const G4double  A,
const G4double  Z 
)
static

Definition at line 246 of file G4LFission.cc.

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 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
Definition: G4Pow.hh:56
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
int G4int
Definition: G4Types.hh:78
G4double A23(G4double A) const
Definition: G4Pow.hh:160
double A(double temperature)
G4double GetPDGMass() const
G4double A13(G4double A) const
Definition: G4Pow.hh:132
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
static G4Deuteron * DeuteronDefinition()
Definition: G4Deuteron.cc:89
static G4Alpha * AlphaDefinition()
Definition: G4Alpha.cc:84
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:99

Here is the call graph for this function:

Here is the caller graph for this function:

const std::pair< G4double, G4double > G4LFission::GetFatalEnergyCheckLevels ( ) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 283 of file G4LFission.cc.

284 {
285  // max energy non-conservation is mass of heavy nucleus
286  return std::pair<G4double, G4double>(5*perCent,250*GeV);
287 }
static constexpr double perCent
Definition: G4SIunits.hh:332
static constexpr double GeV
Definition: G4SIunits.hh:217
void G4LFission::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 68 of file G4LFission.cc.

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 }

The documentation for this class was generated from the following files: