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: G4GEMProbability.cc,v 1.15 2010/11/05 14:43:27 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-04-ref-00 $ |
---|
28 | // |
---|
29 | //--------------------------------------------------------------------- |
---|
30 | // |
---|
31 | // Geant4 class G4GEMProbability |
---|
32 | // |
---|
33 | // |
---|
34 | // Hadronic Process: Nuclear De-excitations |
---|
35 | // by V. Lara (Sept 2001) |
---|
36 | // |
---|
37 | // |
---|
38 | // Hadronic Process: Nuclear De-excitations |
---|
39 | // by V. Lara (Sept 2001) |
---|
40 | // |
---|
41 | // J. M. Quesada : several fixes in total GEM width |
---|
42 | // J. M. Quesada 14/07/2009 bug fixed in total emission width (hbarc) |
---|
43 | // J. M. Quesada (September 2009) several fixes: |
---|
44 | // -level density parameter of residual calculated at its right excitation energy. |
---|
45 | // -InitialLeveldensity calculated according to the right conditions of the |
---|
46 | // initial excited nucleus. |
---|
47 | // J. M. Quesada 19/04/2010 fix in emission probability calculation. |
---|
48 | // V.Ivanchenko 20/04/2010 added usage of G4Pow and use more safe computation |
---|
49 | // V.Ivanchenko 18/05/2010 trying to speedup the most slow method |
---|
50 | // by usage of G4Pow, integer Z and A; moved constructor, |
---|
51 | // destructor and virtual functions to source |
---|
52 | |
---|
53 | #include "G4GEMProbability.hh" |
---|
54 | #include "G4PairingCorrection.hh" |
---|
55 | #include "G4Pow.hh" |
---|
56 | #include "G4IonTable.hh" |
---|
57 | |
---|
58 | G4GEMProbability:: G4GEMProbability() : |
---|
59 | theA(0), theZ(0), Spin(0.0), theCoulombBarrierPtr(0), |
---|
60 | ExcitationEnergies(0), ExcitationSpins(0), ExcitationLifetimes(0), Normalization(1.0) |
---|
61 | { |
---|
62 | theEvapLDPptr = new G4EvaporationLevelDensityParameter; |
---|
63 | fG4pow = G4Pow::GetInstance(); |
---|
64 | fPairCorr = G4PairingCorrection::GetInstance(); |
---|
65 | } |
---|
66 | |
---|
67 | G4GEMProbability:: G4GEMProbability(G4int anA, G4int aZ, G4double aSpin) : |
---|
68 | theA(anA), theZ(aZ), Spin(aSpin), theCoulombBarrierPtr(0), |
---|
69 | ExcitationEnergies(0), ExcitationSpins(0), ExcitationLifetimes(0), Normalization(1.0) |
---|
70 | { |
---|
71 | theEvapLDPptr = new G4EvaporationLevelDensityParameter; |
---|
72 | fG4pow = G4Pow::GetInstance(); |
---|
73 | fPairCorr = G4PairingCorrection::GetInstance(); |
---|
74 | } |
---|
75 | |
---|
76 | G4GEMProbability::~G4GEMProbability() |
---|
77 | { |
---|
78 | delete theEvapLDPptr; |
---|
79 | } |
---|
80 | |
---|
81 | G4double G4GEMProbability::CalcAlphaParam(const G4Fragment & ) const |
---|
82 | { |
---|
83 | return 1.0; |
---|
84 | } |
---|
85 | |
---|
86 | G4double G4GEMProbability::CalcBetaParam(const G4Fragment & ) const |
---|
87 | { |
---|
88 | return 1.0; |
---|
89 | } |
---|
90 | |
---|
91 | G4double G4GEMProbability::CCoeficient(const G4double ) const |
---|
92 | { |
---|
93 | return 0.0; |
---|
94 | } |
---|
95 | |
---|
96 | G4double G4GEMProbability::EmissionProbability(const G4Fragment & fragment, |
---|
97 | G4double MaximalKineticEnergy) |
---|
98 | { |
---|
99 | G4double EmissionProbability = 0.0; |
---|
100 | |
---|
101 | if (MaximalKineticEnergy > 0.0 && fragment.GetExcitationEnergy() > 0.0) { |
---|
102 | G4double CoulombBarrier = GetCoulombBarrier(fragment); |
---|
103 | |
---|
104 | EmissionProbability = CalcProbability(fragment,MaximalKineticEnergy,CoulombBarrier); |
---|
105 | Normalization = EmissionProbability; |
---|
106 | // Next there is a loop over excited states for this channel summing probabilities |
---|
107 | if (ExcitationEnergies && ExcitationSpins && ExcitationLifetimes) { |
---|
108 | G4double SavedSpin = Spin; |
---|
109 | for (size_t i = 0; i < ExcitationEnergies->size(); ++i) { |
---|
110 | Spin = ExcitationSpins->operator[](i); |
---|
111 | // substract excitation energies |
---|
112 | G4double Tmax = MaximalKineticEnergy - ExcitationEnergies->operator[](i); |
---|
113 | if (Tmax > 0.0) { |
---|
114 | G4double width = CalcProbability(fragment,Tmax,CoulombBarrier); |
---|
115 | //JMQ April 2010 added condition to prevent reported crash |
---|
116 | // update probability |
---|
117 | if (width > 0. && |
---|
118 | hbar_Planck*fG4pow->logZ(2) < width*ExcitationLifetimes->operator[](i)) { |
---|
119 | EmissionProbability += width; |
---|
120 | } |
---|
121 | } |
---|
122 | } |
---|
123 | // Restore Spin |
---|
124 | Spin = SavedSpin; |
---|
125 | } |
---|
126 | } |
---|
127 | return EmissionProbability; |
---|
128 | } |
---|
129 | |
---|
130 | |
---|
131 | G4double G4GEMProbability::CalcProbability(const G4Fragment & fragment, |
---|
132 | G4double MaximalKineticEnergy, |
---|
133 | G4double V) |
---|
134 | |
---|
135 | // Calculate integrated probability (width) for evaporation channel |
---|
136 | { |
---|
137 | G4int A = fragment.GetA_asInt(); |
---|
138 | G4int Z = fragment.GetZ_asInt(); |
---|
139 | |
---|
140 | G4int ResidualA = A - theA; |
---|
141 | G4int ResidualZ = Z - theZ; |
---|
142 | G4double U = fragment.GetExcitationEnergy(); |
---|
143 | |
---|
144 | G4double NuclearMass = fragment.ComputeGroundStateMass(theZ, theA); |
---|
145 | |
---|
146 | G4double Alpha = CalcAlphaParam(fragment); |
---|
147 | G4double Beta = CalcBetaParam(fragment); |
---|
148 | |
---|
149 | // ***RESIDUAL*** |
---|
150 | //JMQ (September 2009) the following quantities refer to the RESIDUAL: |
---|
151 | |
---|
152 | G4double delta0 = fPairCorr->GetPairingCorrection(ResidualA, ResidualZ); |
---|
153 | |
---|
154 | G4double a = theEvapLDPptr->LevelDensityParameter(ResidualA, |
---|
155 | ResidualZ,MaximalKineticEnergy+V-delta0); |
---|
156 | G4double Ux = (2.5 + 150.0/G4double(ResidualA))*MeV; |
---|
157 | G4double Ex = Ux + delta0; |
---|
158 | G4double T = 1.0/(std::sqrt(a/Ux) - 1.5/Ux); |
---|
159 | //JMQ fixed bug in units |
---|
160 | G4double E0 = Ex - T*(std::log(T/MeV) - std::log(a*MeV)/4.0 |
---|
161 | - 1.25*std::log(Ux/MeV) + 2.0*std::sqrt(a*Ux)); |
---|
162 | // ***end RESIDUAL *** |
---|
163 | |
---|
164 | // ***PARENT*** |
---|
165 | //JMQ (September 2009) the following quantities refer to the PARENT: |
---|
166 | |
---|
167 | G4double deltaCN = fPairCorr->GetPairingCorrection(A, Z); |
---|
168 | G4double aCN = theEvapLDPptr->LevelDensityParameter(A, Z, U-deltaCN); |
---|
169 | G4double UxCN = (2.5 + 150.0/G4double(A))*MeV; |
---|
170 | G4double ExCN = UxCN + deltaCN; |
---|
171 | G4double TCN = 1.0/(std::sqrt(aCN/UxCN) - 1.5/UxCN); |
---|
172 | // ***end PARENT*** |
---|
173 | |
---|
174 | G4double Width; |
---|
175 | G4double InitialLevelDensity; |
---|
176 | G4double t = MaximalKineticEnergy/T; |
---|
177 | if ( MaximalKineticEnergy < Ex ) { |
---|
178 | //JMQ 190709 bug in I1 fixed (T was missing) |
---|
179 | Width = (I1(t,t)*T + (Beta+V)*I0(t))/std::exp(E0/T); |
---|
180 | //JMQ 160909 fix: InitialLevelDensity has been taken away |
---|
181 | //(different conditions for initial CN..) |
---|
182 | } else { |
---|
183 | |
---|
184 | //VI minor speedup |
---|
185 | G4double expE0T = std::exp(E0/T); |
---|
186 | const G4double sqrt2 = std::sqrt(2.0); |
---|
187 | |
---|
188 | G4double tx = Ex/T; |
---|
189 | G4double s = 2.0*std::sqrt(a*(MaximalKineticEnergy-delta0)); |
---|
190 | G4double sx = 2.0*std::sqrt(a*(Ex-delta0)); |
---|
191 | Width = I1(t,tx)*T/expE0T + I3(s,sx)*std::exp(s)/(sqrt2*a); |
---|
192 | // For charged particles (Beta+V) = 0 beacuse Beta = -V |
---|
193 | if (theZ == 0) { |
---|
194 | Width += (Beta+V)*(I0(tx)/expE0T + 2.0*sqrt2*I2(s,sx)*std::exp(s)); |
---|
195 | } |
---|
196 | } |
---|
197 | |
---|
198 | //JMQ 14/07/2009 BIG BUG : NuclearMass is in MeV => hbarc instead of hbar_planck must be used |
---|
199 | // G4double g = (2.0*Spin+1.0)*NuclearMass/(pi2* hbar_Planck*hbar_Planck); |
---|
200 | G4double g = (2.0*Spin+1.0)*NuclearMass/(pi2* hbarc*hbarc); |
---|
201 | |
---|
202 | //JMQ 190709 fix on Rb and geometrical cross sections according to Furihata's paper |
---|
203 | // (JAERI-Data/Code 2001-105, p6) |
---|
204 | // G4double RN = 0.0; |
---|
205 | G4double Rb = 0.0; |
---|
206 | if (theA > 4) |
---|
207 | { |
---|
208 | // G4double R1 = std::pow(ResidualA,1.0/3.0); |
---|
209 | // G4double R2 = std::pow(G4double(theA),1.0/3.0); |
---|
210 | G4double Ad = fG4pow->Z13(ResidualA); |
---|
211 | G4double Aj = fG4pow->Z13(theA); |
---|
212 | // RN = 1.12*(R1 + R2) - 0.86*((R1+R2)/(R1*R2)); |
---|
213 | Rb = 1.12*(Aj + Ad) - 0.86*((Aj+Ad)/(Aj*Ad))+2.85; |
---|
214 | Rb *= fermi; |
---|
215 | } |
---|
216 | else if (theA>1) |
---|
217 | { |
---|
218 | G4double Ad = fG4pow->Z13(ResidualA); |
---|
219 | G4double Aj = fG4pow->Z13(theA); |
---|
220 | Rb=1.5*(Aj+Ad)*fermi; |
---|
221 | } |
---|
222 | else |
---|
223 | { |
---|
224 | G4double Ad = fG4pow->Z13(ResidualA); |
---|
225 | Rb = 1.5*Ad*fermi; |
---|
226 | } |
---|
227 | // G4double GeometricalXS = pi*RN*RN*std::pow(ResidualA,2./3.); |
---|
228 | G4double GeometricalXS = pi*Rb*Rb; |
---|
229 | //end of JMQ fix on Rb by 190709 |
---|
230 | |
---|
231 | //JMQ 160909 fix: initial level density must be calculated according to the |
---|
232 | // conditions at the initial compound nucleus |
---|
233 | // (it has been removed from previous "if" for the residual) |
---|
234 | |
---|
235 | if ( U < ExCN ) |
---|
236 | { |
---|
237 | //JMQ fixed bug in units |
---|
238 | //VI moved the computation here |
---|
239 | G4double E0CN = ExCN - TCN*(std::log(TCN/MeV) - std::log(aCN*MeV)/4.0 |
---|
240 | - 1.25*std::log(UxCN/MeV) + 2.0*std::sqrt(aCN*UxCN)); |
---|
241 | InitialLevelDensity = (pi/12.0)*std::exp((U-E0CN)/TCN)/TCN; |
---|
242 | } |
---|
243 | else |
---|
244 | { |
---|
245 | //InitialLevelDensity = (pi/12.0)*std::exp(2*std::sqrt(aCN*(U-deltaCN)))/std::pow(aCN*std::pow(U-deltaCN,5.0),1.0/4.0); |
---|
246 | //VI speedup |
---|
247 | G4double x = U-deltaCN; |
---|
248 | G4double x1 = std::sqrt(aCN*x); |
---|
249 | InitialLevelDensity = (pi/12.0)*std::exp(2*x1)/(x*std::sqrt(x1)); |
---|
250 | } |
---|
251 | |
---|
252 | //JMQ 190709 BUG : pi instead of sqrt(pi) must be here according to Furihata's report: |
---|
253 | // Width *= std::sqrt(pi)*g*GeometricalXS*Alpha/(12.0*InitialLevelDensity); |
---|
254 | Width *= pi*g*GeometricalXS*Alpha/(12.0*InitialLevelDensity); |
---|
255 | |
---|
256 | return Width; |
---|
257 | } |
---|
258 | |
---|
259 | G4double G4GEMProbability::I3(G4double s, G4double sx) |
---|
260 | { |
---|
261 | G4double s2 = s*s; |
---|
262 | G4double sx2 = sx*sx; |
---|
263 | G4double S = 1.0/std::sqrt(s); |
---|
264 | G4double S2 = S*S; |
---|
265 | G4double Sx = 1.0/std::sqrt(sx); |
---|
266 | G4double Sx2 = Sx*Sx; |
---|
267 | |
---|
268 | G4double p1 = S *(2.0 + S2 *( 4.0 + S2 *( 13.5 + S2 *( 60.0 + S2 * 325.125 )))); |
---|
269 | G4double p2 = Sx*Sx2 *( |
---|
270 | (s2-sx2) + Sx2 *( |
---|
271 | (1.5*s2+0.5*sx2) + Sx2 *( |
---|
272 | (3.75*s2+0.25*sx2) + Sx2 *( |
---|
273 | (12.875*s2+0.625*sx2) + Sx2 *( |
---|
274 | (59.0625*s2+0.9375*sx2) + Sx2 *(324.8*s2+3.28*sx2)))))); |
---|
275 | |
---|
276 | p2 *= std::exp(sx-s); |
---|
277 | |
---|
278 | return p1-p2; |
---|
279 | } |
---|