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