Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4RPGAntiSigmaMinusInelastic.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: G4RPGAntiSigmaMinusInelastic.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  if (originalIncident->GetKineticEnergy()<= 0.1*MeV)
41  {
45  return &theParticleChange;
46  }
47 
48  // Choose the target particle
49 
50  G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
51 
52  if( verboseLevel > 1 )
53  {
54  const G4Material *targetMaterial = aTrack.GetMaterial();
55  G4cout << "G4RPGAntiSigmaMinusInelastic::ApplyYourself called" << G4endl;
56  G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
57  G4cout << "target material = " << targetMaterial->GetName() << ", ";
58  G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
59  << G4endl;
60  }
61 
62  // Fermi motion and evaporation
63  // As of Geant3, the Fermi energy calculation had not been Done
64 
65  G4double ek = originalIncident->GetKineticEnergy()/MeV;
66  G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
67  G4ReactionProduct modifiedOriginal;
68  modifiedOriginal = *originalIncident;
69 
70  G4double tkin = targetNucleus.Cinema( ek );
71  ek += tkin;
72  modifiedOriginal.SetKineticEnergy( ek*MeV );
73  G4double et = ek + amas;
74  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
75  G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
76  if( pp > 0.0 )
77  {
78  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
79  modifiedOriginal.SetMomentum( momentum * (p/pp) );
80  }
81  //
82  // calculate black track energies
83  //
84  tkin = targetNucleus.EvaporationEffects( ek );
85  ek -= tkin;
86  modifiedOriginal.SetKineticEnergy( ek*MeV );
87  et = ek + amas;
88  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
89  pp = modifiedOriginal.GetMomentum().mag()/MeV;
90  if( pp > 0.0 )
91  {
92  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
93  modifiedOriginal.SetMomentum( momentum * (p/pp) );
94  }
95  G4ReactionProduct currentParticle = modifiedOriginal;
96  G4ReactionProduct targetParticle;
97  targetParticle = *originalTarget;
98  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
99  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
100  G4bool incidentHasChanged = false;
101  G4bool targetHasChanged = false;
102  G4bool quasiElastic = false;
103  G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
104  G4int vecLen = 0;
105  vec.Initialize( 0 );
106 
107  const G4double cutOff = 0.1;
108  const G4double anni = std::min( 1.3*currentParticle.GetTotalMomentum()/GeV, 0.4 );
109  if( (currentParticle.GetKineticEnergy()/MeV > cutOff) || (G4UniformRand() > anni) )
110  Cascade( vec, vecLen,
111  originalIncident, currentParticle, targetParticle,
112  incidentHasChanged, targetHasChanged, quasiElastic );
113 
114  CalculateMomenta( vec, vecLen,
115  originalIncident, originalTarget, modifiedOriginal,
116  targetNucleus, currentParticle, targetParticle,
117  incidentHasChanged, targetHasChanged, quasiElastic );
118 
119  SetUpChange( vec, vecLen,
120  currentParticle, targetParticle,
121  incidentHasChanged );
122 
123  delete originalTarget;
124  return &theParticleChange;
125 }
126 
127 
130  G4int& vecLen,
131  const G4HadProjectile *originalIncident,
132  G4ReactionProduct &currentParticle,
133  G4ReactionProduct &targetParticle,
134  G4bool &incidentHasChanged,
135  G4bool &targetHasChanged,
136  G4bool &quasiElastic )
137 {
138  // Derived from H. Fesefeldt's original FORTRAN code CASASM
139  // AntiSigmaMinus undergoes interaction with nucleon within a nucleus. Check if it is
140  // energetically possible to produce pions/kaons. In not, assume nuclear excitation
141  // occurs and input particle is degraded in energy. No other particles are produced.
142  // If reaction is possible, find the correct number of pions/protons/neutrons
143  // produced using an interpolation to multiplicity data. Replace some pions or
144  // protons/neutrons by kaons or strange baryons according to the average
145  // multiplicity per Inelastic reaction.
146 
147  const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
148  const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
149  const G4double pOriginal = originalIncident->GetTotalMomentum()/MeV;
150  const G4double targetMass = targetParticle.GetMass()/MeV;
151  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
152  targetMass*targetMass +
153  2.0*targetMass*etOriginal );
154  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
155 
156  static G4ThreadLocal G4bool first = true;
157  const G4int numMul = 1200;
158  const G4int numMulA = 400;
159  const G4int numSec = 60;
160  static G4ThreadLocal G4double protmul[numMul], protnorm[numSec]; // proton constants
161  static G4ThreadLocal G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
162  static G4ThreadLocal G4double protmulA[numMulA], protnormA[numSec]; // proton constants
163  static G4ThreadLocal G4double neutmulA[numMulA], neutnormA[numSec]; // neutron constants
164  // np = number of pi+, nneg = number of pi-, nz = number of pi0
165  G4int counter, nt=0, np=0, nneg=0, nz=0;
166  G4double test;
167  const G4double c = 1.25;
168  const G4double b[2] = { 0.7, 0.7 };
169  if( first ) // compute normalization constants, this will only be Done once
170  {
171  first = false;
172  G4int i;
173  for( i=0; i<numMul; ++i )protmul[i] = 0.0;
174  for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
175  counter = -1;
176  for( np=0; np<(numSec/3); ++np )
177  {
178  for( nneg=std::max(0,np-2); nneg<=np; ++nneg )
179  {
180  for( nz=0; nz<numSec/3; ++nz )
181  {
182  if( ++counter < numMul )
183  {
184  nt = np+nneg+nz;
185  if( nt>0 && nt<=numSec )
186  {
187  protmul[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
188  protnorm[nt-1] += protmul[counter];
189  }
190  }
191  }
192  }
193  }
194  for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
195  for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
196  counter = -1;
197  for( np=0; np<numSec/3; ++np )
198  {
199  for( nneg=std::max(0,np-1); nneg<=(np+1); ++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  for( i=0; i<numSec; ++i )
216  {
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  for( np=2; np<(numSec/3); ++np )
227  {
228  nneg = np-2;
229  for( nz=0; nz<numSec/3; ++nz )
230  {
231  if( ++counter < numMulA )
232  {
233  nt = np+nneg+nz;
234  if( nt>1 && nt<=numSec )
235  {
236  protmulA[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
237  protnormA[nt-1] += protmulA[counter];
238  }
239  }
240  }
241  }
242  for( i=0; i<numMulA; ++i )neutmulA[i] = 0.0;
243  for( i=0; i<numSec; ++i )neutnormA[i] = 0.0;
244  counter = -1;
245  for( np=1; np<numSec/3; ++np )
246  {
247  nneg = np-1;
248  for( nz=0; nz<numSec/3; ++nz )
249  {
250  if( ++counter < numMulA )
251  {
252  nt = np+nneg+nz;
253  if( nt>1 && nt<=numSec )
254  {
255  neutmulA[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
256  neutnormA[nt-1] += neutmulA[counter];
257  }
258  }
259  }
260  }
261  for( i=0; i<numSec; ++i )
262  {
263  if( protnormA[i] > 0.0 )protnormA[i] = 1.0/protnormA[i];
264  if( neutnormA[i] > 0.0 )neutnormA[i] = 1.0/neutnormA[i];
265  }
266  } // end of initialization
267  const G4double expxu = 82.; // upper bound for arg. of exp
268  const G4double expxl = -expxu; // lower bound for arg. of exp
277  const G4double anhl[] = {1.00,1.00,1.00,1.00,1.00,1.00,1.00,1.00,0.97,0.88,
278  0.85,0.81,0.75,0.64,0.64,0.55,0.55,0.45,0.47,0.40,
279  0.39,0.36,0.33,0.10,0.01};
280  G4int iplab = G4int( pOriginal/GeV*10.0 );
281  if( iplab > 9 )iplab = G4int( (pOriginal/GeV- 1.0)*5.0 ) + 10;
282  if( iplab > 14 )iplab = G4int( pOriginal/GeV- 2.0 ) + 15;
283  if( iplab > 23 )iplab = G4int( (pOriginal/GeV-10.0)/10.0 ) + 23;
284  if( iplab > 24 )iplab = 24;
285  if( G4UniformRand() > anhl[iplab] )
286  {
287  if( availableEnergy <= aPiPlus->GetPDGMass()/MeV )
288  {
289  quasiElastic = true;
290  return;
291  }
292  G4double n, anpn;
293  GetNormalizationConstant( availableEnergy, n, anpn );
294  G4double ran = G4UniformRand();
295  G4double dum, excs = 0.0;
296  if( targetParticle.GetDefinition() == aProton )
297  {
298  counter = -1;
299  for( np=0; np<numSec/3 && ran>=excs; ++np )
300  {
301  for( nneg=std::max(0,np-2); nneg<=np && ran>=excs; ++nneg )
302  {
303  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
304  {
305  if( ++counter < numMul )
306  {
307  nt = np+nneg+nz;
308  if( nt>0 && nt<=numSec )
309  {
310  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
311  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
312  if( std::fabs(dum) < 1.0 )
313  {
314  if( test >= 1.0e-10 )excs += dum*test;
315  }
316  else
317  excs += dum*test;
318  }
319  }
320  }
321  }
322  }
323  if( ran >= excs ) // 3 previous loops continued to the end
324  {
325  quasiElastic = true;
326  return;
327  }
328  np--; nneg--; nz--;
329  G4int ncht = std::min( 3, std::max( 1, np-nneg+1 ) );
330  switch( ncht )
331  {
332  case 1:
333  break;
334  case 2:
335  if( G4UniformRand() < 0.5 )
336  {
337  targetParticle.SetDefinitionAndUpdateE( aNeutron );
338  targetHasChanged = true;
339  }
340  else
341  {
342  if( G4UniformRand() < 0.5 )
343  currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
344  else
345  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
346  incidentHasChanged = true;
347  }
348  break;
349  case 3:
350  if( G4UniformRand() < 0.5 )
351  currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
352  else
353  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
354  incidentHasChanged = true;
355  targetParticle.SetDefinitionAndUpdateE( aNeutron );
356  targetHasChanged = true;
357  break;
358  }
359  }
360  else // target must be a neutron
361  {
362  counter = -1;
363  for( np=0; np<numSec/3 && ran>=excs; ++np )
364  {
365  for( nneg=std::max(0,np-1); nneg<=(np+1) && ran>=excs; ++nneg )
366  {
367  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
368  {
369  if( ++counter < numMul )
370  {
371  nt = np+nneg+nz;
372  if( nt>0 && nt<=numSec )
373  {
374  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
375  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
376  if( std::fabs(dum) < 1.0 )
377  {
378  if( test >= 1.0e-10 )excs += dum*test;
379  }
380  else
381  excs += dum*test;
382  }
383  }
384  }
385  }
386  }
387  if( ran >= excs ) // 3 previous loops continued to the end
388  {
389  quasiElastic = true;
390  return;
391  }
392  np--; nneg--; nz--;
393  G4int ncht = std::min( 3, std::max( 1, np-nneg+2 ) );
394  switch( ncht )
395  {
396  case 1:
397  {
398  targetParticle.SetDefinitionAndUpdateE( aProton );
399  targetHasChanged = true;
400  }
401  break;
402  case 2:
403  if( G4UniformRand() < 0.5 )
404  {
405  if( G4UniformRand() < 0.5 )
406  {
407  currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
408  incidentHasChanged = true;
409  targetParticle.SetDefinitionAndUpdateE( aProton );
410  targetHasChanged = true;
411  }
412  }
413  else
414  {
415  if( G4UniformRand() < 0.5 )
416  {
417  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
418  incidentHasChanged = true;
419  targetParticle.SetDefinitionAndUpdateE( aProton );
420  targetHasChanged = true;
421  }
422  }
423  break;
424  case 3:
425  if( G4UniformRand() < 0.5 )
426  currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
427  else
428  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
429  incidentHasChanged = true;
430  break;
431  }
432  }
433  }
434  else // random number <= anhl[iplab]
435  {
436  if( centerofmassEnergy <= aPiPlus->GetPDGMass()/MeV+aKaonPlus->GetPDGMass()/MeV )
437  {
438  quasiElastic = true;
439  return;
440  }
441  G4double n, anpn;
442  GetNormalizationConstant( -centerofmassEnergy, n, anpn );
443  G4double ran = G4UniformRand();
444  G4double dum, excs = 0.0;
445  if( targetParticle.GetDefinition() == aProton )
446  {
447  counter = -1;
448  for( np=2; np<numSec/3 && ran>=excs; ++np )
449  {
450  nneg=np-2;
451  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
452  {
453  if( ++counter < numMulA )
454  {
455  nt = np+nneg+nz;
456  if( nt>1 && nt<=numSec )
457  {
458  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
459  dum = (pi/anpn)*nt*protmulA[counter]*protnormA[nt-1]/(2.0*n*n);
460  if( std::fabs(dum) < 1.0 )
461  {
462  if( test >= 1.0e-10 )excs += dum*test;
463  }
464  else
465  excs += dum*test;
466  }
467  }
468  }
469  }
470  if( ran >= excs ) // 3 previous loops continued to the end
471  {
472  quasiElastic = true;
473  return;
474  }
475  np--; nz--;
476  }
477  else // target must be a neutron
478  {
479  counter = -1;
480  for( np=1; np<numSec/3 && ran>=excs; ++np )
481  {
482  nneg = np-1;
483  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
484  {
485  if( ++counter < numMulA )
486  {
487  nt = np+nneg+nz;
488  if( nt>1 && nt<=numSec )
489  {
490  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
491  dum = (pi/anpn)*nt*neutmulA[counter]*neutnormA[nt-1]/(2.0*n*n);
492  if( std::fabs(dum) < 1.0 )
493  {
494  if( test >= 1.0e-10 )excs += dum*test;
495  }
496  else
497  excs += dum*test;
498  }
499  }
500  }
501  }
502  if( ran >= excs ) // 3 previous loops continued to the end
503  {
504  quasiElastic = true;
505  return;
506  }
507  np--; nz--;
508  }
509  if( nz > 0 )
510  {
511  if( nneg > 0 )
512  {
513  if( G4UniformRand() < 0.5 )
514  {
515  vec.Initialize( 1 );
517  p->SetDefinition( aKaonMinus );
518  (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
519  vec.SetElement( vecLen++, p );
520  --nneg;
521  }
522  else // random number >= 0.5
523  {
524  vec.Initialize( 1 );
526  p->SetDefinition( aKaonZL );
527  (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
528  vec.SetElement( vecLen++, p );
529  --nz;
530  }
531  }
532  else // nneg == 0
533  {
534  vec.Initialize( 1 );
536  p->SetDefinition( aKaonZL );
537  (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
538  vec.SetElement( vecLen++, p );
539  --nz;
540  }
541  }
542  else // nz == 0
543  {
544  if( nneg > 0 )
545  {
546  vec.Initialize( 1 );
548  p->SetDefinition( aKaonMinus );
549  (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
550  vec.SetElement( vecLen++, p );
551  --nneg;
552  }
553  }
554  currentParticle.SetMass( 0.0 );
555  targetParticle.SetMass( 0.0 );
556  }
557 
558  SetUpPions( np, nneg, nz, vec, vecLen );
559  return;
560 }
561 
562  /* end of file */
563 
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
static constexpr double MeV
Definition: G4SIunits.hh:214
G4double Cinema(G4double kineticEnergy)
Definition: G4Nucleus.cc:382
void SetMomentumChange(const G4ThreeVector &aV)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
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
const G4ParticleDefinition * GetDefinition() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double GetTotalEnergy() const
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
#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
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:241
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)
static G4KaonZeroLong * KaonZeroLong()
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
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
#define G4UniformRand()
Definition: Randomize.hh:53
G4double GetTotalMomentum() const
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)
static G4AntiSigmaZero * AntiSigmaZero()
const G4LorentzVector & Get4Momentum() const
Hep3Vector unit() const
G4double GetKineticEnergy() const
void Cascade(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool &quasiElastic)
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
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
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)
static G4AntiLambda * AntiLambda()