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