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