Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4RPGAntiSigmaPlusInelastic.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: G4RPGAntiSigmaPlusInelastic.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 << "G4RPGAntiSigmaPlusInelastic::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 
128 void G4RPGAntiSigmaPlusInelastic::Cascade(
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 CASASP
139  // AntiSigmaPlus 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[] = { 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-1); nneg<=(np+1); ++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=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  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=1; np<(numSec/3); ++np )
227  {
228  nneg = np;
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=0; 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 
268  const G4double expxu = 82.; // upper bound for arg. of exp
269  const G4double expxl = -expxu; // lower bound for arg. of exp
278  const G4double anhl[] = {1.00,1.00,1.00,1.00,1.00,1.00,1.00,1.00,0.97,0.88,
279  0.85,0.81,0.75,0.64,0.64,0.55,0.55,0.45,0.47,0.40,
280  0.39,0.36,0.33,0.10,0.01};
281  G4int iplab = G4int( pOriginal/GeV*10.0 );
282  if( iplab > 9 )iplab = G4int( (pOriginal/GeV- 1.0)*5.0 ) + 10;
283  if( iplab > 14 )iplab = G4int( pOriginal/GeV- 2.0 ) + 15;
284  if( iplab > 22 )iplab = G4int( (pOriginal/GeV-10.0)/10.0 ) + 23;
285  if( iplab > 24 )iplab = 24;
286  if( G4UniformRand() > anhl[iplab] )
287  {
288  if( availableEnergy <= aPiPlus->GetPDGMass()/MeV )
289  {
290  quasiElastic = true;
291  return;
292  }
293  G4double n, anpn;
294  GetNormalizationConstant( availableEnergy, n, anpn );
295  G4double ran = G4UniformRand();
296  G4double dum, excs = 0.0;
297  if( targetParticle.GetDefinition() == aProton )
298  {
299  counter = -1;
300  for( np=0; np<numSec/3 && ran>=excs; ++np )
301  {
302  for( nneg=std::max(0,np-1); nneg<=(np+1) && ran>=excs; ++nneg )
303  {
304  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
305  {
306  if( ++counter < numMul )
307  {
308  nt = np+nneg+nz;
309  if( (nt>0) && (nt<=numSec) )
310  {
311  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
312  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
313  if( std::fabs(dum) < 1.0 )
314  {
315  if( test >= 1.0e-10 )excs += dum*test;
316  }
317  else
318  excs += dum*test;
319  }
320  }
321  }
322  }
323  }
324  if( ran >= excs ) // 3 previous loops continued to the end
325  {
326  quasiElastic = true;
327  return;
328  }
329  np--; nneg--; nz--;
330  G4int ncht = std::min( 3, std::max( 1, np-nneg+2 ) );
331  switch( ncht )
332  {
333  case 1:
334  if( G4UniformRand() < 0.5 )
335  currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
336  else
337  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
338  incidentHasChanged = true;
339  break;
340  case 2:
341  if( G4UniformRand() >= 0.5 )
342  {
343  if( G4UniformRand() < 0.5 )
344  currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
345  else
346  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
347  incidentHasChanged = true;
348  }
349  targetParticle.SetDefinitionAndUpdateE( aNeutron );
350  targetHasChanged = true;
351  break;
352  case 3:
353  targetParticle.SetDefinitionAndUpdateE( aNeutron );
354  targetHasChanged = true;
355  break;
356  }
357  }
358  else // target must be a neutron
359  {
360  counter = -1;
361  for( np=0; np<numSec/3 && ran>=excs; ++np )
362  {
363  for( nneg=np; nneg<=(np+2) && ran>=excs; ++nneg )
364  {
365  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
366  {
367  if( ++counter < numMul )
368  {
369  nt = np+nneg+nz;
370  if( (nt>0) && (nt<=numSec) )
371  {
372  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
373  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
374  if( std::fabs(dum) < 1.0 )
375  {
376  if( test >= 1.0e-10 )excs += dum*test;
377  }
378  else
379  excs += dum*test;
380  }
381  }
382  }
383  }
384  }
385  if( ran >= excs ) // 3 previous loops continued to the end
386  {
387  quasiElastic = true;
388  return;
389  }
390  np--; nneg--; nz--;
391  G4int ncht = std::min( 3, std::max( 1, np-nneg+3 ) );
392  switch( ncht )
393  {
394  case 1:
395  if( G4UniformRand() < 0.5 )
396  currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
397  else
398  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
399  incidentHasChanged = true;
400  targetParticle.SetDefinitionAndUpdateE( aProton );
401  targetHasChanged = true;
402  break;
403  case 2:
404  if( G4UniformRand() < 0.5 )
405  {
406  if( G4UniformRand() < 0.5 )
407  {
408  currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
409  incidentHasChanged = true;
410  }
411  else
412  {
413  targetParticle.SetDefinitionAndUpdateE( aProton );
414  targetHasChanged = true;
415  }
416  }
417  else
418  {
419  if( G4UniformRand() < 0.5 )
420  {
421  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
422  incidentHasChanged = true;
423  }
424  else
425  {
426  targetParticle.SetDefinitionAndUpdateE( aProton );
427  targetHasChanged = true;
428  }
429  }
430  break;
431  case 3:
432  break;
433  }
434  }
435  }
436  else // random number <= anhl[iplab]
437  {
438  if( centerofmassEnergy <= aPiPlus->GetPDGMass()/MeV+aKaonPlus->GetPDGMass()/MeV )
439  {
440  quasiElastic = true;
441  return;
442  }
443  G4double n, anpn;
444  GetNormalizationConstant( -centerofmassEnergy, n, anpn );
445  G4double ran = G4UniformRand();
446  G4double dum, excs = 0.0;
447  if( targetParticle.GetDefinition() == aProton )
448  {
449  counter = -1;
450  for( np=1; np<numSec/3 && ran>=excs; ++np )
451  {
452  nneg = np;
453  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
454  {
455  if( ++counter < numMulA )
456  {
457  nt = np+nneg+nz;
458  if( nt>1 && nt<=numSec )
459  {
460  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
461  dum = (pi/anpn)*nt*protmulA[counter]*protnormA[nt-1]/(2.0*n*n);
462  if( std::fabs(dum) < 1.0 )
463  {
464  if( test >= 1.0e-10 )excs += dum*test;
465  }
466  else
467  excs += dum*test;
468  }
469  }
470  }
471  }
472  if( ran >= excs ) // 3 previous loops continued to the end
473  {
474  quasiElastic = true;
475  return;
476  }
477  np--; nz--;
478  }
479  else // target must be a neutron
480  {
481  counter = -1;
482  for( np=0; np<numSec/3 && ran>=excs; ++np )
483  {
484  nneg = np+1;
485  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
486  {
487  if( ++counter < numMulA )
488  {
489  nt = np+nneg+nz;
490  if( nt>1 && nt<=numSec )
491  {
492  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
493  dum = (pi/anpn)*nt*neutmulA[counter]*neutnormA[nt-1]/(2.0*n*n);
494  if( std::fabs(dum) < 1.0 )
495  {
496  if( test >= 1.0e-10 )excs += dum*test;
497  }
498  else
499  excs += dum*test;
500  }
501  }
502  }
503  }
504  if( ran >= excs ) // 3 previous loops continued to the end
505  {
506  quasiElastic = true;
507  return;
508  }
509  np--; nz--;
510  }
511  if( nz > 0 )
512  {
513  if( nneg > 0 )
514  {
515  if( G4UniformRand() < 0.5 )
516  {
517  vec.Initialize( 1 );
519  p->SetDefinition( aKaonMinus );
520  (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
521  vec.SetElement( vecLen++, p );
522  --nneg;
523  }
524  else
525  {
526  vec.Initialize( 1 );
528  p->SetDefinition( aKaonZL );
529  (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
530  vec.SetElement( vecLen++, p );
531  --nz;
532  }
533  }
534  else // nneg == 0
535  {
536  vec.Initialize( 1 );
538  p->SetDefinition( aKaonZL );
539  (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
540  vec.SetElement( vecLen++, p );
541  --nz;
542  }
543  }
544  else // nz == 0
545  {
546  if( nneg > 0 )
547  {
548  vec.Initialize( 1 );
550  p->SetDefinition( aKaonMinus );
551  (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
552  vec.SetElement( vecLen++, p );
553  --nneg;
554  }
555  }
556  currentParticle.SetMass( 0.0 );
557  targetParticle.SetMass( 0.0 );
558  }
559 
560  SetUpPions( np, nneg, nz, vec, vecLen );
561  return;
562 }
563 
564  /* end of file */
565 
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:278
G4double GetTotalMomentum() const
void SetUpChange(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
void SetKineticEnergy(const G4double en)
void SetMomentum(const G4double x, const G4double y, const G4double z)
const char * p
Definition: xmltok.h:285
const G4String & GetName() const
Definition: G4Material.hh:178
static G4KaonZeroLong * KaonZeroLong()
void SetSide(const G4int sid)
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)
G4ParticleDefinition * GetDefinition() const
#define G4ThreadLocal
Definition: tls.hh:89
void Initialize(G4int items)
Definition: G4FastVector.hh:63
int G4int
Definition: G4Types.hh:78
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:241
const G4String & GetParticleName() const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
void SetStatusChange(G4HadFinalStateStatus aS)
G4double Pmltpc(G4int np, G4int nm, G4int nz, G4int n, G4double b, G4double c)
const G4ParticleDefinition * GetDefinition() const
void SetMass(const G4double mas)
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
const G4ParticleDefinition * GetDefinition() const
bool G4bool
Definition: G4Types.hh:79
G4double GetKineticEnergy() const
G4double ek
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
void SetDefinitionAndUpdateE(const G4ParticleDefinition *aParticleDefinition)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
const G4LorentzVector & Get4Momentum() const
G4double GetKineticEnergy() const
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
void SetEnergyChange(G4double anEnergy)
G4double GetPDGMass() const
static G4AntiLambda * AntiLambda()
T max(const T t1, const T t2)
brief Return the largest of the two arguments
Hep3Vector unit() const
G4double Cinema(G4double kineticEnergy)
Definition: G4Nucleus.cc:382
static G4AntiSigmaZero * AntiSigmaZero()
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static constexpr double GeV
Definition: G4SIunits.hh:217
G4ThreeVector GetMomentum() const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
const G4Material * GetMaterial() const
static constexpr double pi
Definition: G4SIunits.hh:75
void GetNormalizationConstant(const G4double availableEnergy, G4double &n, G4double &anpn)
double G4double
Definition: G4Types.hh:76
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
void SetUpPions(const G4int np, const G4int nm, const G4int nz, G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen)
double mag() const
void SetMomentumChange(const G4ThreeVector &aV)
G4double GetMass() const
G4double GetTotalMomentum() const
G4double GetTotalEnergy() const