Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4RPGAntiLambdaInelastic.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: G4RPGAntiLambdaInelastic.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  // Choose the target particle
42 
43  G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
44 
45  if( verboseLevel > 1 )
46  {
47  const G4Material *targetMaterial = aTrack.GetMaterial();
48  G4cout << "G4RPGAntiLambdaInelastic::ApplyYourself called" << G4endl;
49  G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
50  G4cout << "target material = " << targetMaterial->GetName() << ", ";
51  G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
52  << G4endl;
53  }
54 
55  // Fermi motion and evaporation
56  // As of Geant3, the Fermi energy calculation had not been Done
57 
58  G4double ek = originalIncident->GetKineticEnergy()/MeV;
59  G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
60  G4ReactionProduct modifiedOriginal;
61  modifiedOriginal = *originalIncident;
62 
63  G4double tkin = targetNucleus.Cinema( ek );
64  ek += tkin;
65  modifiedOriginal.SetKineticEnergy( ek*MeV );
66  G4double et = ek + amas;
67  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
68  G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
69  if( pp > 0.0 )
70  {
71  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
72  modifiedOriginal.SetMomentum( momentum * (p/pp) );
73  }
74  //
75  // calculate black track energies
76  //
77  tkin = targetNucleus.EvaporationEffects( ek );
78  ek -= tkin;
79  modifiedOriginal.SetKineticEnergy( ek*MeV );
80  et = ek + amas;
81  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
82  pp = modifiedOriginal.GetMomentum().mag()/MeV;
83  if( pp > 0.0 )
84  {
85  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
86  modifiedOriginal.SetMomentum( momentum * (p/pp) );
87  }
88 
89  G4ReactionProduct currentParticle = modifiedOriginal;
90  G4ReactionProduct targetParticle;
91  targetParticle = *originalTarget;
92  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
93  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
94  G4bool incidentHasChanged = false;
95  G4bool targetHasChanged = false;
96  G4bool quasiElastic = false;
97  G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
98  G4int vecLen = 0;
99  vec.Initialize( 0 );
100 
101  const G4double cutOff = 0.1;
102  const G4double anni = std::min( 1.3*currentParticle.GetTotalMomentum()/GeV, 0.4 );
103  if( (originalIncident->GetKineticEnergy()/MeV > cutOff) || (G4UniformRand() > anni) )
104  Cascade( vec, vecLen,
105  originalIncident, currentParticle, targetParticle,
106  incidentHasChanged, targetHasChanged, quasiElastic );
107 
108  CalculateMomenta( vec, vecLen,
109  originalIncident, originalTarget, modifiedOriginal,
110  targetNucleus, currentParticle, targetParticle,
111  incidentHasChanged, targetHasChanged, quasiElastic );
112 
113  SetUpChange( vec, vecLen,
114  currentParticle, targetParticle,
115  incidentHasChanged );
116 
117  delete originalTarget;
118  return &theParticleChange;
119 }
120 
121 
122 void G4RPGAntiLambdaInelastic::Cascade(
124  G4int &vecLen,
125  const G4HadProjectile *originalIncident,
126  G4ReactionProduct &currentParticle,
127  G4ReactionProduct &targetParticle,
128  G4bool &incidentHasChanged,
129  G4bool &targetHasChanged,
130  G4bool &quasiElastic )
131 {
132  // Derived from H. Fesefeldt's original FORTRAN code CASAL0
133  // AntiLambda undergoes interaction with nucleon within a nucleus. Check if it is
134  // energetically possible to produce pions/kaons. In not, assume nuclear excitation
135  // occurs and input particle is degraded in energy. No other particles are produced.
136  // If reaction is possible, find the correct number of pions/protons/neutrons
137  // produced using an interpolation to multiplicity data. Replace some pions or
138  // protons/neutrons by kaons or strange baryons according to the average
139  // multiplicity per Inelastic reaction.
140 
141  const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
142  const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
143  const G4double targetMass = targetParticle.GetMass()/MeV;
144  const G4double pOriginal = originalIncident->GetTotalMomentum()/GeV;
145  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
146  targetMass*targetMass +
147  2.0*targetMass*etOriginal );
148  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
149 
150  static G4ThreadLocal G4bool first = true;
151  const G4int numMul = 1200;
152  const G4int numMulA = 400;
153  const G4int numSec = 60;
154  static G4ThreadLocal G4double protmul[numMul], protnorm[numSec]; // proton constants
155  static G4ThreadLocal G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
156  static G4ThreadLocal G4double protmulA[numMulA], protnormA[numSec]; // proton constants
157  static G4ThreadLocal G4double neutmulA[numMulA], neutnormA[numSec]; // neutron constants
158  // np = number of pi+, nneg = number of pi-, nz = number of pi0
159  G4int nt=0, np=0, nneg=0, nz=0;
160  G4double test;
161  const G4double c = 1.25;
162  const G4double b[] = { 0.7, 0.7 };
163  if( first ) // compute normalization constants, this will only be Done once
164  {
165  first = false;
166  G4int i;
167  for( i=0; i<numMul; ++i )protmul[i] = 0.0;
168  for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
169  G4int counter = -1;
170  for( np=0; np<(numSec/3); ++np )
171  {
172  for( nneg=std::max(0,np-2); nneg<=(np+1); ++nneg )
173  {
174  for( nz=0; nz<numSec/3; ++nz )
175  {
176  if( ++counter < numMul )
177  {
178  nt = np+nneg+nz;
179  if( nt>0 && nt<=numSec )
180  {
181  protmul[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
182  protnorm[nt-1] += protmul[counter];
183  }
184  }
185  }
186  }
187  }
188  for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
189  for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
190  counter = -1;
191  for( np=0; np<numSec/3; ++np )
192  {
193  for( nneg=std::max(0,np-1); nneg<=(np+2); ++nneg )
194  {
195  for( nz=0; nz<numSec/3; ++nz )
196  {
197  if( ++counter < numMul )
198  {
199  nt = np+nneg+nz;
200  if( nt>0 && nt<=numSec )
201  {
202  neutmul[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
203  neutnorm[nt-1] += neutmul[counter];
204  }
205  }
206  }
207  }
208  }
209  for( i=0; i<numSec; ++i )
210  {
211  if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
212  if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
213  }
214  //
215  // do the same for annihilation channels
216  //
217  for( i=0; i<numMulA; ++i )protmulA[i] = 0.0;
218  for( i=0; i<numSec; ++i )protnormA[i] = 0.0;
219  counter = -1;
220  for( np=1; np<(numSec/3); ++np )
221  {
222  nneg = np-1;
223  for( nz=0; nz<numSec/3; ++nz )
224  {
225  if( ++counter < numMulA )
226  {
227  nt = np+nneg+nz;
228  if( nt>1 && nt<=numSec )
229  {
230  protmulA[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
231  protnormA[nt-1] += protmulA[counter];
232  }
233  }
234  }
235  }
236  for( i=0; i<numMulA; ++i )neutmulA[i] = 0.0;
237  for( i=0; i<numSec; ++i )neutnormA[i] = 0.0;
238  counter = -1;
239  for( np=0; np<numSec/3; ++np )
240  {
241  nneg = np;
242  for( nz=0; nz<numSec/3; ++nz )
243  {
244  if( ++counter < numMulA )
245  {
246  nt = np+nneg+nz;
247  if( nt>1 && nt<=numSec )
248  {
249  neutmulA[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
250  neutnormA[nt-1] += neutmulA[counter];
251  }
252  }
253  }
254  }
255  for( i=0; i<numSec; ++i )
256  {
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  const G4double expxu = 82.; // upper bound for arg. of exp
262  const G4double expxl = -expxu; // lower bound for arg. of exp
263 
275  const G4double anhl[] = {1.00,1.00,1.00,1.00,1.00,1.00,1.00,1.00,0.97,0.88,
276  0.85,0.81,0.75,0.64,0.64,0.55,0.55,0.45,0.47,0.40,
277  0.39,0.36,0.33,0.10,0.01};
278  G4int iplab = G4int( pOriginal*10.0 );
279  if( iplab > 9 )iplab = G4int( (pOriginal- 1.0)*5.0 ) + 10;
280  if( iplab > 14 )iplab = G4int( pOriginal- 2.0 ) + 15;
281  if( iplab > 22 )iplab = G4int( (pOriginal-10.0)/10.0 ) + 23;
282  if( iplab > 24 )iplab = 24;
283  if( G4UniformRand() > anhl[iplab] )
284  {
285  if( availableEnergy <= aPiPlus->GetPDGMass()/MeV )
286  { // not energetically possible to produce pion(s)
287  quasiElastic = true;
288  return;
289  }
290  G4double n, anpn;
291  GetNormalizationConstant( availableEnergy, n, anpn );
292  G4double ran = G4UniformRand();
293  G4double dum, excs = 0.0;
294  if( targetParticle.GetDefinition() == aProton )
295  {
296  G4int counter = -1;
297  for( np=0; np<numSec/3 && ran>=excs; ++np )
298  {
299  for( nneg=std::max(0,np-2); nneg<=(np+1) && ran>=excs; ++nneg )
300  {
301  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
302  {
303  if( ++counter < numMul )
304  {
305  nt = np+nneg+nz;
306  if( nt>0 && nt<=numSec )
307  {
308  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
309  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
310  if( std::fabs(dum) < 1.0 )
311  {
312  if( test >= 1.0e-10 )excs += dum*test;
313  }
314  else
315  excs += dum*test;
316  }
317  }
318  }
319  }
320  }
321  if( ran >= excs ) // 3 previous loops continued to the end
322  {
323  quasiElastic = true;
324  return;
325  }
326  np--; nneg--; nz--;
327  G4int ncht = std::min( 4, std::max( 1, np-nneg+2 ) );
328  switch( ncht )
329  {
330  case 1:
331  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaMinus );
332  incidentHasChanged = true;
333  break;
334  case 2:
335  if( G4UniformRand() < 0.5 )
336  {
337  if( G4UniformRand() < 0.5 )
338  {
339  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
340  incidentHasChanged = true;
341  }
342  else
343  {
344  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaMinus );
345  incidentHasChanged = true;
346  targetParticle.SetDefinitionAndUpdateE( aNeutron );
347  targetHasChanged = true;
348  }
349  }
350  else
351  {
352  if( G4UniformRand() >= 0.5 )
353  {
354  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaMinus );
355  incidentHasChanged = true;
356  targetParticle.SetDefinitionAndUpdateE( aNeutron );
357  targetHasChanged = true;
358  }
359  }
360  break;
361  case 3:
362  if( G4UniformRand() < 0.5 )
363  {
364  if( G4UniformRand() < 0.5 )
365  {
366  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
367  incidentHasChanged = true;
368  targetParticle.SetDefinitionAndUpdateE( aNeutron );
369  targetHasChanged = true;
370  }
371  else
372  {
373  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaPlus );
374  incidentHasChanged = true;
375  }
376  }
377  else
378  {
379  if( G4UniformRand() < 0.5 )
380  {
381  targetParticle.SetDefinitionAndUpdateE( aNeutron );
382  targetHasChanged = true;
383  }
384  else
385  {
386  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaPlus );
387  incidentHasChanged = true;
388  }
389  }
390  break;
391  case 4:
392  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaPlus );
393  incidentHasChanged = true;
394  targetParticle.SetDefinitionAndUpdateE( aNeutron );
395  targetHasChanged = true;
396  break;
397  }
398  }
399  else // target must be a neutron
400  {
401  G4int counter = -1;
402  for( np=0; np<numSec/3 && ran>=excs; ++np )
403  {
404  for( nneg=std::max(0,np-1); nneg<=(np+2) && ran>=excs; ++nneg )
405  {
406  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
407  {
408  if( ++counter < numMul )
409  {
410  nt = np+nneg+nz;
411  if( nt>0 && nt<=numSec )
412  {
413  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
414  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
415  if( std::fabs(dum) < 1.0 )
416  {
417  if( test >= 1.0e-10 )excs += dum*test;
418  }
419  else
420  excs += dum*test;
421  }
422  }
423  }
424  }
425  }
426  if( ran >= excs ) // 3 previous loops continued to the end
427  {
428  quasiElastic = true;
429  return;
430  }
431  np--; nneg--; nz--;
432  G4int ncht = std::min( 4, std::max( 1, np-nneg+3 ) );
433  switch( ncht )
434  {
435  case 1:
436  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaMinus );
437  incidentHasChanged = true;
438  targetParticle.SetDefinitionAndUpdateE( aProton );
439  targetHasChanged = true;
440  break;
441  case 2:
442  if( G4UniformRand() < 0.5 )
443  {
444  if( G4UniformRand() < 0.5 )
445  {
446  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
447  incidentHasChanged = true;
448  targetParticle.SetDefinitionAndUpdateE( aProton );
449  targetHasChanged = true;
450  }
451  else
452  {
453  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaMinus );
454  incidentHasChanged = true;
455  }
456  }
457  else
458  {
459  if( G4UniformRand() < 0.5 )
460  {
461  targetParticle.SetDefinitionAndUpdateE( aProton );
462  targetHasChanged = true;
463  }
464  else
465  {
466  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaMinus );
467  incidentHasChanged = true;
468  }
469  }
470  break;
471  case 3:
472  if( G4UniformRand() < 0.5 )
473  {
474  if( G4UniformRand() < 0.5 )
475  {
476  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
477  incidentHasChanged = true;
478  }
479  else
480  {
481  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaPlus );
482  incidentHasChanged = true;
483  targetParticle.SetDefinitionAndUpdateE( aProton );
484  targetHasChanged = true;
485  }
486  }
487  else
488  {
489  if( G4UniformRand() >= 0.5 )
490  {
491  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaPlus );
492  incidentHasChanged = true;
493  targetParticle.SetDefinitionAndUpdateE( aProton );
494  targetHasChanged = true;
495  }
496  }
497  break;
498  default:
499  currentParticle.SetDefinitionAndUpdateE( anAntiSigmaPlus );
500  incidentHasChanged = true;
501  break;
502  }
503  }
504  }
505  else // random number <= anhl[iplab]
506  {
507  if( centerofmassEnergy <= aPiPlus->GetPDGMass()/MeV+aKaonPlus->GetPDGMass()/MeV )
508  {
509  quasiElastic = true;
510  return;
511  }
512  //
513  // annihilation channels
514  //
515  G4double n, anpn;
516  GetNormalizationConstant( -centerofmassEnergy, n, anpn );
517  G4double ran = G4UniformRand();
518  G4double dum, excs = 0.0;
519  if( targetParticle.GetDefinition() == aProton )
520  {
521  G4int counter = -1;
522  for( np=1; np<numSec/3 && ran>=excs; ++np )
523  {
524  nneg = np-1;
525  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
526  {
527  if( ++counter < numMulA )
528  {
529  nt = np+nneg+nz;
530  if( nt>1 && nt<=numSec )
531  {
532  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
533  dum = (pi/anpn)*nt*protmulA[counter]*protnormA[nt-1]/(2.0*n*n);
534  if( std::fabs(dum) < 1.0 )
535  {
536  if( test >= 1.0e-10 )excs += dum*test;
537  }
538  else
539  excs += dum*test;
540  }
541  }
542  }
543  }
544  }
545  else // target must be a neutron
546  {
547  G4int counter = -1;
548  for( np=0; np<numSec/3 && ran>=excs; ++np )
549  {
550  nneg = np;
551  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
552  {
553  if( ++counter < numMulA )
554  {
555  nt = np+nneg+nz;
556  if( nt>1 && nt<=numSec )
557  {
558  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
559  dum = (pi/anpn)*nt*neutmulA[counter]*neutnormA[nt-1]/(2.0*n*n);
560  if( std::fabs(dum) < 1.0 )
561  {
562  if( test >= 1.0e-10 )excs += dum*test;
563  }
564  else
565  excs += dum*test;
566  }
567  }
568  }
569  }
570  }
571  if( ran >= excs ) // 3 previous loops continued to the end
572  {
573  quasiElastic = true;
574  return;
575  }
576  np--; nz--;
577  currentParticle.SetMass( 0.0 );
578  targetParticle.SetMass( 0.0 );
579  }
580  SetUpPions( np, nneg, nz, vec, vecLen );
581  if( currentParticle.GetMass() == 0.0 )
582  {
583  if( nz == 0 )
584  {
585  if( nneg > 0 )
586  {
587  for( G4int i=0; i<vecLen; ++i )
588  {
589  if( vec[i]->GetDefinition() == aPiMinus )
590  {
591  vec[i]->SetDefinitionAndUpdateE( aKaonMinus );
592  break;
593  }
594  }
595  }
596  }
597  else // nz > 0
598  {
599  if( nneg == 0 )
600  {
601  for( G4int i=0; i<vecLen; ++i )
602  {
603  if( vec[i]->GetDefinition() == aPiZero )
604  {
605  vec[i]->SetDefinitionAndUpdateE( aKaonZL );
606  break;
607  }
608  }
609  }
610  else // nneg > 0
611  {
612  if( G4UniformRand() < 0.5 )
613  {
614  if( nneg > 0 )
615  {
616  for( G4int i=0; i<vecLen; ++i )
617  {
618  if( vec[i]->GetDefinition() == aPiMinus )
619  {
620  vec[i]->SetDefinitionAndUpdateE( aKaonMinus );
621  break;
622  }
623  }
624  }
625  }
626  else // random number >= 0.5
627  {
628  for( G4int i=0; i<vecLen; ++i )
629  {
630  if( vec[i]->GetDefinition() == aPiZero )
631  {
632  vec[i]->SetDefinitionAndUpdateE( aKaonZL );
633  break;
634  }
635  }
636  }
637  }
638  }
639  }
640  return;
641 }
642 
643  /* end of file */
644 
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:278
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
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
static G4AntiSigmaPlus * AntiSigmaPlus()
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
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
G4double Pmltpc(G4int np, G4int nm, G4int nz, G4int n, G4double b, G4double c)
static G4AntiSigmaMinus * AntiSigmaMinus()
const G4ParticleDefinition * GetDefinition() const
void SetMass(const G4double mas)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
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
static G4PionZero * PionZero()
Definition: G4PionZero.cc:108
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double Cinema(G4double kineticEnergy)
Definition: G4Nucleus.cc:382
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
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
G4double GetMass() const
G4double GetTotalMomentum() const
G4double GetTotalEnergy() const