Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LFission Class Reference

#include <G4LFission.hh>

+ Inheritance 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, G4doubleGetFatalEnergyCheckLevels () 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, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void InitialiseModel ()
 
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
 
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 

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 58 of file G4LFission.hh.

Constructor & Destructor Documentation

◆ G4LFission()

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

Definition at line 53 of file G4LFission.cc.

54 : G4HadronicInteraction(name), secID(-1)
55{
56 init();
57 SetMinEnergy(0.0*GeV);
60}
void SetMinEnergy(G4double anEnergy)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
const G4String & GetModelName() const
void SetMaxEnergy(const G4double anEnergy)
static G4int GetModelID(const G4int modelIndex)
#define DBL_MAX
Definition templates.hh:62

◆ ~G4LFission()

G4LFission::~G4LFission ( )

Definition at line 63 of file G4LFission.cc.

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 98 of file G4LFission.cc.

100{
102 const G4HadProjectile* aParticle = &aTrack;
103
104 G4double N = targetNucleus.GetA_asInt();
105 G4double Z = targetNucleus.GetZ_asInt();
107
108 G4double P = aParticle->GetTotalMomentum()/MeV;
109 G4double Px = aParticle->Get4Momentum().vect().x();
110 G4double Py = aParticle->Get4Momentum().vect().y();
111 G4double Pz = aParticle->Get4Momentum().vect().z();
112 G4double E = aParticle->GetTotalEnergy()/MeV;
113 G4double E0 = aParticle->GetDefinition()->GetPDGMass()/MeV;
114 G4double Q = aParticle->GetDefinition()->GetPDGCharge();
115 if (verboseLevel > 1) {
116 G4cout << "G4LFission:ApplyYourself: incident particle:" << G4endl;
117 G4cout << "P " << P << " MeV/c" << G4endl;
118 G4cout << "Px " << Px << " MeV/c" << G4endl;
119 G4cout << "Py " << Py << " MeV/c" << G4endl;
120 G4cout << "Pz " << Pz << " MeV/c" << G4endl;
121 G4cout << "E " << E << " MeV" << G4endl;
122 G4cout << "mass " << E0 << " MeV" << G4endl;
123 G4cout << "charge " << Q << G4endl;
124 }
125 // GHEISHA ADD operation to get total energy, mass, charge:
126 if (verboseLevel > 1) {
127 G4cout << "G4LFission:ApplyYourself: material:" << G4endl;
128 G4cout << "A " << N << G4endl;
129 G4cout << "Z " << Z << G4endl;
130 G4cout << "atomic mass " <<
131 Atomas(N, Z) << "MeV" << G4endl;
132 }
133 E = E + Atomas(N, Z);
134 G4double E02 = E*E - P*P;
135 E0 = std::sqrt(std::abs(E02));
136 if (E02 < 0) E0 = -E0;
137 Q = Q + Z;
138 if (verboseLevel > 1) {
139 G4cout << "G4LFission:ApplyYourself: total:" << G4endl;
140 G4cout << "E " << E << " MeV" << G4endl;
141 G4cout << "mass " << E0 << " MeV" << G4endl;
142 G4cout << "charge " << Q << G4endl;
143 }
144 Px = -Px;
145 Py = -Py;
146 Pz = -Pz;
147
148 G4double e1 = aParticle->GetKineticEnergy()/MeV;
149 if (e1 < 1.) e1 = 1.;
150
151// Average number of neutrons
152 G4double avern = 2.569 + 0.559*G4Log(e1);
153 G4bool photofission = 0; // For now
154// Take the following value if photofission is not included
155 if (!photofission) avern = 2.569 + 0.900*G4Log(e1);
156
157// Average number of gammas
158 G4double averg = 9.500 + 0.600*G4Log(e1);
159
160 G4double ran = G4RandGauss::shoot();
161// Number of neutrons
162 G4int nn = static_cast<G4int>(avern + ran*1.23 + 0.5);
163 ran = G4RandGauss::shoot();
164// Number of gammas
165 G4int ng = static_cast<G4int>(averg + ran*3. + 0.5);
166 if (nn < 1) nn = 1;
167 if (ng < 1) ng = 1;
168 G4double exn = 0.;
169 G4double exg = 0.;
170
171// Make secondary neutrons and distribute kinetic energy
172 G4DynamicParticle* aNeutron;
173 G4int i;
174 for (i = 1; i <= nn; i++) {
175 ran = G4UniformRand();
176 G4int j;
177 for (j = 1; j <= 10; j++) {
178 if (ran < spneut[j-1]) goto label12;
179 }
180 j = 10;
181 label12:
182 ran = G4UniformRand();
183 G4double ekin = (j - 1)*1. + ran;
184 exn = exn + ekin;
186 G4ParticleMomentum(1.,0.,0.),
187 ekin*MeV);
188 theParticleChange.AddSecondary(aNeutron, secID);
189 }
190
191// Make secondary gammas and distribute kinetic energy
192 G4DynamicParticle* aGamma;
193 for (i = 1; i <= ng; i++) {
194 ran = G4UniformRand();
195 G4double ekin = -0.87*G4Log(ran);
196 exg = exg + ekin;
198 G4ParticleMomentum(1.,0.,0.),
199 ekin*MeV);
200 theParticleChange.AddSecondary(aGamma, secID);
201 }
202
203// Distribute momentum vectors and do Lorentz transformation
204
205 G4HadSecondary* theSecondary;
206
207 for (i = 1; i <= nn + ng; i++) {
208 G4double ran1 = G4UniformRand();
209 G4double ran2 = G4UniformRand();
210 G4double cost = -1. + 2.*ran1;
211 G4double sint = std::sqrt(std::abs(1. - cost*cost));
212 G4double phi = ran2*twopi;
213 // G4cout << ran1 << " " << ran2 << G4endl;
214 // G4cout << cost << " " << sint << " " << phi << G4endl;
215 theSecondary = theParticleChange.GetSecondary(i - 1);
216 G4double pp = theSecondary->GetParticle()->GetTotalMomentum()/MeV;
217 G4double px = pp*sint*std::sin(phi);
218 G4double py = pp*sint*std::cos(phi);
219 G4double pz = pp*cost;
220 // G4cout << pp << G4endl;
221 // G4cout << px << " " << py << " " << pz << G4endl;
222 G4double e = theSecondary->GetParticle()->GetTotalEnergy()/MeV;
223 G4double e0 = theSecondary->GetParticle()->GetDefinition()->GetPDGMass()/MeV;
224
225 G4double a = px*Px + py*Py + pz*Pz;
226 a = (a/(E + E0) - e)/E0;
227
228 px = px + a*Px;
229 py = py + a*Py;
230 pz = pz + a*Pz;
231 G4double p2 = px*px + py*py + pz*pz;
232 pp = std::sqrt(p2);
233 e = std::sqrt(e0*e0 + p2);
234 G4double ekin = e - theSecondary->GetParticle()->GetDefinition()->GetPDGMass()/MeV;
236 py/pp,
237 pz/pp));
238 theSecondary->GetParticle()->SetKineticEnergy(ekin*MeV);
239 }
240
241 return &theParticleChange;
242}
@ stopAndKill
G4double G4Log(G4double x)
Definition G4Log.hh:227
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:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
double z() const
double x() const
double y() const
Hep3Vector vect() const
void SetMomentumDirection(const G4ThreeVector &aDirection)
G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
G4double GetTotalMomentum() const
void SetKineticEnergy(G4double aEnergy)
static G4Gamma * GammaDefinition()
Definition G4Gamma.cc:76
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()
static G4double Atomas(const G4double A, const G4double Z)
static G4Neutron * NeutronDefinition()
Definition G4Neutron.cc:96
G4int GetA_asInt() const
Definition G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition G4Nucleus.hh:105
#define N
Definition crc32.c:57

