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

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

geant4 tag 9.4

File size: 18.8 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: G4HENeutronInelastic.cc,v 1.17 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// First work done by J.L.Chuma and F.W.Jones, TRIUMF, June 96.
39// H. Fesefeldt, RWTH-Aachen, 23-October-1996
40// Last modified: 29-July-1998
41
42#include "G4HENeutronInelastic.hh"
43
44G4HadFinalState*
45G4HENeutronInelastic::ApplyYourself(const G4HadProjectile& aTrack,
46 G4Nucleus& targetNucleus)
47{
48 G4HEVector* pv = new G4HEVector[MAXPART];
49 const G4HadProjectile* aParticle = &aTrack;
50 const G4double A = targetNucleus.GetN();
51 const G4double Z = targetNucleus.GetZ();
52 G4HEVector incidentParticle(aParticle);
53
54 G4double atomicNumber = Z;
55 G4double atomicWeight = A;
56
57 G4int incidentCode = incidentParticle.getCode();
58 G4double incidentMass = incidentParticle.getMass();
59 G4double incidentTotalEnergy = incidentParticle.getEnergy();
60 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
61 G4double incidentKineticEnergy = incidentTotalEnergy - incidentMass;
62
63 if (incidentKineticEnergy < 1.)
64 G4cout << "GHENeutronInelastic: incident energy < 1 GeV" << G4endl;
65
66 if (verboseLevel > 1) {
67 G4cout << "G4HENeutronInelastic::ApplyYourself" << G4endl;
68 G4cout << "incident particle " << incidentParticle.getName()
69 << "mass " << incidentMass
70 << "kinetic energy " << incidentKineticEnergy
71 << G4endl;
72 G4cout << "target material with (A,Z) = ("
73 << atomicWeight << "," << atomicNumber << ")" << G4endl;
74 }
75
76 G4double inelasticity = NuclearInelasticity(incidentKineticEnergy,
77 atomicWeight, atomicNumber);
78 if (verboseLevel > 1)
79 G4cout << "nuclear inelasticity = " << inelasticity << G4endl;
80
81 incidentKineticEnergy -= inelasticity;
82
83 G4double excitationEnergyGNP = 0.;
84 G4double excitationEnergyDTA = 0.;
85
86 G4double excitation = NuclearExcitation(incidentKineticEnergy,
87 atomicWeight, atomicNumber,
88 excitationEnergyGNP,
89 excitationEnergyDTA);
90 if (verboseLevel > 1)
91 G4cout << "nuclear excitation = " << excitation << excitationEnergyGNP
92 << excitationEnergyDTA << G4endl;
93
94 incidentKineticEnergy -= excitation;
95 incidentTotalEnergy = incidentKineticEnergy + incidentMass;
96 incidentTotalMomentum = std::sqrt( (incidentTotalEnergy-incidentMass)
97 *(incidentTotalEnergy+incidentMass));
98
99 G4HEVector targetParticle;
100 if (G4UniformRand() < atomicNumber/atomicWeight) {
101 targetParticle.setDefinition("Proton");
102 } else {
103 targetParticle.setDefinition("Neutron");
104 }
105
106 G4double targetMass = targetParticle.getMass();
107 G4double centerOfMassEnergy = std::sqrt(incidentMass*incidentMass
108 + targetMass*targetMass
109 + 2.0*targetMass*incidentTotalEnergy);
110 G4double availableEnergy = centerOfMassEnergy - targetMass - incidentMass;
111
112 // In the original Gheisha code the inElastic flag was set as follows
113 // G4bool inElastic = InElasticCrossSectionInFirstInt
114 // (availableEnergy, incidentCode, incidentTotalMomentum);
115 // For unknown reasons, it has been replaced by the following code in Geant???
116 //
117 G4bool inElastic = true;
118 // if (G4UniformRand() < elasticCrossSection/totalCrossSection) inElastic = false;
119
120 vecLength = 0;
121
122 if (verboseLevel > 1)
123 G4cout << "ApplyYourself: CallFirstIntInCascade for particle "
124 << incidentCode << G4endl;
125
126 G4bool successful = false;
127
128 FirstIntInCasNeutron(inElastic, availableEnergy, pv, vecLength,
129 incidentParticle, targetParticle, atomicWeight);
130
131 if (verboseLevel > 1)
132 G4cout << "ApplyYourself::StrangeParticlePairProduction" << G4endl;
133
134 if ((vecLength > 0) && (availableEnergy > 1.))
135 StrangeParticlePairProduction(availableEnergy, centerOfMassEnergy,
136 pv, vecLength,
137 incidentParticle, targetParticle);
138
139 HighEnergyCascading(successful, pv, vecLength,
140 excitationEnergyGNP, excitationEnergyDTA,
141 incidentParticle, targetParticle,
142 atomicWeight, atomicNumber);
143 if (!successful)
144 HighEnergyClusterProduction(successful, pv, vecLength,
145 excitationEnergyGNP, excitationEnergyDTA,
146 incidentParticle, targetParticle,
147 atomicWeight, atomicNumber);
148 if (!successful)
149 MediumEnergyCascading(successful, pv, vecLength,
150 excitationEnergyGNP, excitationEnergyDTA,
151 incidentParticle, targetParticle,
152 atomicWeight, atomicNumber);
153
154 if (!successful)
155 MediumEnergyClusterProduction(successful, pv, vecLength,
156 excitationEnergyGNP, excitationEnergyDTA,
157 incidentParticle, targetParticle,
158 atomicWeight, atomicNumber);
159 if (!successful)
160 QuasiElasticScattering(successful, pv, vecLength,
161 excitationEnergyGNP, excitationEnergyDTA,
162 incidentParticle, targetParticle,
163 atomicWeight, atomicNumber);
164 if (!successful)
165 ElasticScattering(successful, pv, vecLength,
166 incidentParticle,
167 atomicWeight, atomicNumber);
168
169 if (!successful)
170 G4cout << "GHEInelasticInteraction::ApplyYourself fails to produce final state particles"
171 << G4endl;
172
173 FillParticleChange(pv, vecLength);
174 delete [] pv;
175 theParticleChange.SetStatusChange(stopAndKill);
176 return &theParticleChange;
177}
178
179
180void
181G4HENeutronInelastic::FirstIntInCasNeutron(G4bool& inElastic,
182 const G4double availableEnergy,
183 G4HEVector pv[],
184 G4int& vecLen,
185 const G4HEVector& incidentParticle,
186 const G4HEVector& targetParticle,
187 const G4double atomicWeight)
188
189// Neutron undergoes interaction with nucleon within a nucleus. Check if it is
190// energetically possible to produce pions/kaons. In not, assume nuclear excitation
191// occurs and input particle is degraded in energy. No other particles are produced.
192// If reaction is possible, find the correct number of pions/protons/neutrons
193// produced using an interpolation to multiplicity data. Replace some pions or
194// protons/neutrons by kaons or strange baryons according to the average
195// multiplicity per inelastic reaction.
196{
197 static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp
198 static const G4double expxl = -expxu; // lower bound for arg. of exp
199
200 static const G4double protb = 0.35;
201 static const G4double neutb = 0.35;
202 static const G4double c = 1.25;
203
204 static const G4int numMul = 1200;
205 static const G4int numSec = 60;
206
207 G4int neutronCode = Neutron.getCode();
208 G4int protonCode = Proton.getCode();
209
210 G4int targetCode = targetParticle.getCode();
211 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
212
213 static G4bool first = true;
214 static G4double protmul[numMul], protnorm[numSec]; // proton constants
215 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
216
217 // misc. local variables
218 // np = number of pi+, nm = number of pi-, nz = number of pi0
219
220 G4int i, counter, nt, np, nm, nz;
221
222 if( first )
223 { // compute normalization constants, this will only be done once
224 first = false;
225 for( i=0; i<numMul; i++ )protmul[i] = 0.0;
226 for( i=0; i<numSec; i++ )protnorm[i] = 0.0;
227 counter = -1;
228 for( np=0; np<(numSec/3); np++ )
229 {
230 for( nm=std::max(0,np-1); nm<=(np+1); nm++ )
231 {
232 for( nz=0; nz<numSec/3; nz++ )
233 {
234 if( ++counter < numMul )
235 {
236 nt = np+nm+nz;
237 if( (nt>0) && (nt<=numSec) )
238 {
239 protmul[counter] = pmltpc(np,nm,nz,nt,protb,c)
240 /(Factorial(1-np+nm)*Factorial(1+np-nm)) ;
241 protnorm[nt-1] += protmul[counter];
242 }
243 }
244 }
245 }
246 }
247 for( i=0; i<numMul; i++ )neutmul[i] = 0.0;
248 for( i=0; i<numSec; i++ )neutnorm[i] = 0.0;
249 counter = -1;
250 for( np=0; np<numSec/3; np++ )
251 {
252 for( nm=np; nm<=(np+2); nm++ )
253 {
254 for( nz=0; nz<numSec/3; nz++ )
255 {
256 if( ++counter < numMul )
257 {
258 nt = np+nm+nz;
259 if( (nt>0) && (nt<=numSec) )
260 {
261 neutmul[counter] = pmltpc(np,nm,nz,nt,neutb,c)
262 /(Factorial(-np+nm)*Factorial(2+np-nm));
263 neutnorm[nt-1] += neutmul[counter];
264 }
265 }
266 }
267 }
268 }
269 for( i=0; i<numSec; i++ )
270 {
271 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
272 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
273 }
274 } // end of initialization
275
276
277 // initialize the first two places
278 // the same as beam and target
279 pv[0] = incidentParticle;
280 pv[1] = targetParticle;
281 vecLen = 2;
282
283 if( !inElastic )
284 { // quasi-elastic scattering, no pions produced
285 if( targetCode == protonCode )
286 {
287 G4double cech[] = {0.50, 0.45, 0.40, 0.35, 0.30, 0.25, 0.06, 0.04, 0.005, 0.};
288 G4int iplab = G4int( std::min( 9.0, incidentTotalMomentum*2.5) );
289 if( G4UniformRand() < cech[iplab]/std::pow(atomicWeight,0.42) )
290 { // charge exchange n p -> p n
291 pv[0] = Proton;
292 pv[1] = Neutron;
293 }
294 }
295 return;
296 }
297 else if (availableEnergy <= PionPlus.getMass())
298 return;
299
300 // inelastic scattering
301
302 np = 0, nm = 0, nz = 0;
303 G4double eab = availableEnergy;
304 G4int ieab = G4int( eab*5.0 );
305
306 G4double supp[] = {0., 0.4, 0.55, 0.65, 0.75, 0.82, 0.86, 0.90, 0.94, 0.98};
307 if( (ieab <= 9) && (G4UniformRand() >= supp[ieab]) )
308 {
309 // suppress high multiplicity events at low momentum
310 // only one additional pion will be produced
311 G4double w0, wp, wm, wt, ran;
312 if( targetCode == neutronCode ) // target is a neutron
313 {
314 w0 = - sqr(1.+neutb)/(2.*c*c);
315 wm = w0 = std::exp(w0);
316 w0 = w0/2.;
317 if( G4UniformRand() < w0/(w0+wm) )
318 { np = 0; nm = 0; nz = 1; }
319 else
320 { np = 0; nm = 1; nz = 0; }
321 }
322 else
323 { // target is a proton
324 w0 = -sqr(1.+protb)/(2.*c*c);
325 w0 = std::exp(w0);
326 wp = w0/2.;
327 wm = -sqr(-1.+protb)/(2.*c*c);
328 wm = std::exp(wm)/2.;
329 wt = w0+wp+wm;
330 wp = w0+wp;
331 ran = G4UniformRand();
332 if( ran < w0/wt)
333 { np = 0; nm = 0; nz = 1; }
334 else if( ran < wp/wt)
335 { np = 1; nm = 0; nz = 0; }
336 else
337 { np = 0; nm = 1; nz = 0; }
338 }
339 }
340 else
341 {
342 // number of total particles vs. centre of mass Energy - 2*proton mass
343
344 G4double aleab = std::log(availableEnergy);
345 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
346 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
347
348 // normalization constant for kno-distribution.
349 // calculate first the sum of all constants, check for numerical problems.
350 G4double test, dum, anpn = 0.0;
351
352 for (nt=1; nt<=numSec; nt++) {
353 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
354 dum = pi*nt/(2.0*n*n);
355 if (std::fabs(dum) < 1.0) {
356 if( test >= 1.0e-10 )anpn += dum*test;
357 } else {
358 anpn += dum*test;
359 }
360 }
361
362 G4double ran = G4UniformRand();
363 G4double excs = 0.0;
364 if( targetCode == protonCode )
365 {
366 counter = -1;
367 for (np=0; np<numSec/3; np++) {
368 for (nm=std::max(0,np-1); nm<=(np+1); nm++) {
369 for (nz=0; nz<numSec/3; nz++) {
370 if (++counter < numMul) {
371 nt = np+nm+nz;
372 if ( (nt>0) && (nt<=numSec) ) {
373 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
374 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
375 if (std::fabs(dum) < 1.0) {
376 if( test >= 1.0e-10 )excs += dum*test;
377 } else {
378 excs += dum*test;
379 }
380 if (ran < excs) goto outOfLoop; //----------------------->
381 }
382 }
383 }
384 }
385 }
386
387 // 3 previous loops continued to the end
388 inElastic = false; // quasi-elastic scattering
389 return;
390 }
391 else
392 { // target must be a neutron
393 counter = -1;
394 for (np=0; np<numSec/3; np++) {
395 for (nm=np; nm<=(np+2); nm++) {
396 for (nz=0; nz<numSec/3; nz++) {
397 if (++counter < numMul) {
398 nt = np+nm+nz;
399 if ( (nt>=1) && (nt<=numSec) ) {
400 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
401 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
402 if (std::fabs(dum) < 1.0) {
403 if( test >= 1.0e-10 )excs += dum*test;
404 } else {
405 excs += dum*test;
406 }
407 if (ran < excs) goto outOfLoop; // -------------------------->
408 }
409 }
410 }
411 }
412 }
413 // 3 previous loops continued to the end
414 inElastic = false; // quasi-elastic scattering.
415 return;
416 }
417 }
418 outOfLoop: // <----------------------------------------------
419
420 if( targetCode == neutronCode)
421 {
422 if( np == nm)
423 {
424 }
425 else if (np == (nm-1))
426 {
427 if( G4UniformRand() < 0.5)
428 {
429 pv[1] = Proton;
430 }
431 else
432 {
433 pv[0] = Proton;
434 }
435 }
436 else
437 {
438 pv[0] = Proton;
439 pv[1] = Proton;
440 }
441 }
442 else
443 {
444 if( np == nm)
445 {
446 if( G4UniformRand() < 0.25)
447 {
448 pv[0] = Proton;
449 pv[1] = Neutron;
450 }
451 else
452 {
453 }
454 }
455 else if ( np == (1+nm))
456 {
457 pv[1] = Neutron;
458 }
459 else
460 {
461 pv[0] = Proton;
462 }
463 }
464
465
466 nt = np + nm + nz;
467 while ( nt > 0)
468 {
469 G4double ran = G4UniformRand();
470 if ( ran < (G4double)np/nt)
471 {
472 if( np > 0 )
473 { pv[vecLen++] = PionPlus;
474 np--;
475 }
476 }
477 else if ( ran < (G4double)(np+nm)/nt)
478 {
479 if( nm > 0 )
480 {
481 pv[vecLen++] = PionMinus;
482 nm--;
483 }
484 }
485 else
486 {
487 if( nz > 0 )
488 {
489 pv[vecLen++] = PionZero;
490 nz--;
491 }
492 }
493 nt = np + nm + nz;
494 }
495 if (verboseLevel > 1)
496 {
497 G4cout << "Particles produced: " ;
498 G4cout << pv[0].getName() << " " ;
499 G4cout << pv[1].getName() << " " ;
500 for (i=2; i < vecLen; i++)
501 {
502 G4cout << pv[i].getName() << " " ;
503 }
504 G4cout << G4endl;
505 }
506 return;
507 }
508
509
510
511
512
513
514
515
516
Note: See TracBrowser for help on using the repository browser.