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