◆ Atomas()

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

Definition at line 247 of file G4LFission.cc.

248{
254
255 G4int ia = static_cast<G4int>(A + 0.5);
256 if (ia < 1) return 0;
257 G4int iz = static_cast<G4int>(Z + 0.5);
258 if (iz < 0) return 0;
259 if (iz > ia) return 0;
260
261 if (ia == 1) {
262 if (iz == 0) return rmn; //neutron
263 if (iz == 1) return rmp + rmel; //Hydrogen
264 }
265 else if (ia == 2 && iz == 1) {
266 return rmd; //Deuteron
267 }
268 else if (ia == 4 && iz == 2) {
269 return rma; //Alpha
270 }
271
273 G4double mass = (A - Z)*rmn + Z*rmp + Z*rmel - 15.67*A
274 + 17.23*Pow->A23(A)
275 + 93.15*(A/2. - Z)*(A/2. - Z)/A
276 + 0.6984523*Z*Z/Pow->A13(A);
277 G4int ipp = (ia - iz)%2;
278 G4int izz = iz%2;
279 if (ipp == izz) mass = mass + (ipp + izz -1)*12.*Pow->powA(A, -0.5);
280
281 return mass;
282}
const G4double A[17]
static G4Alpha * AlphaDefinition()
Definition G4Alpha.cc:78
static G4Deuteron * DeuteronDefinition()
Definition G4Deuteron.cc:85
static G4Electron * ElectronDefinition()
Definition G4Electron.cc:86
Definition G4Pow.hh:49
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double A13(G4double A) const
Definition G4Pow.cc:116
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:85

Referenced by ApplyYourself().

◆ GetFatalEnergyCheckLevels()

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

Reimplemented from G4HadronicInteraction.

Definition at line 284 of file G4LFission.cc.

285{
286 // max energy non-conservation is mass of heavy nucleus
287 return std::pair<G4double, G4double>(10.0*perCent, 350.0*CLHEP::GeV);
288}

◆ ModelDescription()

void G4LFission::ModelDescription ( std::ostream & outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 69 of file G4LFission.cc.

70{
71 outFile << "G4LFission is one of the Low Energy Parameterized\n"
72 << "(LEP) models used to implement neutron-induced fission of\n"
73 << "nuclei. It is a re-engineered version of the GHEISHA code\n"
74 << "of H. Fesefeldt which emits neutrons and gammas but no\n"
75 << "nuclear fragments. The model is applicable to all incident\n"
76 << "neutron energies.\n";
77}

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