Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CompetitiveFission.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// Hadronic Process: Nuclear De-excitations
28// by V. Lara (Oct 1998)
29//
30// J. M. Quesada (March 2009). Bugs fixed:
31// - Full relativistic calculation (Lorentz boosts)
32// - Fission pairing energy is included in fragment excitation energies
33// Now Energy and momentum are conserved in fission
34
37#include "G4ParticleMomentum.hh"
38#include "G4NuclearLevelData.hh"
39#include "G4VFissionBarrier.hh"
40#include "G4FissionBarrier.hh"
44#include "G4Pow.hh"
45#include "Randomize.hh"
46#include "G4RandomDirection.hh"
49
51{
52 theFissionBarrierPtr = new G4FissionBarrier;
53 myOwnFissionBarrier = true;
54
55 theFissionProbabilityPtr = new G4FissionProbability;
56 myOwnFissionProbability = true;
57
58 theLevelDensityPtr = new G4FissionLevelDensityParameter;
59 myOwnLevelDensity = true;
60
61 maxKineticEnergy = fissionBarrier = fissionProbability = 0.0;
63
64 theSecID = G4PhysicsModelCatalog::GetModelID("model_G4CompetitiveFission");
65}
66
68{
69 if (myOwnFissionBarrier) delete theFissionBarrierPtr;
70 if (myOwnFissionProbability) delete theFissionProbabilityPtr;
71 if (myOwnLevelDensity) delete theLevelDensityPtr;
72}
73
75{
76 G4int Z = fragment->GetZ_asInt();
77 G4int A = fragment->GetA_asInt();
78 fissionProbability = 0.0;
79 // Saddle point excitation energy ---> A = 65
80 if (A >= 65 && Z > 16) {
81 G4double exEnergy = fragment->GetExcitationEnergy() -
82 pairingCorrection->GetFissionPairingCorrection(A, Z);
83
84 if (exEnergy > 0.0) {
85 fissionBarrier = theFissionBarrierPtr->FissionBarrier(A, Z, exEnergy);
86 maxKineticEnergy = exEnergy - fissionBarrier;
87 fissionProbability =
88 theFissionProbabilityPtr->EmissionProbability(*fragment,
89 maxKineticEnergy);
90 }
91 }
92 return fissionProbability;
93}
94
96{
97 G4Fragment * Fragment1 = nullptr;
98 // Nucleus data
99 // Atomic number of nucleus
100 G4int A = theNucleus->GetA_asInt();
101 // Charge of nucleus
102 G4int Z = theNucleus->GetZ_asInt();
103 // Excitation energy (in MeV)
104 G4double U = theNucleus->GetExcitationEnergy();
105 G4double pcorr = pairingCorrection->GetFissionPairingCorrection(A,Z);
106 if (U <= pcorr) { return Fragment1; }
107
108 // Atomic Mass of Nucleus (in MeV)
109 G4double M = theNucleus->GetGroundStateMass();
110
111 // Nucleus Momentum
112 G4LorentzVector theNucleusMomentum = theNucleus->GetMomentum();
113
114 // Calculate fission parameters
115 theParam.DefineParameters(A, Z, U-pcorr, fissionBarrier);
116
117 // First fragment
118 G4int A1 = 0;
119 G4int Z1 = 0;
120 G4double M1 = 0.0;
121
122 // Second fragment
123 G4int A2 = 0;
124 G4int Z2 = 0;
125 G4double M2 = 0.0;
126
127 G4double FragmentsExcitationEnergy = 0.0;
128 G4double FragmentsKineticEnergy = 0.0;
129
130 G4int Trials = 0;
131 do {
132
133 // First fragment
134 A1 = FissionAtomicNumber(A);
135 Z1 = FissionCharge(A, Z, A1);
137
138 // Second Fragment
139 A2 = A - A1;
140 Z2 = Z - Z1;
141 if (A2 < 1 || Z2 < 0 || Z2 > A2) {
142 FragmentsExcitationEnergy = -1.0;
143 continue;
144 }
146 // Maximal Kinetic Energy (available energy for fragments)
147 G4double Tmax = M + U - M1 - M2 - pcorr;
148
149 // Check that fragment masses are less or equal than total energy
150 if (Tmax < 0.0) {
151 FragmentsExcitationEnergy = -1.0;
152 continue;
153 }
154
155 FragmentsKineticEnergy = FissionKineticEnergy( A , Z,
156 A1, Z1,
157 A2, Z2,
158 U , Tmax);
159
160 // Excitation Energy
161 // FragmentsExcitationEnergy = Tmax - FragmentsKineticEnergy;
162 // JMQ 04/03/09 BUG FIXED: in order to fulfill energy conservation the
163 // fragments carry the fission pairing energy in form of
164 // excitation energy
165
166 FragmentsExcitationEnergy =
167 Tmax - FragmentsKineticEnergy + pcorr;
168
169 // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
170 } while (FragmentsExcitationEnergy < 0.0 && ++Trials < 100);
171
172 if (FragmentsExcitationEnergy <= 0.0) {
173 throw G4HadronicException(__FILE__, __LINE__,
174 "G4CompetitiveFission::BreakItUp: Excitation energy for fragments < 0.0!");
175 }
176
177 // Fragment 1
178 M1 += FragmentsExcitationEnergy * A1/static_cast<G4double>(A);
179 // Fragment 2
180 M2 += FragmentsExcitationEnergy * A2/static_cast<G4double>(A);
181 // primary
182 M += U;
183
184 G4double etot1 = ((M - M2)*(M + M2) + M1*M1)/(2*M);
185 G4ParticleMomentum Momentum1 =
186 std::sqrt((etot1 - M1)*(etot1+M1))*G4RandomDirection();
187 G4LorentzVector FourMomentum1(Momentum1, etot1);
188 FourMomentum1.boost(theNucleusMomentum.boostVector());
189
190 // Create Fragments
191 Fragment1 = new G4Fragment( A1, Z1, FourMomentum1);
192 if (Fragment1 != nullptr) { Fragment1->SetCreatorModelID(theSecID); }
193 theNucleusMomentum -= FourMomentum1;
194 theNucleus->SetZandA_asInt(Z2, A2);
195 theNucleus->SetMomentum(theNucleusMomentum);
196 theNucleus->SetCreatorModelID(theSecID);
197 return Fragment1;
198}
199
200G4int
201G4CompetitiveFission::FissionAtomicNumber(G4int A)
202 // Calculates the atomic number of a fission product
203{
204
205 // For Simplicity reading code
206 G4int A1 = theParam.GetA1();
207 G4int A2 = theParam.GetA2();
208 G4double As = theParam.GetAs();
209 G4double Sigma2 = theParam.GetSigma2();
210 G4double SigmaS = theParam.GetSigmaS();
211 G4double w = theParam.GetW();
212
213 G4double C2A = A2 + 3.72*Sigma2;
214 G4double C2S = As + 3.72*SigmaS;
215
216 G4double C2 = 0.0;
217 if (w > 1000.0 ) { C2 = C2S; }
218 else if (w < 0.001) { C2 = C2A; }
219 else { C2 = std::max(C2A,C2S); }
220
221 G4double C1 = A-C2;
222 if (C1 < 30.0) {
223 C2 = A-30.0;
224 C1 = 30.0;
225 }
226
227 G4double Am1 = (As + A1)*0.5;
228 G4double Am2 = (A1 + A2)*0.5;
229
230 // Get Mass distributions as sum of symmetric and asymmetric Gasussians
231 G4double Mass1 = MassDistribution(As,A);
232 G4double Mass2 = MassDistribution(Am1,A);
233 G4double Mass3 = MassDistribution(G4double(A1),A);
234 G4double Mass4 = MassDistribution(Am2,A);
235 G4double Mass5 = MassDistribution(G4double(A2),A);
236 // get maximal value among Mass1,...,Mass5
237 G4double MassMax = Mass1;
238 if (Mass2 > MassMax) { MassMax = Mass2; }
239 if (Mass3 > MassMax) { MassMax = Mass3; }
240 if (Mass4 > MassMax) { MassMax = Mass4; }
241 if (Mass5 > MassMax) { MassMax = Mass5; }
242
243 // Sample a fragment mass number, which lies between C1 and C2
244 G4double xm;
245 G4double Pm;
246 do {
247 xm = C1+G4UniformRand()*(C2-C1);
248 Pm = MassDistribution(xm,A);
249 // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
250 } while (MassMax*G4UniformRand() > Pm);
251 G4int ires = G4lrint(xm);
252
253 return ires;
254}
255
257G4CompetitiveFission::MassDistribution(G4double x, G4int A)
258 // This method gives mass distribution F(x) = F_{asym}(x)+w*F_{sym}(x)
259 // which consist of symmetric and asymmetric sum of gaussians components.
260{
261 G4double y0 = (x-theParam.GetAs())/theParam.GetSigmaS();
262 G4double Xsym = LocalExp(y0);
263
264 G4double y1 = (x - theParam.GetA1())/theParam.GetSigma1();
265 G4double y2 = (x - theParam.GetA2())/theParam.GetSigma2();
266 G4double z1 = (x - A + theParam.GetA1())/theParam.GetSigma1();
267 G4double z2 = (x - A + theParam.GetA2())/theParam.GetSigma2();
268 G4double Xasym = LocalExp(y1) + LocalExp(y2)
269 + 0.5*(LocalExp(z1) + LocalExp(z2));
270
271 G4double res;
272 G4double w = theParam.GetW();
273 if (w > 1000) { res = Xsym; }
274 else if (w < 0.001) { res = Xasym; }
275 else { res = w*Xsym+Xasym; }
276 return res;
277}
278
279G4int G4CompetitiveFission::FissionCharge(G4int A, G4int Z, G4double Af)
280 // Calculates the charge of a fission product for a given atomic number Af
281{
282 static const G4double sigma = 0.6;
283 G4double DeltaZ = 0.0;
284 if (Af >= 134.0) { DeltaZ = -0.45; }
285 else if (Af <= (A-134.0)) { DeltaZ = 0.45; }
286 else { DeltaZ = -0.45*(Af-A*0.5)/(134.0-A*0.5); }
287
288 G4double Zmean = (Af/A)*Z + DeltaZ;
289
290 G4double theZ;
291 do {
292 theZ = G4RandGauss::shoot(Zmean,sigma);
293 // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
294 } while (theZ < 1.0 || theZ > (Z-1.0) || theZ > Af);
295
296 return G4lrint(theZ);
297}
298
300G4CompetitiveFission::FissionKineticEnergy(G4int A, G4int Z,
301 G4int Af1, G4int /*Zf1*/,
302 G4int Af2, G4int /*Zf2*/,
303 G4double /*U*/, G4double Tmax)
304 // Gives the kinetic energy of fission products
305{
306 // Find maximal value of A for fragments
307 G4int AfMax = std::max(Af1,Af2);
308
309 // Weights for symmetric and asymmetric components
310 G4double Pas = 0.0;
311 if (theParam.GetW() <= 1000) {
312 G4double x1 = (AfMax-theParam.GetA1())/theParam.GetSigma1();
313 G4double x2 = (AfMax-theParam.GetA2())/theParam.GetSigma2();
314 Pas = 0.5*LocalExp(x1) + LocalExp(x2);
315 }
316
317 G4double Ps = 0.0;
318 if (theParam.GetW() >= 0.001) {
319 G4double xs = (AfMax-theParam.GetAs())/theParam.GetSigmaS();
320 Ps = theParam.GetW()*LocalExp(xs);
321 }
322 G4double Psy = (Pas + Ps > 0.0) ? Ps/(Pas+Ps) : 0.5;
323
324 // Fission fractions Xsy and Xas formed in symmetric and asymmetric modes
325 G4double PPas = theParam.GetSigma1() + 2.0 * theParam.GetSigma2();
326 G4double PPsy = theParam.GetW() * theParam.GetSigmaS();
327 G4double Xas = (PPas + PPsy > 0.0) ? PPas/(PPas+PPsy) : 0.5;
328 G4double Xsy = 1.0 - Xas;
329
330 // Average kinetic energy for symmetric and asymmetric components
331 G4double Eaverage = (0.1071*(Z*Z)/G4Pow::GetInstance()->Z13(A) + 22.2)*CLHEP::MeV;
332
333 // Compute maximal average kinetic energy of fragments and Energy Dispersion
334 G4double TaverageAfMax;
335 G4double ESigma = 10*CLHEP::MeV;
336 // Select randomly fission mode (symmetric or asymmetric)
337 if (G4UniformRand() > Psy) { // Asymmetric Mode
338 G4double A11 = theParam.GetA1()-0.7979*theParam.GetSigma1();
339 G4double A12 = theParam.GetA1()+0.7979*theParam.GetSigma1();
340 G4double A21 = theParam.GetA2()-0.7979*theParam.GetSigma2();
341 G4double A22 = theParam.GetA2()+0.7979*theParam.GetSigma2();
342 // scale factor
343 G4double ScaleFactor = 0.5*theParam.GetSigma1()*
344 (AsymmetricRatio(A,A11)+AsymmetricRatio(A,A12))+
345 theParam.GetSigma2()*(AsymmetricRatio(A,A21)+AsymmetricRatio(A,A22));
346 // Compute average kinetic energy for fragment with AfMax
347 TaverageAfMax = (Eaverage + 12.5 * Xsy) * (PPas/ScaleFactor) *
348 AsymmetricRatio(A,G4double(AfMax));
349
350 } else { // Symmetric Mode
351 G4double As0 = theParam.GetAs() + 0.7979*theParam.GetSigmaS();
352 // Compute average kinetic energy for fragment with AfMax
353 TaverageAfMax = (Eaverage - 12.5*CLHEP::MeV*Xas)
354 *SymmetricRatio(A, G4double(AfMax))/SymmetricRatio(A, As0);
355 ESigma = 8.0*CLHEP::MeV;
356 }
357
358 // Select randomly, in accordance with Gaussian distribution,
359 // fragment kinetic energy
361 G4int i = 0;
362 do {
363 KineticEnergy = G4RandGauss::shoot(TaverageAfMax, ESigma);
364 if (++i > 100) return Eaverage;
365 // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
366 } while (KineticEnergy < Eaverage-3.72*ESigma ||
367 KineticEnergy > Eaverage+3.72*ESigma ||
368 KineticEnergy > Tmax);
369
370 return KineticEnergy;
371}
372
374{
375 if (myOwnFissionBarrier) delete theFissionBarrierPtr;
376 theFissionBarrierPtr = aBarrier;
377 myOwnFissionBarrier = false;
378}
379
380void
382{
383 if (myOwnFissionProbability) delete theFissionProbabilityPtr;
384 theFissionProbabilityPtr = aFissionProb;
385 myOwnFissionProbability = false;
386}
387
388void
390{
391 if (myOwnLevelDensity) delete theLevelDensityPtr;
392 theLevelDensityPtr = aLevelDensity;
393 myOwnLevelDensity = false;
394}
395
#define A22
#define A12
#define A21
#define A11
#define M(row, col)
G4ThreeVector G4RandomDirection()
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
const double C2
#define C1
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
void SetEmissionStrategy(G4VEmissionProbability *aFissionProb)
void SetLevelDensityParameter(G4VLevelDensityParameter *aLevelDensity)
G4double GetEmissionProbability(G4Fragment *theNucleus) override
void SetFissionBarrier(G4VFissionBarrier *aBarrier)
G4Fragment * EmittedFragment(G4Fragment *theNucleus) override
G4double GetAs(void) const
G4double GetSigma1(void) const
G4double GetW(void) const
G4int GetA1(void) const
G4double GetSigmaS(void) const
void DefineParameters(G4int A, G4int Z, G4double ExEnergy, G4double FissionBarrier)
G4int GetA2(void) const
G4double GetSigma2(void) const
G4double GetGroundStateMass() const
void SetZandA_asInt(G4int Znew, G4int Anew, G4int Lnew=0)
G4double GetExcitationEnergy() const
const G4LorentzVector & GetMomentum() const
void SetCreatorModelID(G4int value)
G4int GetZ_asInt() const
void SetMomentum(const G4LorentzVector &value)
G4int GetA_asInt() const
G4PairingCorrection * GetPairingCorrection()
static G4NuclearLevelData * GetInstance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetFissionPairingCorrection(G4int A, G4int Z) const
static G4int GetModelID(const G4int modelIndex)
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double Z13(G4int Z) const
Definition G4Pow.hh:123
virtual G4double EmissionProbability(const G4Fragment &fragment, G4double anEnergy)
virtual G4double FissionBarrier(G4int A, G4int Z, G4double U) const =0
int G4lrint(double ad)
Definition templates.hh:134