Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4RPGAntiProtonInelastic.cc
이 파일의 문서화 페이지로 가기
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: G4RPGAntiProtonInelastic.cc 94214 2015-11-09 08:18:05Z gcosmo $
27 //
28 
30 #include "G4Exp.hh"
31 #include "G4PhysicalConstants.hh"
32 #include "G4SystemOfUnits.hh"
33 #include "Randomize.hh"
34 
37  G4Nucleus &targetNucleus )
38 {
39  const G4HadProjectile* originalIncident = &aTrack;
40 
41  if (originalIncident->GetKineticEnergy() <= 0.1*MeV) {
45  return &theParticleChange;
46  }
47 
48  //
49  // create the target particle
50  //
51  G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
52 
53  if( verboseLevel > 1 )
54  {
55  const G4Material *targetMaterial = aTrack.GetMaterial();
56  G4cout << "G4RPGAntiProtonInelastic::ApplyYourself called" << G4endl;
57  G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
58  G4cout << "target material = " << targetMaterial->GetName() << ", ";
59  G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
60  << G4endl;
61  }
62  //
63  // Fermi motion and evaporation
64  // As of Geant3, the Fermi energy calculation had not been Done
65  //
66  G4double ek = originalIncident->GetKineticEnergy()/MeV;
67  G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
68  G4ReactionProduct modifiedOriginal;
69  modifiedOriginal = *originalIncident;
70 
71  G4double tkin = targetNucleus.Cinema( ek );
72  ek += tkin;
73  modifiedOriginal.SetKineticEnergy( ek*MeV );
74  G4double et = ek + amas;
75  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
76  G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
77  if( pp > 0.0 )
78  {
79  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
80  modifiedOriginal.SetMomentum( momentum * (p/pp) );
81  }
82  //
83  // calculate black track energies
84  //
85  tkin = targetNucleus.EvaporationEffects( ek );
86  ek -= tkin;
87  modifiedOriginal.SetKineticEnergy( ek*MeV );
88  et = ek + amas;
89  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
90  pp = modifiedOriginal.GetMomentum().mag()/MeV;
91  if( pp > 0.0 )
92  {
93  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
94  modifiedOriginal.SetMomentum( momentum * (p/pp) );
95  }
96  G4ReactionProduct currentParticle = modifiedOriginal;
97  G4ReactionProduct targetParticle;
98  targetParticle = *originalTarget;
99  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
100  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
101  G4bool incidentHasChanged = false;
102  G4bool targetHasChanged = false;
103  G4bool quasiElastic = false;
104  G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
105  G4int vecLen = 0;
106  vec.Initialize( 0 );
107 
108  const G4double cutOff = 0.1;
109  const G4double anni = std::min( 1.3*originalIncident->GetTotalMomentum()/GeV, 0.4 );
110 
111  if( (currentParticle.GetKineticEnergy()/MeV > cutOff) ||
112  (G4UniformRand() > anni) )
113  Cascade( vec, vecLen,
114  originalIncident, currentParticle, targetParticle,
115  incidentHasChanged, targetHasChanged, quasiElastic );
116  else
117  quasiElastic = true;
118 
119  CalculateMomenta( vec, vecLen,
120  originalIncident, originalTarget, modifiedOriginal,
121  targetNucleus, currentParticle, targetParticle,
122  incidentHasChanged, targetHasChanged, quasiElastic );
123 
124  SetUpChange( vec, vecLen,
125  currentParticle, targetParticle,
126  incidentHasChanged );
127 
128  delete originalTarget;
129  return &theParticleChange;
130 }
131 
132 
135  G4int& vecLen,
136  const G4HadProjectile *originalIncident,
137  G4ReactionProduct &currentParticle,
138  G4ReactionProduct &targetParticle,
139  G4bool &incidentHasChanged,
140  G4bool &targetHasChanged,
141  G4bool &quasiElastic )
142 {
143  // Derived from H. Fesefeldt's original FORTRAN code CASPB
144  // AntiProton undergoes interaction with nucleon within a nucleus. Check if it is
145  // energetically possible to produce pions/kaons. In not, assume nuclear excitation
146  // occurs and input particle is degraded in energy. No other particles are produced.
147  // If reaction is possible, find the correct number of pions/protons/neutrons
148  // produced using an interpolation to multiplicity data. Replace some pions or
149  // protons/neutrons by kaons or strange baryons according to the average
150  // multiplicity per Inelastic reaction.
151 
152  const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
153  const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
154  const G4double pOriginal = originalIncident->GetTotalMomentum()/MeV;
155  const G4double targetMass = targetParticle.GetMass()/MeV;
156  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
157  targetMass*targetMass +
158  2.0*targetMass*etOriginal );
159  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
160 
161  static G4ThreadLocal G4bool first = true;
162  const G4int numMul = 1200;
163  const G4int numMulA = 400;
164  const G4int numSec = 60;
165  static G4ThreadLocal G4double protmul[numMul], protnorm[numSec]; // proton constants
166  static G4ThreadLocal G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
167  static G4ThreadLocal G4double protmulA[numMulA], protnormA[numSec]; // proton constants
168  static G4ThreadLocal G4double neutmulA[numMulA], neutnormA[numSec]; // neutron constants
169  // np = number of pi+, nneg = number of pi-, nz = number of pi0
170  G4int counter, nt=0, np=0, nneg=0, nz=0;
171  G4double test;
172  const G4double c = 1.25;
173  const G4double b[] = { 0.7, 0.7 };
174 
175  if (first) { // compute normalization constants, this will only be Done once
176  first = false;
177  G4int i;
178  for( i=0; i<numMul; ++i )protmul[i] = 0.0;
179  for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
180  counter = -1;
181  for (np=0; np<(numSec/3); ++np ) {
182  for (nneg=std::max(0,np-1); nneg<=(np+1); ++nneg ) {
183  for (nz=0; nz<numSec/3; ++nz ) {
184  if (++counter < numMul ) {
185  nt = np+nneg+nz;
186  if (nt>0 && nt<=numSec ) {
187  protmul[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
188  protnorm[nt-1] += protmul[counter];
189  }
190  }
191  }
192  }
193  }
194 
195  for (i=0; i<numMul; ++i )neutmul[i] = 0.0;
196  for (i=0; i<numSec; ++i )neutnorm[i] = 0.0;
197  counter = -1;
198  for (np=0; np<numSec/3; ++np ) {
199  for( nneg=np; nneg<=(np+2); ++nneg )
200  {
201  for( nz=0; nz<numSec/3; ++nz )
202  {
203  if( ++counter < numMul )
204  {
205  nt = np+nneg+nz;
206  if( (nt>0) && (nt<=numSec) )
207  {
208  neutmul[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
209  neutnorm[nt-1] += neutmul[counter];
210  }
211  }
212  }
213  }
214  }
215 
216  for (i=0; i<numSec; ++i ) {
217  if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
218  if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
219  }
220 
221  // do the same for annihilation channels
222 
223  for (i=0; i<numMulA; ++i )protmulA[i] = 0.0;
224  for (i=0; i<numSec; ++i )protnormA[i] = 0.0;
225  counter = -1;
226 
227  for (np=0; np<(numSec/3); ++np ) {
228  nneg = np;
229  for (nz=0; nz<numSec/3; ++nz ) {
230  if ( ++counter < numMulA ) {
231  nt = np+nneg+nz;
232  if ( nt>1 && nt<=numSec ) {
233  protmulA[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
234  protnormA[nt-1] += protmulA[counter];
235  }
236  }
237  }
238  }
239 
240  for( i=0; i<numMulA; ++i )neutmulA[i] = 0.0;
241  for( i=0; i<numSec; ++i )neutnormA[i] = 0.0;
242  counter = -1;
243 
244  for (np=0; np<numSec/3; ++np ) {
245  nneg = np+1;
246  for (nz=0; nz<numSec/3; ++nz ) {
247  if (++counter < numMulA ) {
248  nt = np+nneg+nz;
249  if ( (nt>1) && (nt<=numSec) ) {
250  neutmulA[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
251  neutnormA[nt-1] += neutmulA[counter];
252  }
253  }
254  }
255  }
256  for (i=0; i<numSec; ++i) {
257  if( protnormA[i] > 0.0 )protnormA[i] = 1.0/protnormA[i];
258  if( neutnormA[i] > 0.0 )neutnormA[i] = 1.0/neutnormA[i];
259  }
260  } // end of initialization
261 
262  // initialization is OK
263 
264  const G4double expxu = 82.; // upper bound for arg. of exp
265  const G4double expxl = -expxu; // lower bound for arg. of exp
266 
271 
272  const G4double anhl[] = {1.0,1.0,1.0,1.0,1.0,1.0,1.0,1.0,1.0,0.90,
273  0.6,0.52,0.47,0.44,0.41,0.39,0.37,0.35,0.34,0.24,
274  0.19,0.15,0.12,0.10,0.09,0.07,0.06,0.05,0.0};
275 
276  G4int iplab = G4int( pOriginal/GeV*10.0 );
277  if( iplab > 9 )iplab = G4int( pOriginal/GeV ) + 9;
278  if( iplab > 18 )iplab = G4int( pOriginal/GeV/10.0 ) + 18;
279  if( iplab > 27 )iplab = 28;
280  if (G4UniformRand() > anhl[iplab] ) {
281  if (availableEnergy <= aPiPlus->GetPDGMass()/MeV ) {
282  quasiElastic = true;
283  return;
284  }
285  G4int ieab = static_cast<G4int>(availableEnergy*5.0/GeV);
286  const G4double supp[] = {0.,0.4,0.55,0.65,0.75,0.82,0.86,0.90,0.94,0.98};
287  G4double w0, wp, wt, wm;
288  if ( (availableEnergy < 2.0*GeV) && (G4UniformRand() >= supp[ieab]) ) {
289  //
290  // suppress high multiplicity events at low momentum
291  // only one pion will be produced
292  //
293  np = nneg = nz = 0;
294  if ( targetParticle.GetDefinition() == aProton ) {
295  test = G4Exp( std::min( expxu, std::max( expxl, -(1.0+b[1])*(1.0+b[1])/(2.0*c*c) ) ) );
296  w0 = test;
297  wp = test;
298  test = G4Exp( std::min( expxu, std::max( expxl, -(-1.0+b[1])*(-1.0+b[1])/(2.0*c*c) ) ) );
299  wm = test;
300  wt = w0+wp+wm;
301  wp += w0;
302  G4double ran = G4UniformRand();
303  if( ran < w0/wt )
304  nz = 1;
305  else if( ran < wp/wt )
306  np = 1;
307  else
308  nneg = 1;
309  } else { // target is a neutron
310  test = G4Exp( std::min( expxu, std::max( expxl, -(1.0+b[0])*(1.0+b[0])/(2.0*c*c) ) ) );
311  w0 = test;
312  test = G4Exp( std::min( expxu, std::max( expxl, -(-1.0+b[0])*(-1.0+b[0])/(2.0*c*c) ) ) );
313  wm = test;
314  G4double ran = G4UniformRand();
315  if( ran < w0/(w0+wm) )
316  nz = 1;
317  else
318  nneg = 1;
319  }
320  } else { // (availableEnergy >= 2.0*GeV) || (random number < supp[ieab])
321  G4double n, anpn;
322  GetNormalizationConstant( availableEnergy, n, anpn );
323  G4double ran = G4UniformRand();
324  G4double dum, excs = 0.0;
325  if (targetParticle.GetDefinition() == aProton ) {
326  counter = -1;
327  for (np=0; np<numSec/3 && ran>=excs; ++np ) {
328  for (nneg=std::max(0,np-1); nneg<=(np+1) && ran>=excs; ++nneg ) {
329  for (nz=0; nz<numSec/3 && ran>=excs; ++nz ) {
330  if (++counter < numMul ) {
331  nt = np+nneg+nz;
332  if ((nt>0) && (nt<=numSec) ) {
333  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
334  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
335  if (std::fabs(dum) < 1.0 ) {
336  if( test >= 1.0e-10 )excs += dum*test;
337  }
338  else
339  excs += dum*test;
340  }
341  }
342  }
343  }
344  }
345  if (ran >= excs ) { // 3 previous loops continued to the end
346  quasiElastic = true;
347  return;
348  }
349  } else { // target must be a neutron
350  counter = -1;
351  for ( np=0; np<numSec/3 && ran>=excs; ++np ) {
352  for ( nneg=np; nneg<=(np+2) && ran>=excs; ++nneg ) {
353  for ( nz=0; nz<numSec/3 && ran>=excs; ++nz ) {
354  if ( ++counter < numMul ) {
355  nt = np+nneg+nz;
356  if ( (nt>0) && (nt<=numSec) ) {
357  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
358  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
359  if ( std::fabs(dum) < 1.0 )
360  {
361  if( test >= 1.0e-10 )excs += dum*test;
362  }
363  else
364  excs += dum*test;
365  }
366  }
367  }
368  }
369  }
370  if( ran >= excs ) { // 3 previous loops continued to the end
371  quasiElastic = true;
372  return;
373  }
374  }
375  np--; nneg--; nz--;
376 
377  } // end if availableEnergy
378 
379  if (targetParticle.GetDefinition() == aProton ) {
380  switch( np-nneg )
381  {
382  case 0:
383  if( G4UniformRand() < 0.33 )
384  {
385  currentParticle.SetDefinitionAndUpdateE( anAntiNeutron );
386  targetParticle.SetDefinitionAndUpdateE( aNeutron );
387  incidentHasChanged = true;
388  targetHasChanged = true;
389  }
390  break;
391  case 1:
392  targetParticle.SetDefinitionAndUpdateE( aNeutron );
393  targetHasChanged = true;
394  break;
395  default:
396  currentParticle.SetDefinitionAndUpdateE( anAntiNeutron );
397  incidentHasChanged = true;
398  break;
399  }
400 
401  } else { // target must be a neutron
402  switch( np-nneg )
403  {
404  case -1:
405  if( G4UniformRand() < 0.5 )
406  {
407  targetParticle.SetDefinitionAndUpdateE( aProton );
408  targetHasChanged = true;
409  }
410  else
411  {
412  currentParticle.SetDefinitionAndUpdateE( anAntiNeutron );
413  incidentHasChanged = true;
414  }
415  break;
416  case 0:
417  break;
418  default:
419  currentParticle.SetDefinitionAndUpdateE( anAntiNeutron );
420  targetParticle.SetDefinitionAndUpdateE( aProton );
421  incidentHasChanged = true;
422  targetHasChanged = true;
423  break;
424  }
425  }
426 
427  } else { // random number <= anhl[iplab]
428  if (centerofmassEnergy <= 2*aPiPlus->GetPDGMass()/MeV ) {
429  quasiElastic = true;
430  return;
431  }
432 
433  // annihilation channels
434 
435  G4double n, anpn;
436  GetNormalizationConstant( -centerofmassEnergy, n, anpn );
437  G4double ran = G4UniformRand();
438  G4double dum, excs = 0.0;
439  if (targetParticle.GetDefinition() == aProton ) {
440  counter = -1;
441  for( np=0; (np<numSec/3) && (ran>=excs); ++np )
442  {
443  nneg = np;
444  for( nz=0; (nz<numSec/3) && (ran>=excs); ++nz )
445  {
446  if( ++counter < numMulA )
447  {
448  nt = np+nneg+nz;
449  if( (nt>1) && (nt<=numSec) )
450  {
451  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
452  dum = (pi/anpn)*nt*protmulA[counter]*protnormA[nt-1]/(2.0*n*n);
453  if( std::abs(dum) < 1.0 )
454  {
455  if( test >= 1.0e-10 )excs += dum*test;
456  }
457  else
458  excs += dum*test;
459  }
460  }
461  }
462  }
463  } else { // target must be a neutron
464  counter = -1;
465  for( np=0; (np<numSec/3) && (ran>=excs); ++np )
466  {
467  nneg = np+1;
468  for( nz=0; (nz<numSec/3) && (ran>=excs); ++nz )
469  {
470  if( ++counter < numMulA )
471  {
472  nt = np+nneg+nz;
473  if( (nt>1) && (nt<=numSec) )
474  {
475  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
476  dum = (pi/anpn)*nt*neutmulA[counter]*neutnormA[nt-1]/(2.0*n*n);
477  if( std::fabs(dum) < 1.0 )
478  {
479  if( test >= 1.0e-10 )excs += dum*test;
480  }
481  else
482  excs += dum*test;
483  }
484  }
485  }
486  }
487  }
488  if ( ran >= excs ) { // 3 previous loops continued to the end
489  quasiElastic = true;
490  return;
491  }
492  np--; nz--;
493  currentParticle.SetMass( 0.0 );
494  targetParticle.SetMass( 0.0 );
495  }
496 
497  G4int loop = 0;
499  ed << " While count exceeded " << G4endl;
500  while(np + nneg + nz < 3) { /* Loop checking, 01.09.2015, D.Wright */
501  nz++;
502  loop++;
503  if (loop > 1000) {
504  G4Exception("G4RPGAntiProtonInelastic::Cascade()", "HAD_RPG_100", JustWarning, ed);
505  break;
506  }
507  }
508 
509  SetUpPions( np, nneg, nz, vec, vecLen );
510  return;
511 }
512 
513  /* end of file */
514 
void Initialize(G4int items)
Definition: G4FastVector.hh:63
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void SetUpChange(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
G4double ek
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
static constexpr double MeV
Definition: G4SIunits.hh:214
static G4AntiNeutron * AntiNeutron()
G4double Cinema(G4double kineticEnergy)
Definition: G4Nucleus.cc:382
void SetMomentumChange(const G4ThreeVector &aV)
void SetMomentum(const G4double x, const G4double y, const G4double z)
#define G4endl
Definition: G4ios.hh:61
const char * p
Definition: xmltok.h:285
void SetDefinitionAndUpdateE(const G4ParticleDefinition *aParticleDefinition)
const G4String & GetParticleName() const
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
void GetNormalizationConstant(const G4double availableEnergy, G4double &n, G4double &anpn)
G4double GetTotalMomentum() const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
const G4ParticleDefinition * GetDefinition() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double GetTotalEnergy() const
#define G4ThreadLocal
Definition: tls.hh:69
void CalculateMomenta(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, const G4DynamicParticle *originalTarget, G4ReactionProduct &modifiedOriginal, G4Nucleus &targetNucleus, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool quasiElastic)
G4double GetPDGMass() const
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:241
static unsigned wp
Definition: csz_inflate.cc:294
void SetMass(const G4double mas)
const G4String & GetName() const
Definition: G4Material.hh:179
void SetEnergyChange(G4double anEnergy)
G4double Pmltpc(G4int np, G4int nm, G4int nz, G4int n, G4double b, G4double c)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
void SetKineticEnergy(const G4double en)
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
const G4Material * GetMaterial() const
#define G4UniformRand()
Definition: Randomize.hh:53
void SetUpPions(const G4int np, const G4int nm, const G4int nz, G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen)
G4ThreeVector GetMomentum() const
void SetSide(const G4int sid)
const G4LorentzVector & Get4Momentum() const
Hep3Vector unit() const
G4double GetKineticEnergy() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
G4double GetMass() const
const G4ParticleDefinition * GetDefinition() const
double mag() const
int G4int
Definition: G4Types.hh:78
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4GLOB_DLL std::ostream G4cout
void Cascade(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool &quasiElastic)
Hep3Vector vect() const
Char_t n[5]
static constexpr double pi
Definition: G4SIunits.hh:75
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:278
static constexpr double GeV
Definition: G4SIunits.hh:217
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
void SetStatusChange(G4HadFinalStateStatus aS)