source: trunk/source/processes/hadronic/models/high_energy/src/G4HEPionPlusInelastic.cc@ 1036

Last change on this file since 1036 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

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