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: G4AblaFissionSimfis18.cc,v 1.5 2010/11/17 20:19:09 kaitanie Exp $ |
---|
27 | // Translation of INCL4.2/ABLA V3 |
---|
28 | // Pekka Kaitaniemi, HIP (translation) |
---|
29 | // Christelle Schmidt, IPNL (fission code) |
---|
30 | // Alain Boudard, CEA (contact person INCL/ABLA) |
---|
31 | // Aatos Heikkinen, HIP (project coordination) |
---|
32 | |
---|
33 | #include "G4AblaFissionSimfis18.hh" |
---|
34 | #include <time.h> |
---|
35 | |
---|
36 | G4AblaFissionSimfis18::G4AblaFissionSimfis18() |
---|
37 | { |
---|
38 | hazard = 0; |
---|
39 | randomGenerator = 0; |
---|
40 | } |
---|
41 | |
---|
42 | G4AblaFissionSimfis18::G4AblaFissionSimfis18(G4Hazard *hzr, G4InclRandomInterface *rndm) |
---|
43 | { |
---|
44 | hazard = hzr; |
---|
45 | randomGenerator = rndm; |
---|
46 | setAboutString("Fission model: Based on ABLA with SimFis18"); |
---|
47 | } |
---|
48 | |
---|
49 | G4AblaFissionSimfis18::~G4AblaFissionSimfis18() |
---|
50 | { |
---|
51 | |
---|
52 | } |
---|
53 | |
---|
54 | void G4AblaFissionSimfis18::doFission(G4double &A, G4double &Z, G4double &E, |
---|
55 | G4double &A1, G4double &Z1, G4double &E1, G4double &K1, |
---|
56 | G4double &A2, G4double &Z2, G4double &E2, G4double &K2) |
---|
57 | { |
---|
58 | fissionDistri(A,Z,E,A1,Z1,E1,K1,A2,Z2,E2,K2); |
---|
59 | } |
---|
60 | |
---|
61 | void G4AblaFissionSimfis18::even_odd(G4double r_origin,G4double r_even_odd,G4int &i_out) |
---|
62 | { |
---|
63 | // Procedure to calculate I_OUT from R_IN in a way that |
---|
64 | // on the average a flat distribution in R_IN results in a |
---|
65 | // fluctuating distribution in I_OUT with an even-odd effect as |
---|
66 | // given by R_EVEN_ODD |
---|
67 | |
---|
68 | // /* ------------------------------------------------------------ */ |
---|
69 | // /* EXAMPLES : */ |
---|
70 | // /* ------------------------------------------------------------ */ |
---|
71 | // /* If R_EVEN_ODD = 0 : */ |
---|
72 | // /* CEIL(R_IN) ---- */ |
---|
73 | // /* */ |
---|
74 | // /* R_IN -> */ |
---|
75 | // /* (somewhere in between CEIL(R_IN) and FLOOR(R_IN)) */ */ |
---|
76 | // /* */ |
---|
77 | // /* FLOOR(R_IN) ---- --> I_OUT */ |
---|
78 | // /* ------------------------------------------------------------ */ |
---|
79 | // /* If R_EVEN_ODD > 0 : */ |
---|
80 | // /* The interval for the above treatment is */ |
---|
81 | // /* larger for FLOOR(R_IN) = even and */ |
---|
82 | // /* smaller for FLOOR(R_IN) = odd */ |
---|
83 | // /* For R_EVEN_ODD < 0 : just opposite treatment */ |
---|
84 | // /* ------------------------------------------------------------ */ |
---|
85 | |
---|
86 | // /* ------------------------------------------------------------ */ |
---|
87 | // /* On input: R_ORIGIN nuclear charge (real number) */ |
---|
88 | // /* R_EVEN_ODD requested even-odd effect */ |
---|
89 | // /* Intermediate quantity: R_IN = R_ORIGIN + 0.5 */ |
---|
90 | // /* On output: I_OUT nuclear charge (integer) */ |
---|
91 | // /* ------------------------------------------------------------ */ |
---|
92 | |
---|
93 | // G4double R_ORIGIN,R_IN,R_EVEN_ODD,R_REST,R_HELP; |
---|
94 | G4double r_in = 0.0, r_rest = 0.0, r_help = 0.0; |
---|
95 | G4double r_floor = 0.0; |
---|
96 | G4double r_middle = 0.0; |
---|
97 | // G4int I_OUT,N_FLOOR; |
---|
98 | G4int n_floor = 0; |
---|
99 | |
---|
100 | r_in = r_origin + 0.5; |
---|
101 | r_floor = (float)((int)(r_in)); |
---|
102 | if (r_even_odd < 0.001) { |
---|
103 | i_out = (int)(r_floor); |
---|
104 | } |
---|
105 | else { |
---|
106 | r_rest = r_in - r_floor; |
---|
107 | r_middle = r_floor + 0.5; |
---|
108 | n_floor = (int)(r_floor); |
---|
109 | if (n_floor%2 == 0) { |
---|
110 | // even before modif. |
---|
111 | r_help = r_middle + (r_rest - 0.5) * (1.0 - r_even_odd); |
---|
112 | } |
---|
113 | else { |
---|
114 | // odd before modification |
---|
115 | r_help = r_middle + (r_rest - 0.5) * (1.0 + r_even_odd); |
---|
116 | } |
---|
117 | i_out = (int)(r_help); |
---|
118 | } |
---|
119 | } |
---|
120 | |
---|
121 | G4double G4AblaFissionSimfis18::umass(G4double z,G4double n,G4double beta) |
---|
122 | { |
---|
123 | // liquid-drop mass, Myers & Swiatecki, Lysekil, 1967 |
---|
124 | // pure liquid drop, without pairing and shell effects |
---|
125 | |
---|
126 | // On input: Z nuclear charge of nucleus |
---|
127 | // N number of neutrons in nucleus |
---|
128 | // beta deformation of nucleus |
---|
129 | // On output: binding energy of nucleus |
---|
130 | |
---|
131 | G4double a = 0.0, umass = 0.0; |
---|
132 | G4double alpha = 0.0; |
---|
133 | G4double xcom = 0.0, xvs = 0.0, xe = 0.0; |
---|
134 | const G4double pi = 3.1416; |
---|
135 | |
---|
136 | a = n + z; |
---|
137 | alpha = ( std::sqrt(5.0/(4.0*pi)) ) * beta; |
---|
138 | |
---|
139 | xcom = 1.0 - 1.7826 * ((a - 2.0*z)/a)*((a - 2.0*z)/a); |
---|
140 | // factor for asymmetry dependence of surface and volume term |
---|
141 | xvs = - xcom * ( 15.4941 * a - |
---|
142 | 17.9439 * std::pow(a,0.66667) * (1.0+0.4*alpha*alpha) ); |
---|
143 | // sum of volume and surface energy |
---|
144 | xe = z*z * (0.7053/(std::pow(a,0.33333)) * (1.0-0.2*alpha*alpha) - 1.1529/a); |
---|
145 | umass = xvs + xe; |
---|
146 | |
---|
147 | return umass; |
---|
148 | } |
---|
149 | |
---|
150 | G4double G4AblaFissionSimfis18::ecoul(G4double z1,G4double n1,G4double beta1,G4double z2,G4double n2,G4double beta2,G4double d) |
---|
151 | { |
---|
152 | // Coulomb potential between two nuclei |
---|
153 | // surfaces are in a distance of d |
---|
154 | // in a tip to tip configuration |
---|
155 | |
---|
156 | // approximate formulation |
---|
157 | // On input: Z1 nuclear charge of first nucleus |
---|
158 | // N1 number of neutrons in first nucleus |
---|
159 | // beta1 deformation of first nucleus |
---|
160 | // Z2 nuclear charge of second nucleus |
---|
161 | // N2 number of neutrons in second nucleus |
---|
162 | // beta2 deformation of second nucleus |
---|
163 | // d distance of surfaces of the nuclei |
---|
164 | |
---|
165 | // G4double Z1,N1,beta1,Z2,N2,beta2,d,ecoul; |
---|
166 | G4double ecoul = 0; |
---|
167 | G4double dtot = 0; |
---|
168 | const G4double r0 = 1.16; |
---|
169 | |
---|
170 | dtot = r0 * ( std::pow((z1+n1),0.33333) * (1.0+(2.0/3.0)*beta1) |
---|
171 | + std::pow((z2+n2),0.33333) * (1.0+(2.0/3.0)*beta2) ) + d; |
---|
172 | ecoul = z1 * z2 * 1.44 / dtot; |
---|
173 | |
---|
174 | return ecoul; |
---|
175 | } |
---|
176 | |
---|
177 | void G4AblaFissionSimfis18::fissionDistri(G4double &a,G4double &z,G4double &e, |
---|
178 | G4double &a1,G4double &z1,G4double &e1,G4double &v1, |
---|
179 | G4double &a2,G4double &z2,G4double &e2,G4double &v2) |
---|
180 | { |
---|
181 | // On input: A, Z, E (mass, atomic number and exc. energy of compound nucleus |
---|
182 | // before fission) |
---|
183 | // On output: Ai, Zi, Ei (mass, atomic number and exc. energy of fragment 1 and 2 |
---|
184 | // after fission) |
---|
185 | |
---|
186 | // Additionally calculated but not put in the parameter list: |
---|
187 | // Kinetic energy of prefragments EkinR1, EkinR2 |
---|
188 | |
---|
189 | // Translation of SIMFIS18.PLI (KHS, 2.1.2001) |
---|
190 | |
---|
191 | // This program calculates isotopic distributions of fission fragments |
---|
192 | // with a semiempirical model |
---|
193 | // Copy from SIMFIS3, KHS, 8. February 1995 |
---|
194 | // Modifications made by Jose Benlliure and KHS in August 1996 |
---|
195 | // Energy counted from lowest barrier (J. Benlliure, KHS 1997) |
---|
196 | // Some bugs corrected (J. Benlliure, KHS 1997) |
---|
197 | // Version used for thesis S. Steinhaueser (August 1997) |
---|
198 | // (Curvature of LD potential increased by factor of 2!) |
---|
199 | |
---|
200 | // Weiter veraendert mit der Absicht, eine Version zu erhalten, die |
---|
201 | // derjenigen entspricht, die von J. Benlliure et al. |
---|
202 | // in Nucl. Phys. A 628 (1998) 458 verwendet wurde, |
---|
203 | // allerdings ohne volle Neutronenabdampfung. |
---|
204 | |
---|
205 | // The excitation energy was calculate now for each fission channel |
---|
206 | // separately. The dissipation from saddle to scission was taken from |
---|
207 | // systematics, the deformation energy at scission considers the shell |
---|
208 | // effects in a simplified way, and the fluctuation is included. |
---|
209 | // KHS, April 1999 |
---|
210 | |
---|
211 | // The width in N/Z was carefully adapted to values given by Lang et al. |
---|
212 | |
---|
213 | // The width and eventually a shift in N/Z (polarization) follows the |
---|
214 | // following rules: |
---|
215 | |
---|
216 | // The line N/Z following UCD has an angle of std::atan(Zcn/Ncn) |
---|
217 | // to the horizontal axis on a chart of nuclides. |
---|
218 | // (For 238U the angle is 32.2 deg.) |
---|
219 | |
---|
220 | // The following relations hold: (from Armbruster) |
---|
221 | // |
---|
222 | // sigma(N) (A=const) = sigma(Z) (A=const) |
---|
223 | // sigma(A) (N=const) = sigma(Z) (N=const) |
---|
224 | // sigma(A) (Z=const) = sigma(N) (Z=const) |
---|
225 | // |
---|
226 | // From this we get: |
---|
227 | // sigma(Z) (N=const) * N = sigma(N) (Z=const) * Z |
---|
228 | // sigma(A) (Z=const) = sigma(Z) (A=const) * A/Z |
---|
229 | // sigma(N) (Z=const) = sigma(Z) (A=const) * A/Z |
---|
230 | // Z*sigma(N) (Z=const) = N*sigma(Z) (N=const) = A*sigma(Z) (A=const) |
---|
231 | |
---|
232 | // Excitation energy now calculated above the lowest potential point |
---|
233 | // Inclusion of a distribution of excitation energies |
---|
234 | |
---|
235 | // Several modifications, starting from SIMFIS12: KHS November 2000 |
---|
236 | // This version seems to work quite well for 238U. |
---|
237 | // The transition from symmetric to asymmetric fission around 226Th |
---|
238 | // is reasonably well reproduced, although St. I is too strong and St. II |
---|
239 | // is too weak. St. I and St. II are also weakly seen for 208Pb. |
---|
240 | |
---|
241 | // Extensions for an event generator of fission events (21.11.2000,KHS) |
---|
242 | |
---|
243 | // Defalt parameters (IPARS) rather carefully adjusted to |
---|
244 | // pre-neutron mass distributions of Vives et al. (238U + n) |
---|
245 | // Die Parameter Fgamma1 und Fgamma2 sind kleiner als die resultierenden |
---|
246 | // Breiten der Massenverteilungen!!! |
---|
247 | // Fgamma1 und Fgamma2 wurden angepaᅵ, so daᅵ |
---|
248 | // Sigma-A(ST-I) = 3.3, Sigma-A(St-II) = 5.8 (nach Vives) |
---|
249 | |
---|
250 | // Parameters of the model carefully adjusted by KHS (2.2.2001) to |
---|
251 | // 238U + 208Pb, 1000 A MeV, Timo Enqvist et al. |
---|
252 | |
---|
253 | |
---|
254 | G4double n = 0.0; |
---|
255 | G4double nlight1 = 0.0, nlight2 = 0.0; |
---|
256 | G4double aheavy1 = 0.0,alight1 = 0.0, aheavy2 = 0.0, alight2 = 0.0; |
---|
257 | G4double eheavy1 = 0.0, elight1 = 0.0, eheavy2 = 0.0, elight2 = 0.0; |
---|
258 | G4double zheavy1_shell = 0.0, zheavy2_shell = 0.0; |
---|
259 | G4double zlight1 = 0.0, zlight2 = 0.0; |
---|
260 | G4double masscurv = 0.0; |
---|
261 | G4double sasymm1 = 0.0, sasymm2 = 0.0, ssymm = 0.0, ysum = 0.0, yasymm = 0.0; |
---|
262 | G4double ssymm_mode1 = 0.0, ssymm_mode2 = 0.0; |
---|
263 | G4double cz_asymm1_saddle = 0.0, cz_asymm2_saddle = 0.0; |
---|
264 | // Curvature at saddle, modified by ld-potential |
---|
265 | G4double wzasymm1_saddle, wzasymm2_saddle, wzsymm_saddle = 0.0; |
---|
266 | G4double wzasymm1_scission = 0.0, wzasymm2_scission = 0.0, wzsymm_scission = 0.0; |
---|
267 | G4double wzasymm1 = 0.0, wzasymm2 = 0.0, wzsymm = 0.0; |
---|
268 | G4double nlight1_eff = 0.0, nlight2_eff = 0.0; |
---|
269 | G4int imode = 0; |
---|
270 | G4double rmode = 0.0; |
---|
271 | G4double z1mean = 0.0, z2mean = 0.0, z1width = 0.0, za1width = 0.0; |
---|
272 | // G4double Z1,Z2,N1R,N2R,A1R,A2R,N1,N2,A1,A2; |
---|
273 | G4double n1r = 0.0, n2r = 0.0, a1r = 0.0, a2r = 0.0, n1 = 0.0, n2 = 0.0; |
---|
274 | |
---|
275 | G4double zsymm = 0.0, nsymm = 0.0, asymm = 0.0; |
---|
276 | G4double n1mean = 0.0, n2mean, n1width; |
---|
277 | G4double dueff = 0.0; |
---|
278 | // effective shell effect at lowest barrier |
---|
279 | G4double eld = 0.0; |
---|
280 | // Excitation energy with respect to ld barrier |
---|
281 | G4double re1 = 0.0, re2 = 0.0, re3 = 0.0; |
---|
282 | G4double eps1 = 0.0, eps2 = 0.0; |
---|
283 | G4double n1ucd = 0.0, n2ucd = 0.0, z1ucd = 0.0, z2ucd = 0.0; |
---|
284 | G4double beta = 0.0, beta1 = 0.0, beta2 = 0.0; |
---|
285 | |
---|
286 | G4double dn1_pol = 0.0; |
---|
287 | // shift of most probable neutron number for given Z, |
---|
288 | // according to polarization |
---|
289 | G4int i_help = 0; |
---|
290 | |
---|
291 | // /* Parameters of the semiempirical fission model */ |
---|
292 | G4double a_levdens = 0.0; |
---|
293 | // /* level-density parameter */ |
---|
294 | G4double a_levdens_light1 = 0.0, a_levdens_light2 = 0.0; |
---|
295 | G4double a_levdens_heavy1 = 0.0, a_levdens_heavy2 = 0.0; |
---|
296 | const G4double r_null = 1.16; |
---|
297 | // /* radius parameter */ |
---|
298 | G4double epsilon_1_saddle = 0.0, epsilon0_1_saddle = 0.0; |
---|
299 | G4double epsilon_2_saddle = 0.0, epsilon0_2_saddle = 0.0, epsilon_symm_saddle = 0.0; |
---|
300 | G4double epsilon_1_scission = 0.0, epsilon0_1_scission = 0.0; |
---|
301 | G4double epsilon_2_scission = 0.0, epsilon0_2_scission = 0.0; |
---|
302 | G4double epsilon_symm_scission = 0.0; |
---|
303 | // /* modified energy */ |
---|
304 | G4double e_eff1_saddle = 0.0, e_eff2_saddle = 0.0; |
---|
305 | G4double epot0_mode1_saddle = 0.0, epot0_mode2_saddle = 0.0, epot0_symm_saddle = 0.0; |
---|
306 | G4double epot_mode1_saddle = 0.0, epot_mode2_saddle = 0.0, epot_symm_saddle = 0.0; |
---|
307 | G4double e_defo = 0.0, e_defo1 = 0.0, e_defo2 = 0.0, e_scission = 0.0, e_asym = 0.0; |
---|
308 | G4double e1exc = 0.0, e2exc = 0.0; |
---|
309 | G4double e1exc_sigma = 0.0, e2exc_sigma = 0.0; |
---|
310 | G4double e1final = 0.0, e2final = 0.0; |
---|
311 | |
---|
312 | const G4double r0 = 1.16; |
---|
313 | G4double tker = 0.0; |
---|
314 | G4double ekin1 = 0.0, ekin2 = 0.0; |
---|
315 | // G4double EkinR1,EkinR2,E1,E2,V1,V2; |
---|
316 | G4double ekinr1 = 0.0, ekinr2 = 0.0; |
---|
317 | G4int icz = 0, k = 0; |
---|
318 | |
---|
319 | // Input parameters: |
---|
320 | //OMMENT(Nuclear charge number); |
---|
321 | // G4double Z; |
---|
322 | //OMMENT(Nuclear mass number); |
---|
323 | // G4double A; |
---|
324 | //OMMENT(Excitation energy above fission barrier); |
---|
325 | // G4double E; |
---|
326 | |
---|
327 | // Model parameters: |
---|
328 | //OMMENT(position of heavy peak valley 1); |
---|
329 | const G4double nheavy1 = 83.0; |
---|
330 | //OMMENT(position of heavy peak valley 2); |
---|
331 | const G4double nheavy2 = 90.0; |
---|
332 | //OMMENT(Shell effect for valley 1); |
---|
333 | const G4double delta_u1_shell = -2.65; |
---|
334 | // Parameter (Delta_U1_shell = -2) |
---|
335 | //OMMENT(Shell effect for valley 2); |
---|
336 | const G4double delta_u2_shell = -3.8; |
---|
337 | // Parameter (Delta_U2_shell = -3.2) |
---|
338 | //OMMENT(I: used shell effect); |
---|
339 | G4double delta_u1 = 0.0; |
---|
340 | //omment(I: used shell effect); |
---|
341 | G4double delta_u2 = 0.0; |
---|
342 | //OMMENT(Curvature of asymmetric valley 1); |
---|
343 | const G4double cz_asymm1_shell = 0.7; |
---|
344 | //OMMENT(Curvature of asymmetric valley 2); |
---|
345 | const G4double cz_asymm2_shell = 0.15; |
---|
346 | //OMMENT(Factor for width of distr. valley 1); |
---|
347 | const G4double fwidth_asymm1 = 0.63; |
---|
348 | //OMMENT(Factor for width of distr. valley 2); |
---|
349 | const G4double fwidth_asymm2 = 0.97; |
---|
350 | // Parameter (CZ_asymm2_scission = 0.12) |
---|
351 | //OMMENT(Parameter x: a = A/x); |
---|
352 | const G4double xlevdens = 12.0; |
---|
353 | //OMMENT(Factor to gamma_heavy1); |
---|
354 | const G4double fgamma1 = 2.0; |
---|
355 | //OMMENT(I: fading of shells (general)); |
---|
356 | G4double gamma = 0.0; |
---|
357 | //OMMENT(I: fading of shell 1); |
---|
358 | G4double gamma_heavy1 = 0.0; |
---|
359 | //OMMENT(I: fading of shell 2); |
---|
360 | G4double gamma_heavy2 = 0.0; |
---|
361 | //OMMENT(Zero-point energy at saddle); |
---|
362 | const G4double e_zero_point = 0.5; |
---|
363 | //OMMENT(I: friction from saddle to scission); |
---|
364 | G4double e_saddle_scission = 0.0; |
---|
365 | //OMMENT(Friction factor); |
---|
366 | const G4double friction_factor = 1.0; |
---|
367 | //OMMENT(I: Internal counter for different modes); INIT(0,0,0) |
---|
368 | // Integer*4 I_MODE(3) |
---|
369 | //OMMENT(I: Yield of symmetric mode); |
---|
370 | G4double ysymm = 0.0; |
---|
371 | //OMMENT(I: Yield of asymmetric mode 1); |
---|
372 | G4double yasymm1 = 0.0; |
---|
373 | //OMMENT(I: Yield of asymmetric mode 2); |
---|
374 | G4double yasymm2 = 0.0; |
---|
375 | //OMMENT(I: Effective position of valley 1); |
---|
376 | G4double nheavy1_eff = 0.0; |
---|
377 | //OMMENT(I: position of heavy peak valley 1); |
---|
378 | G4double zheavy1 = 0.0; |
---|
379 | //omment(I: Effective position of valley 2); |
---|
380 | G4double nheavy2_eff = 0.0; |
---|
381 | //OMMENT(I: position of heavy peak valley 2); |
---|
382 | G4double zheavy2 = 0.0; |
---|
383 | //omment(I: Excitation energy above saddle 1); |
---|
384 | G4double eexc1_saddle = 0.0; |
---|
385 | //omment(I: Excitation energy above saddle 2); |
---|
386 | G4double eexc2_saddle = 0.0; |
---|
387 | //omment(I: Excitation energy above lowest saddle); |
---|
388 | G4double eexc_max = 0.0; |
---|
389 | //omment(I: Effective mass mode 1); |
---|
390 | G4double aheavy1_mean = 0.0; |
---|
391 | //omment(I: Effective mass mode 2); |
---|
392 | G4double aheavy2_mean = 0.0; |
---|
393 | //omment(I: Width of symmetric mode); |
---|
394 | G4double wasymm_saddle = 0.0; |
---|
395 | //OMMENT(I: Width of asymmetric mode 1); |
---|
396 | G4double waheavy1_saddle = 0.0; |
---|
397 | //OMMENT(I: Width of asymmetric mode 2); |
---|
398 | G4double waheavy2_saddle = 0.0; |
---|
399 | //omment(I: Width of symmetric mode); |
---|
400 | G4double wasymm = 0.0; |
---|
401 | //OMMENT(I: Width of asymmetric mode 1); |
---|
402 | G4double waheavy1 = 0.0; |
---|
403 | //OMMENT(I: Width of asymmetric mode 2); |
---|
404 | G4double waheavy2 = 0.0; |
---|
405 | //OMMENT(I: Even-odd effect in Z); |
---|
406 | G4double r_e_o = 0.0, r_e_o_exp = 0.0; |
---|
407 | //OMMENT(I: Curveture of symmetric valley); |
---|
408 | G4double cz_symm = 0.0; |
---|
409 | //OMMENT(I: Curvature of mass distribution for fixed Z); |
---|
410 | G4double cn = 0.0; |
---|
411 | //OMMENT(I: Curvature of Z distribution for fixed A); |
---|
412 | G4double cz = 0.0; |
---|
413 | //OMMENT(Minimum neutron width for constant Z); |
---|
414 | const G4double sigzmin = 1.16; |
---|
415 | //OMMENT(Surface distance of scission configuration); |
---|
416 | const G4double d = 2.0; |
---|
417 | |
---|
418 | // /* Charge polarisation from Wagemanns p. 397: */ |
---|
419 | //OMMENT(Charge polarisation standard I); |
---|
420 | const G4double cpol1 = 0.65; |
---|
421 | //OMMENT(Charge polarisation standard II); |
---|
422 | const G4double cpol2 = 0.55; |
---|
423 | //OMMENT(=1: Polarisation simult. in N and Z); |
---|
424 | const G4int nzpol = 1; |
---|
425 | //OMMENT(=1: test output, =0: no test output); |
---|
426 | const G4int itest = 0; |
---|
427 | |
---|
428 | // G4double UMASS, ECOUL, reps1, reps2, rn1_pol; |
---|
429 | G4double reps1 = 0.0, reps2 = 0.0, rn1_pol = 0.0; |
---|
430 | // Float_t HAZ,GAUSSHAZ; |
---|
431 | G4int kkk = 0; |
---|
432 | // G4int kkk = 10; // PK |
---|
433 | |
---|
434 | // I_MODE = 0; |
---|
435 | |
---|
436 | if(itest == 1) { |
---|
437 | G4cout << " cn mass " << a << G4endl; |
---|
438 | G4cout << " cn charge " << z << G4endl; |
---|
439 | G4cout << " cn energy " << e << G4endl; |
---|
440 | } |
---|
441 | |
---|
442 | // /* average Z of asymmetric and symmetric components: */ |
---|
443 | n = a - z; /* neutron number of the fissioning nucleus */ |
---|
444 | |
---|
445 | k = 0; |
---|
446 | icz = 0; |
---|
447 | if ( (std::pow(z,2)/a < 25.0) || (n < nheavy2) || (e > 500.0) ) { |
---|
448 | icz = -1; |
---|
449 | // GOTO 1002; |
---|
450 | goto milledeux; |
---|
451 | } |
---|
452 | |
---|
453 | nlight1 = n - nheavy1; |
---|
454 | nlight2 = n - nheavy2; |
---|
455 | |
---|
456 | // /* Polarisation assumed for standard I and standard II: |
---|
457 | // Z - Zucd = cpol (for A = const); |
---|
458 | // from this we get (see Armbruster) |
---|
459 | // Z - Zucd = Acn/Ncn * cpol (for N = const) */ |
---|
460 | |
---|
461 | zheavy1_shell = ((nheavy1/n) * z) - ((a/n) * cpol1); |
---|
462 | zheavy2_shell = ((nheavy2/n) * z) - ((a/n) * cpol2); |
---|
463 | |
---|
464 | e_saddle_scission = |
---|
465 | (-24.0 + 0.02227 * (std::pow(z,2))/(std::pow(a,0.33333)) ) * friction_factor; |
---|
466 | |
---|
467 | // /* Energy dissipated from saddle to scission */ |
---|
468 | // /* F. Rejmund et al., Nucl. Phys. A 678 (2000) 215, fig. 4 b */ |
---|
469 | // E_saddle_scission = DMAX1(0.,E_saddle_scission); |
---|
470 | if (e_saddle_scission > 0.) { |
---|
471 | e_saddle_scission = e_saddle_scission; |
---|
472 | } |
---|
473 | else { |
---|
474 | e_saddle_scission = 0.; |
---|
475 | } |
---|
476 | // /* Semiempirical fission model: */ |
---|
477 | |
---|
478 | // /* Fit to experimental result on curvature of potential at saddle */ |
---|
479 | // /* reference: */ |
---|
480 | // /* IF Z**2/A < 33.15E0 THEN |
---|
481 | // MassCurv = 30.5438538E0 - 4.00212049E0 * Z**2/A |
---|
482 | // + 0.11983384E0 * Z**4 / (A**2) ; |
---|
483 | // ELSE |
---|
484 | // MassCurv = 10.E0 ** (7.16993332E0 - 0.26602401E0 * Z**2/A |
---|
485 | // + 0.00283802E0 * Z**4 / (A**2)) ; */ |
---|
486 | // /* New parametrization of T. Enqvist according to Mulgin et al. 1998 */ |
---|
487 | if ( (std::pow(z,2))/a < 34.0) { |
---|
488 | masscurv = std::pow( 10.0,(-1.093364 + 0.082933 * (std::pow(z,2)/a) |
---|
489 | - 0.0002602 * (std::pow(z,4)/std::pow(a,2))) ); |
---|
490 | } else { |
---|
491 | masscurv = std::pow( 10.0,(3.053536 - 0.056477 * (std::pow(z,2)/a) |
---|
492 | + 0.0002454 * (std::pow(z,4)/std::pow(a,2))) ); |
---|
493 | } |
---|
494 | |
---|
495 | cz_symm = (8.0/std::pow(z,2)) * masscurv; |
---|
496 | |
---|
497 | if(itest == 1) { |
---|
498 | G4cout << "cz_symmetry= " << cz_symm << G4endl; |
---|
499 | } |
---|
500 | |
---|
501 | if (cz_symm < 0) { |
---|
502 | icz = -1; |
---|
503 | // GOTO 1002; |
---|
504 | goto milledeux; |
---|
505 | } |
---|
506 | |
---|
507 | // /* proton number in symmetric fission (centre) */ |
---|
508 | zsymm = z/2.0; |
---|
509 | nsymm = n/2.0; |
---|
510 | asymm = nsymm + zsymm; |
---|
511 | |
---|
512 | zheavy1 = (cz_symm*zsymm + cz_asymm1_shell*zheavy1_shell)/(cz_symm + cz_asymm1_shell); |
---|
513 | zheavy2 = (cz_symm*zsymm + cz_asymm2_shell*zheavy2_shell)/(cz_symm + cz_asymm2_shell); |
---|
514 | // /* position of valley due to influence of liquid-drop potential */ |
---|
515 | nheavy1_eff = (zheavy1 + (a/n * cpol1))*(n/z); |
---|
516 | nheavy2_eff = (zheavy2 + (a/n * cpol2))*(n/z); |
---|
517 | nlight1_eff = n - nheavy1_eff; |
---|
518 | nlight2_eff = n - nheavy2_eff; |
---|
519 | // /* proton number of light fragments (centre) */ |
---|
520 | zlight1 = z - zheavy1; |
---|
521 | // /* proton number of light fragments (centre) */ |
---|
522 | zlight2 = z - zheavy2; |
---|
523 | aheavy1 = nheavy1_eff + zheavy1; |
---|
524 | aheavy2 = nheavy2_eff + zheavy2; |
---|
525 | aheavy1_mean = aheavy1; |
---|
526 | aheavy2_mean = aheavy2; |
---|
527 | alight1 = nlight1_eff + zlight1; |
---|
528 | alight2 = nlight2_eff + zlight2; |
---|
529 | |
---|
530 | a_levdens = a / xlevdens; |
---|
531 | a_levdens_heavy1 = aheavy1 / xlevdens; |
---|
532 | a_levdens_heavy2 = aheavy2 / xlevdens; |
---|
533 | a_levdens_light1 = alight1 / xlevdens; |
---|
534 | a_levdens_light2 = alight2 / xlevdens; |
---|
535 | gamma = a_levdens / (0.4 * (std::pow(a,1.3333)) ); |
---|
536 | gamma_heavy1 = ( a_levdens_heavy1 / (0.4 * (std::pow(aheavy1,1.3333)) ) ) * fgamma1; |
---|
537 | gamma_heavy2 = a_levdens_heavy2 / (0.4 * (std::pow(aheavy2,1.3333)) ); |
---|
538 | |
---|
539 | cz_asymm1_saddle = cz_asymm1_shell + cz_symm; |
---|
540 | cz_asymm2_saddle = cz_asymm2_shell + cz_symm; |
---|
541 | |
---|
542 | // Up to here: Ok! Checked CS 10/10/05 |
---|
543 | |
---|
544 | cn = umass(zsymm,(nsymm+1.),0.0) + umass(zsymm,(nsymm-1.),0.0) |
---|
545 | + 1.44 * (std::pow(zsymm,2))/ |
---|
546 | ( (std::pow(r_null,2)) * |
---|
547 | ( std::pow((asymm+1.0),0.33333) + std::pow((asymm-1.0),0.33333) ) * |
---|
548 | ( std::pow((asymm+1.0),0.33333) + std::pow((asymm-1.0),0.33333) ) ) |
---|
549 | - 2.0 * umass(zsymm,nsymm,0.0) |
---|
550 | - 1.44 * (std::pow(zsymm,2))/ |
---|
551 | ( ( 2.0 * r_null * (std::pow(asymm,0.33333)) ) * |
---|
552 | ( 2.0 * r_null * (std::pow(asymm,0.33333)) ) ); |
---|
553 | |
---|
554 | // /* shell effect in valley of mode 1 */ |
---|
555 | delta_u1 = delta_u1_shell + (std::pow((zheavy1_shell-zheavy1),2))*cz_asymm1_shell; |
---|
556 | // /* shell effect in valley of mode 2 */ |
---|
557 | delta_u2 = delta_u2_shell + (std::pow((zheavy2_shell-zheavy2),2))*cz_asymm2_shell; |
---|
558 | |
---|
559 | // /* liquid drop energies |
---|
560 | // at the centres of the different shell effects |
---|
561 | // with respect to liquid drop at symmetry: */ |
---|
562 | epot0_mode1_saddle = (std::pow((zheavy1-zsymm),2)) * cz_symm; |
---|
563 | epot0_mode2_saddle = (std::pow((zheavy2-zsymm),2)) * cz_symm; |
---|
564 | epot0_symm_saddle = 0.0; |
---|
565 | |
---|
566 | if (itest == 1) { |
---|
567 | G4cout << "check zheavy1 = " << zheavy1 << G4endl; |
---|
568 | G4cout << "check zheavy2 = " << zheavy2 << G4endl; |
---|
569 | G4cout << "check zsymm = " << zsymm << G4endl; |
---|
570 | G4cout << "check czsymm = " << cz_symm << G4endl; |
---|
571 | G4cout << "check epot0_mode1_saddle = " << epot0_mode1_saddle << G4endl; |
---|
572 | G4cout << "check epot0_mode2_saddle = " << epot0_mode2_saddle << G4endl; |
---|
573 | G4cout << "check epot0_symm_saddle = " << epot0_symm_saddle << G4endl; |
---|
574 | G4cout << "delta_u1 = " << delta_u1 << G4endl; |
---|
575 | G4cout << "delta_u2 = " << delta_u2 << G4endl; |
---|
576 | } |
---|
577 | |
---|
578 | // /* energies including shell effects |
---|
579 | // at the centres of the different shell effects |
---|
580 | // with respect to liquid drop at symmetry: */ |
---|
581 | epot_mode1_saddle = epot0_mode1_saddle + delta_u1; |
---|
582 | epot_mode2_saddle = epot0_mode2_saddle + delta_u2; |
---|
583 | epot_symm_saddle = epot0_symm_saddle; |
---|
584 | if (itest == 1) { |
---|
585 | G4cout << "check epot_mode1_saddle = " << epot_mode1_saddle << G4endl; |
---|
586 | G4cout << "check epot_mode2_saddle = " << epot_mode2_saddle << G4endl; |
---|
587 | G4cout << "check epot_symm_saddle = " << epot_symm_saddle << G4endl; |
---|
588 | } |
---|
589 | |
---|
590 | // /* Minimum of potential with respect to ld potential at symmetry */ |
---|
591 | dueff = min(epot_mode1_saddle,epot_mode2_saddle); |
---|
592 | dueff = min(dueff,epot_symm_saddle); |
---|
593 | dueff = dueff - epot_symm_saddle; |
---|
594 | |
---|
595 | eld = e + dueff + e_zero_point; |
---|
596 | |
---|
597 | if (itest == 1) { |
---|
598 | G4cout << "check dueff = " << dueff << G4endl; |
---|
599 | G4cout << "check e = " << e << G4endl; |
---|
600 | G4cout << "check e_zero_point = " << e_zero_point << G4endl; |
---|
601 | G4cout << "check eld = " << eld << G4endl; |
---|
602 | } |
---|
603 | // Up to here: Ok! Checked CS 10/10/05 |
---|
604 | |
---|
605 | // /* E = energy above lowest effective barrier */ |
---|
606 | // /* Eld = energy above liquid-drop barrier */ |
---|
607 | |
---|
608 | // /* Due to this treatment the energy E on input means the excitation */ |
---|
609 | // /* energy above the lowest saddle. */ |
---|
610 | |
---|
611 | // /* These energies are not used */ |
---|
612 | eheavy1 = e * aheavy1 / a; |
---|
613 | eheavy2 = e * aheavy2 / a; |
---|
614 | elight1 = e * alight1 / a; |
---|
615 | elight2 = e * alight2 / a; |
---|
616 | |
---|
617 | epsilon0_1_saddle = eld - e_zero_point - epot0_mode1_saddle; |
---|
618 | // /* excitation energy at saddle mode 1 without shell effect */ |
---|
619 | epsilon0_2_saddle = eld - e_zero_point - epot0_mode2_saddle; |
---|
620 | // /* excitation energy at saddle mode 2 without shell effect */ |
---|
621 | |
---|
622 | epsilon_1_saddle = eld - e_zero_point - epot_mode1_saddle; |
---|
623 | // /* excitation energy at saddle mode 1 with shell effect */ |
---|
624 | epsilon_2_saddle = eld - e_zero_point - epot_mode2_saddle; |
---|
625 | // /* excitation energy at saddle mode 2 with shell effect */ |
---|
626 | epsilon_symm_saddle = eld - e_zero_point - epot_symm_saddle; |
---|
627 | |
---|
628 | // /* global parameters */ |
---|
629 | eexc1_saddle = epsilon_1_saddle; |
---|
630 | eexc2_saddle = epsilon_2_saddle; |
---|
631 | eexc_max = max(eexc1_saddle,eexc2_saddle); |
---|
632 | eexc_max = max(eexc_max,eld); |
---|
633 | |
---|
634 | // /* EEXC_MAX is energy above the lowest saddle */ |
---|
635 | |
---|
636 | |
---|
637 | epsilon0_1_scission = eld + e_saddle_scission - epot0_mode1_saddle; |
---|
638 | // /* excitation energy without shell effect */ |
---|
639 | epsilon0_2_scission = eld + e_saddle_scission - epot0_mode2_saddle; |
---|
640 | // /* excitation energy without shell effect */ |
---|
641 | |
---|
642 | epsilon_1_scission = eld + e_saddle_scission - epot_mode1_saddle; |
---|
643 | // /* excitation energy at scission */ |
---|
644 | epsilon_2_scission = eld+ e_saddle_scission - epot_mode2_saddle; |
---|
645 | // /* excitation energy at scission */ |
---|
646 | epsilon_symm_scission = eld + e_saddle_scission - epot_symm_saddle; |
---|
647 | // /* excitation energy of symmetric fragment at scission */ |
---|
648 | |
---|
649 | // /* Calculate widhts at the saddle: */ |
---|
650 | |
---|
651 | e_eff1_saddle = epsilon0_1_saddle - delta_u1 * (std::exp((-epsilon_1_saddle*gamma))); |
---|
652 | |
---|
653 | if (e_eff1_saddle > 0.0) { |
---|
654 | wzasymm1_saddle = std::sqrt( (0.5 * |
---|
655 | (std::sqrt(1.0/a_levdens*e_eff1_saddle)) / |
---|
656 | (cz_asymm1_shell * std::exp((-epsilon_1_saddle*gamma)) + cz_symm) ) ); |
---|
657 | } |
---|
658 | else { |
---|
659 | wzasymm1_saddle = 1.0; |
---|
660 | } |
---|
661 | |
---|
662 | e_eff2_saddle = epsilon0_2_saddle - delta_u2 * (std::exp((-epsilon_2_saddle*gamma))); |
---|
663 | if (e_eff2_saddle > 0.0) { |
---|
664 | wzasymm2_saddle = std::sqrt( (0.5 * |
---|
665 | (std::sqrt(1.0/a_levdens*e_eff2_saddle)) / |
---|
666 | (cz_asymm2_shell * std::exp((-epsilon_2_saddle*gamma)) + cz_symm) ) ); |
---|
667 | } |
---|
668 | else { |
---|
669 | wzasymm2_saddle = 1.0; |
---|
670 | } |
---|
671 | |
---|
672 | if (eld > e_zero_point) { |
---|
673 | if ( (eld + epsilon_symm_saddle) < 0.0) { |
---|
674 | G4cout << "<e> eld + epsilon_symm_saddle < 0" << G4endl; |
---|
675 | } |
---|
676 | wzsymm_saddle = std::sqrt( (0.5 * |
---|
677 | (std::sqrt(1.0/a_levdens*(eld+epsilon_symm_saddle))) / cz_symm ) ); |
---|
678 | } else { |
---|
679 | wzsymm_saddle = 1.0; |
---|
680 | } |
---|
681 | |
---|
682 | if (itest == 1) { |
---|
683 | G4cout << "wz1(saddle) = " << wzasymm1_saddle << G4endl; |
---|
684 | G4cout << "wz2(saddle) = " << wzasymm2_saddle << G4endl; |
---|
685 | G4cout << "wzsymm(saddle) = " << wzsymm_saddle << G4endl; |
---|
686 | } |
---|
687 | |
---|
688 | // /* Calculate widhts at the scission point: */ |
---|
689 | // /* fits of ref. Beizin 1991 (Plots brought to GSI by Sergei Zhdanov) */ |
---|
690 | |
---|
691 | wzsymm_scission = wzsymm_saddle; |
---|
692 | |
---|
693 | if (e_saddle_scission == 0.0) { |
---|
694 | |
---|
695 | wzasymm1_scission = wzasymm1_saddle; |
---|
696 | wzasymm2_scission = wzasymm2_saddle; |
---|
697 | |
---|
698 | } |
---|
699 | else { |
---|
700 | |
---|
701 | if (nheavy1_eff > 75.0) { |
---|
702 | wzasymm1_scission = (std::sqrt(21.0)) * z/a; |
---|
703 | wzasymm2_scission = (std::sqrt (max( (70.0-28.0)/3.0*(z*z/a-35.0)+28.,0.0 )) ) * z/a; |
---|
704 | } |
---|
705 | else { |
---|
706 | wzasymm1_scission = wzasymm1_saddle; |
---|
707 | wzasymm2_scission = wzasymm2_saddle; |
---|
708 | } |
---|
709 | |
---|
710 | } |
---|
711 | |
---|
712 | wzasymm1_scission = max(wzasymm1_scission,wzasymm1_saddle); |
---|
713 | wzasymm2_scission = max(wzasymm2_scission,wzasymm2_saddle); |
---|
714 | |
---|
715 | wzasymm1 = wzasymm1_scission * fwidth_asymm1; |
---|
716 | wzasymm2 = wzasymm2_scission * fwidth_asymm2; |
---|
717 | wzsymm = wzsymm_scission; |
---|
718 | |
---|
719 | /* if (ITEST == 1) { |
---|
720 | G4cout << "WZ1(scission) = " << WZasymm1_scission << G4endl; |
---|
721 | G4cout << "WZ2(scission) = " << WZasymm2_scission << G4endl; |
---|
722 | G4cout << "WZsymm(scission) = " << WZsymm_scission << G4endl; |
---|
723 | } |
---|
724 | if (ITEST == 1) { |
---|
725 | G4cout << "WZ1(scission) final= " << WZasymm1 << G4endl; |
---|
726 | G4cout << "WZ2(scission) final= " << WZasymm2 << G4endl; |
---|
727 | G4cout << "WZsymm(scission) final= " << WZsymm << G4endl; |
---|
728 | } */ |
---|
729 | |
---|
730 | wasymm = wzsymm * a/z; |
---|
731 | waheavy1 = wzasymm1 * a/z; |
---|
732 | waheavy2 = wzasymm2 * a/z; |
---|
733 | |
---|
734 | wasymm_saddle = wzsymm_saddle * a/z; |
---|
735 | waheavy1_saddle = wzasymm1_saddle * a/z; |
---|
736 | waheavy2_saddle = wzasymm2_saddle * a/z; |
---|
737 | |
---|
738 | if (itest == 1) { |
---|
739 | G4cout << "wasymm = " << wzsymm << G4endl; |
---|
740 | G4cout << "waheavy1 = " << waheavy1 << G4endl; |
---|
741 | G4cout << "waheavy2 = " << waheavy2 << G4endl; |
---|
742 | } |
---|
743 | // Up to here: Ok! Checked CS 11/10/05 |
---|
744 | |
---|
745 | if ( (epsilon0_1_saddle - delta_u1*std::exp((-epsilon_1_saddle*gamma_heavy1))) < 0.0) { |
---|
746 | sasymm1 = -10.0; |
---|
747 | } |
---|
748 | else { |
---|
749 | sasymm1 = 2.0 * std::sqrt( a_levdens * (epsilon0_1_saddle - |
---|
750 | delta_u1*(std::exp((-epsilon_1_saddle*gamma_heavy1))) ) ); |
---|
751 | } |
---|
752 | |
---|
753 | if ( (epsilon0_2_saddle - delta_u2*std::exp((-epsilon_2_saddle*gamma_heavy2))) < 0.0) { |
---|
754 | sasymm2 = -10.0; |
---|
755 | } |
---|
756 | else { |
---|
757 | sasymm2 = 2.0 * std::sqrt( a_levdens * (epsilon0_2_saddle - |
---|
758 | delta_u2*(std::exp((-epsilon_2_saddle*gamma_heavy2))) ) ); |
---|
759 | } |
---|
760 | |
---|
761 | if (epsilon_symm_saddle > 0.0) { |
---|
762 | ssymm = 2.0 * std::sqrt( a_levdens*(epsilon_symm_saddle) ); |
---|
763 | } |
---|
764 | else { |
---|
765 | ssymm = -10.0; |
---|
766 | } |
---|
767 | |
---|
768 | if (ssymm > -10.0) { |
---|
769 | ysymm = 1.0; |
---|
770 | |
---|
771 | if (epsilon0_1_saddle < 0.0) { |
---|
772 | // /* low energy */ |
---|
773 | yasymm1 = std::exp((sasymm1-ssymm)) * wzasymm1_saddle / wzsymm_saddle * 2.0; |
---|
774 | // /* factor of 2 for symmetry classes */ |
---|
775 | } |
---|
776 | else { |
---|
777 | // /* high energy */ |
---|
778 | ssymm_mode1 = 2.0 * std::sqrt( a_levdens*(epsilon0_1_saddle) ); |
---|
779 | yasymm1 = ( std::exp((sasymm1-ssymm)) - std::exp((ssymm_mode1 - ssymm)) ) |
---|
780 | * wzasymm1_saddle / wzsymm_saddle * 2.0; |
---|
781 | } |
---|
782 | |
---|
783 | if (epsilon0_2_saddle < 0.0) { |
---|
784 | // /* low energy */ |
---|
785 | yasymm2 = std::exp((sasymm2-ssymm)) * wzasymm2_saddle / wzsymm_saddle * 2.0; |
---|
786 | // /* factor of 2 for symmetry classes */ |
---|
787 | } |
---|
788 | else { |
---|
789 | // /* high energy */ |
---|
790 | ssymm_mode2 = 2.0 * std::sqrt( a_levdens*(epsilon0_2_saddle) ); |
---|
791 | yasymm2 = ( std::exp((sasymm2-ssymm)) - std::exp((ssymm_mode2 - ssymm)) ) |
---|
792 | * wzasymm2_saddle / wzsymm_saddle * 2.0; |
---|
793 | } |
---|
794 | // /* difference in the exponent in order */ |
---|
795 | // /* to avoid numerical overflow */ |
---|
796 | |
---|
797 | } |
---|
798 | else { |
---|
799 | if ( (sasymm1 > -10.0) && (sasymm2 > -10.0) ) { |
---|
800 | ysymm = 0.0; |
---|
801 | yasymm1 = std::exp(sasymm1) * wzasymm1_saddle * 2.0; |
---|
802 | yasymm2 = std::exp(sasymm2) * wzasymm2_saddle * 2.0; |
---|
803 | } |
---|
804 | } |
---|
805 | |
---|
806 | // /* normalize */ |
---|
807 | ysum = ysymm + yasymm1 + yasymm2; |
---|
808 | if (ysum > 0.0) { |
---|
809 | ysymm = ysymm / ysum; |
---|
810 | yasymm1 = yasymm1 / ysum; |
---|
811 | yasymm2 = yasymm2 / ysum; |
---|
812 | yasymm = yasymm1 + yasymm2; |
---|
813 | } |
---|
814 | else { |
---|
815 | ysymm = 0.0; |
---|
816 | yasymm1 = 0.0; |
---|
817 | yasymm2 = 0.0; |
---|
818 | // /* search minimum threshold and attribute all events to this mode */ |
---|
819 | if ( (epsilon_symm_saddle < epsilon_1_saddle) && (epsilon_symm_saddle < epsilon_2_saddle) ) { |
---|
820 | ysymm = 1.0; |
---|
821 | } |
---|
822 | else { |
---|
823 | if (epsilon_1_saddle < epsilon_2_saddle) { |
---|
824 | yasymm1 = 1.0; |
---|
825 | } |
---|
826 | else { |
---|
827 | yasymm2 = 1.0; |
---|
828 | } |
---|
829 | } |
---|
830 | } |
---|
831 | |
---|
832 | if (itest == 1) { |
---|
833 | G4cout << "ysymm normalized= " << ysymm << G4endl; |
---|
834 | G4cout << "yasymm1 normalized= " << yasymm1 << G4endl; |
---|
835 | G4cout << "yasymm2 normalized= " << yasymm2 << G4endl; |
---|
836 | } |
---|
837 | // Up to here: Ok! Ckecked CS 11/10/05 |
---|
838 | |
---|
839 | // /* even-odd effect */ |
---|
840 | // /* simple parametrization KHS, Nov. 2000. From Rejmund et al. */ |
---|
841 | if ((int)(z) % 2 == 0) { |
---|
842 | r_e_o_exp = -0.017 * (e_saddle_scission + eld) * (e_saddle_scission + eld); |
---|
843 | if ( r_e_o_exp < -307.0) { |
---|
844 | r_e_o_exp = -307.0; |
---|
845 | r_e_o = std::pow(10.0,r_e_o_exp); |
---|
846 | } |
---|
847 | else { |
---|
848 | r_e_o = std::pow(10.0,r_e_o_exp); |
---|
849 | } |
---|
850 | } |
---|
851 | else { |
---|
852 | r_e_o = 0.0; |
---|
853 | } |
---|
854 | |
---|
855 | // $LOOP; /* event loop */ |
---|
856 | // I_COUNT = I_COUNT + 1; |
---|
857 | |
---|
858 | // /* random decision: symmetric or asymmetric */ |
---|
859 | // /* IMODE = 1 means asymmetric fission, mode 1, |
---|
860 | // IMODE = 2 means asymmetric fission, mode 2, |
---|
861 | // IMODE = 3 means symmetric */ |
---|
862 | // RMODE = dble(HAZ(k)); |
---|
863 | // rmode = rnd.rndm(); |
---|
864 | |
---|
865 | // Safety check added to make sure we always select well defined |
---|
866 | // fission mode. |
---|
867 | do { |
---|
868 | rmode = haz(k); |
---|
869 | // Cast for test CS 11/10/05 |
---|
870 | // RMODE = 0.54; |
---|
871 | // rmode = 0.54; |
---|
872 | if (rmode < yasymm1) { |
---|
873 | imode = 1; |
---|
874 | } else if ( (rmode > yasymm1) && (rmode < (yasymm1+yasymm2)) ) { |
---|
875 | imode = 2; |
---|
876 | } else if ( (rmode > yasymm1) && (rmode > (yasymm1+yasymm2)) ) { |
---|
877 | imode = 3; |
---|
878 | } |
---|
879 | } while(imode == 0); |
---|
880 | |
---|
881 | // /* determine parameters of the Z distribution */ |
---|
882 | // force imode (for testing, PK) |
---|
883 | // imode = 3; |
---|
884 | if (imode == 1) { |
---|
885 | z1mean = zheavy1; |
---|
886 | z1width = wzasymm1; |
---|
887 | } |
---|
888 | if (imode == 2) { |
---|
889 | z1mean = zheavy2; |
---|
890 | z1width = wzasymm2; |
---|
891 | } |
---|
892 | if (imode == 3) { |
---|
893 | z1mean = zsymm; |
---|
894 | z1width = wzsymm; |
---|
895 | } |
---|
896 | |
---|
897 | if (itest == 1) { |
---|
898 | G4cout << "nbre aleatoire tire " << rmode << G4endl; |
---|
899 | G4cout << "fission mode " << imode << G4endl; |
---|
900 | G4cout << "z1mean= " << z1mean << G4endl; |
---|
901 | G4cout << "z1width= " << z1width << G4endl; |
---|
902 | } |
---|
903 | |
---|
904 | // /* random decision: Z1 and Z2 at scission: */ |
---|
905 | z1 = 1.0; |
---|
906 | z2 = 1.0; |
---|
907 | while ( (z1<5.0) || (z2<5.0) ) { |
---|
908 | // Z1 = dble(GAUSSHAZ(K,sngl(Z1mean),sngl(Z1width))); |
---|
909 | // z1 = rnd.gaus(z1mean,z1width); |
---|
910 | z1 = gausshaz(k, z1mean, z1width); |
---|
911 | z2 = z - z1; |
---|
912 | } |
---|
913 | if (itest == 1) { |
---|
914 | G4cout << "ff charge sample " << G4endl; |
---|
915 | G4cout << "z1 = " << z1 << G4endl; |
---|
916 | G4cout << "z2 = " << z2 << G4endl; |
---|
917 | } |
---|
918 | |
---|
919 | // CALL EVEN_ODD(Z1,R_E_O,I_HELP); |
---|
920 | // /* Integer proton number with even-odd effect */ |
---|
921 | // Z1 = REAL(I_HELP) |
---|
922 | // /* Z1 = INT(Z1+0.5E0); */ |
---|
923 | z2 = z - z1; |
---|
924 | |
---|
925 | // /* average N of both fragments: */ |
---|
926 | if (imode == 1) { |
---|
927 | n1mean = (z1 + cpol1 * a/n) * n/z; |
---|
928 | } |
---|
929 | if (imode == 2) { |
---|
930 | n1mean = (z1 + cpol2 * a/n) * n/z; |
---|
931 | } |
---|
932 | /* CASE(99) ! only for testing; |
---|
933 | N1UCD = Z1 * N/Z; |
---|
934 | N2UCD = Z2 * N/Z; |
---|
935 | re1 = UMASS(Z1,N1UCD,0.6) +; |
---|
936 | & UMASS(Z2,N2UCD,0.6) +; |
---|
937 | & ECOUL(Z1,N1UCD,0.6,Z2,N2UCD,0.6,d); |
---|
938 | re2 = UMASS(Z1,N1UCD+1.,0.6) +; |
---|
939 | & UMASS(Z2,N2UCD-1.,0.6) +; |
---|
940 | & ECOUL(Z1,N1UCD+1.,0.6,Z2,N2UCD-1.,0.6,d); |
---|
941 | re3 = UMASS(Z1,N1UCD+2.,0.6) +; |
---|
942 | & UMASS(Z2,N2UCD-2.,0.6) +; |
---|
943 | & ECOUL(Z1,N1UCD+2.,0.6,Z2,N2UCD-2.,0.6,d); |
---|
944 | eps2 = (re1-2.0*re2+re3) / 2.0; |
---|
945 | eps1 = re2 - re1 - eps2; |
---|
946 | DN1_POL = - eps1 / (2.0 * eps2); |
---|
947 | N1mean = N1UCD + DN1_POL; */ |
---|
948 | if (imode == 3) { |
---|
949 | n1ucd = z1 * n/z; |
---|
950 | n2ucd = z2 * n/z; |
---|
951 | re1 = umass(z1,n1ucd,0.6) + umass(z2,n2ucd,0.6) + ecoul(z1,n1ucd,0.6,z2,n2ucd,0.6,d); |
---|
952 | re2 = umass(z1,n1ucd+1.,0.6) + umass(z2,n2ucd-1.,0.6) + ecoul(z1,n1ucd+1.,0.6,z2,n2ucd-1.,0.6,d); |
---|
953 | re3 = umass(z1,n1ucd+2.,0.6) + umass(z2,n2ucd-2.,0.6) + ecoul(z1,n1ucd+2.,0.6,z2,n2ucd-2.,0.6,d); |
---|
954 | eps2 = (re1-2.0*re2+re3) / 2.0; |
---|
955 | eps1 = re2 - re1 - eps2; |
---|
956 | dn1_pol = - eps1 / (2.0 * eps2); |
---|
957 | n1mean = n1ucd + dn1_pol; |
---|
958 | } |
---|
959 | // all fission modes features have been checked CS 11/10/05 |
---|
960 | n2mean = n - n1mean; |
---|
961 | z2mean = z - z1mean; |
---|
962 | |
---|
963 | // /* Excitation energies */ |
---|
964 | // /* formulated in energies in close consistency with the fission model */ |
---|
965 | |
---|
966 | // /* E_defo = UMASS(Z*0.5E0,N*0.5E0,0.6E0) - |
---|
967 | // UMASS(Z*0.5E0,N*0.5E0,0); */ |
---|
968 | // /* calculates the deformation energy of the liquid drop for |
---|
969 | // deformation beta = 0.6 which is most probable at scission */ |
---|
970 | |
---|
971 | // /* N1R and N2R provisionaly taken without fluctuations in |
---|
972 | // polarisation: */ |
---|
973 | n1r = n1mean; |
---|
974 | n2r = n2mean; |
---|
975 | a1r = n1r + z1; |
---|
976 | a2r = n2r + z2; |
---|
977 | |
---|
978 | if (imode == 1) { /* N = 82 */; |
---|
979 | //! /* Eexc at scission */ |
---|
980 | e_scission = max(epsilon_1_scission,1.0); |
---|
981 | if (n1mean > (n * 0.5) ) { |
---|
982 | //! /* 1. fragment is spherical */ |
---|
983 | beta1 = 0.0; |
---|
984 | beta2 = 0.6; |
---|
985 | e1exc = epsilon_1_scission * a1r / a; |
---|
986 | e_defo = umass(z2,n2r,beta2) - umass(z2,n2r,0.0); |
---|
987 | e2exc = epsilon_1_scission * a2r / a + e_defo; |
---|
988 | } |
---|
989 | else { |
---|
990 | //! /* 2. fragment is spherical */ |
---|
991 | beta1 = 0.6; |
---|
992 | beta2 = 0.0; |
---|
993 | e_defo = umass(z1,n1r,beta1) - umass(z1,n1r,0.0); |
---|
994 | e1exc = epsilon_1_scission * a1r / a + e_defo; |
---|
995 | e2exc = epsilon_1_scission * a2r / a; |
---|
996 | } |
---|
997 | } |
---|
998 | |
---|
999 | if (imode == 2) { |
---|
1000 | //! /* N appr. 86 */ |
---|
1001 | e_scission = max(epsilon_2_scission,1.0); |
---|
1002 | if (n1mean > (n * 0.5) ) { |
---|
1003 | //! /* 2. fragment is spherical */ |
---|
1004 | beta1 = (n1r - nheavy2) * 0.034 + 0.3; |
---|
1005 | e_defo = umass(z1,n1r,beta1) - umass(z1,n1r,0.0); |
---|
1006 | e1exc = epsilon_2_scission * a1r / a + e_defo; |
---|
1007 | beta2 = 0.6 - beta1; |
---|
1008 | e_defo = umass(z2,n2r,beta2) - umass(z2,n2r,0.0); |
---|
1009 | e2exc = epsilon_2_scission * a2r / a + e_defo; |
---|
1010 | } |
---|
1011 | else { |
---|
1012 | //! /* 1. fragment is spherical */ |
---|
1013 | beta2 = (n2r - nheavy2) * 0.034 + 0.3; |
---|
1014 | e_defo = umass(z2,n2r,beta2) - umass(z2,n2r,0.0); |
---|
1015 | e2exc = epsilon_2_scission * a2r / a + e_defo; |
---|
1016 | beta1 = 0.6 - beta2; |
---|
1017 | e_defo = umass(z1,n1r,beta1) - umass(z1,n1r,0.0); |
---|
1018 | e1exc = epsilon_2_scission * a1r / a + e_defo; |
---|
1019 | } |
---|
1020 | } |
---|
1021 | |
---|
1022 | if (imode == 3) { |
---|
1023 | // ! /* Symmetric fission channel */ |
---|
1024 | |
---|
1025 | // /* the fit function for beta is the deformation for |
---|
1026 | // optimum energy at the scission point, d = 2 */ |
---|
1027 | // /* beta : deformation of symmetric fragments */ |
---|
1028 | // /* beta1 : deformation of first fragment */ |
---|
1029 | // /* beta2 : deformation of second fragment */ |
---|
1030 | beta = 0.177963 + 0.0153241 * zsymm - 0.000162037 * zsymm*zsymm; |
---|
1031 | beta1 = 0.177963 + 0.0153241 * z1 - 0.000162037 * z1*z1; |
---|
1032 | // beta1 = 0.6 |
---|
1033 | e_defo1 = umass(z1,n1r,beta1) - umass(z1,n1r,0.0); |
---|
1034 | beta2 = 0.177963 + 0.0153241 * z2 - 0.000162037 * z2*z2; |
---|
1035 | // beta2 = 0.6 |
---|
1036 | e_defo2 = umass(z2,n2r,beta2) - umass(z2,n2r,0.0); |
---|
1037 | e_asym = umass(z1 , n1r, beta1) + umass(z2, n2r ,beta2) |
---|
1038 | + ecoul(z1,n1r,beta1,z2,n2r,beta2,2.0) |
---|
1039 | - 2.0 * umass(zsymm,nsymm,beta) |
---|
1040 | - ecoul(zsymm,nsymm,beta,zsymm,nsymm,beta,2.0); |
---|
1041 | // E_asym = CZ_symm * (Z1 - Zsymm)**2 |
---|
1042 | e_scission = max((epsilon_symm_scission - e_asym),1.0); |
---|
1043 | // /* $LIST(Z1,N1R,Z2,N2R,E_asym,E_scission); */ |
---|
1044 | e1exc = e_scission * a1r / a + e_defo1; |
---|
1045 | e2exc = e_scission * a2r / a + e_defo2; |
---|
1046 | } |
---|
1047 | // Energies checked for all the modes CS 11/10/05 |
---|
1048 | |
---|
1049 | // /* random decision: N1R and N2R at scission, before evaporation: */ |
---|
1050 | // /* CN = UMASS(Zsymm , Nsymm + 1.E0,0) + |
---|
1051 | // UMASS(Zsymm, Nsymm - 1.E0,0) |
---|
1052 | // + 1.44E0 * (Zsymm)**2 / |
---|
1053 | // (r_null**2 * ((Asymm+1)**1/3 + (Asymm-1)**1/3)**2 ) |
---|
1054 | // - 2.E0 * UMASS(Zsymm,Nsymm,0) |
---|
1055 | // - 1.44E0 * (Zsymm)**2 / (r_null * 2.E0 * (Asymm)**1/3)**2; */ |
---|
1056 | |
---|
1057 | |
---|
1058 | // /* N1width = std::sqrt(0.5E0 * std::sqrt(1.E0/A_levdens*(Eld+E_saddle_scission)) / CN); */ |
---|
1059 | // /* 8. 9. 1998: KHS (see also consideration in the first comment block) |
---|
1060 | // sigma_N(Z=const) = A/Z * sigma_Z(A=const) |
---|
1061 | // sigma_Z(A=const) = 0.4 to 0.5 (from Lang paper Nucl Phys. A345 (1980) 34) |
---|
1062 | // sigma_N(Z=const) = 0.45 * A/Z (= 1.16 for 238U) |
---|
1063 | // therefore: SIGZMIN = 1.16 */ |
---|
1064 | |
---|
1065 | if ( (imode == 1) || (imode == 2) ) { |
---|
1066 | cn=(umass(z1,n1mean+1.,beta1) + umass(z1,n1mean-1.,beta1) |
---|
1067 | + umass(z2,n2mean+1.,beta2) + umass(z2,n2mean-1.,beta2) |
---|
1068 | + ecoul(z1,n1mean+1.,beta1,z2,n2mean-1.,beta2,2.0) |
---|
1069 | + ecoul(z1,n1mean-1.,beta1,z2,n2mean+1.,beta2,2.0) |
---|
1070 | - 2.0 * ecoul(z1,n1mean,beta1,z2,n2mean,beta2,2.0) |
---|
1071 | - 2.0 * umass(z1, n1mean, beta1) |
---|
1072 | - 2.0 * umass(z2, n2mean, beta2) ) * 0.5; |
---|
1073 | // /* Coulomb energy neglected for the moment! */ |
---|
1074 | // IF (E_scission.lt.0.) Then |
---|
1075 | // write(6,*)'<E> E_scission < 0, MODE 1,2' |
---|
1076 | // ENDIF |
---|
1077 | // IF (CN.lt.0.) Then |
---|
1078 | // write(6,*)'CN < 0, MODE 1,2' |
---|
1079 | // ENDIF |
---|
1080 | n1width=std::sqrt( (0.5 * (std::sqrt(1.0/a_levdens*(e_scission)))/cn) ); |
---|
1081 | n1width=max(n1width, sigzmin); |
---|
1082 | |
---|
1083 | // /* random decision: N1R and N2R at scission, before evaporation: */ |
---|
1084 | n1r = 1.0; |
---|
1085 | n2r = 1.0; |
---|
1086 | while ( (n1r<5.0) || (n2r<5.0) ) { |
---|
1087 | // n1r = dble(gausshaz(k,sngl(n1mean),sngl(n1width))); |
---|
1088 | // n1r = rnd.gaus(n1mean,n1width); |
---|
1089 | n1r = gausshaz(k, n1mean, n1width); |
---|
1090 | n2r = n - n1r; |
---|
1091 | } |
---|
1092 | // N1R = GAUSSHAZ(K,N1mean,N1width) |
---|
1093 | if (itest == 1) { |
---|
1094 | G4cout << "after neutron sample " << n1r << G4endl; |
---|
1095 | } |
---|
1096 | n1r = (float)( (int)((n1r+0.5)) ); |
---|
1097 | n2r = n - n1r; |
---|
1098 | |
---|
1099 | even_odd(z1,r_e_o,i_help); |
---|
1100 | // /* proton number with even-odd effect */ |
---|
1101 | z1 = (float)(i_help); |
---|
1102 | z2 = z - z1; |
---|
1103 | |
---|
1104 | a1r = z1 + n1r; |
---|
1105 | a2r = z2 + n2r; |
---|
1106 | } |
---|
1107 | |
---|
1108 | if (imode == 3) { |
---|
1109 | //! /* When(3) */ |
---|
1110 | if (nzpol > 0.0) { |
---|
1111 | // /* We treat a simultaneous split in Z and N to determine polarisation */ |
---|
1112 | cz = ( umass(z1-1., n1mean+1.,beta1) |
---|
1113 | + umass(z2+1., n2mean-1.,beta1) |
---|
1114 | + umass(z1+1., n1mean-1.,beta2) |
---|
1115 | + umass(z2 - 1., n2mean + 1.,beta2) |
---|
1116 | + ecoul(z1-1.,n1mean+1.,beta1,z2+1.,n2mean-1.,beta2,2.0) |
---|
1117 | + ecoul(z1+1.,n1mean-1.,beta1,z2-1.,n2mean+1.,beta2,2.0) |
---|
1118 | - 2.0 * ecoul(z1,n1mean,beta1,z2,n2mean,beta2,2.0) |
---|
1119 | - 2.0 * umass(z1, n1mean,beta1) |
---|
1120 | - 2.0 * umass(z2, n2mean,beta2) ) * 0.5; |
---|
1121 | // IF (E_scission.lt.0.) Then |
---|
1122 | // write(6,*) '<E> E_scission < 0, MODE 1,2' |
---|
1123 | // ENDIF |
---|
1124 | // IF (CZ.lt.0.) Then |
---|
1125 | // write(6,*) 'CZ < 0, MODE 1,2' |
---|
1126 | // ENDIF |
---|
1127 | za1width=std::sqrt( (0.5 * std::sqrt(1.0/a_levdens*(e_scission)) / cz) ); |
---|
1128 | za1width=std::sqrt( (max((za1width*za1width-(1.0/12.0)),0.1)) ); |
---|
1129 | // /* Check the value of 0.1 ! */ |
---|
1130 | // /* Shephard correction */ |
---|
1131 | a1r = z1 + n1mean; |
---|
1132 | a1r = (float)((int)((a1r+0.5))); |
---|
1133 | a2r = a - a1r; |
---|
1134 | // /* A1R and A2R are integer numbers now */ |
---|
1135 | // /* $LIST(A1R,A2R,ZA1WIDTH); */ |
---|
1136 | |
---|
1137 | n1ucd = n/a * a1r; |
---|
1138 | n2ucd = n/a * a2r; |
---|
1139 | z1ucd = z/a * a1r; |
---|
1140 | z2ucd = z/a * a2r; |
---|
1141 | |
---|
1142 | re1 = umass(z1ucd-1.,n1ucd+1.,beta1) + umass(z2ucd+1.,n2ucd-1.,beta2) |
---|
1143 | + ecoul(z1ucd-1.,n1ucd+1.,beta1,z2ucd+1.,n2ucd-1.,beta2,d); |
---|
1144 | re2 = umass(z1ucd,n1ucd,beta1) + umass(z2ucd,n2ucd,beta2) |
---|
1145 | + ecoul(z1ucd,n1ucd,beta1,z2ucd,n2ucd,beta2,d); |
---|
1146 | re3 = umass(z1ucd+1.,n1ucd-1.,beta1) + umass(z2ucd-1.,n2ucd+1.,beta2) + |
---|
1147 | + ecoul(z1ucd+1.,n1ucd-1.,beta1,z2ucd-1.,n2ucd+1.,beta2,d); |
---|
1148 | |
---|
1149 | eps2 = (re1-2.0*re2+re3) / 2.0; |
---|
1150 | eps1 = (re3 - re1)/2.0; |
---|
1151 | dn1_pol = - eps1 / (2.0 * eps2); |
---|
1152 | z1 = z1ucd + dn1_pol; |
---|
1153 | if (itest == 1) { |
---|
1154 | G4cout << "before proton sample " << z1 << G4endl; |
---|
1155 | } |
---|
1156 | // Z1 = dble(GAUSSHAZ(k,sngl(Z1),sngl(ZA1width))); |
---|
1157 | // z1 = rnd.gaus(z1,za1width); |
---|
1158 | z1 = gausshaz(k, z1, za1width); |
---|
1159 | if (itest == 1) { |
---|
1160 | G4cout << "after proton sample " << z1 << G4endl; |
---|
1161 | } |
---|
1162 | even_odd(z1,r_e_o,i_help); |
---|
1163 | // /* proton number with even-odd effect */ |
---|
1164 | z1 = (float)(i_help); |
---|
1165 | z2 = (float)((int)( (z - z1 + 0.5)) ); |
---|
1166 | |
---|
1167 | n1r = a1r - z1; |
---|
1168 | n2r = n - n1r; |
---|
1169 | } |
---|
1170 | else { |
---|
1171 | // /* First division of protons, then adjustment of neutrons */ |
---|
1172 | cn = ( umass(z1, n1mean+1.,beta1) + umass(z1, n1mean-1., beta1) |
---|
1173 | + umass(z2, n2mean+1.,beta2) + umass(z2, n2mean-1., beta2) |
---|
1174 | + ecoul(z1,n1mean+1.,beta1,z2,n2mean-1.,beta2,2.0) |
---|
1175 | + ecoul(z1,n1mean-1.,beta1,z2,n2mean+1.,beta2,2.0) |
---|
1176 | - 2.0 * ecoul(z1,n1mean,beta1,z2,n2mean,beta2,2.0) |
---|
1177 | - 2.0 * umass(z1, n1mean, 0.6) |
---|
1178 | - 2.0 * umass(z2, n2mean, 0.6) ) * 0.5; |
---|
1179 | // /* Coulomb energy neglected for the moment! */ |
---|
1180 | // IF (E_scission.lt.0.) Then |
---|
1181 | // write(6,*) '<E> E_scission < 0, MODE 1,2' |
---|
1182 | // Endif |
---|
1183 | // IF (CN.lt.0.) Then |
---|
1184 | // write(6,*) 'CN < 0, MODE 1,2' |
---|
1185 | // Endif |
---|
1186 | n1width=std::sqrt( (0.5 * std::sqrt(1.0/a_levdens*(e_scission)) / cn) ); |
---|
1187 | n1width=max(n1width, sigzmin); |
---|
1188 | |
---|
1189 | // /* random decision: N1R and N2R at scission, before evaporation: */ |
---|
1190 | // N1R = dble(GAUSSHAZ(k,sngl(N1mean),sngl(N1width))); |
---|
1191 | // n1r = rnd.gaus(n1mean,n1width); |
---|
1192 | n1r = gausshaz(k, n1mean, n1width); |
---|
1193 | n1r = (float)( (int)((n1r+0.5)) ); |
---|
1194 | n2r = n - n1r; |
---|
1195 | |
---|
1196 | even_odd(z1,r_e_o,i_help); |
---|
1197 | // /* Integer proton number with even-odd effect */ |
---|
1198 | z1 = (float)(i_help); |
---|
1199 | z2 = z - z1; |
---|
1200 | |
---|
1201 | a1r = z1 + n1r; |
---|
1202 | a2r = z2 + n2r; |
---|
1203 | |
---|
1204 | } |
---|
1205 | } |
---|
1206 | |
---|
1207 | if (itest == 1) { |
---|
1208 | G4cout << "remid imode = " << imode << G4endl; |
---|
1209 | G4cout << "n1width = " << n1width << G4endl; |
---|
1210 | G4cout << "n1r = " << n1r << G4endl; |
---|
1211 | G4cout << "a1r = " << a1r << G4endl; |
---|
1212 | G4cout << "n2r = " << n2r << G4endl; |
---|
1213 | G4cout << "a2r = " << a2r << G4endl; |
---|
1214 | } |
---|
1215 | // Up to here: checked CS 11/10/05 |
---|
1216 | |
---|
1217 | // /* Extracted from Lang et al. Nucl. Phys. A 345 (1980) 34 */ |
---|
1218 | e1exc_sigma = 5.5; |
---|
1219 | e2exc_sigma = 5.5; |
---|
1220 | |
---|
1221 | neufcentquatrevingtsept: |
---|
1222 | // E1final = dble(Gausshaz(k,sngl(E1exc),sngl(E1exc_sigma))); |
---|
1223 | // E2final = dble(Gausshaz(k,sngl(E2exc),sngl(E2exc_sigma))); |
---|
1224 | // e1final = rnd.gaus(e1exc,e1exc_sigma); |
---|
1225 | // e2final = rnd.gaus(e2exc,e2exc_sigma); |
---|
1226 | e1final = gausshaz(k, e1exc, e1exc_sigma); |
---|
1227 | e2final = gausshaz(k, e2exc, e2exc_sigma); |
---|
1228 | if ( (e1final < 0.0) || (e2final < 0.0) ) goto neufcentquatrevingtsept; |
---|
1229 | if (itest == 1) { |
---|
1230 | G4cout << "sampled exc 1 " << e1final << G4endl; |
---|
1231 | G4cout << "sampled exc 2 " << e2final << G4endl; |
---|
1232 | } |
---|
1233 | |
---|
1234 | // /* OUTPUT QUANTITIES OF THE EVENT GENERATOR: */ |
---|
1235 | |
---|
1236 | // /* Quantities before neutron evaporation */ |
---|
1237 | |
---|
1238 | // /* Neutron number of prefragments: N1R and N2R */ |
---|
1239 | // /* Atomic number of fragments: Z1 and Z2 */ |
---|
1240 | // /* Kinetic energy of fragments: EkinR1, EkinR2 *7 |
---|
1241 | |
---|
1242 | // /* Quantities after neutron evaporation: */ |
---|
1243 | |
---|
1244 | // /* Neutron number of fragments: N1 and N2 */ |
---|
1245 | // /* Mass number of fragments: A1 and A2 */ |
---|
1246 | // /* Atomic number of fragments: Z1 and Z2 */ |
---|
1247 | // /* Number of evaporated neutrons: N1R-N1 and N2R-N2 */ |
---|
1248 | // /* Kinetic energy of fragments: EkinR1*A1/A1R and |
---|
1249 | // EkinR2*A2/A2R */ |
---|
1250 | |
---|
1251 | n1 = n1r; |
---|
1252 | n2 = n2r; |
---|
1253 | a1 = n1 + z1; |
---|
1254 | a2 = n2 + z2; |
---|
1255 | e1 = e1final; |
---|
1256 | e2 = e2final; |
---|
1257 | |
---|
1258 | // /* Pre-neutron-emission total kinetic energy: */ |
---|
1259 | tker = (z1 * z2 * 1.44) / |
---|
1260 | ( r0 * std::pow(a1,0.33333) * (1.0 + 2.0/3.0 * beta1) + |
---|
1261 | r0 * std::pow(a2,0.33333) * (1.0 + 2.0/3.0 * beta2) + 2.0 ); |
---|
1262 | // /* Pre-neutron-emission kinetic energy of 1. fragment: */ |
---|
1263 | ekinr1 = tker * a2 / a; |
---|
1264 | // /* Pre-neutron-emission kinetic energy of 2. fragment: */ |
---|
1265 | ekinr2 = tker * a1 / a; |
---|
1266 | |
---|
1267 | v1 = std::sqrt( (ekinr1/a1) ) * 1.3887; |
---|
1268 | v2 = std::sqrt( (ekinr2/a2) ) * 1.3887; |
---|
1269 | |
---|
1270 | if (itest == 1) { |
---|
1271 | G4cout << "ekinr1 " << ekinr1 << G4endl; |
---|
1272 | G4cout << "ekinr2 " << ekinr2 << G4endl; |
---|
1273 | } |
---|
1274 | |
---|
1275 | milledeux: |
---|
1276 | //************************** |
---|
1277 | //*** only symmetric fission |
---|
1278 | //************************** |
---|
1279 | // Symmetric fission: Ok! Checked CS 10/10/05 |
---|
1280 | if ( (icz == -1) || (a1 < 0.0) || (a2 < 0.0) ) { |
---|
1281 | // IF (z.eq.92) THEN |
---|
1282 | // write(6,*)'symmetric fission' |
---|
1283 | // write(6,*)'Z,A,E,A1,A2,icz,Atot',Z,A,E,A1,A2,icz,Atot |
---|
1284 | // END IF |
---|
1285 | |
---|
1286 | if (itest == 1) { |
---|
1287 | G4cout << "milledeux: liquid-drop option " << G4endl; |
---|
1288 | } |
---|
1289 | |
---|
1290 | n = a-z; |
---|
1291 | // proton number in symmetric fission (centre) * |
---|
1292 | zsymm = z / 2.0; |
---|
1293 | nsymm = n / 2.0; |
---|
1294 | asymm = nsymm + zsymm; |
---|
1295 | |
---|
1296 | a_levdens = a / xlevdens; |
---|
1297 | |
---|
1298 | masscurv = 2.0; |
---|
1299 | cz_symm = 8.0 / std::pow(z,2) * masscurv; |
---|
1300 | |
---|
1301 | wzsymm = std::sqrt( (0.5 * std::sqrt(1.0/a_levdens*e) / cz_symm) ) ; |
---|
1302 | |
---|
1303 | if (itest == 1) { |
---|
1304 | G4cout << " symmetric high energy fission " << G4endl; |
---|
1305 | G4cout << "wzsymm " << wzsymm << G4endl; |
---|
1306 | } |
---|
1307 | |
---|
1308 | z1mean = zsymm; |
---|
1309 | z1width = wzsymm; |
---|
1310 | |
---|
1311 | // random decision: Z1 and Z2 at scission: */ |
---|
1312 | z1 = 1.0; |
---|
1313 | z2 = 1.0; |
---|
1314 | while ( (z1 < 5.0) || (z2 < 5.0) ) { |
---|
1315 | // z1 = dble(gausshaz(kkk,sngl(z1mean),sngl(z1width))); |
---|
1316 | // z1 = rnd.gaus(z1mean,z1width); |
---|
1317 | z1 = gausshaz(kkk, z1mean, z1width); |
---|
1318 | z2 = z - z1; |
---|
1319 | } |
---|
1320 | |
---|
1321 | if (itest == 1) { |
---|
1322 | G4cout << " z1 " << z1 << G4endl; |
---|
1323 | G4cout << " z2 " << z2 << G4endl; |
---|
1324 | } |
---|
1325 | if (itest == 1) { |
---|
1326 | G4cout << " zsymm " << zsymm << G4endl; |
---|
1327 | G4cout << " nsymm " << nsymm << G4endl; |
---|
1328 | G4cout << " asymm " << asymm << G4endl; |
---|
1329 | } |
---|
1330 | // CN = UMASS(Zsymm , Nsymm + 1.E0) + UMASS(Zsymm, Nsymm - 1.E0) |
---|
1331 | // # + 1.44E0 * (Zsymm)**2 / |
---|
1332 | // # (r_null**2 * ((Asymm+1)**(1./3.) + |
---|
1333 | // # (Asymm-1)**(1./3.))**2 ) |
---|
1334 | // # - 2.E0 * UMASS(Zsymm,Nsymm) |
---|
1335 | // # - 1.44E0 * (Zsymm)**2 / |
---|
1336 | // # (r_null * 2.E0 * (Asymm)**(1./3.))**2 |
---|
1337 | |
---|
1338 | n1ucd = z1 * n/z; |
---|
1339 | n2ucd = z2 * n/z; |
---|
1340 | re1 = umass(z1,n1ucd,0.6) + umass(z2,n2ucd,0.6) + |
---|
1341 | ecoul(z1,n1ucd,0.6,z2,n2ucd,0.6,2.0); |
---|
1342 | re2 = umass(z1,n1ucd+1.,0.6) + umass(z2,n2ucd-1.,0.6) + |
---|
1343 | ecoul(z1,n1ucd+1.,0.6,z2,n2ucd-1.,0.6,2.0); |
---|
1344 | re3 = umass(z1,n1ucd+2.,0.6) + umass(z2,n2ucd-2.,0.6) + |
---|
1345 | ecoul(z1,n1ucd+2.,0.6,z2,n2ucd-2.,0.6,2.0); |
---|
1346 | reps2 = (re1-2.0*re2+re3)/2.0; |
---|
1347 | reps1 = re2 - re1 -reps2; |
---|
1348 | rn1_pol = -reps1/(2.0*reps2); |
---|
1349 | n1mean = n1ucd + rn1_pol; |
---|
1350 | n2mean = n - n1mean; |
---|
1351 | |
---|
1352 | if (itest == 1) { |
---|
1353 | G4cout << " n1mean " << n1mean << G4endl; |
---|
1354 | G4cout << " n2mean " << n2mean << G4endl; |
---|
1355 | } |
---|
1356 | |
---|
1357 | cn = (umass(z1,n1mean+1.,0.0) + umass(z1,n1mean-1.,0.0) + |
---|
1358 | + umass(z2,n2mean+1.,0.0) + umass(z2,n2mean-1.,0.0) |
---|
1359 | - 2.0 * umass(z1,n1mean,0.0) + |
---|
1360 | - 2.0 * umass(z2,n2mean,0.0) ) * 0.5; |
---|
1361 | // This is an approximation! Coulomb energy is neglected. |
---|
1362 | |
---|
1363 | n1width = std::sqrt( (0.5 * std::sqrt(1.0/a_levdens*e) / cn) ); |
---|
1364 | |
---|
1365 | if (itest == 1) { |
---|
1366 | G4cout << " cn " << cn << G4endl; |
---|
1367 | G4cout << " n1width " << n1width << G4endl; |
---|
1368 | } |
---|
1369 | |
---|
1370 | // random decision: N1R and N2R at scission, before evaporation: */ |
---|
1371 | // N1R = dfloat(NINT(GAUSSHAZ(KKK,sngl(N1mean),sngl(N1width)))); |
---|
1372 | // n1r = (float)( (int)(rnd.gaus(n1mean,n1width)) ); |
---|
1373 | n1r = (float)( (int)(gausshaz(k, n1mean,n1width)) ); |
---|
1374 | n2r = n - n1r; |
---|
1375 | // Mass of first and second fragment */ |
---|
1376 | a1 = z1 + n1r; |
---|
1377 | a2 = z2 + n2r; |
---|
1378 | |
---|
1379 | e1 = e*a1/(a1+a2); |
---|
1380 | e2 = e - e*a1/(a1+a2); |
---|
1381 | if (itest == 1) { |
---|
1382 | G4cout << " n1r " << n1r << G4endl; |
---|
1383 | G4cout << " n2r " << n2r << G4endl; |
---|
1384 | } |
---|
1385 | |
---|
1386 | } |
---|
1387 | |
---|
1388 | if (itest == 1) { |
---|
1389 | G4cout << " a1 " << a1 << G4endl; |
---|
1390 | G4cout << " z1 " << z1 << G4endl; |
---|
1391 | G4cout << " a2 " << a2 << G4endl; |
---|
1392 | G4cout << " z2 " << z2 << G4endl; |
---|
1393 | G4cout << " e1 " << e1 << G4endl; |
---|
1394 | G4cout << " e2 " << e << G4endl; |
---|
1395 | } |
---|
1396 | |
---|
1397 | // /* Pre-neutron-emission total kinetic energy: */ |
---|
1398 | tker = (z1 * z2 * 1.44) / |
---|
1399 | ( r0 * std::pow(a1,0.33333) * (1.0 + 2.0/3.0 * beta1) + |
---|
1400 | r0 * std::pow(a2,0.33333) * (1.0 + 2.0/3.0 * beta2) + 2.0 ); |
---|
1401 | // /* Pre-neutron-emission kinetic energy of 1. fragment: */ |
---|
1402 | ekin1 = tker * a2 / a; |
---|
1403 | // /* Pre-neutron-emission kinetic energy of 2. fragment: */ |
---|
1404 | ekin2 = tker * a1 / a; |
---|
1405 | |
---|
1406 | v1 = std::sqrt( (ekin1/a1) ) * 1.3887; |
---|
1407 | v2 = std::sqrt( (ekin2/a2) ) * 1.3887; |
---|
1408 | |
---|
1409 | if (itest == 1) { |
---|
1410 | G4cout << " kinetic energies " << G4endl; |
---|
1411 | G4cout << " ekin1 " << ekin1 << G4endl; |
---|
1412 | G4cout << " ekin2 " << ekin2 << G4endl; |
---|
1413 | } |
---|
1414 | } |
---|
1415 | |
---|
1416 | void G4AblaFissionSimfis18::standardRandom(G4double *rndm, G4long*) |
---|
1417 | { |
---|
1418 | // Use Geant4 G4UniformRand |
---|
1419 | (*rndm) = randomGenerator->getRandom(); |
---|
1420 | } |
---|
1421 | |
---|
1422 | G4double G4AblaFissionSimfis18::haz(G4int k) |
---|
1423 | { |
---|
1424 | const G4int pSize = 110; |
---|
1425 | static G4double p[pSize]; |
---|
1426 | static G4long ix = 0, i = 0; |
---|
1427 | static G4double x = 0.0, y = 0.0, a = 0.0, haz = 0.0; |
---|
1428 | // k =< -1 on initialise |
---|
1429 | // k = -1 c'est reproductible |
---|
1430 | // k < -1 || k > -1 ce n'est pas reproductible |
---|
1431 | |
---|
1432 | // Zero is invalid random seed. Set proper value from our random seed collection: |
---|
1433 | if(ix == 0) { |
---|
1434 | ix = hazard->ial; |
---|
1435 | } |
---|
1436 | |
---|
1437 | if (k <= -1) { //then |
---|
1438 | if(k == -1) { //then |
---|
1439 | ix = 0; |
---|
1440 | } |
---|
1441 | else { |
---|
1442 | x = 0.0; |
---|
1443 | y = secnds(int(x)); |
---|
1444 | ix = int(y * 100 + 43543000); |
---|
1445 | if(mod(ix,2) == 0) { |
---|
1446 | ix = ix + 1; |
---|
1447 | } |
---|
1448 | } |
---|
1449 | |
---|
1450 | // Here we are using random number generator copied from INCL code |
---|
1451 | // instead of the CERNLIB one! This causes difficulties for |
---|
1452 | // automatic testing since the random number generators, and thus |
---|
1453 | // the behavior of the routines in C++ and FORTRAN versions is no |
---|
1454 | // longer exactly the same! |
---|
1455 | x = randomGenerator->getRandom(); |
---|
1456 | // standardRandom(&x, &ix); |
---|
1457 | for(G4int i = 0; i < pSize; i++) { //do i=1,110 |
---|
1458 | p[i] = randomGenerator->getRandom(); |
---|
1459 | // standardRandom(&(p[i]), &ix); |
---|
1460 | } |
---|
1461 | a = randomGenerator->getRandom(); |
---|
1462 | standardRandom(&a, &ix); |
---|
1463 | k = 0; |
---|
1464 | } |
---|
1465 | |
---|
1466 | i = nint(100*a)+1; |
---|
1467 | haz = p[i]; |
---|
1468 | a = randomGenerator->getRandom(); |
---|
1469 | // standardRandom(&a, &ix); |
---|
1470 | p[i] = a; |
---|
1471 | |
---|
1472 | hazard->ial = ix; |
---|
1473 | return haz; |
---|
1474 | } |
---|
1475 | |
---|
1476 | |
---|
1477 | G4double G4AblaFissionSimfis18::gausshaz(int k, double xmoy, double sig) |
---|
1478 | { |
---|
1479 | // Gaussian random numbers: |
---|
1480 | |
---|
1481 | // 1005 C*** TIRAGE ALEATOIRE DANS UNE GAUSSIENNE DE LARGEUR SIG ET MOYENNE XMOY |
---|
1482 | static G4int iset = 0; |
---|
1483 | static G4double v1,v2,r,fac,gset,gausshaz; |
---|
1484 | |
---|
1485 | if(iset == 0) { //then |
---|
1486 | do { |
---|
1487 | v1 = 2.0*haz(k) - 1.0; |
---|
1488 | v2 = 2.0*haz(k) - 1.0; |
---|
1489 | r = std::pow(v1,2) + std::pow(v2,2); |
---|
1490 | } while(r >= 1); |
---|
1491 | |
---|
1492 | fac = std::sqrt(-2.*std::log(r)/r); |
---|
1493 | gset = v1*fac; |
---|
1494 | gausshaz = v2*fac*sig+xmoy; |
---|
1495 | iset = 1; |
---|
1496 | } |
---|
1497 | else { |
---|
1498 | gausshaz=gset*sig+xmoy; |
---|
1499 | iset=0; |
---|
1500 | } |
---|
1501 | return gausshaz; |
---|
1502 | } |
---|
1503 | |
---|
1504 | |
---|
1505 | // Utilities |
---|
1506 | |
---|
1507 | G4double G4AblaFissionSimfis18::min(G4double a, G4double b) |
---|
1508 | { |
---|
1509 | if(a < b) { |
---|
1510 | return a; |
---|
1511 | } |
---|
1512 | else { |
---|
1513 | return b; |
---|
1514 | } |
---|
1515 | } |
---|
1516 | |
---|
1517 | G4int G4AblaFissionSimfis18::min(G4int a, G4int b) |
---|
1518 | { |
---|
1519 | if(a < b) { |
---|
1520 | return a; |
---|
1521 | } |
---|
1522 | else { |
---|
1523 | return b; |
---|
1524 | } |
---|
1525 | } |
---|
1526 | |
---|
1527 | G4double G4AblaFissionSimfis18::max(G4double a, G4double b) |
---|
1528 | { |
---|
1529 | if(a > b) { |
---|
1530 | return a; |
---|
1531 | } |
---|
1532 | else { |
---|
1533 | return b; |
---|
1534 | } |
---|
1535 | } |
---|
1536 | |
---|
1537 | G4int G4AblaFissionSimfis18::max(G4int a, G4int b) |
---|
1538 | { |
---|
1539 | if(a > b) { |
---|
1540 | return a; |
---|
1541 | } |
---|
1542 | else { |
---|
1543 | return b; |
---|
1544 | } |
---|
1545 | } |
---|
1546 | |
---|
1547 | G4int G4AblaFissionSimfis18::nint(G4double number) |
---|
1548 | { |
---|
1549 | G4double intpart = 0.0; |
---|
1550 | G4double fractpart = 0.0; |
---|
1551 | fractpart = std::modf(number, &intpart); |
---|
1552 | if(number == 0) { |
---|
1553 | return 0; |
---|
1554 | } |
---|
1555 | if(number > 0) { |
---|
1556 | if(fractpart < 0.5) { |
---|
1557 | return int(std::floor(number)); |
---|
1558 | } |
---|
1559 | else { |
---|
1560 | return int(std::ceil(number)); |
---|
1561 | } |
---|
1562 | } |
---|
1563 | if(number < 0) { |
---|
1564 | if(fractpart < -0.5) { |
---|
1565 | return int(std::floor(number)); |
---|
1566 | } |
---|
1567 | else { |
---|
1568 | return int(std::ceil(number)); |
---|
1569 | } |
---|
1570 | } |
---|
1571 | |
---|
1572 | return int(std::floor(number)); |
---|
1573 | } |
---|
1574 | |
---|
1575 | G4int G4AblaFissionSimfis18::secnds(G4int x) |
---|
1576 | { |
---|
1577 | time_t mytime; |
---|
1578 | tm *mylocaltime; |
---|
1579 | |
---|
1580 | time(&mytime); |
---|
1581 | mylocaltime = localtime(&mytime); |
---|
1582 | |
---|
1583 | if(x == 0) { |
---|
1584 | return(mylocaltime->tm_hour*60*60 + mylocaltime->tm_min*60 + mylocaltime->tm_sec); |
---|
1585 | } |
---|
1586 | else { |
---|
1587 | return(mytime - x); |
---|
1588 | } |
---|
1589 | } |
---|
1590 | |
---|
1591 | G4int G4AblaFissionSimfis18::mod(G4int a, G4int b) |
---|
1592 | { |
---|
1593 | if(b != 0) { |
---|
1594 | return (a - (a/b)*b); |
---|
1595 | } |
---|
1596 | else { |
---|
1597 | return 0; |
---|
1598 | } |
---|
1599 | } |
---|
1600 | |
---|
1601 | G4double G4AblaFissionSimfis18::dmod(G4double a, G4double b) |
---|
1602 | { |
---|
1603 | if(b != 0) { |
---|
1604 | return (a - (a/b)*b); |
---|
1605 | } |
---|
1606 | else { |
---|
1607 | return 0.0; |
---|
1608 | } |
---|
1609 | } |
---|
1610 | |
---|
1611 | G4double G4AblaFissionSimfis18::dint(G4double a) |
---|
1612 | { |
---|
1613 | G4double value = 0.0; |
---|
1614 | |
---|
1615 | if(a < 0.0) { |
---|
1616 | value = double(std::ceil(a)); |
---|
1617 | } |
---|
1618 | else { |
---|
1619 | value = double(std::floor(a)); |
---|
1620 | } |
---|
1621 | |
---|
1622 | return value; |
---|
1623 | } |
---|
1624 | |
---|
1625 | G4int G4AblaFissionSimfis18::idint(G4double a) |
---|
1626 | { |
---|
1627 | G4int value = 0; |
---|
1628 | |
---|
1629 | if(a < 0) { |
---|
1630 | value = int(std::ceil(a)); |
---|
1631 | } |
---|
1632 | else { |
---|
1633 | value = int(std::floor(a)); |
---|
1634 | } |
---|
1635 | |
---|
1636 | return value; |
---|
1637 | } |
---|
1638 | |
---|
1639 | G4int G4AblaFissionSimfis18::idnint(G4double value) |
---|
1640 | { |
---|
1641 | G4double valueCeil = int(std::ceil(value)); |
---|
1642 | G4double valueFloor = int(std::floor(value)); |
---|
1643 | |
---|
1644 | if(std::fabs(value - valueCeil) < std::fabs(value - valueFloor)) { |
---|
1645 | return int(valueCeil); |
---|
1646 | } |
---|
1647 | else { |
---|
1648 | return int(valueFloor); |
---|
1649 | } |
---|
1650 | } |
---|
1651 | |
---|
1652 | G4double G4AblaFissionSimfis18::dmin1(G4double a, G4double b, G4double c) |
---|
1653 | { |
---|
1654 | if(a < b && a < c) { |
---|
1655 | return a; |
---|
1656 | } |
---|
1657 | if(b < a && b < c) { |
---|
1658 | return b; |
---|
1659 | } |
---|
1660 | if(c < a && c < b) { |
---|
1661 | return c; |
---|
1662 | } |
---|
1663 | return a; |
---|
1664 | } |
---|
1665 | |
---|
1666 | G4double G4AblaFissionSimfis18::utilabs(G4double a) |
---|
1667 | { |
---|
1668 | if(a > 0) { |
---|
1669 | return a; |
---|
1670 | } |
---|
1671 | if(a < 0) { |
---|
1672 | return (-1*a); |
---|
1673 | } |
---|
1674 | if(a == 0) { |
---|
1675 | return a; |
---|
1676 | } |
---|
1677 | |
---|
1678 | return a; |
---|
1679 | } |
---|