Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 // $Id: G4CompetitiveFission.cc 98739 2016-08-09 12:56:55Z gcosmo $
28 //
29 // Hadronic Process: Nuclear De-excitations
30 // by V. Lara (Oct 1998)
31 //
32 // J. M. Quesada (March 2009). Bugs fixed:
33 // - Full relativistic calculation (Lorentz boosts)
34 // - Fission pairing energy is included in fragment excitation energies
35 // Now Energy and momentum are conserved in fission
36 
37 #include "G4CompetitiveFission.hh"
38 #include "G4PairingCorrection.hh"
39 #include "G4ParticleMomentum.hh"
40 #include "G4Pow.hh"
41 #include "G4Exp.hh"
42 #include "G4PhysicalConstants.hh"
43 
45 {
46  theFissionBarrierPtr = new G4FissionBarrier;
47  MyOwnFissionBarrier = true;
48 
49  theFissionProbabilityPtr = new G4FissionProbability;
50  MyOwnFissionProbability = true;
51 
52  theLevelDensityPtr = new G4FissionLevelDensityParameter;
53  MyOwnLevelDensity = true;
54 
55  MaximalKineticEnergy = -1000.0;
56  FissionBarrier = 0.0;
57  FissionProbability = 0.0;
58  LevelDensityParameter = 0.0;
59  pairingCorrection = G4PairingCorrection::GetInstance();
60 }
61 
63 {
64  if (MyOwnFissionBarrier) delete theFissionBarrierPtr;
65  if (MyOwnFissionProbability) delete theFissionProbabilityPtr;
66  if (MyOwnLevelDensity) delete theLevelDensityPtr;
67 }
68 
70 {
71  G4int anA = fragment->GetA_asInt();
72  G4int aZ = fragment->GetZ_asInt();
73  G4double ExEnergy = fragment->GetExcitationEnergy() -
74  pairingCorrection->GetFissionPairingCorrection(anA,aZ);
75 
76  // Saddle point excitation energy ---> A = 65
77  // Fission is excluded for A < 65
78  if (anA >= 65 && ExEnergy > 0.0) {
79  FissionBarrier = theFissionBarrierPtr->FissionBarrier(anA,aZ,ExEnergy);
80  MaximalKineticEnergy = ExEnergy - FissionBarrier;
81  LevelDensityParameter =
82  theLevelDensityPtr->LevelDensityParameter(anA,aZ,ExEnergy);
83  FissionProbability =
84  theFissionProbabilityPtr->EmissionProbability(*fragment,
85  MaximalKineticEnergy);
86  }
87  else {
88  MaximalKineticEnergy = -1000.0;
89  LevelDensityParameter = 0.0;
90  FissionProbability = 0.0;
91  }
92  return FissionProbability;
93 }
94 
96 {
97  G4Fragment * Fragment1 = 0;
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;
168  Tmax - FragmentsKineticEnergy + pcorr;
169 
170  // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
171  } while (FragmentsExcitationEnergy < 0.0
172  && ++Trials < 100);
173 
174  if (FragmentsExcitationEnergy <= 0.0) {
175  throw G4HadronicException(__FILE__, __LINE__,
176  "G4CompetitiveFission::BreakItUp: Excitation energy for fragments < 0.0!");
177  }
178 
179  // Fragment 1
180  M1 += FragmentsExcitationEnergy * A1/static_cast<G4double>(A);
181  // Fragment 2
182  M2 += FragmentsExcitationEnergy * A2/static_cast<G4double>(A);
183  // primary
184  M += U;
185 
186  G4double etot1 = ((M - M2)*(M + M2) + M1*M1)/(2*M);
187  G4ParticleMomentum Momentum1(IsotropicVector(std::sqrt((etot1 - M1)*(etot1+M1))));
188  G4LorentzVector FourMomentum1(Momentum1, etot1);
189  FourMomentum1.boost(theNucleusMomentum.boostVector());
190 
191  // Create Fragments
192  Fragment1 = new G4Fragment( A1, Z1, FourMomentum1);
193  theNucleusMomentum -= FourMomentum1;
194  theNucleus->SetZandA_asInt(Z2, A2);
195  theNucleus->SetMomentum(theNucleusMomentum);
196  return Fragment1;
197 }
198 
199 G4int
200 G4CompetitiveFission::FissionAtomicNumber(G4int A)
201  // Calculates the atomic number of a fission product
202 {
203 
204  // For Simplicity reading code
205  G4int A1 = theParam.GetA1();
206  G4int A2 = theParam.GetA2();
207  G4double As = theParam.GetAs();
208  G4double Sigma2 = theParam.GetSigma2();
209  G4double SigmaS = theParam.GetSigmaS();
210  G4double w = theParam.GetW();
211 
212  G4double C2A = A2 + 3.72*Sigma2;
213  G4double C2S = As + 3.72*SigmaS;
214 
215  G4double C2 = 0.0;
216  if (w > 1000.0 ) { C2 = C2S; }
217  else if (w < 0.001) { C2 = C2A; }
218  else { C2 = std::max(C2A,C2S); }
219 
220  G4double C1 = A-C2;
221  if (C1 < 30.0) {
222  C2 = A-30.0;
223  C1 = 30.0;
224  }
225 
226  G4double Am1 = (As + A1)*0.5;
227  G4double Am2 = (A1 + A2)*0.5;
228 
229  // Get Mass distributions as sum of symmetric and asymmetric Gasussians
230  G4double Mass1 = MassDistribution(As,A);
231  G4double Mass2 = MassDistribution(Am1,A);
232  G4double Mass3 = MassDistribution(G4double(A1),A);
233  G4double Mass4 = MassDistribution(Am2,A);
234  G4double Mass5 = MassDistribution(G4double(A2),A);
235  // get maximal value among Mass1,...,Mass5
236  G4double MassMax = Mass1;
237  if (Mass2 > MassMax) { MassMax = Mass2; }
238  if (Mass3 > MassMax) { MassMax = Mass3; }
239  if (Mass4 > MassMax) { MassMax = Mass4; }
240  if (Mass5 > MassMax) { MassMax = Mass5; }
241 
242  // Sample a fragment mass number, which lies between C1 and C2
243  G4double xm;
244  G4double Pm;
245  do {
246  xm = C1+G4UniformRand()*(C2-C1);
247  Pm = MassDistribution(xm,A);
248  // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
249  } while (MassMax*G4UniformRand() > Pm);
250  G4int ires = G4lrint(xm);
251 
252  return ires;
253 }
254 
255 G4double
256 G4CompetitiveFission::MassDistribution(G4double x, G4int A)
257  // This method gives mass distribution F(x) = F_{asym}(x)+w*F_{sym}(x)
258  // which consist of symmetric and asymmetric sum of gaussians components.
259 {
260  G4double y0 = (x-theParam.GetAs())/theParam.GetSigmaS();
261  G4double Xsym = G4Exp(-0.5*y0*y0);
262 
263  G4double y1 = (x - theParam.GetA1())/theParam.GetSigma1();
264  G4double y2 = (x - theParam.GetA2())/theParam.GetSigma2();
265  G4double z1 = (x - A + theParam.GetA1())/theParam.GetSigma1();
266  G4double z2 = (x - A + theParam.GetA2())/theParam.GetSigma2();
267  G4double Xasym = G4Exp(-0.5*y1*y1) + G4Exp(-0.5*y2*y2)
268  + 0.5*( G4Exp(-0.5*z1*z1) + G4Exp(-0.5*z2*z2));
269 
270  G4double res;
271  G4double w = theParam.GetW();
272  if (w > 1000) { res = Xsym; }
273  else if (w < 0.001) { res = Xasym; }
274  else { res = w*Xsym+Xasym; }
275  return res;
276 }
277 
278 G4int G4CompetitiveFission::FissionCharge(G4int A, G4int Z, G4double Af)
279  // Calculates the charge of a fission product for a given atomic number Af
280 {
281  static const G4double sigma = 0.6;
282  G4double DeltaZ = 0.0;
283  if (Af >= 134.0) { DeltaZ = -0.45; }
284  else if (Af <= (A-134.0)) { DeltaZ = 0.45; }
285  else { DeltaZ = -0.45*(Af-A*0.5)/(134.0-A*0.5); }
286 
287  G4double Zmean = (Af/A)*Z + DeltaZ;
288 
289  G4double theZ;
290  do {
291  theZ = G4RandGauss::shoot(Zmean,sigma);
292  // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
293  } while (theZ < 1.0 || theZ > (Z-1.0) || theZ > Af);
294 
295  return G4lrint(theZ);
296 }
297 
298 G4double
299 G4CompetitiveFission::FissionKineticEnergy(G4int A, G4int Z,
300  G4int Af1, G4int /*Zf1*/,
301  G4int Af2, G4int /*Zf2*/,
302  G4double /*U*/, G4double Tmax)
303  // Gives the kinetic energy of fission products
304 {
305  // Find maximal value of A for fragments
306  G4int AfMax = std::max(Af1,Af2);
307 
308  // Weights for symmetric and asymmetric components
309  G4double Pas = 0.0;
310  if (theParam.GetW() <= 1000) {
311  G4double x1 = (AfMax-theParam.GetA1())/theParam.GetSigma1();
312  G4double x2 = (AfMax-theParam.GetA2())/theParam.GetSigma2();
313  Pas = 0.5*G4Exp(-0.5*x1*x1) + G4Exp(-0.5*x2*x2);
314  }
315 
316  G4double Ps = 0.0;
317  if (theParam.GetW() >= 0.001) {
318  G4double xs = (AfMax-theParam.GetAs())/theParam.GetSigmaS();
319  Ps = theParam.GetW()*G4Exp(-0.5*xs*xs);
320  }
321  G4double Psy = Ps/(Pas+Ps);
322 
323  // Fission fractions Xsy and Xas formed in symmetric and asymmetric modes
324  G4double PPas = theParam.GetSigma1() + 2.0 * theParam.GetSigma2();
325  G4double PPsy = theParam.GetW() * theParam.GetSigmaS();
326  G4double Xas = PPas / (PPas+PPsy);
327  G4double Xsy = PPsy / (PPas+PPsy);
328 
329  // Average kinetic energy for symmetric and asymmetric components
330  G4double Eaverage = (0.1071*(Z*Z)/G4Pow::GetInstance()->Z13(A) + 22.2)*CLHEP::MeV;
331 
332  // Compute maximal average kinetic energy of fragments and Energy Dispersion
333  G4double TaverageAfMax;
334  G4double ESigma = 10*CLHEP::MeV;
335  // Select randomly fission mode (symmetric or asymmetric)
336  if (G4UniformRand() > Psy) { // Asymmetric Mode
337  G4double A11 = theParam.GetA1()-0.7979*theParam.GetSigma1();
338  G4double A12 = theParam.GetA1()+0.7979*theParam.GetSigma1();
339  G4double A21 = theParam.GetA2()-0.7979*theParam.GetSigma2();
340  G4double A22 = theParam.GetA2()+0.7979*theParam.GetSigma2();
341  // scale factor
342  G4double ScaleFactor = 0.5*theParam.GetSigma1()*
343  (AsymmetricRatio(A,A11)+AsymmetricRatio(A,A12))+
344  theParam.GetSigma2()*(AsymmetricRatio(A,A21)+AsymmetricRatio(A,A22));
345  // Compute average kinetic energy for fragment with AfMax
346  TaverageAfMax = (Eaverage + 12.5 * Xsy) * (PPas/ScaleFactor) *
347  AsymmetricRatio(A,G4double(AfMax));
348 
349  } else { // Symmetric Mode
350  G4double As0 = theParam.GetAs() + 0.7979*theParam.GetSigmaS();
351  // Compute average kinetic energy for fragment with AfMax
352  TaverageAfMax = (Eaverage - 12.5*CLHEP::MeV*Xas)
353  *SymmetricRatio(A, G4double(AfMax))/SymmetricRatio(A, As0);
354  ESigma = 8.0*CLHEP::MeV;
355  }
356 
357  // Select randomly, in accordance with Gaussian distribution,
358  // fragment kinetic energy
359  G4double KineticEnergy;
360  G4int i = 0;
361  do {
362  KineticEnergy = G4RandGauss::shoot(TaverageAfMax, ESigma);
363  if (++i > 100) return Eaverage;
364  // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
365  } while (KineticEnergy < Eaverage-3.72*ESigma ||
366  KineticEnergy > Eaverage+3.72*ESigma ||
367  KineticEnergy > Tmax);
368 
369  return KineticEnergy;
370 }
371 
372 
#define A21
const double C2
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
ThreeVector shoot(const G4int Ap, const G4int Af)
Hep3Vector boostVector() const
virtual G4double FissionBarrier(G4int A, G4int Z, G4double U) const =0
G4double GetSigma2(void) const
static G4double GetNuclearMass(const G4double A, const G4double Z)
const double C1
G4double GetSigma1(void) const
#define A22
virtual G4double GetEmissionProbability(G4Fragment *theNucleus)
tuple x
Definition: test.py:50
int G4int
Definition: G4Types.hh:78
#define A12
G4int GetA2(void) const
#define G4UniformRand()
Definition: Randomize.hh:97
G4double Z13(G4int Z) const
Definition: G4Pow.hh:127
double A(double temperature)
G4int GetA_asInt() const
Definition: G4Fragment.hh:266
G4double GetSigmaS(void) const
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:307
HepLorentzVector & boost(double, double, double)
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:312
virtual G4Fragment * EmittedFragment(G4Fragment *theNucleus)
static constexpr double MeV
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:288
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static G4PairingCorrection * GetInstance()
G4int GetA1(void) const
G4double GetFissionPairingCorrection(G4int A, G4int Z) const
G4double GetAs(void) const
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void SetZandA_asInt(G4int Znew, G4int Anew)
Definition: G4Fragment.hh:276
virtual G4double LevelDensityParameter(G4int A, G4int Z, G4double U) const =0
#define A11
G4int GetZ_asInt() const
Definition: G4Fragment.hh:271
double G4double
Definition: G4Types.hh:76
G4double GetW(void) const
virtual G4double EmissionProbability(const G4Fragment &fragment, const G4double anEnergy)=0
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:283
void DefineParameters(G4int A, G4int Z, G4double ExEnergy, G4double FissionBarrier)