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 | // $Id: G4BigBanger.cc,v 1.29 2010/05/21 17:56:34 mkelsey Exp $ |
---|
26 | // Geant4 tag: $Name: geant4-09-04-beta-cand-01 $ |
---|
27 | // |
---|
28 | // 20100114 M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly |
---|
29 | // 20100301 M. Kelsey -- In generateBangInSCM(), restore old G4CascMom calcs. |
---|
30 | // for (N-1)th outgoing nucleon. |
---|
31 | // 20100319 M. Kelsey -- Use new generateWithRandomAngles for theta,phi stuff |
---|
32 | // 20100407 M. Kelsey -- Replace std::vector<> returns with data members. |
---|
33 | // 20100413 M. Kelsey -- Pass G4CollisionOutput by ref to ::collide() |
---|
34 | // 20100517 M. Kelsey -- Inherit from common base class, clean up code |
---|
35 | |
---|
36 | #include "G4BigBanger.hh" |
---|
37 | #include "G4CollisionOutput.hh" |
---|
38 | #include "G4InuclNuclei.hh" |
---|
39 | #include "G4InuclElementaryParticle.hh" |
---|
40 | #include "G4InuclSpecialFunctions.hh" |
---|
41 | #include "G4ParticleLargerEkin.hh" |
---|
42 | #include "G4LorentzConvertor.hh" |
---|
43 | #include "G4HadTmpUtil.hh" |
---|
44 | #include "G4NucleiProperties.hh" |
---|
45 | #include <algorithm> |
---|
46 | |
---|
47 | using namespace G4InuclSpecialFunctions; |
---|
48 | |
---|
49 | typedef std::vector<G4InuclElementaryParticle>::iterator particleIterator; |
---|
50 | |
---|
51 | G4BigBanger::G4BigBanger() : G4VCascadeCollider("G4BigBanger") {} |
---|
52 | |
---|
53 | void |
---|
54 | G4BigBanger::collide(G4InuclParticle* /*bullet*/, G4InuclParticle* target, |
---|
55 | G4CollisionOutput& output) { |
---|
56 | |
---|
57 | if (verboseLevel > 3) { |
---|
58 | G4cout << " >>> G4BigBanger::collide" << G4endl; |
---|
59 | } |
---|
60 | |
---|
61 | // primitive explosion model A -> nucleons to prevent too exotic evaporation |
---|
62 | |
---|
63 | const G4double small_ekin = 1.0e-6; |
---|
64 | |
---|
65 | G4LorentzVector totscm; |
---|
66 | G4LorentzVector totlab; |
---|
67 | |
---|
68 | G4InuclNuclei* nuclei_target = dynamic_cast<G4InuclNuclei*>(target); |
---|
69 | if (!nuclei_target) { |
---|
70 | G4cerr << " BigBanger -> try to bang not nuclei " << G4endl; |
---|
71 | return; |
---|
72 | } |
---|
73 | |
---|
74 | G4double A = nuclei_target->getA(); |
---|
75 | G4double Z = nuclei_target->getZ(); |
---|
76 | G4LorentzVector PEX = nuclei_target->getMomentum(); |
---|
77 | G4double EEXS = nuclei_target->getExitationEnergy(); |
---|
78 | G4InuclElementaryParticle dummy(small_ekin, 1); |
---|
79 | G4LorentzConvertor toTheNucleiSystemRestFrame; |
---|
80 | |
---|
81 | toTheNucleiSystemRestFrame.setBullet(dummy); |
---|
82 | toTheNucleiSystemRestFrame.setTarget(nuclei_target); |
---|
83 | toTheNucleiSystemRestFrame.toTheTargetRestFrame(); |
---|
84 | |
---|
85 | G4double etot = (EEXS - G4NucleiProperties::GetBindingEnergy(G4lrint(A), G4lrint(Z) ) ) * MeV/GeV; // To Bertini units |
---|
86 | if (etot < 0.0) etot = 0.0; |
---|
87 | |
---|
88 | if (verboseLevel > 2) { |
---|
89 | G4cout << " BigBanger: target " << G4endl; |
---|
90 | nuclei_target->printParticle(); |
---|
91 | G4cout << " BigBanger: a " << A << " z " << Z << " eexs " << EEXS << " etot " << |
---|
92 | etot << " nm " << nuclei_target->getMass() << G4endl; |
---|
93 | } |
---|
94 | |
---|
95 | generateBangInSCM(etot, A, Z); |
---|
96 | |
---|
97 | if (verboseLevel > 2) { |
---|
98 | G4cout << " particles " << particles.size() << G4endl; |
---|
99 | for(G4int i = 0; i < G4int(particles.size()); i++) |
---|
100 | particles[i].printParticle(); |
---|
101 | } |
---|
102 | |
---|
103 | if(!particles.empty()) { // convert back to Lab |
---|
104 | particleIterator ipart; |
---|
105 | |
---|
106 | for(ipart = particles.begin(); ipart != particles.end(); ipart++) { |
---|
107 | if (verboseLevel > 2) { |
---|
108 | totscm += ipart->getMomentum(); |
---|
109 | } |
---|
110 | G4LorentzVector mom = |
---|
111 | toTheNucleiSystemRestFrame.backToTheLab(ipart->getMomentum()); |
---|
112 | ipart->setMomentum(mom); |
---|
113 | |
---|
114 | if (verboseLevel > 2) { |
---|
115 | mom = ipart->getMomentum(); |
---|
116 | totlab += mom; |
---|
117 | } |
---|
118 | } |
---|
119 | |
---|
120 | std::sort(particles.begin(), particles.end(), G4ParticleLargerEkin()); |
---|
121 | |
---|
122 | if (verboseLevel > 2) { |
---|
123 | G4cout << " In SCM: total outgoing momentum " << G4endl |
---|
124 | << " E " << totscm.e() << " px " << totscm.x() |
---|
125 | << " py " << totscm.y() << " pz " << totscm.z() << G4endl; |
---|
126 | G4cout << " In Lab: mom cons " << G4endl |
---|
127 | << " E " << PEX.e() + 0.001 * EEXS - totlab.e() |
---|
128 | << " px " << PEX.x() - totlab.x() |
---|
129 | << " py " << PEX.y() - totlab.y() |
---|
130 | << " pz " << PEX.z() - totlab.z() << G4endl; |
---|
131 | } |
---|
132 | } |
---|
133 | |
---|
134 | output.addOutgoingParticles(particles); |
---|
135 | |
---|
136 | return; |
---|
137 | } |
---|
138 | |
---|
139 | void G4BigBanger::generateBangInSCM(G4double etot, G4double a, G4double z) { |
---|
140 | if (verboseLevel > 3) { |
---|
141 | G4cout << " >>> G4BigBanger::generateBangInSCM" << G4endl; |
---|
142 | } |
---|
143 | |
---|
144 | // Proton and neutron masses |
---|
145 | const G4double mp = G4InuclElementaryParticle::getParticleMass(1); |
---|
146 | const G4double mn = G4InuclElementaryParticle::getParticleMass(2); |
---|
147 | |
---|
148 | const G4double ang_cut = 0.9999; |
---|
149 | const G4int itry_max = 1000; |
---|
150 | |
---|
151 | G4int ia = int(a + 0.1); |
---|
152 | G4int iz = int(z + 0.1); |
---|
153 | |
---|
154 | if (verboseLevel > 2) { |
---|
155 | G4cout << " ia " << ia << " iz " << iz << G4endl; |
---|
156 | } |
---|
157 | |
---|
158 | particles.clear(); // Reset output vector before filling |
---|
159 | |
---|
160 | if(ia == 1) { |
---|
161 | // abnormal situation |
---|
162 | G4double m = iz > 0 ? mp : mn; |
---|
163 | G4double pmod = std::sqrt((etot + 2.0 * m) * etot); |
---|
164 | G4LorentzVector mom = generateWithRandomAngles(pmod, m); |
---|
165 | |
---|
166 | G4int knd = iz > 0 ? 1 : 2; |
---|
167 | |
---|
168 | particles.push_back(G4InuclElementaryParticle(mom, knd, 8)); // modelId included |
---|
169 | |
---|
170 | return; |
---|
171 | }; |
---|
172 | |
---|
173 | generateMomentumModules(etot, a, z); |
---|
174 | G4bool bad = true; |
---|
175 | G4int itry = 0; |
---|
176 | |
---|
177 | while(bad && itry < itry_max) { |
---|
178 | itry++; |
---|
179 | std::vector<G4LorentzVector> scm_momentums; |
---|
180 | G4LorentzVector tot_mom; |
---|
181 | |
---|
182 | if(ia == 2) { |
---|
183 | // FIXME: This isn't actually a correct four-vector, wrong mass! |
---|
184 | G4LorentzVector mom = generateWithRandomAngles(momModules[0]); |
---|
185 | |
---|
186 | tot_mom += mom; |
---|
187 | |
---|
188 | scm_momentums.push_back(mom); |
---|
189 | |
---|
190 | G4LorentzVector mom1 = -mom; |
---|
191 | |
---|
192 | scm_momentums.push_back(mom1); |
---|
193 | bad = false; |
---|
194 | } else { |
---|
195 | for(G4int i = 0; i < ia - 2; i++) { |
---|
196 | // FIXME: This isn't actually a correct four-vector, wrong mass! |
---|
197 | G4LorentzVector mom = generateWithRandomAngles(momModules[i]); |
---|
198 | |
---|
199 | tot_mom += mom; |
---|
200 | |
---|
201 | scm_momentums.push_back(mom); |
---|
202 | }; |
---|
203 | |
---|
204 | // handle last two |
---|
205 | G4double tot_mod = tot_mom.rho(); |
---|
206 | G4double ct = -0.5 * (tot_mod * tot_mod + momModules[ia - 2] * momModules[ia - 2] - |
---|
207 | momModules[ia - 1] * momModules[ia - 1]) / tot_mod / momModules[ia - 2]; |
---|
208 | |
---|
209 | if (verboseLevel > 2) { |
---|
210 | G4cout << " ct last " << ct << G4endl; |
---|
211 | } |
---|
212 | |
---|
213 | if(std::fabs(ct) < ang_cut) { |
---|
214 | // FIXME: These are not actually four-vectors, just three-momenta |
---|
215 | G4LorentzVector mom2 = generateWithFixedTheta(ct, momModules[ia - 2]); |
---|
216 | // rotate to the normal system |
---|
217 | G4LorentzVector apr = tot_mom/tot_mod; |
---|
218 | G4double a_tr = std::sqrt(apr.x()*apr.x() + apr.y()*apr.y()); |
---|
219 | G4LorentzVector mom; |
---|
220 | mom.setX(mom2.z()*apr.x() + ( mom2.x()*apr.y() + mom2.y()*apr.z()*apr.x())/a_tr); |
---|
221 | mom.setY(mom2.z()*apr.y() + (-mom2.x()*apr.x() + mom2.y()*apr.z()*apr.y())/a_tr); |
---|
222 | mom.setZ(mom2.z()*apr.z() - mom2.y()*a_tr); |
---|
223 | |
---|
224 | scm_momentums.push_back(mom); |
---|
225 | // and the last one |
---|
226 | G4LorentzVector mom1= - mom - tot_mom; |
---|
227 | scm_momentums.push_back(mom1); |
---|
228 | bad = false; |
---|
229 | }; |
---|
230 | }; |
---|
231 | if(!bad) { |
---|
232 | for(G4int i = 0; i < ia; i++) { |
---|
233 | G4int knd = i < iz ? 1 : 2; |
---|
234 | |
---|
235 | particles.push_back(G4InuclElementaryParticle(scm_momentums[i], knd, 8)); |
---|
236 | }; |
---|
237 | }; |
---|
238 | }; |
---|
239 | if (verboseLevel > 2) { |
---|
240 | if(itry == itry_max) G4cout << " BigBanger -> can not generate bang " << G4endl; |
---|
241 | } |
---|
242 | |
---|
243 | return; |
---|
244 | } |
---|
245 | |
---|
246 | void G4BigBanger::generateMomentumModules(G4double etot, G4double a, G4double z) { |
---|
247 | if (verboseLevel > 3) { |
---|
248 | G4cout << " >>> G4BigBanger::generateMomentumModules" << G4endl; |
---|
249 | } |
---|
250 | |
---|
251 | // Proton and neutron masses |
---|
252 | const G4double mp = G4InuclElementaryParticle::getParticleMass(1); |
---|
253 | const G4double mn = G4InuclElementaryParticle::getParticleMass(2); |
---|
254 | |
---|
255 | momModules.clear(); // Reset buffer for filling |
---|
256 | |
---|
257 | G4int ia = G4int(a + 0.1); |
---|
258 | G4int iz = G4int(z + 0.1); |
---|
259 | |
---|
260 | G4double xtot = 0.0; |
---|
261 | G4double promax = maxProbability(a); |
---|
262 | |
---|
263 | G4int i; |
---|
264 | for(i = 0; i < ia; i++) { |
---|
265 | G4double x = generateX(ia, a, promax); |
---|
266 | |
---|
267 | if (verboseLevel > 2) { |
---|
268 | G4cout << " i " << i << " x " << x << G4endl; |
---|
269 | } |
---|
270 | momModules.push_back(x); |
---|
271 | xtot += x; |
---|
272 | }; |
---|
273 | for(i = 0; i < ia; i++) { |
---|
274 | G4double m = i < iz ? mp : mn; |
---|
275 | |
---|
276 | momModules[i] = momModules[i] * etot / xtot; |
---|
277 | momModules[i] = std::sqrt(momModules[i] * (momModules[i] + 2.0 * m)); |
---|
278 | |
---|
279 | if (verboseLevel > 2) { |
---|
280 | G4cout << " i " << i << " pmod " << momModules[i] << G4endl; |
---|
281 | } |
---|
282 | }; |
---|
283 | |
---|
284 | return; |
---|
285 | } |
---|
286 | |
---|
287 | G4double G4BigBanger::xProbability(G4double x, G4int ia) const { |
---|
288 | |
---|
289 | |
---|
290 | if (verboseLevel > 3) { |
---|
291 | G4cout << " >>> G4BigBanger::xProbability" << G4endl; |
---|
292 | } |
---|
293 | |
---|
294 | G4int ihalf = ia / 2; |
---|
295 | G4double ekpr = 0.0; |
---|
296 | |
---|
297 | if(x < 1.0 || x > 0.0) { |
---|
298 | ekpr = x * x; |
---|
299 | |
---|
300 | if(2 * ihalf == ia) { // even A |
---|
301 | ekpr *= std::sqrt(1.0 - x) * std::pow((1.0 - x), G4int(G4double(3 * ia - 6) / 2.0)); |
---|
302 | } |
---|
303 | else { |
---|
304 | ekpr *= std::pow((1.0 - x), G4int(G4double(3 * ia - 5) / 2.0)); |
---|
305 | }; |
---|
306 | }; |
---|
307 | |
---|
308 | return ekpr; |
---|
309 | } |
---|
310 | |
---|
311 | G4double G4BigBanger::maxProbability(G4double a) const { |
---|
312 | |
---|
313 | if (verboseLevel > 3) { |
---|
314 | G4cout << " >>> G4BigBanger::maxProbability" << G4endl; |
---|
315 | } |
---|
316 | |
---|
317 | return xProbability(1.0 / (a - 1.0) / 1.5, G4int(a + 0.1)); |
---|
318 | } |
---|
319 | |
---|
320 | G4double G4BigBanger::generateX(G4int ia, |
---|
321 | G4double a, |
---|
322 | G4double promax) const { |
---|
323 | |
---|
324 | if (verboseLevel > 3) { |
---|
325 | G4cout << " >>> G4BigBanger::generateX" << G4endl; |
---|
326 | } |
---|
327 | |
---|
328 | const G4int itry_max = 1000; |
---|
329 | G4int itry = 0; |
---|
330 | G4double x; |
---|
331 | |
---|
332 | while(itry < itry_max) { |
---|
333 | itry++; |
---|
334 | x = inuclRndm(); |
---|
335 | |
---|
336 | if(xProbability(x, ia) >= promax * inuclRndm()) return x; |
---|
337 | }; |
---|
338 | if (verboseLevel > 2) { |
---|
339 | G4cout << " BigBanger -> can not generate x " << G4endl; |
---|
340 | } |
---|
341 | |
---|
342 | return maxProbability(a); |
---|
343 | } |
---|