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: G4HEKaonZeroLongInelastic.cc,v 1.13 2010/11/29 05:44:44 dennis Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-04-ref-00 $ |
---|
28 | // |
---|
29 | |
---|
30 | #include "globals.hh" |
---|
31 | #include "G4ios.hh" |
---|
32 | |
---|
33 | // G4 Process: Gheisha High Energy Collision model. |
---|
34 | // This includes the high energy cascading model, the two-body-resonance model |
---|
35 | // and the low energy two-body model. Not included are the low energy stuff |
---|
36 | // like nuclear reactions, nuclear fission without any cascading and all |
---|
37 | // processes for particles at rest. |
---|
38 | // |
---|
39 | // New version by D.H. Wright (SLAC) to fix seg fault in old version |
---|
40 | // 26 January 2010 |
---|
41 | |
---|
42 | |
---|
43 | #include "G4HEKaonZeroLongInelastic.hh" |
---|
44 | |
---|
45 | G4HadFinalState* |
---|
46 | G4HEKaonZeroLongInelastic::ApplyYourself(const G4HadProjectile& aTrack, |
---|
47 | G4Nucleus& targetNucleus) |
---|
48 | { |
---|
49 | G4HEVector* pv = new G4HEVector[MAXPART]; |
---|
50 | const G4HadProjectile* aParticle = &aTrack; |
---|
51 | const G4double atomicWeight = targetNucleus.GetN(); |
---|
52 | const G4double atomicNumber = targetNucleus.GetZ(); |
---|
53 | G4HEVector incidentParticle(aParticle); |
---|
54 | |
---|
55 | G4int incidentCode = incidentParticle.getCode(); |
---|
56 | G4double incidentMass = incidentParticle.getMass(); |
---|
57 | G4double incidentTotalEnergy = incidentParticle.getEnergy(); |
---|
58 | G4double incidentTotalMomentum = incidentParticle.getTotalMomentum(); |
---|
59 | G4double incidentKineticEnergy = incidentTotalEnergy - incidentMass; |
---|
60 | |
---|
61 | if(incidentKineticEnergy < 1) |
---|
62 | G4cout << "GHEKaonZeroLongInelastic: incident energy < 1 GeV " << G4endl; |
---|
63 | |
---|
64 | if(verboseLevel > 1) { |
---|
65 | G4cout << "G4HEKaonZeroLongInelastic::ApplyYourself" << G4endl; |
---|
66 | G4cout << "incident particle " << incidentParticle.getName() |
---|
67 | << "mass " << incidentMass |
---|
68 | << "kinetic energy " << incidentKineticEnergy |
---|
69 | << G4endl; |
---|
70 | G4cout << "target material with (A,Z) = (" |
---|
71 | << atomicWeight << "," << atomicNumber << ")" << G4endl; |
---|
72 | } |
---|
73 | |
---|
74 | G4double inelasticity = NuclearInelasticity(incidentKineticEnergy, |
---|
75 | atomicWeight, atomicNumber); |
---|
76 | if(verboseLevel > 1) |
---|
77 | G4cout << "nuclear inelasticity = " << inelasticity << G4endl; |
---|
78 | |
---|
79 | incidentKineticEnergy -= inelasticity; |
---|
80 | |
---|
81 | G4double excitationEnergyGNP = 0.; |
---|
82 | G4double excitationEnergyDTA = 0.; |
---|
83 | |
---|
84 | G4double excitation = NuclearExcitation(incidentKineticEnergy, |
---|
85 | atomicWeight, atomicNumber, |
---|
86 | excitationEnergyGNP, |
---|
87 | excitationEnergyDTA); |
---|
88 | if(verboseLevel > 1) |
---|
89 | G4cout << "nuclear excitation = " << excitation << excitationEnergyGNP |
---|
90 | << excitationEnergyDTA << G4endl; |
---|
91 | |
---|
92 | incidentKineticEnergy -= excitation; |
---|
93 | incidentTotalEnergy = incidentKineticEnergy + incidentMass; |
---|
94 | incidentTotalMomentum = std::sqrt( (incidentTotalEnergy-incidentMass) |
---|
95 | *(incidentTotalEnergy+incidentMass)); |
---|
96 | |
---|
97 | G4HEVector targetParticle; |
---|
98 | if (G4UniformRand() < atomicNumber/atomicWeight) { |
---|
99 | targetParticle.setDefinition("Proton"); |
---|
100 | } else { |
---|
101 | targetParticle.setDefinition("Neutron"); |
---|
102 | } |
---|
103 | |
---|
104 | G4double targetMass = targetParticle.getMass(); |
---|
105 | G4double centerOfMassEnergy = std::sqrt(incidentMass*incidentMass |
---|
106 | + targetMass*targetMass |
---|
107 | + 2.0*targetMass*incidentTotalEnergy); |
---|
108 | G4double availableEnergy = centerOfMassEnergy - targetMass - incidentMass; |
---|
109 | |
---|
110 | G4bool inElastic = true; |
---|
111 | vecLength = 0; |
---|
112 | |
---|
113 | if(verboseLevel > 1) |
---|
114 | G4cout << "ApplyYourself: CallFirstIntInCascade for particle " |
---|
115 | << incidentCode << G4endl; |
---|
116 | |
---|
117 | G4bool successful = false; |
---|
118 | |
---|
119 | // Split K0L into K0 and K0bar |
---|
120 | if (G4UniformRand() < 0.5) |
---|
121 | FirstIntInCasAntiKaonZero(inElastic, availableEnergy, pv, vecLength, |
---|
122 | incidentParticle, targetParticle); |
---|
123 | else |
---|
124 | FirstIntInCasKaonZero(inElastic, availableEnergy, pv, vecLength, |
---|
125 | incidentParticle, targetParticle, atomicWeight); |
---|
126 | |
---|
127 | // Do nuclear interaction with either K0 or K0bar |
---|
128 | if ((vecLength > 0) && (availableEnergy > 1.)) |
---|
129 | StrangeParticlePairProduction(availableEnergy, centerOfMassEnergy, |
---|
130 | pv, vecLength, |
---|
131 | incidentParticle, targetParticle); |
---|
132 | |
---|
133 | HighEnergyCascading(successful, pv, vecLength, |
---|
134 | excitationEnergyGNP, excitationEnergyDTA, |
---|
135 | incidentParticle, targetParticle, |
---|
136 | atomicWeight, atomicNumber); |
---|
137 | if (!successful) |
---|
138 | HighEnergyClusterProduction(successful, pv, vecLength, |
---|
139 | excitationEnergyGNP, excitationEnergyDTA, |
---|
140 | incidentParticle, targetParticle, |
---|
141 | atomicWeight, atomicNumber); |
---|
142 | if (!successful) |
---|
143 | MediumEnergyCascading(successful, pv, vecLength, |
---|
144 | excitationEnergyGNP, excitationEnergyDTA, |
---|
145 | incidentParticle, targetParticle, |
---|
146 | atomicWeight, atomicNumber); |
---|
147 | |
---|
148 | if (!successful) |
---|
149 | MediumEnergyClusterProduction(successful, pv, vecLength, |
---|
150 | excitationEnergyGNP, excitationEnergyDTA, |
---|
151 | incidentParticle, targetParticle, |
---|
152 | atomicWeight, atomicNumber); |
---|
153 | if (!successful) |
---|
154 | QuasiElasticScattering(successful, pv, vecLength, |
---|
155 | excitationEnergyGNP, excitationEnergyDTA, |
---|
156 | incidentParticle, targetParticle, |
---|
157 | atomicWeight, atomicNumber); |
---|
158 | |
---|
159 | if (!successful) |
---|
160 | ElasticScattering(successful, pv, vecLength, |
---|
161 | incidentParticle, |
---|
162 | atomicWeight, atomicNumber); |
---|
163 | |
---|
164 | if (!successful) |
---|
165 | G4cout << "GHEInelasticInteraction::ApplyYourself fails to produce final state particles" |
---|
166 | << G4endl; |
---|
167 | |
---|
168 | // Check for K0, K0bar and change particle types to K0L, K0S if necessary |
---|
169 | G4int kcode; |
---|
170 | for (G4int i = 0; i < vecLength; i++) { |
---|
171 | kcode = pv[i].getCode(); |
---|
172 | if (kcode == KaonZero.getCode() || kcode == AntiKaonZero.getCode()) { |
---|
173 | if (G4UniformRand() < 0.5) |
---|
174 | pv[i] = KaonZeroShort; |
---|
175 | else |
---|
176 | pv[i] = KaonZeroLong; |
---|
177 | } |
---|
178 | } |
---|
179 | |
---|
180 | // ................ |
---|
181 | |
---|
182 | FillParticleChange(pv, vecLength); |
---|
183 | delete [] pv; |
---|
184 | theParticleChange.SetStatusChange(stopAndKill); |
---|
185 | return &theParticleChange; |
---|
186 | } |
---|
187 | |
---|
188 | |
---|
189 | void |
---|
190 | G4HEKaonZeroLongInelastic::FirstIntInCasKaonZero(G4bool& inElastic, |
---|
191 | const G4double availableEnergy, |
---|
192 | G4HEVector pv[], |
---|
193 | G4int& vecLen, |
---|
194 | const G4HEVector& incidentParticle, |
---|
195 | const G4HEVector& targetParticle, |
---|
196 | const G4double atomicWeight) |
---|
197 | |
---|
198 | // Kaon0 undergoes interaction with nucleon within a nucleus. Check if it is |
---|
199 | // energetically possible to produce pions/kaons. In not, assume nuclear excitation |
---|
200 | // occurs and input particle is degraded in energy. No other particles are produced. |
---|
201 | // If reaction is possible, find the correct number of pions/protons/neutrons |
---|
202 | // produced using an interpolation to multiplicity data. Replace some pions or |
---|
203 | // protons/neutrons by kaons or strange baryons according to the average |
---|
204 | // multiplicity per inelastic reaction. |
---|
205 | { |
---|
206 | static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp |
---|
207 | static const G4double expxl = -expxu; // lower bound for arg. of exp |
---|
208 | |
---|
209 | static const G4double protb = 0.7; |
---|
210 | static const G4double neutb = 0.7; |
---|
211 | static const G4double c = 1.25; |
---|
212 | |
---|
213 | static const G4int numMul = 1200; |
---|
214 | static const G4int numSec = 60; |
---|
215 | |
---|
216 | G4int neutronCode = Neutron.getCode(); |
---|
217 | G4int protonCode = Proton.getCode(); |
---|
218 | |
---|
219 | G4int targetCode = targetParticle.getCode(); |
---|
220 | G4double incidentTotalMomentum = incidentParticle.getTotalMomentum(); |
---|
221 | |
---|
222 | static G4bool first = true; |
---|
223 | static G4double protmul[numMul], protnorm[numSec]; // proton constants |
---|
224 | static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants |
---|
225 | |
---|
226 | // misc. local variables |
---|
227 | // np = number of pi+, nm = number of pi-, nz = number of pi0 |
---|
228 | |
---|
229 | G4int i, counter, nt, np, nm, nz; |
---|
230 | |
---|
231 | if (first) { |
---|
232 | // compute normalization constants, this will only be done once |
---|
233 | first = false; |
---|
234 | for( i=0; i<numMul; i++ )protmul[i] = 0.0; |
---|
235 | for( i=0; i<numSec; i++ )protnorm[i] = 0.0; |
---|
236 | counter = -1; |
---|
237 | for (np=0; np<(numSec/3); np++) { |
---|
238 | for (nm=std::max(0,np-1); nm<=(np+1); nm++) { |
---|
239 | for (nz=0; nz<numSec/3; nz++) { |
---|
240 | if (++counter < numMul) { |
---|
241 | nt = np+nm+nz; |
---|
242 | if( (nt>0) && (nt<=numSec) ) { |
---|
243 | protmul[counter] = pmltpc(np,nm,nz,nt,protb,c) ; |
---|
244 | protnorm[nt-1] += protmul[counter]; |
---|
245 | } |
---|
246 | } |
---|
247 | } |
---|
248 | } |
---|
249 | } |
---|
250 | |
---|
251 | for( i=0; i<numMul; i++ )neutmul[i] = 0.0; |
---|
252 | for( i=0; i<numSec; i++ )neutnorm[i] = 0.0; |
---|
253 | counter = -1; |
---|
254 | for (np=0; np<numSec/3; np++) { |
---|
255 | for (nm=np; nm<=(np+2); nm++) { |
---|
256 | for (nz=0; nz<numSec/3; nz++) { |
---|
257 | if (++counter < numMul) { |
---|
258 | nt = np+nm+nz; |
---|
259 | if( (nt>0) && (nt<=numSec) ) { |
---|
260 | neutmul[counter] = pmltpc(np,nm,nz,nt,neutb,c); |
---|
261 | neutnorm[nt-1] += neutmul[counter]; |
---|
262 | } |
---|
263 | } |
---|
264 | } |
---|
265 | } |
---|
266 | } |
---|
267 | |
---|
268 | for (i=0; i<numSec; i++) { |
---|
269 | if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i]; |
---|
270 | if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i]; |
---|
271 | } |
---|
272 | } // end of initialization |
---|
273 | |
---|
274 | |
---|
275 | // Initialize the first two particles |
---|
276 | // the same as beam and target |
---|
277 | pv[0] = incidentParticle; |
---|
278 | pv[1] = targetParticle; |
---|
279 | vecLen = 2; |
---|
280 | |
---|
281 | if( !inElastic ) { |
---|
282 | // quasi-elastic scattering, no pions produced |
---|
283 | if( targetCode == protonCode) { |
---|
284 | G4double cech[] = {0.33,0.27,0.29,0.31,0.27,0.18,0.13,0.10,0.09,0.07}; |
---|
285 | G4int iplab = G4int( std::min( 9.0, incidentTotalMomentum*5. ) ); |
---|
286 | if( G4UniformRand() < cech[iplab]/std::pow(atomicWeight,0.42)) { |
---|
287 | // charge exchange K+ n -> K0 p |
---|
288 | pv[0] = KaonPlus; |
---|
289 | pv[1] = Neutron; |
---|
290 | } |
---|
291 | } |
---|
292 | return; |
---|
293 | } else if (availableEnergy <= PionPlus.getMass()) { |
---|
294 | return; |
---|
295 | } |
---|
296 | |
---|
297 | // Inelastic scattering |
---|
298 | |
---|
299 | np = 0, nm = 0, nz = 0; |
---|
300 | G4double eab = availableEnergy; |
---|
301 | G4int ieab = G4int( eab*5.0 ); |
---|
302 | |
---|
303 | G4double supp[] = {0., 0.4, 0.55, 0.65, 0.75, 0.82, 0.86, 0.90, 0.94, 0.98}; |
---|
304 | if( (ieab <= 9) && (G4UniformRand() >= supp[ieab])) { |
---|
305 | // Suppress high multiplicity events at low momentum |
---|
306 | // only one additional pion will be produced |
---|
307 | G4double w0, wp, wm, wt, ran; |
---|
308 | if (targetCode == neutronCode) { |
---|
309 | // target is a neutron |
---|
310 | w0 = - sqr(1.+protb)/(2.*c*c); |
---|
311 | w0 = std::exp(w0); |
---|
312 | wm = - sqr(-1.+protb)/(2.*c*c); |
---|
313 | wm = std::exp(wm); |
---|
314 | w0 = w0/2.; |
---|
315 | wm = wm*1.5; |
---|
316 | if (G4UniformRand() < w0/(w0+wm) ) { |
---|
317 | np = 0; |
---|
318 | nm = 0; |
---|
319 | nz = 1; |
---|
320 | } else { |
---|
321 | np = 0; |
---|
322 | nm = 1; |
---|
323 | nz = 0; |
---|
324 | } |
---|
325 | |
---|
326 | } else { |
---|
327 | // target is a proton |
---|
328 | w0 = -sqr(1.+neutb)/(2.*c*c); |
---|
329 | wp = w0 = std::exp(w0); |
---|
330 | wm = -sqr(-1.+neutb)/(2.*c*c); |
---|
331 | wm = std::exp(wm); |
---|
332 | wt = w0+wp+wm; |
---|
333 | wp = w0+wp; |
---|
334 | ran = G4UniformRand(); |
---|
335 | if ( ran < w0/wt) { |
---|
336 | np = 0; |
---|
337 | nm = 0; |
---|
338 | nz = 1; |
---|
339 | } else if (ran < wp/wt) { |
---|
340 | np = 1; |
---|
341 | nm = 0; |
---|
342 | nz = 0; |
---|
343 | } else { |
---|
344 | np = 0; |
---|
345 | nm = 1; |
---|
346 | nz = 0; |
---|
347 | } |
---|
348 | } |
---|
349 | } else { |
---|
350 | // number of total particles vs. centre of mass Energy - 2*proton mass |
---|
351 | |
---|
352 | G4double aleab = std::log(availableEnergy); |
---|
353 | G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514 |
---|
354 | + aleab*(0.117712+0.0136912*aleab))) - 2.0; |
---|
355 | |
---|
356 | // Normalization constant for kno-distribution. |
---|
357 | // Calculate first the sum of all constants, check for numerical problems. |
---|
358 | G4double test, dum, anpn = 0.0; |
---|
359 | |
---|
360 | for (nt=1; nt<=numSec; nt++) { |
---|
361 | test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) ); |
---|
362 | dum = pi*nt/(2.0*n*n); |
---|
363 | if (std::fabs(dum) < 1.0) { |
---|
364 | if( test >= 1.0e-10 )anpn += dum*test; |
---|
365 | } else { |
---|
366 | anpn += dum*test; |
---|
367 | } |
---|
368 | } |
---|
369 | |
---|
370 | G4double ran = G4UniformRand(); |
---|
371 | G4double excs = 0.0; |
---|
372 | if( targetCode == protonCode ) |
---|
373 | { |
---|
374 | counter = -1; |
---|
375 | for( np=0; np<numSec/3; np++ ) |
---|
376 | { |
---|
377 | for( nm=std::max(0,np-1); nm<=(np+1); nm++ ) |
---|
378 | { |
---|
379 | for (nz=0; nz<numSec/3; nz++) { |
---|
380 | if (++counter < numMul) { |
---|
381 | nt = np+nm+nz; |
---|
382 | if ( (nt>0) && (nt<=numSec) ) { |
---|
383 | test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) ); |
---|
384 | dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n); |
---|
385 | if (std::fabs(dum) < 1.0) { |
---|
386 | if( test >= 1.0e-10 )excs += dum*test; |
---|
387 | } else { |
---|
388 | excs += dum*test; |
---|
389 | } |
---|
390 | if (ran < excs) goto outOfLoop; //-----------------------> |
---|
391 | } |
---|
392 | } |
---|
393 | } |
---|
394 | } |
---|
395 | } |
---|
396 | |
---|
397 | // 3 previous loops continued to the end |
---|
398 | inElastic = false; // quasi-elastic scattering |
---|
399 | return; |
---|
400 | } |
---|
401 | else |
---|
402 | { // target must be a neutron |
---|
403 | counter = -1; |
---|
404 | for( np=0; np<numSec/3; np++ ) |
---|
405 | { |
---|
406 | for( nm=np; nm<=(np+2); nm++ ) |
---|
407 | { |
---|
408 | for (nz=0; nz<numSec/3; nz++) { |
---|
409 | if (++counter < numMul) { |
---|
410 | nt = np+nm+nz; |
---|
411 | if ( (nt>=1) && (nt<=numSec) ) { |
---|
412 | test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) ); |
---|
413 | dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n); |
---|
414 | if (std::fabs(dum) < 1.0) { |
---|
415 | if( test >= 1.0e-10 )excs += dum*test; |
---|
416 | } else { |
---|
417 | excs += dum*test; |
---|
418 | } |
---|
419 | if (ran < excs) goto outOfLoop; // --------------------------> |
---|
420 | } |
---|
421 | } |
---|
422 | } |
---|
423 | } |
---|
424 | } |
---|
425 | // 3 previous loops continued to the end |
---|
426 | inElastic = false; // quasi-elastic scattering. |
---|
427 | return; |
---|
428 | } |
---|
429 | } |
---|
430 | outOfLoop: // <----------------------------------------------- |
---|
431 | |
---|
432 | if( targetCode == neutronCode) |
---|
433 | { |
---|
434 | if( np == nm) |
---|
435 | { |
---|
436 | } |
---|
437 | else if (np == (nm-1)) |
---|
438 | { |
---|
439 | if( G4UniformRand() < 0.5) |
---|
440 | { |
---|
441 | pv[0] = KaonPlus; |
---|
442 | } |
---|
443 | else |
---|
444 | { |
---|
445 | pv[1] = Proton; |
---|
446 | } |
---|
447 | } |
---|
448 | else |
---|
449 | { |
---|
450 | pv[0] = KaonPlus; |
---|
451 | pv[1] = Proton; |
---|
452 | } |
---|
453 | } |
---|
454 | else |
---|
455 | { |
---|
456 | if( np == nm ) |
---|
457 | { |
---|
458 | if( G4UniformRand() < 0.25) |
---|
459 | { |
---|
460 | pv[0] = KaonPlus; |
---|
461 | pv[1] = Neutron; |
---|
462 | } |
---|
463 | else |
---|
464 | { |
---|
465 | } |
---|
466 | } |
---|
467 | else if ( np == (nm+1)) |
---|
468 | { |
---|
469 | pv[1] = Neutron; |
---|
470 | } |
---|
471 | else |
---|
472 | { |
---|
473 | pv[0] = KaonPlus; |
---|
474 | } |
---|
475 | } |
---|
476 | |
---|
477 | nt = np + nm + nz; |
---|
478 | while (nt > 0) { |
---|
479 | G4double ran = G4UniformRand(); |
---|
480 | if (ran < (G4double)np/nt) { |
---|
481 | if (np > 0) { |
---|
482 | pv[vecLen++] = PionPlus; |
---|
483 | np--; |
---|
484 | } |
---|
485 | } else if ( ran < (G4double)(np+nm)/nt) { |
---|
486 | if (nm > 0) { |
---|
487 | pv[vecLen++] = PionMinus; |
---|
488 | nm--; |
---|
489 | } |
---|
490 | } else { |
---|
491 | if (nz > 0) { |
---|
492 | pv[vecLen++] = PionZero; |
---|
493 | nz--; |
---|
494 | } |
---|
495 | } |
---|
496 | nt = np + nm + nz; |
---|
497 | } |
---|
498 | |
---|
499 | if (verboseLevel > 1) { |
---|
500 | G4cout << "Particles produced: " ; |
---|
501 | G4cout << pv[0].getName() << " " ; |
---|
502 | G4cout << pv[1].getName() << " " ; |
---|
503 | for (i=2; i < vecLen; i++) G4cout << pv[i].getName() << " " ; |
---|
504 | G4cout << G4endl; |
---|
505 | } |
---|
506 | |
---|
507 | return; |
---|
508 | } |
---|
509 | |
---|
510 | |
---|
511 | void |
---|
512 | G4HEKaonZeroLongInelastic::FirstIntInCasAntiKaonZero(G4bool& inElastic, |
---|
513 | const G4double availableEnergy, |
---|
514 | G4HEVector pv[], |
---|
515 | G4int& vecLen, |
---|
516 | const G4HEVector& incidentParticle, |
---|
517 | const G4HEVector& targetParticle) |
---|
518 | |
---|
519 | // AntiKaon0 undergoes interaction with nucleon within a nucleus. Check if it is |
---|
520 | // energetically possible to produce pions/kaons. In not, assume nuclear excitation |
---|
521 | // occurs and input particle is degraded in energy. No other particles are produced. |
---|
522 | // If reaction is possible, find the correct number of pions/protons/neutrons |
---|
523 | // produced using an interpolation to multiplicity data. Replace some pions or |
---|
524 | // protons/neutrons by kaons or strange baryons according to the average |
---|
525 | // multiplicity per inelastic reaction. |
---|
526 | |
---|
527 | { |
---|
528 | static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp |
---|
529 | static const G4double expxl = -expxu; // lower bound for arg. of exp |
---|
530 | |
---|
531 | static const G4double protb = 0.7; |
---|
532 | static const G4double neutb = 0.7; |
---|
533 | static const G4double c = 1.25; |
---|
534 | |
---|
535 | static const G4int numMul = 1200; |
---|
536 | static const G4int numSec = 60; |
---|
537 | |
---|
538 | G4int neutronCode = Neutron.getCode(); |
---|
539 | G4int protonCode = Proton.getCode(); |
---|
540 | G4int kaonMinusCode = KaonMinus.getCode(); |
---|
541 | G4int kaonZeroCode = KaonZero.getCode(); |
---|
542 | G4int antiKaonZeroCode = AntiKaonZero.getCode(); |
---|
543 | |
---|
544 | G4int targetCode = targetParticle.getCode(); |
---|
545 | G4double incidentTotalMomentum = incidentParticle.getTotalMomentum(); |
---|
546 | |
---|
547 | static G4bool first = true; |
---|
548 | static G4double protmul[numMul], protnorm[numSec]; // proton constants |
---|
549 | static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants |
---|
550 | |
---|
551 | // misc. local variables |
---|
552 | // np = number of pi+, nm = number of pi-, nz = number of pi0 |
---|
553 | |
---|
554 | G4int i, counter, nt, np, nm, nz; |
---|
555 | |
---|
556 | if (first) { |
---|
557 | // compute normalization constants, this will only be done once |
---|
558 | first = false; |
---|
559 | for( i=0; i<numMul; i++ )protmul[i] = 0.0; |
---|
560 | for( i=0; i<numSec; i++ )protnorm[i] = 0.0; |
---|
561 | counter = -1; |
---|
562 | for(np=0; np<(numSec/3); np++) { |
---|
563 | for(nm=std::max(0,np-2); nm<=np; nm++) { |
---|
564 | for(nz=0; nz<numSec/3; nz++) { |
---|
565 | if(++counter < numMul) { |
---|
566 | nt = np+nm+nz; |
---|
567 | if( (nt>0) && (nt<=numSec) ) { |
---|
568 | protmul[counter] = pmltpc(np,nm,nz,nt,protb,c) ; |
---|
569 | protnorm[nt-1] += protmul[counter]; |
---|
570 | } |
---|
571 | } |
---|
572 | } |
---|
573 | } |
---|
574 | } |
---|
575 | |
---|
576 | for( i=0; i<numMul; i++ )neutmul[i] = 0.0; |
---|
577 | for( i=0; i<numSec; i++ )neutnorm[i] = 0.0; |
---|
578 | counter = -1; |
---|
579 | for(np=0; np<numSec/3; np++) { |
---|
580 | for(nm=std::max(0,np-1); nm<=(np+1); nm++) { |
---|
581 | for(nz=0; nz<numSec/3; nz++) { |
---|
582 | if(++counter < numMul) { |
---|
583 | nt = np+nm+nz; |
---|
584 | if( (nt>0) && (nt<=numSec) ) { |
---|
585 | neutmul[counter] = pmltpc(np,nm,nz,nt,neutb,c); |
---|
586 | neutnorm[nt-1] += neutmul[counter]; |
---|
587 | } |
---|
588 | } |
---|
589 | } |
---|
590 | } |
---|
591 | } |
---|
592 | |
---|
593 | for(i=0; i<numSec; i++) { |
---|
594 | if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i]; |
---|
595 | if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i]; |
---|
596 | } |
---|
597 | } // end of initialization |
---|
598 | |
---|
599 | // initialize the first two particles |
---|
600 | // the same as beam and target |
---|
601 | pv[0] = incidentParticle; |
---|
602 | pv[1] = targetParticle; |
---|
603 | vecLen = 2; |
---|
604 | |
---|
605 | if (!inElastic || (availableEnergy <= PionPlus.getMass())) |
---|
606 | return; |
---|
607 | |
---|
608 | // Inelastic scattering |
---|
609 | |
---|
610 | np = 0, nm = 0, nz = 0; |
---|
611 | G4double cech[] = { 1., 1., 1., 0.70, 0.60, 0.55, 0.35, 0.25, 0.18, 0.15}; |
---|
612 | G4int iplab = G4int( incidentTotalMomentum*5.); |
---|
613 | if( (iplab < 10) && (G4UniformRand() < cech[iplab]) ) { |
---|
614 | G4int iplab = std::min(19, G4int( incidentTotalMomentum*5.)); |
---|
615 | G4double cnk0[] = {0.17, 0.18, 0.17, 0.24, 0.26, 0.20, 0.22, 0.21, 0.34, 0.45, |
---|
616 | 0.58, 0.55, 0.36, 0.29, 0.29, 0.32, 0.32, 0.33, 0.33, 0.33}; |
---|
617 | if(G4UniformRand() < cnk0[iplab]) { |
---|
618 | if(targetCode == protonCode) { |
---|
619 | return; |
---|
620 | } else { |
---|
621 | pv[0] = KaonMinus; |
---|
622 | pv[1] = Proton; |
---|
623 | return; |
---|
624 | } |
---|
625 | } |
---|
626 | |
---|
627 | G4double ran = G4UniformRand(); |
---|
628 | if(targetCode == protonCode) { |
---|
629 | |
---|
630 | // target is a proton |
---|
631 | if( ran < 0.25 ) { |
---|
632 | ; |
---|
633 | } else if (ran < 0.50) { |
---|
634 | pv[0] = PionPlus; |
---|
635 | pv[1] = SigmaZero; |
---|
636 | } else if (ran < 0.75) { |
---|
637 | ; |
---|
638 | } else { |
---|
639 | pv[0] = PionPlus; |
---|
640 | pv[1] = Lambda; |
---|
641 | } |
---|
642 | } else { |
---|
643 | |
---|
644 | // target is a neutron |
---|
645 | if( ran < 0.25 ) { |
---|
646 | pv[0] = PionMinus; |
---|
647 | pv[1] = SigmaPlus; |
---|
648 | } else if (ran < 0.50) { |
---|
649 | pv[0] = PionZero; |
---|
650 | pv[1] = SigmaZero; |
---|
651 | } else if (ran < 0.75) { |
---|
652 | pv[0] = PionPlus; |
---|
653 | pv[1] = SigmaMinus; |
---|
654 | } else { |
---|
655 | pv[0] = PionZero; |
---|
656 | pv[1] = Lambda; |
---|
657 | } |
---|
658 | } |
---|
659 | return; |
---|
660 | |
---|
661 | } else { |
---|
662 | // number of total particles vs. centre of mass Energy - 2*proton mass |
---|
663 | |
---|
664 | G4double aleab = std::log(availableEnergy); |
---|
665 | G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514 |
---|
666 | + aleab*(0.117712+0.0136912*aleab))) - 2.0; |
---|
667 | |
---|
668 | // Normalization constant for kno-distribution. |
---|
669 | // Calculate first the sum of all constants, check for numerical problems. |
---|
670 | G4double test, dum, anpn = 0.0; |
---|
671 | |
---|
672 | for (nt=1; nt<=numSec; nt++) { |
---|
673 | test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) ); |
---|
674 | dum = pi*nt/(2.0*n*n); |
---|
675 | if (std::fabs(dum) < 1.0) { |
---|
676 | if( test >= 1.0e-10 )anpn += dum*test; |
---|
677 | } else { |
---|
678 | anpn += dum*test; |
---|
679 | } |
---|
680 | } |
---|
681 | |
---|
682 | G4double ran = G4UniformRand(); |
---|
683 | G4double excs = 0.0; |
---|
684 | if (targetCode == protonCode) { |
---|
685 | counter = -1; |
---|
686 | for (np=0; np<numSec/3; np++) { |
---|
687 | for (nm=std::max(0,np-2); nm<=np; nm++) { |
---|
688 | for (nz=0; nz<numSec/3; nz++) { |
---|
689 | if (++counter < numMul) { |
---|
690 | nt = np+nm+nz; |
---|
691 | if( (nt>0) && (nt<=numSec) ) { |
---|
692 | test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) ); |
---|
693 | dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n); |
---|
694 | |
---|
695 | if (std::fabs(dum) < 1.0) { |
---|
696 | if( test >= 1.0e-10 )excs += dum*test; |
---|
697 | } else { |
---|
698 | excs += dum*test; |
---|
699 | } |
---|
700 | |
---|
701 | if (ran < excs) goto outOfLoop; //-----------------------> |
---|
702 | } |
---|
703 | } |
---|
704 | } |
---|
705 | } |
---|
706 | } |
---|
707 | // 3 previous loops continued to the end |
---|
708 | inElastic = false; // quasi-elastic scattering |
---|
709 | return; |
---|
710 | |
---|
711 | } else { // target must be a neutron |
---|
712 | counter = -1; |
---|
713 | for (np=0; np<numSec/3; np++) { |
---|
714 | for (nm=std::max(0,np-1); nm<=(np+1); nm++) { |
---|
715 | for (nz=0; nz<numSec/3; nz++) { |
---|
716 | if (++counter < numMul) { |
---|
717 | nt = np+nm+nz; |
---|
718 | if( (nt>=1) && (nt<=numSec) ) { |
---|
719 | test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) ); |
---|
720 | dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n); |
---|
721 | |
---|
722 | if (std::fabs(dum) < 1.0) { |
---|
723 | if( test >= 1.0e-10 )excs += dum*test; |
---|
724 | } else { |
---|
725 | excs += dum*test; |
---|
726 | } |
---|
727 | |
---|
728 | if (ran < excs) goto outOfLoop; // --------------------------> |
---|
729 | } |
---|
730 | } |
---|
731 | } |
---|
732 | } |
---|
733 | } |
---|
734 | // 3 previous loops continued to the end |
---|
735 | inElastic = false; // quasi-elastic scattering. |
---|
736 | return; |
---|
737 | } |
---|
738 | } |
---|
739 | outOfLoop: // <------------------------------------------------------------------------ |
---|
740 | |
---|
741 | if( targetCode == protonCode) |
---|
742 | { |
---|
743 | if( np == nm) |
---|
744 | { |
---|
745 | } |
---|
746 | else if (np == (1+nm)) |
---|
747 | { |
---|
748 | if( G4UniformRand() < 0.5) |
---|
749 | { |
---|
750 | pv[0] = KaonMinus; |
---|
751 | } |
---|
752 | else |
---|
753 | { |
---|
754 | pv[1] = Neutron; |
---|
755 | } |
---|
756 | } |
---|
757 | else |
---|
758 | { |
---|
759 | pv[0] = KaonMinus; |
---|
760 | pv[1] = Neutron; |
---|
761 | } |
---|
762 | } |
---|
763 | else |
---|
764 | { |
---|
765 | if( np == nm) |
---|
766 | { |
---|
767 | if( G4UniformRand() < 0.75) |
---|
768 | { |
---|
769 | } |
---|
770 | else |
---|
771 | { |
---|
772 | pv[0] = KaonMinus; |
---|
773 | pv[1] = Proton; |
---|
774 | } |
---|
775 | } |
---|
776 | else if ( np == (1+nm)) |
---|
777 | { |
---|
778 | pv[0] = KaonMinus; |
---|
779 | } |
---|
780 | else |
---|
781 | { |
---|
782 | pv[1] = Proton; |
---|
783 | } |
---|
784 | } |
---|
785 | |
---|
786 | |
---|
787 | if( G4UniformRand() < 0.5 ) |
---|
788 | { |
---|
789 | if( ( (pv[0].getCode() == kaonMinusCode) |
---|
790 | && (pv[1].getCode() == neutronCode) ) |
---|
791 | || ( (pv[0].getCode() == kaonZeroCode) |
---|
792 | && (pv[1].getCode() == protonCode) ) |
---|
793 | || ( (pv[0].getCode() == antiKaonZeroCode) |
---|
794 | && (pv[1].getCode() == protonCode) ) ) |
---|
795 | { |
---|
796 | G4double ran = G4UniformRand(); |
---|
797 | if( pv[1].getCode() == protonCode) |
---|
798 | { |
---|
799 | if(ran < 0.68) |
---|
800 | { |
---|
801 | pv[0] = PionPlus; |
---|
802 | pv[1] = Lambda; |
---|
803 | } |
---|
804 | else if (ran < 0.84) |
---|
805 | { |
---|
806 | pv[0] = PionZero; |
---|
807 | pv[1] = SigmaPlus; |
---|
808 | } |
---|
809 | else |
---|
810 | { |
---|
811 | pv[0] = PionPlus; |
---|
812 | pv[1] = SigmaZero; |
---|
813 | } |
---|
814 | } |
---|
815 | else |
---|
816 | { |
---|
817 | if(ran < 0.68) |
---|
818 | { |
---|
819 | pv[0] = PionMinus; |
---|
820 | pv[1] = Lambda; |
---|
821 | } |
---|
822 | else if (ran < 0.84) |
---|
823 | { |
---|
824 | pv[0] = PionMinus; |
---|
825 | pv[1] = SigmaZero; |
---|
826 | } |
---|
827 | else |
---|
828 | { |
---|
829 | pv[0] = PionZero; |
---|
830 | pv[1] = SigmaMinus; |
---|
831 | } |
---|
832 | } |
---|
833 | } |
---|
834 | else |
---|
835 | { |
---|
836 | G4double ran = G4UniformRand(); |
---|
837 | if (ran < 0.67) |
---|
838 | { |
---|
839 | pv[0] = PionZero; |
---|
840 | pv[1] = Lambda; |
---|
841 | } |
---|
842 | else if (ran < 0.78) |
---|
843 | { |
---|
844 | pv[0] = PionMinus; |
---|
845 | pv[1] = SigmaPlus; |
---|
846 | } |
---|
847 | else if (ran < 0.89) |
---|
848 | { |
---|
849 | pv[0] = PionZero; |
---|
850 | pv[1] = SigmaZero; |
---|
851 | } |
---|
852 | else |
---|
853 | { |
---|
854 | pv[0] = PionPlus; |
---|
855 | pv[1] = SigmaMinus; |
---|
856 | } |
---|
857 | } |
---|
858 | } |
---|
859 | |
---|
860 | nt = np + nm + nz; |
---|
861 | while ( nt > 0) { |
---|
862 | G4double ran = G4UniformRand(); |
---|
863 | if ( ran < (G4double)np/nt) { |
---|
864 | if( np > 0 ) { |
---|
865 | pv[vecLen++] = PionPlus; |
---|
866 | np--; |
---|
867 | } |
---|
868 | } else if (ran < (G4double)(np+nm)/nt) { |
---|
869 | if( nm > 0 ) { |
---|
870 | pv[vecLen++] = PionMinus; |
---|
871 | nm--; |
---|
872 | } |
---|
873 | } else { |
---|
874 | if( nz > 0 ) { |
---|
875 | pv[vecLen++] = PionZero; |
---|
876 | nz--; |
---|
877 | } |
---|
878 | } |
---|
879 | nt = np + nm + nz; |
---|
880 | } |
---|
881 | |
---|
882 | if (verboseLevel > 1) { |
---|
883 | G4cout << "Particles produced: " ; |
---|
884 | G4cout << pv[0].getName() << " " ; |
---|
885 | G4cout << pv[1].getName() << " " ; |
---|
886 | for (i=2; i < vecLen; i++) G4cout << pv[i].getName() << " " ; |
---|
887 | G4cout << G4endl; |
---|
888 | } |
---|
889 | |
---|
890 | return; |
---|
891 | } |
---|