Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EquilibriumEvaporator.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// 20100114 M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly
28// 20100308 M. Kelsey -- Bug fix for setting masses of evaporating nuclei
29// 20100319 M. Kelsey -- Use new generateWithRandomAngles for theta,phi stuff;
30// eliminate some unnecessary std::pow()
31// 20100319 M. Kelsey -- Bug fix in new GetBindingEnergy() use right after
32// goodRemnant() -- Q1 should be outside call.
33// 20100412 M. Kelsey -- Pass G4CollisionOutput by ref to ::collide()
34// 20100413 M. Kelsey -- Pass buffers to paraMaker[Truncated]
35// 20100419 M. Kelsey -- Handle fission output list via const-ref
36// 20100517 M. Kelsey -- Use G4CascadeInterpolator for QFF
37// 20100520 M. Kelsey -- Inherit from common base class, make other colliders
38// simple data members. Rename timeToBigBang() to override
39// base explosion().
40// 20100617 M. Kelsey -- Remove "RUN" preprocessor flag and all "#else" code,
41// pass verbosity to colliders.
42// 20100620 M. Kelsey -- Use local "bindingEnergy()" function to call through.
43// 20100701 M. Kelsey -- Don't need to add excitation to nuclear mass; compute
44// new excitation energies properly (mass differences)
45// 20100702 M. Kelsey -- Simplify if-cascades, indentation
46// 20100712 M. Kelsey -- Add conservation checking
47// 20100714 M. Kelsey -- Move conservation checking to base class. Use
48// _generated_ evaporate energy (S) to adjust EEXS directly,
49// and test for S < EEXS before any particle generation; shift
50// nucleus momentum (PEX) by evaporate momentum directly
51// 20100719 M. Kelsey -- Remove duplicative EESX_new calculation.
52// 20100923 M. Kelsey -- Migrate to integer A and Z
53// 20110214 M. Kelsey -- Follow G4InuclParticle::Model enumerator migration
54// 20110728 M. Kelsey -- Fix Coverity #22951: check for "icase = -1" after
55// loop which is supposed to set it. Otherwise indexing is wrong.
56// 20110801 M. Kelsey -- Move "parms" buffer to data member, allocate in
57// constructor.
58// 20110809 M. Kelsey -- Move "foutput" to data member, get list by reference;
59// create final-state particles within "push_back" to avoid
60// creation of temporaries.
61// 20110922 M. Kelsey -- Follow G4InuclParticle::print(ostream&) migration
62// 20111007 M. Kelsey -- Add G4InuclParticleNames, replace hardcoded numbers
63// 20120608 M. Kelsey -- Fix variable-name "shadowing" compiler warnings.
64// 20121009 M. Kelsey -- Improve debugging output
65// 20130129 M. Kelsey -- Move QF interpolation to global statics for
66// multi-threaded shared memory.
67// 20130621 M. Kelsey -- Remove turning off conservation checks after fission
68// 20130622 Inherit from G4CascadeDeexciteBase, move to deExcite() interface
69// with G4Fragment
70// 20130628 Fissioner produces G4Fragment output directly in G4CollisionOutput
71// 20130808 M. Kelsey -- Use new object-version of paraMaker, for thread safety
72// 20130924 M. Kelsey -- Use G4Log, G4Exp for CPU speedup
73// 20150608 M. Kelsey -- Label all while loops as terminating.
74// 20150622 M. Kelsey -- For new G4cbrt(int), move powers of A outside.
75
77#include "G4SystemOfUnits.hh"
78#include "Randomize.hh"
79#include "G4BigBanger.hh"
81#include "G4CollisionOutput.hh"
82#include "G4Fissioner.hh"
83#include "G4Fragment.hh"
84#include "G4InuclNuclei.hh"
87#include "G4LorentzConvertor.hh"
88#include "G4LorentzVector.hh"
89#include "G4ThreeVector.hh"
90#include "G4Log.hh"
91#include "G4Exp.hh"
92
93using namespace G4InuclParticleNames;
94using namespace G4InuclSpecialFunctions;
95
96namespace { // Interpolation arrays for QF
97 const G4double QFREP[72] = {
98 // TL201 * * * *
99 // 1 2 3 4 5
100 22.5, 22.0, 21.0, 21.0, 20.0,
101 // BI209 BI207 PO210 AT213 * TH234
102 // 6 7 8 9 10 11
103 20.6, 20.6, 18.6, 15.8, 13.5, 6.5,
104 // TH233 TH232 TH231 TH230 TX229 PA233 PA232 PA231 PA230 U240
105 // 12 13 14 15 16 17 18 19 20 21
106 6.65, 6.22, 6.27, 6.5, 6.7, 6.2, 6.25, 5.9, 6.1, 5.75,
107 // U239 U238 U237 U236 U235 U234 U233 U232 U231
108 // 22 23 24 25 26 27 28 29 30
109 6.46, 5.7, 6.28, 5.8, 6.15, 5.6, 5.8, 5.2, 5.8,
110 // NP238 NP237 NP236 NP235 PU245 NP234 PU244 NP233
111 // 31 32 33 34 35 36 37 38
112 6.2 , 5.9 , 5.9, 6.0, 5.8, 5.7, 5.4, 5.4,
113 // PU242 PU241 PU240 PU239 PU238 AM247 PU236 AM245 AM244 AM243
114 // 39 40 41 42 43 44 45 46 47 48
115 5.6, 6.1, 5.57, 6.3, 5.5, 5.8, 4.7, 6.2, 6.4, 6.2,
116 // AM242 AM241 AM240 CM250 AM239 CM249 CM248 CM247 CM246
117 // 49 50 51 52 53 54 55 56 57
118 6.5, 6.2, 6.5, 5.3, 6.4, 5.7, 5.7, 6.2, 5.7,
119 // CM245 CM244 CM243 CM242 CM241 BK250 CM240
120 // 58 59 60 61 62 63 64
121 6.3, 5.8, 6.7, 5.8, 6.6, 6.1, 4.3,
122 // BK249 CF252 CF250 CF248 CF246 ES254 ES253 FM254
123 // 65 66 67 68 69 70 71 72
124 6.2, 3.8, 5.6, 4.0, 4.0, 4.2, 4.2, 3.5 };
125
126 static const G4double XREP[72] = {
127 // 1 2 3 4 5
128 0.6761, 0.677, 0.6788, 0.6803, 0.685,
129 // 6 7 8 9 10 11
130 0.6889, 0.6914, 0.6991, 0.7068, 0.725, 0.7391,
131 // 12 13 14 15 16 17 18 19 20 21
132 0.74, 0.741, 0.742, 0.743, 0.744, 0.7509, 0.752, 0.7531, 0.7543, 0.7548,
133 // 22 23 24
134 0.7557, 0.7566, 0.7576,
135 // 25 26 27 28 29 30 31 32 33 34
136 0.7587, 0.7597, 0.7608, 0.762, 0.7632, 0.7644, 0.7675, 0.7686, 0.7697, 0.7709,
137 // 35 36 37 38 39 40 41
138 0.7714, 0.7721, 0.7723, 0.7733, 0.7743, 0.7753, 0.7764,
139 // 42 43 44 45 46 47 48 49
140 0.7775, 0.7786, 0.7801, 0.781, 0.7821, 0.7831, 0.7842, 0.7852,
141 // 50 51 52 53 54 55 56 57 58
142 0.7864, 0.7875, 0.7880, 0.7887, 0.7889, 0.7899, 0.7909, 0.7919, 0.7930,
143 // 59 60 61 62 63 64
144 0.7941, 0.7953, 0.7965, 0.7977, 0.7987, 0.7989,
145 // 65 66 67 68 69 70 71 72
146 0.7997, 0.8075, 0.8097, 0.8119, 0.8143, 0.8164, 0.8174, 0.8274 };
147}
148
149
150// Constructor and destructor
151
153 : G4CascadeDeexciteBase("G4EquilibriumEvaporator"),
154 theParaMaker(verboseLevel), QFinterp(XREP) {
155 parms.first.resize(6,0.);
156 parms.second.resize(6,0.);
157}
158
160
163 theFissioner.setVerboseLevel(verbose);
164 theBigBanger.setVerboseLevel(verbose);
165}
166
167
168// Main operation
169
171 G4CollisionOutput& output) {
172 if (verboseLevel) {
173 G4cout << " >>> G4EquilibriumEvaporator::deExcite" << G4endl;
174 }
175
176 if (verboseLevel>1) G4cout << " evaporating target: \n" << target << G4endl;
177
178 // Implementation of the equilibium evaporation according to Dostrovsky
179 // (Phys. Rev. 116 (1959) 683.
180 // Note that pairing energy which is of order 1 MeV for odd-even and even-odd
181 // nuclei is not included
182
183 // Element in array: 0 : neutron
184 // 1 : proton
185 // 2 : deuteron
186 // 3 : triton
187 // 4 : 3He
188 // 5 : alpha
189
190 const G4double huge_num = 50.0;
191 const G4double small = -50.0;
192 const G4double prob_cut_off = 1.0e-15;
193 const G4double Q1[6] = { 0.0, 0.0, 2.23, 8.49, 7.72, 28.3 }; // binding energy
194 const G4int AN[6] = { 1, 1, 2, 3, 3, 4 }; // mass number
195 const G4int Q[6] = { 0, 1, 1, 1, 2, 2 }; // charge
196 const G4double G[6] = { 2.0, 2.0, 6.0, 6.0, 6.0, 4.0 };
197 const G4double BE = 0.0063;
198 const G4double fisssion_cut = 1000.0;
199 const G4double cut_off_energy = 0.1;
200
201 const G4double BF = 0.0242;
202 const G4int itry_max = 1000;
203 const G4int itry_global_max = 1000;
204 const G4double small_ekin = 1.0e-6;
205 const G4int itry_gam_max = 100;
206
207 G4double W[8];
208 G4double u[6]; // Level density for each particle type: (Atarg - Aemitted)*parlev
209 G4double V[6]; // Coulomb energy
210 G4double TM[6]; // Maximum possible kinetic energy of emitted particle
211 G4int A1[6]; // A of remnant nucleus
212 G4int Z1[6]; // Z of remnant nucleus
213
214 getTargetData(target);
215 if (verboseLevel > 3) G4cout << " after noeq: eexs " << EEXS << G4endl;
216
217 G4InuclElementaryParticle dummy(small_ekin, proton);
218 G4LorentzConvertor toTheNucleiSystemRestFrame;
219 //*** toTheNucleiSystemRestFrame.setVerbose(verboseLevel);
220 toTheNucleiSystemRestFrame.setBullet(dummy);
221
222 G4LorentzVector ppout;
223
224 // See if fragment should just be dispersed
225 if (explosion(A, Z, EEXS)) {
226 if (verboseLevel > 1) G4cout << " big bang in eql start " << G4endl;
227 theBigBanger.deExcite(target, output);
228
229 validateOutput(target, output); // Check energy conservation
230 return;
231 }
232
233 // If nucleus is in ground state, no evaporation
234 if (EEXS < cut_off_energy) {
235 if (verboseLevel > 1) G4cout << " no energy for evaporation" << G4endl;
236 output.addRecoilFragment(target);
237
238 validateOutput(target, output); // Check energy conservation
239 return;
240 }
241
242 // Initialize evaporation attempts
243 G4double coul_coeff = (A >= 100.0) ? 1.4 : 1.2;
244
245 G4LorentzVector pin = PEX; // Save original target for testing
246
247 G4bool try_again = true;
248 G4bool fission_open = true;
249 G4int itry_global = 0;
250
251 /* Loop checking 08.06.2015 MHK */
252 while (try_again && itry_global < itry_global_max) {
253 itry_global++;
254
255 // Set rest frame of current (recoiling) nucleus
256 toTheNucleiSystemRestFrame.setTarget(PEX);
257 toTheNucleiSystemRestFrame.toTheTargetRestFrame();
258
259 if (verboseLevel > 2) {
261 G4cout << " A " << A << " Z " << Z << " mass " << nuc_mass
262 << " EEXS " << EEXS << G4endl;
263 }
264
265 if (explosion(A, Z, EEXS)) { // big bang
266 if (verboseLevel > 2)
267 G4cout << " big bang in eql step " << itry_global << G4endl;
268
269 theBigBanger.deExcite(makeFragment(PEX,A,Z,EEXS), output);
270
271 validateOutput(target, output); // Check energy conservation
272 return;
273 }
274
275 if (EEXS < cut_off_energy) { // Evaporation not possible
276 if (verboseLevel > 2)
277 G4cout << " no energy for evaporation in eql step " << itry_global
278 << G4endl;
279
280 try_again = false;
281 break;
282 }
283
284 // Normal evaporation chain
285 G4double E0 = getE0(A);
286 G4double parlev = getPARLEVDEN(A, Z);
287 G4double u1 = parlev * A;
288
289 theParaMaker.getParams(Z, parms);
290 const std::vector<G4double>& AK = parms.first;
291 const std::vector<G4double>& CPA = parms.second;
292
293 G4double DM0 = bindingEnergy(A,Z);
294 G4int i(0);
295
296 for (i = 0; i < 6; i++) {
297 A1[i] = A - AN[i];
298 Z1[i] = Z - Q[i];
299 u[i] = parlev * A1[i];
300 V[i] = 0.0;
301 TM[i] = -0.1;
302
303 if (goodRemnant(A1[i], Z1[i])) {
304 G4double QB = DM0 - bindingEnergy(A1[i],Z1[i]) - Q1[i];
305 V[i] = coul_coeff * Z * Q[i] * AK[i] / (1.0 + EEXS / E0) /
306 (G4cbrt(A1[i]) + G4cbrt(AN[i]));
307 TM[i] = EEXS - QB - V[i] * A / A1[i];
308 if (verboseLevel > 3) {
309 G4cout << " i=" << i << " : QB " << QB << " u " << u[i]
310 << " V " << V[i] << " TM " << TM[i] << G4endl;
311 }
312 }
313 }
314
315 G4double ue = 2.0 * std::sqrt(u1 * EEXS);
316 G4double prob_sum = 0.0;
317
318 W[0] = 0.0;
319 if (TM[0] > cut_off_energy) {
320 G4double AL = getAL(A);
321 G4double A13 = G4cbrt(A1[0]);
322 W[0] = BE * A13*A13 * G[0] * AL;
323 G4double TM1 = 2.0 * std::sqrt(u[0] * TM[0]) - ue;
324
325 if (TM1 > huge_num) TM1 = huge_num;
326 else if (TM1 < small) TM1 = small;
327
328 W[0] *= G4Exp(TM1);
329 prob_sum += W[0];
330 }
331
332 for (i = 1; i < 6; i++) {
333 W[i] = 0.0;
334 if (TM[i] > cut_off_energy) {
335 G4double A13 = G4cbrt(A1[i]);
336 W[i] = BE * A13*A13 * G[i] * (1.0 + CPA[i]);
337 G4double TM1 = 2.0 * std::sqrt(u[i] * TM[i]) - ue;
338
339 if (TM1 > huge_num) TM1 = huge_num;
340 else if (TM1 < small) TM1 = small;
341
342 W[i] *= G4Exp(TM1);
343 prob_sum += W[i];
344 }
345 }
346
347 // fisson part
348 W[6] = 0.0;
349 if (A >= 100.0 && fission_open) {
350 G4double X2 = G4double(Z*Z)/G4double(A);
351 G4double X1 = 1.0 - 2.0*G4double(Z)/G4double(A);
352 G4double X = 0.019316 * X2 / (1.0 - 1.79 * X1 * X1);
353 G4double EF = EEXS - getQF(X, X2, A, Z, EEXS);
354
355 if (EF > 0.0) {
356 G4double AF = u1 * getAF(X, A, Z, EEXS);
357 G4double TM1 = 2.0 * std::sqrt(AF * EF) - ue;
358
359 if (TM1 > huge_num) TM1 = huge_num;
360 else if (TM1 < small) TM1 = small;
361
362 W[6] = BF * G4Exp(TM1);
363 if (W[6] > fisssion_cut*W[0]) W[6] = fisssion_cut*W[0];
364
365 prob_sum += W[6];
366 }
367 }
368
369 // again time to decide what next
370 if (verboseLevel > 2) {
371 G4cout << " Evaporation probabilities: sum = " << prob_sum
372 << "\n\t n " << W[0] << " p " << W[1] << " d " << W[2]
373 << " He3 " << W[3] << "\n\t t " << W[4] << " alpha " << W[5]
374 << " fission " << W[6] << G4endl;
375 }
376
377 G4int icase = -1;
378
379 if (prob_sum < prob_cut_off) { // photon emission chain
380 G4double UCR0 = 2.5 + 150.0 / A;
381 G4double T00 = 1.0 / (std::sqrt(u1 / UCR0) - 1.25 / UCR0);
382
383 if (verboseLevel > 3)
384 G4cout << " UCR0 " << UCR0 << " T00 " << T00 << G4endl;
385
386 G4int itry_gam = 0;
387 while (EEXS > cut_off_energy && try_again) {
388 itry_gam++;
389 G4int itry = 0;
390 G4double T04 = 4.0 * T00;
391 G4double FMAX;
392
393 if (T04 < EEXS) {
394 FMAX = (T04*T04*T04*T04) * G4Exp((EEXS - T04) / T00);
395 } else {
396 FMAX = EEXS*EEXS*EEXS*EEXS;
397 };
398
399 if (verboseLevel > 3)
400 G4cout << " T04 " << T04 << " FMAX (EEXS^4) " << FMAX << G4endl;
401
402 G4double S(0), X1(0);
403 while (itry < itry_max) {
404 itry++;
405 S = EEXS*G4UniformRand();
406 X1 = (S*S*S*S) * G4Exp((EEXS - S) / T00);
407
408 if (X1 > FMAX*G4UniformRand() ) break;
409 }
410
411 if (itry == itry_max) { // Maximum attempts exceeded
412 try_again = false;
413 break;
414 }
415
416 if (verboseLevel > 2) G4cout << " photon escape ?" << G4endl;
417
418 if (S < EEXS) { // Valid evaporate
419 S /= GeV; // Convert to Bertini units
420
421 G4double pmod = S;
423
424 // Push evaporated particle into current rest frame
425 mom = toTheNucleiSystemRestFrame.backToTheLab(mom);
426
427 if (verboseLevel > 3) {
428 G4cout << " nucleus px " << PEX.px() << " py " << PEX.py()
429 << " pz " << PEX.pz() << " E " << PEX.e() << G4endl
430 << " evaporate px " << mom.px() << " py " << mom.py()
431 << " pz " << mom.pz() << " E " << mom.e() << G4endl;
432 }
433
434 PEX -= mom; // Remaining four-momentum
435 EEXS -= S*GeV; // New excitation energy (in MeV)
436
437 // NOTE: In-situ construction will be optimized away (no copying)
440
441 if (verboseLevel > 3)
442 G4cout << output.getOutgoingParticles().back() << G4endl;
443
444 ppout += mom;
445 } else {
446 if (itry_gam == itry_gam_max) try_again = false;
447 }
448 } // while (EEXS > cut_off
449 try_again = false;
450 } else { // if (prob_sum < prob_cut_off)
451 G4double SL = prob_sum*G4UniformRand();
452
453 if (verboseLevel > 3) G4cout << " random SL " << SL << G4endl;
454
455 G4double S1 = 0.0;
456 for (i = 0; i < 7; i++) { // Select evaporation scenario
457 S1 += W[i];
458 if (SL <= S1) {
459 icase = i;
460 break;
461 };
462 };
463
464 if (icase < 0) continue; // Failed to choose scenario, try again
465
466 if (icase < 6) { // particle or light nuclei escape
467 if (verboseLevel > 2) {
468 G4cout << " particle/light-ion escape ("
469 << (icase==0 ? "neutron" : icase==1 ? "proton" :
470 icase==2 ? "deuteron" : icase==3 ? "He3" :
471 icase==4 ? "triton" : icase==5 ? "alpha" : "ERROR!")
472 << ")?" << G4endl;
473 }
474
475 G4double Xmax = (std::sqrt(u[icase]*TM[icase] + 0.25) - 0.5)/u[icase];
476 G4int itry1 = 0;
477 G4bool bad = true;
478
479 while (itry1 < itry_max && bad) {
480 itry1++;
481 G4int itry = 0;
482 G4double S = 0.0;
483 G4double X = 0.0;
484 G4double Ptest = 0.0;
485
486 while (itry < itry_max) {
487 itry++;
488
489 // Sampling from eq. 17 of Dostrovsky
490 X = G4UniformRand()*TM[icase];
491 Ptest = (X/Xmax)*G4Exp(-2.*u[icase]*Xmax + 2.*std::sqrt(u[icase]*(TM[icase] - X)));
492 if (G4UniformRand() < Ptest) {
493 S = X + V[icase];
494 break;
495 }
496 }
497
498 if (S > V[icase] && S < EEXS) { // real escape
499
500 if (verboseLevel > 2)
501 G4cout << " escape itry1 " << itry1 << " icase "
502 << icase << " S (MeV) " << S << G4endl;
503
504 S /= GeV; // Convert to Bertini units
505
506 if (icase < 2) { // particle escape
507 G4int ptype = 2 - icase;
508 if (verboseLevel > 2)
509 G4cout << " particle " << ptype << " escape" << G4endl;
510
511 // generate particle momentum
512 G4double mass =
514
515 G4double pmod = std::sqrt((2.0 * mass + S) * S);
517
518 // Push evaporated particle into current rest frame
519 mom = toTheNucleiSystemRestFrame.backToTheLab(mom);
520
521 if (verboseLevel > 2) {
522 G4cout << " nucleus px " << PEX.px() << " py " << PEX.py()
523 << " pz " << PEX.pz() << " E " << PEX.e() << G4endl
524 << " evaporate px " << mom.px() << " py " << mom.py()
525 << " pz " << mom.pz() << " E " << mom.e() << G4endl;
526 }
527
528 // New excitation energy depends on residual nuclear state
529 G4double mass_new =
530 G4InuclNuclei::getNucleiMass(A1[icase],Z1[icase]);
531
532 G4double EEXS_new = ((PEX-mom).m() - mass_new)*GeV;
533 if (EEXS_new < 0.0) continue; // Sanity check for new nucleus
534
535 PEX -= mom; // Remaining four-momentum
536 EEXS = EEXS_new;
537
538 A = A1[icase];
539 Z = Z1[icase];
540
541 // NOTE: In-situ construction optimizes away (no copying)
544
545 if (verboseLevel > 3)
546 G4cout << output.getOutgoingParticles().back() << G4endl;
547
548 ppout += mom;
549 bad = false;
550 } else { // if (icase < 2)
551 if (verboseLevel > 2) {
552 G4cout << " nucleus A " << AN[icase] << " Z " << Q[icase]
553 << " escape icase " << icase << G4endl;
554 }
555
556 G4double mass =
557 G4InuclNuclei::getNucleiMass(AN[icase],Q[icase]);
558
559 // generate particle momentum
560 G4double pmod = std::sqrt((2.0 * mass + S) * S);
562
563 // Push evaporated particle into current rest frame
564 mom = toTheNucleiSystemRestFrame.backToTheLab(mom);
565
566 if (verboseLevel > 2) {
567 G4cout << " nucleus px " << PEX.px() << " py " << PEX.py()
568 << " pz " << PEX.pz() << " E " << PEX.e() << G4endl
569 << " evaporate px " << mom.px() << " py " << mom.py()
570 << " pz " << mom.pz() << " E " << mom.e() << G4endl;
571 }
572
573 // New excitation energy depends on residual nuclear state
574 G4double mass_new =
575 G4InuclNuclei::getNucleiMass(A1[icase],Z1[icase]);
576
577 G4double EEXS_new = ((PEX-mom).m() - mass_new)*GeV;
578 if (EEXS_new < 0.0) continue; // Sanity check for new nucleus
579
580 PEX -= mom; // Remaining four-momentum
581 EEXS = EEXS_new;
582
583 A = A1[icase];
584 Z = Z1[icase];
585
586 // NOTE: In-situ constructor optimizes away (no copying)
588 AN[icase], Q[icase], 0.*GeV,
590
591 if (verboseLevel > 3)
592 G4cout << output.getOutgoingNuclei().back() << G4endl;
593
594 ppout += mom;
595 bad = false;
596 } // if (icase < 2)
597 } // if (EPR > 0.0 ...
598 } // while (itry1 ...
599
600 if (itry1 == itry_max || bad) try_again = false;
601 } else { // if (icase < 6)
602 if (verboseLevel > 2) {
603 G4cout << " fission: A " << A << " Z " << Z << " eexs " << EEXS
604 << " Wn " << W[0] << " Wf " << W[6] << G4endl;
605 }
606
607 // Catch fission output separately for verification
608 fission_output.reset();
609 theFissioner.deExcite(makeFragment(A,Z,EEXS), fission_output);
610
611 if (fission_output.numberOfFragments() == 2) { // fission ok
612 if (verboseLevel > 2) G4cout << " fission done in eql" << G4endl;
613
614 // Move fission fragments to lab frame for processing
615 fission_output.boostToLabFrame(toTheNucleiSystemRestFrame);
616
617 // Now evaporate the fission fragments individually
618 this->deExcite(fission_output.getRecoilFragment(0), output);
619 this->deExcite(fission_output.getRecoilFragment(1), output);
620
621 validateOutput(target, output); // Check energy conservation
622 return;
623 } else { // fission forbidden now
624 fission_open = false;
625 }
626 } // End of fission case
627 } // if (prob_sum < prob_cut_off)
628 } // while (try_again
629
630 // this time it's final nuclei
631
632 if (itry_global == itry_global_max) {
633 if (verboseLevel > 3) {
634 G4cout << " ! itry_global " << itry_global_max << G4endl;
635 }
636 }
637
638 G4LorentzVector pnuc = pin - ppout;
639
640 // NOTE: In-situ constructor will be optimized away (no copying)
641 output.addOutgoingNucleus(G4InuclNuclei(pnuc, A, Z, 0.,
643
644 if (verboseLevel > 3) {
645 G4cout << " remaining nucleus \n" << output.getOutgoingNuclei().back()
646 << G4endl;
647 }
648
649
650 validateOutput(target, output); // Check energy conservation
651 return;
652}
653
654G4bool G4EquilibriumEvaporator::explosion(G4int a,
655 G4int z,
656 G4double e) const {
657 if (verboseLevel > 3) {
658 G4cout << " >>> G4EquilibriumEvaporator::explosion? ";
659 }
660
661 const G4double be_cut = 3.0;
662
663 // Different criteria from base class, since nucleus more "agitated"
664 G4bool bigb = (!(a >= 12 && z >= 0 && z < 3*(a-z)) &&
665 (e >= be_cut * bindingEnergy(a,z))
666 );
667
668 if (verboseLevel > 3) G4cout << bigb << G4endl;
669
670 return bigb;
671}
672
673G4bool G4EquilibriumEvaporator::goodRemnant(G4int a,
674 G4int z) const {
675 if (verboseLevel > 3) {
676 G4cout << " >>> G4EquilibriumEvaporator::goodRemnant(" << a << "," << z
677 << ")? " << (a>1 && z>0 && a>z) << G4endl;
678 }
679
680 return a > 1 && z > 0 && a > z;
681}
682
683G4double G4EquilibriumEvaporator::getQF(G4double x,
684 G4double x2,
685 G4int a,
686 G4int /*z*/,
687 G4double ) const {
688 if (verboseLevel > 3) {
689 G4cout << " >>> G4EquilibriumEvaporator::getQF ";
690 }
691
692 const G4double G0 = 20.4;
693 const G4double XMIN = 0.6761;
694 const G4double XMAX = 0.8274;
695
696 G4double QFF = 0.0;
697
698 if (x < XMIN || x > XMAX) {
699 G4double X1 = 1.0 - 0.02 * x2;
700 G4double FX = (0.73 + (3.33 * X1 - 0.66) * X1) * (X1*X1*X1);
701 G4double A13 = G4cbrt(a);
702 QFF = G0 * FX * A13*A13;
703 } else {
704 QFF = QFinterp.interpolate(x, QFREP);
705 }
706
707 if (QFF < 0.0) QFF = 0.0;
708
709 if (verboseLevel>3) G4cout << " returns " << QFF << G4endl;
710
711 return QFF;
712}
713
714G4double G4EquilibriumEvaporator::getAF(G4double ,
715 G4int /*a*/,
716 G4int /*z*/,
717 G4double e) const {
718
719 if (verboseLevel > 3) {
720 G4cout << " >>> G4EquilibriumEvaporator::getAF" << G4endl;
721 }
722
723 // ugly parameterisation to fit the experimental fission cs for Hg - Bi nuclei
724 G4double AF = 1.285 * (1.0 - e / 1100.0);
725
726 if(AF < 1.06) AF = 1.06;
727 // if(AF < 1.08) AF = 1.08;
728
729 return AF;
730}
731
732G4double G4EquilibriumEvaporator::getPARLEVDEN(G4int /*a*/,
733 G4int /*z*/) const {
734
735 if (verboseLevel > 3) {
736 G4cout << " >>> G4EquilibriumEvaporator::getPARLEVDEN" << G4endl;
737 }
738
739 const G4double par = 0.125;
740
741 return par;
742}
743
744G4double G4EquilibriumEvaporator::getE0(G4int /*a*/) const {
745
746 if (verboseLevel > 3) {
747 G4cout << " >>> G4EquilibriumEvaporator::getE0" << G4endl;
748 }
749
750 const G4double e0 = 200.0;
751
752 return e0;
753}
G4double S(G4double temp)
#define A13
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
CLHEP::HepLorentzVector G4LorentzVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
G4CascadeDeexciteBase(const char *name)
void getTargetData(const G4Fragment &target)
const G4Fragment & makeFragment(G4LorentzVector mom, G4int A, G4int Z, G4double EX=0.)
virtual void setVerboseLevel(G4int verbose=0)
virtual G4bool validateOutput(const G4Fragment &target, G4CollisionOutput &output)
void addRecoilFragment(const G4Fragment *aFragment)
const std::vector< G4InuclNuclei > & getOutgoingNuclei() const
void addOutgoingParticle(const G4InuclElementaryParticle &particle)
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
void addOutgoingNucleus(const G4InuclNuclei &nuclei)
virtual void setVerboseLevel(G4int verbose)
virtual void deExcite(const G4Fragment &target, G4CollisionOutput &output)
static G4double getParticleMass(G4int type)
G4double getNucleiMass() const
void setBullet(const G4InuclParticle *bullet)
G4LorentzVector backToTheLab(const G4LorentzVector &mom) const
void setTarget(const G4InuclParticle *target)
#define W
Definition crc32.c:85
G4double bindingEnergy(G4int A, G4int Z)
G4LorentzVector generateWithRandomAngles(G4double p, G4double mass=0.)