Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LEProtonInelastic.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 // Hadronic Process: Low Energy Proton Inelastic Process
27 // J.L. Chuma, TRIUMF, 19-Nov-1996
28 
29 #include <iostream>
30 
31 #include "G4LEProtonInelastic.hh"
32 #include "G4PhysicalConstants.hh"
33 #include "G4SystemOfUnits.hh"
34 #include "Randomize.hh"
35 
38 {
39  SetMinEnergy(0.0);
40  SetMaxEnergy(55.*GeV);
41  G4cout << "WARNING: model G4LEProtonInelastic is being deprecated and will\n"
42  << "disappear in Geant4 version 10.0" << G4endl;
43 }
44 
45 
47 {
48  outFile << "G4LEProtonInelastic is one of the Low Energy Parameterized\n"
49  << "(LEP) models used to implement inelastic proton scattering\n"
50  << "from nuclei. It is a re-engineered version of the GHEISHA\n"
51  << "code of H. Fesefeldt. It divides the initial collision\n"
52  << "products into backward- and forward-going clusters which are\n"
53  << "then decayed into final state hadrons. The model does not\n"
54  << "conserve energy on an event-by-event basis. It may be\n"
55  << "applied to protons with initial energies between 0 and 25\n"
56  << "GeV.\n";
57 }
58 
59 
62  G4Nucleus& targetNucleus)
63 {
65  const G4HadProjectile *originalIncident = &aTrack;
66  if (originalIncident->GetKineticEnergy()<= 0.1*MeV) {
70  return &theParticleChange;
71  }
72 
73  // Create the target particle
74 
75  G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
76  if (verboseLevel > 1) {
77  const G4Material *targetMaterial = aTrack.GetMaterial();
78  G4cout << "G4LEProtonInelastic::ApplyYourself called" << G4endl;
79  G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
80  G4cout << "target material = " << targetMaterial->GetName() << ", ";
81  G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
82  << G4endl;
83  }
84 
85  if (originalIncident->GetKineticEnergy()/GeV < 0.01+2.*G4UniformRand()/9.) {
86  SlowProton(originalIncident, targetNucleus);
87  if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
88  delete originalTarget;
89  return &theParticleChange;
90  }
91 
92  // Fermi motion and evaporation
93  // As of Geant3, the Fermi energy calculation had not been Done
94 
95  G4double ek = originalIncident->GetKineticEnergy()/MeV;
96  G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
97  G4ReactionProduct modifiedOriginal;
98  modifiedOriginal = *originalIncident;
99 
100  G4double tkin = targetNucleus.Cinema( ek );
101  ek += tkin;
102  modifiedOriginal.SetKineticEnergy( ek*MeV );
103  G4double et = ek + amas;
104  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
105  G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
106  if (pp > 0.0) {
107  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
108  modifiedOriginal.SetMomentum( momentum * (p/pp) );
109  }
110 
111  // calculate black track energies
112 
113  tkin = targetNucleus.EvaporationEffects( ek );
114  ek -= tkin;
115  modifiedOriginal.SetKineticEnergy( ek*MeV );
116  et = ek + amas;
117  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
118  pp = modifiedOriginal.GetMomentum().mag()/MeV;
119  if (pp > 0.0) {
120  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
121  modifiedOriginal.SetMomentum( momentum * (p/pp) );
122  }
123  const G4double cutOff = 0.1;
124  if (modifiedOriginal.GetKineticEnergy()/MeV <= cutOff) {
125  SlowProton( originalIncident, targetNucleus);
126  if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
127  delete originalTarget;
128  return &theParticleChange;
129  }
130  G4ReactionProduct currentParticle = modifiedOriginal;
131  G4ReactionProduct targetParticle;
132  targetParticle = *originalTarget;
133  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
134  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
135  G4bool incidentHasChanged = false;
136  G4bool targetHasChanged = false;
137  G4bool quasiElastic = false;
138  G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the sec. particles
139  G4int vecLen = 0;
140  vec.Initialize( 0 );
141 
142  Cascade(vec, vecLen,
143  originalIncident, currentParticle, targetParticle,
144  incidentHasChanged, targetHasChanged, quasiElastic);
145 
146  CalculateMomenta(vec, vecLen,
147  originalIncident, originalTarget, modifiedOriginal,
148  targetNucleus, currentParticle, targetParticle,
149  incidentHasChanged, targetHasChanged, quasiElastic);
150 
151  SetUpChange(vec, vecLen,
152  currentParticle, targetParticle,
153  incidentHasChanged);
154 
155  if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
156  delete originalTarget;
157  return &theParticleChange;
158 }
159 
160 
161 void
162 G4LEProtonInelastic::SlowProton(const G4HadProjectile* originalIncident,
163  G4Nucleus& targetNucleus)
164 {
165  const G4double A = targetNucleus.GetA_asInt(); // atomic weight
166  const G4double Z = targetNucleus.GetZ_asInt(); // atomic number
167  //
168  // calculate Q-value of reactions
169  //
170  G4double theAtomicMass = targetNucleus.AtomicMass( A, Z );
171  G4double massVec[9];
172  massVec[0] = targetNucleus.AtomicMass( A+1.0, Z+1.0 );
173  massVec[1] = 0.;
174  if (A > Z+1.0)
175  massVec[1] = targetNucleus.AtomicMass( A , Z+1.0 );
176  massVec[2] = theAtomicMass;
177  massVec[3] = 0.;
178  if (A > 1.0 && A-1.0 > Z)
179  massVec[3] = targetNucleus.AtomicMass( A-1.0, Z );
180  massVec[4] = 0.;
181  if (A > 2.0 && A-2.0 > Z)
182  massVec[4] = targetNucleus.AtomicMass( A-2.0, Z );
183  massVec[5] = 0.;
184  if (A > 3.0 && Z > 1.0 && A-3.0 > Z-1.0)
185  massVec[5] = targetNucleus.AtomicMass( A-3.0, Z-1.0 );
186  massVec[6] = 0.;
187  if (A > 1.0 && A-1.0 > Z+1.0)
188  massVec[6] = targetNucleus.AtomicMass( A-1.0, Z+1.0 );
189  massVec[7] = massVec[3];
190  massVec[8] = 0.;
191  if (A > 1.0 && Z > 1.0)
192  massVec[8] = targetNucleus.AtomicMass( A-1.0, Z-1.0 );
193 
194  G4FastVector<G4ReactionProduct,4> vec; // vec will contain the secondary particles
195  G4int vecLen = 0;
196  vec.Initialize( 0 );
197 
198  theReactionDynamics.NuclearReaction( vec, vecLen, originalIncident,
199  targetNucleus, theAtomicMass, massVec );
200 
203 
204  G4DynamicParticle *pd;
205  for( G4int i=0; i<vecLen; ++i )
206  {
207  pd = new G4DynamicParticle();
208  pd->SetDefinition( vec[i]->GetDefinition() );
209  pd->SetMomentum( vec[i]->GetMomentum() );
211  delete vec[i];
212  }
213  }
214 
215  void
216  G4LEProtonInelastic::Cascade(
218  G4int &vecLen,
219  const G4HadProjectile *originalIncident,
220  G4ReactionProduct &currentParticle,
221  G4ReactionProduct &targetParticle,
222  G4bool &incidentHasChanged,
223  G4bool &targetHasChanged,
224  G4bool &quasiElastic )
225  {
226  // derived from original FORTRAN code CASP by H. Fesefeldt (13-Sep-1987)
227  //
228  // Proton undergoes interaction with nucleon within a nucleus. Check if it is
229  // energetically possible to produce pions/kaons. In not, assume nuclear excitation
230  // occurs and input particle is degraded in energy. No other particles are produced.
231  // If reaction is possible, find the correct number of pions/protons/neutrons
232  // produced using an interpolation to multiplicity data. Replace some pions or
233  // protons/neutrons by kaons or strange baryons according to the average
234  // multiplicity per Inelastic reaction.
235  //
236  // the center of mass energy is based on the initial energy, before
237  // Fermi motion and evaporation effects are taken into account
238  //
239  const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
240  const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
241  const G4double targetMass = targetParticle.GetMass()/MeV;
242  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
243  targetMass*targetMass +
244  2.0*targetMass*etOriginal );
245  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
246  if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV )
247  { // not energetically possible to produce pion(s)
248  quasiElastic = true;
249  return;
250  }
251  static G4bool first = true;
252  const G4int numMul = 1200;
253  const G4int numSec = 60;
254  static G4double protmul[numMul], protnorm[numSec]; // proton constants
255  static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
256  // npos = number of pi+, nneg = number of pi-, nzero = number of pi0
257  G4int counter, nt=0, npos=0, nneg=0, nzero=0;
258  const G4double c = 1.25;
259  const G4double b[] = { 0.70, 0.35 };
260  if( first ) // compute normalization constants, this will only be Done once
261  {
262  first = false;
263  G4int i;
264  for( i=0; i<numMul; ++i )protmul[i] = 0.0;
265  for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
266  counter = -1;
267  for( npos=0; npos<(numSec/3); ++npos )
268  {
269  for( nneg=std::max(0,npos-2); nneg<=npos; ++nneg )
270  {
271  for( nzero=0; nzero<numSec/3; ++nzero )
272  {
273  if( ++counter < numMul )
274  {
275  nt = npos+nneg+nzero;
276  if( nt>0 && nt<=numSec )
277  {
278  protmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[0],c) /
279  ( theReactionDynamics.Factorial(2-npos+nneg)*
280  theReactionDynamics.Factorial(npos-nneg) );
281  protnorm[nt-1] += protmul[counter];
282  }
283  }
284  }
285  }
286  }
287  for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
288  for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
289  counter = -1;
290  for( npos=0; npos<numSec/3; ++npos )
291  {
292  for( nneg=std::max(0,npos-1); nneg<=(npos+1); ++nneg )
293  {
294  for( nzero=0; nzero<numSec/3; ++nzero )
295  {
296  if( ++counter < numMul )
297  {
298  nt = npos+nneg+nzero;
299  if( nt>0 && nt<=numSec )
300  {
301  neutmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[1],c) /
302  ( theReactionDynamics.Factorial(1-npos+nneg)*
303  theReactionDynamics.Factorial(1+npos-nneg) );
304  neutnorm[nt-1] += neutmul[counter];
305  }
306  }
307  }
308  }
309  }
310  for( i=0; i<numSec; ++i )
311  {
312  if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
313  if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
314  }
315  } // end of initialization
316 
317  const G4double expxu = 82.; // upper bound for arg. of exp
318  const G4double expxl = -expxu; // lower bound for arg. of exp
321  G4int ieab = static_cast<G4int>(availableEnergy*5.0/GeV);
322  const G4double supp[] = {0.,0.4,0.55,0.65,0.75,0.82,0.86,0.90,0.94,0.98};
323  G4double test, w0, wp, wt, wm;
324  if( (availableEnergy < 2.0*GeV) && (G4UniformRand() >= supp[ieab]) )
325  {
326  // suppress high multiplicity events at low momentum
327  // only one pion will be produced
328 
329  npos = nneg = nzero = 0;
330  if( targetParticle.GetDefinition() == aProton )
331  {
332  test = std::exp( std::min( expxu, std::max(
333  expxl, -(1.0+b[0])*(1.0+b[0])/(2.0*c*c) ) ) );
334  w0 = test/2.0;
335  wp = test;
336  if( G4UniformRand() < w0/(w0+wp) )
337  nzero = 1;
338  else
339  npos = 1;
340  }
341  else // target is a neutron
342  {
343  test = std::exp( std::min( expxu, std::max(
344  expxl, -(1.0+b[1])*(1.0+b[1])/(2.0*c*c) ) ) );
345  w0 = test;
346  wp = test/2.0;
347  test = std::exp( std::min( expxu, std::max(
348  expxl, -(-1.0+b[1])*(-1.0+b[1])/(2.0*c*c) ) ) );
349  wm = test/2.0;
350  wt = w0+wp+wm;
351  wp += w0;
352  G4double ran = G4UniformRand();
353  if( ran < w0/wt )
354  nzero = 1;
355  else if( ran < wp/wt )
356  npos = 1;
357  else
358  nneg = 1;
359  }
360  }
361  else // (availableEnergy >= 2.0*GeV) || (random number < supp[ieab])
362  {
363  G4double n, anpn;
364  GetNormalizationConstant( availableEnergy, n, anpn );
365  G4double ran = G4UniformRand();
366  G4double dum, excs = 0.0;
367  if( targetParticle.GetDefinition() == aProton )
368  {
369  counter = -1;
370  for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
371  {
372  for( nneg=std::max(0,npos-2); nneg<=npos && ran>=excs; ++nneg )
373  {
374  for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
375  {
376  if( ++counter < numMul )
377  {
378  nt = npos+nneg+nzero;
379  if( nt>0 && nt<=numSec )
380  {
381  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
382  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
383  if( std::fabs(dum) < 1.0 )
384  {
385  if( test >= 1.0e-10 )excs += dum*test;
386  } else {
387  excs += dum*test;
388  }
389  }
390  }
391  }
392  }
393  }
394  if( ran >= excs ) // 3 previous loops continued to the end
395  {
396  quasiElastic = true;
397  return;
398  }
399  }
400  else // target must be a neutron
401  {
402  counter = -1;
403  for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
404  {
405  for( nneg=std::max(0,npos-1); nneg<=(npos+1) && ran>=excs; ++nneg )
406  {
407  for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
408  {
409  if( ++counter < numMul )
410  {
411  nt = npos+nneg+nzero;
412  if( nt>0 && nt<=numSec )
413  {
414  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
415  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
416  if( std::fabs(dum) < 1.0 )
417  {
418  if( test >= 1.0e-10 )excs += dum*test;
419  } else {
420  excs += dum*test;
421  }
422  }
423  }
424  }
425  }
426  }
427  if( ran >= excs ) // 3 previous loops continued to the end
428  {
429  quasiElastic = true;
430  return;
431  }
432  }
433  npos--; nneg--; nzero--;
434  }
435  if( targetParticle.GetDefinition() == aProton )
436  {
437  switch( npos-nneg )
438  {
439  case 1:
440  if( G4UniformRand() < 0.5 )
441  {
442  targetParticle.SetDefinitionAndUpdateE( aNeutron );
443  targetHasChanged = true;
444  } else {
445  currentParticle.SetDefinitionAndUpdateE( aNeutron );
446  incidentHasChanged = true;
447  }
448  break;
449  case 2:
450  currentParticle.SetDefinitionAndUpdateE( aNeutron );
451  targetParticle.SetDefinitionAndUpdateE( aNeutron );
452  incidentHasChanged = true;
453  targetHasChanged = true;
454  break;
455  default:
456  break;
457  }
458  }
459  else // target is a neutron
460  {
461  switch( npos-nneg )
462  {
463  case 0:
464  if( G4UniformRand() < 0.333333 )
465  {
466  currentParticle.SetDefinitionAndUpdateE( aNeutron );
467  targetParticle.SetDefinitionAndUpdateE( aProton );
468  incidentHasChanged = true;
469  targetHasChanged = true;
470  }
471  break;
472  case 1:
473  currentParticle.SetDefinitionAndUpdateE( aNeutron );
474  incidentHasChanged = true;
475  break;
476  default:
477  targetParticle.SetDefinitionAndUpdateE( aProton );
478  targetHasChanged = true;
479  break;
480  }
481  }
482  SetUpPions( npos, nneg, nzero, vec, vecLen );
483  return;
484  }
485 
486  /* end of file */
487