source: trunk/source/processes/hadronic/models/high_energy/src/G4HEKaonZeroShortInelastic.cc@ 1350

Last change on this file since 1350 was 1347, checked in by garnier, 15 years ago

geant4 tag 9.4

File size: 28.7 KB
Line 
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: G4HEKaonZeroShortInelastic.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// 21 January 2010
41
42
43#include "G4HEKaonZeroShortInelastic.hh"
44
45G4HadFinalState*
46G4HEKaonZeroShortInelastic::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 << "GHEKaonZeroShortInelastic: incident energy < 1 GeV" << G4endl;
63
64 if(verboseLevel > 1) {
65 G4cout << "G4HEKaonZeroShortInelastic::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
189void
190G4HEKaonZeroShortInelastic::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
511void
512G4HEKaonZeroShortInelastic::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 static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp
528 static const G4double expxl = -expxu; // lower bound for arg. of exp
529
530 static const G4double protb = 0.7;
531 static const G4double neutb = 0.7;
532 static const G4double c = 1.25;
533
534 static const G4int numMul = 1200;
535 static const G4int numSec = 60;
536
537 G4int neutronCode = Neutron.getCode();
538 G4int protonCode = Proton.getCode();
539 G4int kaonMinusCode = KaonMinus.getCode();
540 G4int kaonZeroCode = KaonZero.getCode();
541 G4int antiKaonZeroCode = AntiKaonZero.getCode();
542
543 G4int targetCode = targetParticle.getCode();
544 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
545
546 static G4bool first = true;
547 static G4double protmul[numMul], protnorm[numSec]; // proton constants
548 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
549
550 // misc. local variables
551 // np = number of pi+, nm = number of pi-, nz = number of pi0
552
553 G4int i, counter, nt, np, nm, nz;
554
555 if (first) {
556 // compute normalization constants, this will only be done once
557 first = false;
558 for( i=0; i<numMul; i++ )protmul[i] = 0.0;
559 for( i=0; i<numSec; i++ )protnorm[i] = 0.0;
560 counter = -1;
561 for(np=0; np<(numSec/3); np++) {
562 for(nm=std::max(0,np-2); nm<=np; nm++) {
563 for(nz=0; nz<numSec/3; nz++) {
564 if(++counter < numMul) {
565 nt = np+nm+nz;
566 if( (nt>0) && (nt<=numSec) ) {
567 protmul[counter] = pmltpc(np,nm,nz,nt,protb,c) ;
568 protnorm[nt-1] += protmul[counter];
569 }
570 }
571 }
572 }
573 }
574
575 for( i=0; i<numMul; i++ )neutmul[i] = 0.0;
576 for( i=0; i<numSec; i++ )neutnorm[i] = 0.0;
577 counter = -1;
578 for(np=0; np<numSec/3; np++) {
579 for(nm=std::max(0,np-1); nm<=(np+1); nm++) {
580 for(nz=0; nz<numSec/3; nz++) {
581 if(++counter < numMul) {
582 nt = np+nm+nz;
583 if( (nt>0) && (nt<=numSec) ) {
584 neutmul[counter] = pmltpc(np,nm,nz,nt,neutb,c);
585 neutnorm[nt-1] += neutmul[counter];
586 }
587 }
588 }
589 }
590 }
591
592 for(i=0; i<numSec; i++) {
593 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
594 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
595 }
596 } // end of initialization
597
598 // initialize the first two particles
599 // the same as beam and target
600 pv[0] = incidentParticle;
601 pv[1] = targetParticle;
602 vecLen = 2;
603
604 if (!inElastic || (availableEnergy <= PionPlus.getMass()))
605 return;
606
607 // Inelastic scattering
608
609 np = 0, nm = 0, nz = 0;
610 G4double cech[] = { 1., 1., 1., 0.70, 0.60, 0.55, 0.35, 0.25, 0.18, 0.15};
611 G4int iplab = G4int( incidentTotalMomentum*5.);
612 if( (iplab < 10) && (G4UniformRand() < cech[iplab]) ) {
613 G4int iplab = std::min(19, G4int( incidentTotalMomentum*5.));
614 G4double cnk0[] = {0.17, 0.18, 0.17, 0.24, 0.26, 0.20, 0.22, 0.21, 0.34, 0.45,
615 0.58, 0.55, 0.36, 0.29, 0.29, 0.32, 0.32, 0.33, 0.33, 0.33};
616 if(G4UniformRand() < cnk0[iplab]) {
617 if(targetCode == protonCode) {
618 return;
619 } else {
620 pv[0] = KaonMinus;
621 pv[1] = Proton;
622 return;
623 }
624 }
625
626 G4double ran = G4UniformRand();
627 if(targetCode == protonCode) {
628
629 // target is a proton
630 if( ran < 0.25 ) {
631 ;
632 } else if (ran < 0.50) {
633 pv[0] = PionPlus;
634 pv[1] = SigmaZero;
635 } else if (ran < 0.75) {
636 ;
637 } else {
638 pv[0] = PionPlus;
639 pv[1] = Lambda;
640 }
641 } else {
642
643 // target is a neutron
644 if( ran < 0.25 ) {
645 pv[0] = PionMinus;
646 pv[1] = SigmaPlus;
647 } else if (ran < 0.50) {
648 pv[0] = PionZero;
649 pv[1] = SigmaZero;
650 } else if (ran < 0.75) {
651 pv[0] = PionPlus;
652 pv[1] = SigmaMinus;
653 } else {
654 pv[0] = PionZero;
655 pv[1] = Lambda;
656 }
657 }
658 return;
659
660 } else {
661 // number of total particles vs. centre of mass Energy - 2*proton mass
662
663 G4double aleab = std::log(availableEnergy);
664 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
665 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
666
667 // Normalization constant for kno-distribution.
668 // Calculate first the sum of all constants, check for numerical problems.
669 G4double test, dum, anpn = 0.0;
670
671 for (nt=1; nt<=numSec; nt++) {
672 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
673 dum = pi*nt/(2.0*n*n);
674 if (std::fabs(dum) < 1.0) {
675 if( test >= 1.0e-10 )anpn += dum*test;
676 } else {
677 anpn += dum*test;
678 }
679 }
680
681 G4double ran = G4UniformRand();
682 G4double excs = 0.0;
683 if (targetCode == protonCode) {
684 counter = -1;
685 for (np=0; np<numSec/3; np++) {
686 for (nm=std::max(0,np-2); nm<=np; nm++) {
687 for (nz=0; nz<numSec/3; nz++) {
688 if (++counter < numMul) {
689 nt = np+nm+nz;
690 if( (nt>0) && (nt<=numSec) ) {
691 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
692 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
693
694 if (std::fabs(dum) < 1.0) {
695 if( test >= 1.0e-10 )excs += dum*test;
696 } else {
697 excs += dum*test;
698 }
699
700 if (ran < excs) goto outOfLoop; //----------------------->
701 }
702 }
703 }
704 }
705 }
706 // 3 previous loops continued to the end
707 inElastic = false; // quasi-elastic scattering
708 return;
709
710 } else { // target must be a neutron
711 counter = -1;
712 for (np=0; np<numSec/3; np++) {
713 for (nm=std::max(0,np-1); nm<=(np+1); nm++) {
714 for (nz=0; nz<numSec/3; nz++) {
715 if (++counter < numMul) {
716 nt = np+nm+nz;
717 if( (nt>=1) && (nt<=numSec) ) {
718 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
719 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
720
721 if (std::fabs(dum) < 1.0) {
722 if( test >= 1.0e-10 )excs += dum*test;
723 } else {
724 excs += dum*test;
725 }
726
727 if (ran < excs) goto outOfLoop; // -------------------------->
728 }
729 }
730 }
731 }
732 }
733 // 3 previous loops continued to the end
734 inElastic = false; // quasi-elastic scattering.
735 return;
736 }
737 }
738 outOfLoop: // <------------------------------------------------------------------------
739
740 if( targetCode == protonCode)
741 {
742 if( np == nm)
743 {
744 }
745 else if (np == (1+nm))
746 {
747 if( G4UniformRand() < 0.5)
748 {
749 pv[0] = KaonMinus;
750 }
751 else
752 {
753 pv[1] = Neutron;
754 }
755 }
756 else
757 {
758 pv[0] = KaonMinus;
759 pv[1] = Neutron;
760 }
761 }
762 else
763 {
764 if( np == nm)
765 {
766 if( G4UniformRand() < 0.75)
767 {
768 }
769 else
770 {
771 pv[0] = KaonMinus;
772 pv[1] = Proton;
773 }
774 }
775 else if ( np == (1+nm))
776 {
777 pv[0] = KaonMinus;
778 }
779 else
780 {
781 pv[1] = Proton;
782 }
783 }
784
785
786 if( G4UniformRand() < 0.5 )
787 {
788 if( ( (pv[0].getCode() == kaonMinusCode)
789 && (pv[1].getCode() == neutronCode) )
790 || ( (pv[0].getCode() == kaonZeroCode)
791 && (pv[1].getCode() == protonCode) )
792 || ( (pv[0].getCode() == antiKaonZeroCode)
793 && (pv[1].getCode() == protonCode) ) )
794 {
795 G4double ran = G4UniformRand();
796 if( pv[1].getCode() == protonCode)
797 {
798 if(ran < 0.68)
799 {
800 pv[0] = PionPlus;
801 pv[1] = Lambda;
802 }
803 else if (ran < 0.84)
804 {
805 pv[0] = PionZero;
806 pv[1] = SigmaPlus;
807 }
808 else
809 {
810 pv[0] = PionPlus;
811 pv[1] = SigmaZero;
812 }
813 }
814 else
815 {
816 if(ran < 0.68)
817 {
818 pv[0] = PionMinus;
819 pv[1] = Lambda;
820 }
821 else if (ran < 0.84)
822 {
823 pv[0] = PionMinus;
824 pv[1] = SigmaZero;
825 }
826 else
827 {
828 pv[0] = PionZero;
829 pv[1] = SigmaMinus;
830 }
831 }
832 }
833 else
834 {
835 G4double ran = G4UniformRand();
836 if (ran < 0.67)
837 {
838 pv[0] = PionZero;
839 pv[1] = Lambda;
840 }
841 else if (ran < 0.78)
842 {
843 pv[0] = PionMinus;
844 pv[1] = SigmaPlus;
845 }
846 else if (ran < 0.89)
847 {
848 pv[0] = PionZero;
849 pv[1] = SigmaZero;
850 }
851 else
852 {
853 pv[0] = PionPlus;
854 pv[1] = SigmaMinus;
855 }
856 }
857 }
858
859 nt = np + nm + nz;
860 while ( nt > 0) {
861 G4double ran = G4UniformRand();
862 if ( ran < (G4double)np/nt) {
863 if( np > 0 ) {
864 pv[vecLen++] = PionPlus;
865 np--;
866 }
867 } else if (ran < (G4double)(np+nm)/nt) {
868 if( nm > 0 ) {
869 pv[vecLen++] = PionMinus;
870 nm--;
871 }
872 } else {
873 if( nz > 0 ) {
874 pv[vecLen++] = PionZero;
875 nz--;
876 }
877 }
878 nt = np + nm + nz;
879 }
880
881 if (verboseLevel > 1) {
882 G4cout << "Particles produced: " ;
883 G4cout << pv[0].getName() << " " ;
884 G4cout << pv[1].getName() << " " ;
885 for (i=2; i < vecLen; i++) G4cout << pv[i].getName() << " " ;
886 G4cout << G4endl;
887 }
888
889 return;
890}
Note: See TracBrowser for help on using the repository browser.