source: trunk/source/processes/hadronic/models/high_energy/src/G4HEPionPlusInelastic.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.3 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: G4HEPionPlusInelastic.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 "G4HEPionPlusInelastic.hh"
43
44G4HadFinalState*
45G4HEPionPlusInelastic::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 << "G4HEPionPlusInelastic: incident energy < 1 GeV" << G4endl;
65
66 if (verboseLevel > 1) {
67 G4cout << "G4HEPionPlusInelastic::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 G4bool inElastic = true;
113
114 vecLength = 0;
115
116 if (verboseLevel > 1)
117 G4cout << "ApplyYourself: CallFirstIntInCascade for particle "
118 << incidentCode << G4endl;
119
120 G4bool successful = false;
121
122 FirstIntInCasPionPlus(inElastic, availableEnergy, pv, vecLength,
123 incidentParticle, targetParticle, atomicWeight);
124
125 if (verboseLevel > 1)
126 G4cout << "ApplyYourself::StrangeParticlePairProduction" << G4endl;
127
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 if (!successful)
159 ElasticScattering(successful, pv, vecLength,
160 incidentParticle,
161 atomicWeight, atomicNumber);
162
163 if (!successful)
164 G4cout << "GHEInelasticInteraction::ApplyYourself fails to produce final state particles"
165 << G4endl;
166
167 FillParticleChange(pv, vecLength);
168 delete [] pv;
169 theParticleChange.SetStatusChange(stopAndKill);
170 return &theParticleChange;
171}
172
173
174void
175G4HEPionPlusInelastic::FirstIntInCasPionPlus(G4bool& inElastic,
176 const G4double availableEnergy,
177 G4HEVector pv[],
178 G4int& vecLen,
179 const G4HEVector& incidentParticle,
180 const G4HEVector& targetParticle,
181 const G4double atomicWeight)
182
183// Pion+ undergoes interaction with nucleon within a nucleus. Check if it is
184// energetically possible to produce pions/kaons. In not, assume nuclear excitation
185// occurs and input particle is degraded in energy. No other particles are produced.
186// If reaction is possible, find the correct number of pions/protons/neutrons
187// produced using an interpolation to multiplicity data. Replace some pions or
188// protons/neutrons by kaons or strange baryons according to the average
189// multiplicity per inelastic reaction.
190{
191 static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp
192 static const G4double expxl = -expxu; // lower bound for arg. of exp
193
194 static const G4double protb = 0.7;
195 static const G4double neutb = 0.7;
196 static const G4double c = 1.25;
197
198 static const G4int numMul = 1200;
199 static const G4int numSec = 60;
200
201 G4int neutronCode = Neutron.getCode();
202 G4int protonCode = Proton.getCode();
203 G4double pionMass = PionPlus.getMass();
204
205 G4int targetCode = targetParticle.getCode();
206 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
207
208 static G4bool first = true;
209 static G4double protmul[numMul], protnorm[numSec]; // proton constants
210 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
211
212 // misc. local variables
213 // np = number of pi+, nm = number of pi-, nz = number of pi0
214
215 G4int i, counter, nt, np, nm, nz;
216
217 if( first )
218 { // compute normalization constants, this will only be done once
219 first = false;
220 for( i=0; i<numMul; i++ )protmul[i] = 0.0;
221 for( i=0; i<numSec; i++ )protnorm[i] = 0.0;
222 counter = -1;
223 for( np=0; np<(numSec/3); np++ )
224 {
225 for( nm=Imax(0,np-2); nm<=np; nm++ )
226 {
227 for( nz=0; nz<numSec/3; nz++ )
228 {
229 if( ++counter < numMul )
230 {
231 nt = np+nm+nz;
232 if( (nt>0) && (nt<=numSec) )
233 {
234 protmul[counter] =
235 pmltpc(np,nm,nz,nt,protb,c) ;
236 protnorm[nt-1] += protmul[counter];
237 }
238 }
239 }
240 }
241 }
242 for( i=0; i<numMul; i++ )neutmul[i] = 0.0;
243 for( i=0; i<numSec; i++ )neutnorm[i] = 0.0;
244 counter = -1;
245 for( np=0; np<numSec/3; np++ )
246 {
247 for( nm=Imax(0,np-1); nm<=(np+1); nm++ )
248 {
249 for( nz=0; nz<numSec/3; nz++ )
250 {
251 if( ++counter < numMul )
252 {
253 nt = np+nm+nz;
254 if( (nt>0) && (nt<=numSec) )
255 {
256 neutmul[counter] =
257 pmltpc(np,nm,nz,nt,neutb,c);
258 neutnorm[nt-1] += neutmul[counter];
259 }
260 }
261 }
262 }
263 }
264 for( i=0; i<numSec; i++ )
265 {
266 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
267 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
268 }
269 } // end of initialization
270
271
272 // initialize the first two places
273 // the same as beam and target
274 pv[0] = incidentParticle;
275 pv[1] = targetParticle;
276 vecLen = 2;
277
278 if( !inElastic )
279 { // quasi-elastic scattering, no pions produced
280 if( targetCode == neutronCode )
281 {
282 G4double cech[] = {0.33,0.27,0.29,0.31,0.27,0.18,0.13,0.10,0.09,0.07};
283 G4int iplab = G4int( Amin( 9.0, incidentTotalMomentum*5. ) );
284 if( G4UniformRand() < cech[iplab]/std::pow(atomicWeight,0.42) )
285 { // charge exchange pi+ n -> pi0 p
286 pv[0] = PionZero;
287 pv[1] = Proton;
288 }
289 }
290 return;
291 }
292 else if (availableEnergy <= pionMass)
293 return;
294
295// inelastic scattering
296
297 np = 0, nm = 0, nz = 0;
298 G4double eab = availableEnergy;
299 G4int ieab = G4int( eab*5.0 );
300
301 G4double supp[] = {0., 0.2, 0.45, 0.55, 0.65, 0.75, 0.85, 0.90, 0.94, 0.98};
302 if( (ieab <= 9) && (G4UniformRand() >= supp[ieab]) )
303 {
304// suppress high multiplicity events at low momentum
305// only one additional pion will be produced
306 G4double w0, wp, wm, wt, ran;
307 if( targetCode == protonCode ) // target is a proton
308 {
309 w0 = - sqr(1.+protb)/(2.*c*c);
310 wp = w0 = std::exp(w0);
311 if( G4UniformRand() < w0/(w0+wp) )
312 { np = 0; nm = 0; nz = 1; }
313 else
314 { np = 1; nm = 0; nz = 0; }
315 }
316 else
317 { // target is a neutron
318 w0 = -sqr(1.+neutb)/(2.*c*c);
319 wp = w0 = std::exp(w0);
320 wm = -sqr(-1.+neutb)/(2.*c*c);
321 wm = std::exp(wm);
322 wt = w0+wp+wm;
323 wp = w0+wp;
324 ran = G4UniformRand();
325 if( ran < w0/wt)
326 { np = 0; nm = 0; nz = 1; }
327 else if( ran < wp/wt)
328 { np = 1; nm = 0; nz = 0; }
329 else
330 { np = 0; nm = 1; nz = 0; }
331 }
332 }
333 else
334 {
335 // number of total particles vs. centre of mass Energy - 2*proton mass
336
337 G4double aleab = std::log(availableEnergy);
338 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
339 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
340
341 // normalization constant for kno-distribution.
342 // calculate first the sum of all constants, check for numerical problems.
343 G4double test, dum, anpn = 0.0;
344
345 for (nt=1; nt<=numSec; nt++) {
346 test = std::exp( Amin( expxu, Amax( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
347 dum = pi*nt/(2.0*n*n);
348 if (std::fabs(dum) < 1.0) {
349 if( test >= 1.0e-10 )anpn += dum*test;
350 } else {
351 anpn += dum*test;
352 }
353 }
354
355 G4double ran = G4UniformRand();
356 G4double excs = 0.0;
357 if( targetCode == protonCode )
358 {
359 counter = -1;
360 for (np=0; np<numSec/3; np++) {
361 for (nm=Imax(0,np-2); nm<=np; nm++) {
362 for (nz=0; nz<numSec/3; nz++) {
363 if (++counter < numMul) {
364 nt = np+nm+nz;
365 if ( (nt>0) && (nt<=numSec) ) {
366 test = std::exp( Amin( expxu, Amax( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
367 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
368 if (std::fabs(dum) < 1.0) {
369 if( test >= 1.0e-10 )excs += dum*test;
370 } else {
371 excs += dum*test;
372 }
373 if (ran < excs) goto outOfLoop; //------------------>
374 }
375 }
376 }
377 }
378 }
379
380 // 3 previous loops continued to the end
381 inElastic = false; // quasi-elastic scattering
382 return;
383 }
384 else
385 { // target must be a neutron
386 counter = -1;
387 for (np=0; np<numSec/3; np++) {
388 for (nm=Imax(0,np-1); nm<=(np+1); nm++) {
389 for (nz=0; nz<numSec/3; nz++) {
390 if (++counter < numMul) {
391 nt = np+nm+nz;
392 if ( (nt>=1) && (nt<=numSec) ) {
393 test = std::exp( Amin( expxu, Amax( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
394 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
395 if (std::fabs(dum) < 1.0) {
396 if( test >= 1.0e-10 )excs += dum*test;
397 } else {
398 excs += dum*test;
399 }
400 if (ran < excs) goto outOfLoop; // --------------------->
401 }
402 }
403 }
404 }
405 }
406 // 3 previous loops continued to the end
407 inElastic = false; // quasi-elastic scattering.
408 return;
409 }
410 }
411 outOfLoop: // <--------------------------------------------
412
413 if( targetCode == protonCode)
414 {
415 if( np == nm)
416 {
417 }
418 else if (np == (1+nm))
419 {
420 if( G4UniformRand() < 0.5)
421 {
422 pv[1] = Neutron;
423 }
424 else
425 {
426 pv[0] = PionZero;
427 }
428 }
429 else
430 {
431 pv[0] = PionZero;
432 pv[1] = Neutron;
433 }
434 }
435 else
436 {
437 if( np == nm)
438 {
439 if( G4UniformRand() < 0.25)
440 {
441 pv[0] = PionZero;
442 pv[1] = Proton;
443 }
444 else
445 {
446 }
447 }
448 else if ( np == (1+nm))
449 {
450 pv[0] = PionZero;
451 }
452 else
453 {
454 pv[1] = Proton;
455 }
456 }
457
458
459 nt = np + nm + nz;
460 while ( nt > 0)
461 {
462 G4double ran = G4UniformRand();
463 if ( ran < (G4double)np/nt)
464 {
465 if( np > 0 )
466 { pv[vecLen++] = PionPlus;
467 np--;
468 }
469 }
470 else if ( ran < (G4double)(np+nm)/nt)
471 {
472 if( nm > 0 )
473 {
474 pv[vecLen++] = PionMinus;
475 nm--;
476 }
477 }
478 else
479 {
480 if( nz > 0 )
481 {
482 pv[vecLen++] = PionZero;
483 nz--;
484 }
485 }
486 nt = np + nm + nz;
487 }
488 if (verboseLevel > 1)
489 {
490 G4cout << "Particles produced: " ;
491 G4cout << pv[0].getName() << " " ;
492 G4cout << pv[1].getName() << " " ;
493 for (i=2; i < vecLen; i++)
494 {
495 G4cout << pv[i].getName() << " " ;
496 }
497 G4cout << G4endl;
498 }
499 return;
500 }
501
502
503
504
505
506
507
508
509
Note: See TracBrowser for help on using the repository browser.