Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4RPGKMinusInelastic.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: G4RPGKMinusInelastic.cc 94214 2015-11-09 08:18:05Z gcosmo $
27 //
28 
29 #include "G4RPGKMinusInelastic.hh"
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  // create the target particle
49 
50  G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
51  G4ReactionProduct targetParticle( originalTarget->GetDefinition() );
52 
53  if( verboseLevel > 1 )
54  {
55  const G4Material *targetMaterial = aTrack.GetMaterial();
56  G4cout << "G4RPGKMinusInelastic::ApplyYourself called" << G4endl;
57  G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy() << "MeV, ";
58  G4cout << "target material = " << targetMaterial->GetName() << ", ";
59  G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
60  << G4endl;
61  }
62 
63  G4ReactionProduct currentParticle(originalIncident->GetDefinition() );
64  currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
65  currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
66 
67  // Fermi motion and evaporation
68  // As of Geant3, the Fermi energy calculation had not been Done
69 
70  G4double ek = originalIncident->GetKineticEnergy();
71  G4double amas = originalIncident->GetDefinition()->GetPDGMass();
72 
73  G4double tkin = targetNucleus.Cinema( ek );
74  ek += tkin;
75  currentParticle.SetKineticEnergy( ek );
76  G4double et = ek + amas;
77  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
78  G4double pp = currentParticle.GetMomentum().mag();
79  if( pp > 0.0 )
80  {
81  G4ThreeVector momentum = currentParticle.GetMomentum();
82  currentParticle.SetMomentum( momentum * (p/pp) );
83  }
84 
85  // calculate black track energies
86 
87  tkin = targetNucleus.EvaporationEffects( ek );
88  ek -= tkin;
89  currentParticle.SetKineticEnergy( ek );
90  et = ek + amas;
91  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
92  pp = currentParticle.GetMomentum().mag();
93  if( pp > 0.0 )
94  {
95  G4ThreeVector momentum = currentParticle.GetMomentum();
96  currentParticle.SetMomentum( momentum * (p/pp) );
97  }
98 
99  G4ReactionProduct modifiedOriginal = currentParticle;
100 
101  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
102  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
103  G4bool incidentHasChanged = false;
104  G4bool targetHasChanged = false;
105  G4bool quasiElastic = false;
106  G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
107  G4int vecLen = 0;
108  vec.Initialize( 0 );
109 
110  const G4double cutOff = 0.1*MeV;
111  if( currentParticle.GetKineticEnergy() > cutOff )
112  Cascade( vec, vecLen,
113  originalIncident, currentParticle, targetParticle,
114  incidentHasChanged, targetHasChanged, quasiElastic );
115 
116  CalculateMomenta( vec, vecLen,
117  originalIncident, originalTarget, modifiedOriginal,
118  targetNucleus, currentParticle, targetParticle,
119  incidentHasChanged, targetHasChanged, quasiElastic );
120 
121  SetUpChange( vec, vecLen,
122  currentParticle, targetParticle,
123  incidentHasChanged );
124 
125  delete originalTarget;
126  return &theParticleChange;
127 }
128 
129 
130 void G4RPGKMinusInelastic::Cascade(
132  G4int& vecLen,
133  const G4HadProjectile *originalIncident,
134  G4ReactionProduct &currentParticle,
135  G4ReactionProduct &targetParticle,
136  G4bool &incidentHasChanged,
137  G4bool &targetHasChanged,
138  G4bool &quasiElastic )
139 {
140  // Derived from H. Fesefeldt's original FORTRAN code CASKM
141  //
142  // K- undergoes interaction with nucleon within a nucleus. Check if it is
143  // energetically possible to produce pions/kaons. In not, assume nuclear excitation
144  // occurs and input particle is degraded in energy. No other particles are produced.
145  // If reaction is possible, find the correct number of pions/protons/neutrons
146  // produced using an interpolation to multiplicity data. Replace some pions or
147  // protons/neutrons by kaons or strange baryons according to the average
148  // multiplicity per Inelastic reaction.
149  //
150  const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass();
151  const G4double etOriginal = originalIncident->GetTotalEnergy();
152  const G4double pOriginal = originalIncident->GetTotalMomentum();
153  const G4double targetMass = targetParticle.GetMass();
154  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
155  targetMass*targetMass +
156  2.0*targetMass*etOriginal );
157  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
158 
159  static G4ThreadLocal G4bool first = true;
160  const G4int numMul = 1200;
161  const G4int numSec = 60;
162  static G4ThreadLocal G4double protmul[numMul], protnorm[numSec]; // proton constants
163  static G4ThreadLocal G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
164  // np = number of pi+, nneg = number of pi-, nz = number of pi0
165  G4int nt(0), np(0), nneg(0), nz(0);
166  const G4double c = 1.25;
167  const G4double b[] = { 0.70, 0.70 };
168  if( first ) // compute normalization constants, this will only be Done once
169  {
170  first = false;
171  G4int i;
172  for( i=0; i<numMul; ++i )protmul[i] = 0.0;
173  for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
174  G4int counter = -1;
175  for( np=0; np<(numSec/3); ++np )
176  {
177  for( nneg=std::max(0,np-1); nneg<=(np+1); ++nneg )
178  {
179  for( nz=0; nz<numSec/3; ++nz )
180  {
181  if( ++counter < numMul )
182  {
183  nt = np+nneg+nz;
184  if( (nt>0) && (nt<=numSec) )
185  {
186  protmul[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
187  protnorm[nt-1] += protmul[counter];
188  }
189  }
190  }
191  }
192  }
193  for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
194  for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
195  counter = -1;
196  for( np=0; np<numSec/3; ++np )
197  {
198  for( nneg=np; nneg<=(np+2); ++nneg )
199  {
200  for( nz=0; nz<numSec/3; ++nz )
201  {
202  if( ++counter < numMul )
203  {
204  nt = np+nneg+nz;
205  if( (nt>0) && (nt<=numSec) )
206  {
207  neutmul[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
208  neutnorm[nt-1] += neutmul[counter];
209  }
210  }
211  }
212  }
213  }
214  for( i=0; i<numSec; ++i )
215  {
216  if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
217  if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
218  }
219  } // end of initialization
220 
221  const G4double expxu = 82.; // upper bound for arg. of exp
222  const G4double expxl = -expxu; // lower bound for arg. of exp
235  const G4double cech[] = {1.,1.,1.,0.70,0.60,0.55,0.35,0.25,0.18,0.15};
236  G4int iplab = G4int(std::min( 9.0, pOriginal/GeV*5.0 ));
237  if( (pOriginal <= 2.0*GeV) && (G4UniformRand() < cech[iplab]) )
238  {
239  np = nneg = nz = nt = 0;
240  iplab = G4int(std::min( 19.0, pOriginal/GeV*10.0 ));
241  const G4double cnk0[] = {0.17,0.18,0.17,0.24,0.26,0.20,0.22,0.21,0.34,0.45,
242  0.58,0.55,0.36,0.29,0.29,0.32,0.32,0.33,0.33,0.33};
243  if( G4UniformRand() <= cnk0[iplab] )
244  {
245  quasiElastic = true;
246  if( targetParticle.GetDefinition() == aProton )
247  {
248  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
249  incidentHasChanged = true;
250  targetParticle.SetDefinitionAndUpdateE( aNeutron );
251  targetHasChanged = true;
252  }
253  }
254  else // random number > cnk0[iplab]
255  {
256  G4double ran = G4UniformRand();
257  if( ran < 0.25 ) // k- p --> pi- s+
258  {
259  if( targetParticle.GetDefinition() == aProton )
260  {
261  currentParticle.SetDefinitionAndUpdateE( aPiMinus );
262  targetParticle.SetDefinitionAndUpdateE( aSigmaPlus );
263  incidentHasChanged = true;
264  targetHasChanged = true;
265  }
266  }
267  else if( ran < 0.50 ) // k- p --> pi0 s0 or k- n --> pi- s0
268  {
269  if( targetParticle.GetDefinition() == aNeutron )
270  currentParticle.SetDefinitionAndUpdateE( aPiMinus );
271  else
272  currentParticle.SetDefinitionAndUpdateE( aPiZero );
273  targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
274  incidentHasChanged = true;
275  targetHasChanged = true;
276  }
277  else if( ran < 0.75 ) // k- p --> pi+ s- or k- n --> pi0 s-
278  {
279  if( targetParticle.GetDefinition() == aNeutron )
280  currentParticle.SetDefinitionAndUpdateE( aPiZero );
281  else
282  currentParticle.SetDefinitionAndUpdateE( aPiPlus );
283  targetParticle.SetDefinitionAndUpdateE( aSigmaMinus );
284  incidentHasChanged = true;
285  targetHasChanged = true;
286  }
287  else // k- p --> pi0 L or k- n --> pi- L
288  {
289  if( targetParticle.GetDefinition() == aNeutron )
290  currentParticle.SetDefinitionAndUpdateE( aPiMinus );
291  else
292  currentParticle.SetDefinitionAndUpdateE( aPiZero );
293  targetParticle.SetDefinitionAndUpdateE( aLambda );
294  incidentHasChanged = true;
295  targetHasChanged = true;
296  }
297  }
298  }
299  else // (pOriginal > 2.0*GeV) || (random number >= cech[iplab])
300  {
301  if( availableEnergy < aPiPlus->GetPDGMass() )
302  { // not energetically possible to produce pion(s)
303  quasiElastic = true;
304  return;
305  }
306  G4double n, anpn;
307  GetNormalizationConstant( availableEnergy, n, anpn );
308  G4double ran = G4UniformRand();
309  G4double dum, test, excs = 0.0;
310  if( targetParticle.GetDefinition() == aProton )
311  {
312  G4int counter = -1;
313  for( np=0; np<numSec/3 && ran>=excs; ++np )
314  {
315  for( nneg=std::max(0,np-1); nneg<=(np+1) && ran>=excs; ++nneg )
316  {
317  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
318  {
319  if( ++counter < numMul )
320  {
321  nt = np+nneg+nz;
322  if( nt > 0 )
323  {
324  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
325  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
326  if( std::fabs(dum) < 1.0 )
327  {
328  if( test >= 1.0e-10 )excs += dum*test;
329  }
330  else
331  excs += dum*test;
332  }
333  }
334  }
335  }
336  }
337  if( ran >= excs ) // 3 previous loops continued to the end
338  {
339  quasiElastic = true;
340  return;
341  }
342  np--; nneg--; nz--;
343  if( np == nneg )
344  {
345  if( G4UniformRand() >= 0.75 )
346  {
347  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
348  targetParticle.SetDefinitionAndUpdateE( aNeutron );
349  incidentHasChanged = true;
350  targetHasChanged = true;
351  }
352  }
353  else if( np == nneg+1 )
354  {
355  targetParticle.SetDefinitionAndUpdateE( aNeutron );
356  targetHasChanged = true;
357  }
358  else
359  {
360  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
361  incidentHasChanged = true;
362  }
363  }
364  else // target must be a neutron
365  {
366  G4int counter = -1;
367  for( np=0; np<numSec/3 && ran>=excs; ++np )
368  {
369  for( nneg=np; nneg<=(np+2) && ran>=excs; ++nneg )
370  {
371  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
372  {
373  if( ++counter < numMul )
374  {
375  nt = np+nneg+nz;
376  if( (nt>=1) && (nt<=numSec) )
377  {
378  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
379  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
380  if( std::fabs(dum) < 1.0 )
381  {
382  if( test >= 1.0e-10 )excs += dum*test;
383  }
384  else
385  excs += dum*test;
386  }
387  }
388  }
389  }
390  }
391  if( ran >= excs ) // 3 previous loops continued to the end
392  {
393  quasiElastic = true;
394  return;
395  }
396  np--; nneg--; nz--;
397  if( np == nneg-1 )
398  {
399  if( G4UniformRand() < 0.5 )
400  {
401  targetParticle.SetDefinitionAndUpdateE( aProton );
402  targetHasChanged = true;
403  }
404  else
405  {
406  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
407  incidentHasChanged = true;
408  }
409  }
410  else if( np != nneg )
411  {
412  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
413  incidentHasChanged = true;
414  }
415  }
416  if( G4UniformRand() >= 0.5 )
417  {
418  if( (currentParticle.GetDefinition() == aKaonMinus &&
419  targetParticle.GetDefinition() == aNeutron ) ||
420  (currentParticle.GetDefinition() == aKaonZL &&
421  targetParticle.GetDefinition() == aProton ) )
422  {
423  ran = G4UniformRand();
424  if( ran < 0.68 )
425  {
426  if( targetParticle.GetDefinition() == aProton )
427  {
428  currentParticle.SetDefinitionAndUpdateE( aPiPlus );
429  targetParticle.SetDefinitionAndUpdateE( aLambda );
430  incidentHasChanged = true;
431  targetHasChanged = true;
432  }
433  else
434  {
435  currentParticle.SetDefinitionAndUpdateE( aPiMinus );
436  targetParticle.SetDefinitionAndUpdateE( aLambda );
437  incidentHasChanged = true;
438  targetHasChanged = true;
439  }
440  }
441  else if( ran < 0.84 )
442  {
443  if( targetParticle.GetDefinition() == aProton )
444  {
445  currentParticle.SetDefinitionAndUpdateE( aPiZero );
446  targetParticle.SetDefinitionAndUpdateE( aSigmaPlus );
447  incidentHasChanged = true;
448  targetHasChanged = true;
449  }
450  else
451  {
452  currentParticle.SetDefinitionAndUpdateE( aPiMinus );
453  targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
454  incidentHasChanged = true;
455  targetHasChanged = true;
456  }
457  }
458  else
459  {
460  if( targetParticle.GetDefinition() == aProton )
461  {
462  currentParticle.SetDefinitionAndUpdateE( aPiPlus );
463  targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
464  incidentHasChanged = true;
465  targetHasChanged = true;
466  }
467  else
468  {
469  currentParticle.SetDefinitionAndUpdateE( aPiZero );
470  targetParticle.SetDefinitionAndUpdateE( aSigmaMinus );
471  incidentHasChanged = true;
472  targetHasChanged = true;
473  }
474  }
475  }
476  else // ( current != aKaonMinus || target != aNeutron ) &&
477  // ( current != aKaonZL || target != aProton )
478  {
479  ran = G4UniformRand();
480  if( ran < 0.67 )
481  {
482  currentParticle.SetDefinitionAndUpdateE( aPiZero );
483  targetParticle.SetDefinitionAndUpdateE( aLambda );
484  incidentHasChanged = true;
485  targetHasChanged = true;
486  }
487  else if( ran < 0.78 )
488  {
489  currentParticle.SetDefinitionAndUpdateE( aPiMinus );
490  targetParticle.SetDefinitionAndUpdateE( aSigmaPlus );
491  incidentHasChanged = true;
492  targetHasChanged = true;
493  }
494  else if( ran < 0.89 )
495  {
496  currentParticle.SetDefinitionAndUpdateE( aPiZero );
497  targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
498  incidentHasChanged = true;
499  targetHasChanged = true;
500  }
501  else
502  {
503  currentParticle.SetDefinitionAndUpdateE( aPiPlus );
504  targetParticle.SetDefinitionAndUpdateE( aSigmaMinus );
505  incidentHasChanged = true;
506  targetHasChanged = true;
507  }
508  }
509  }
510  }
511 
512  if (currentParticle.GetDefinition() == aKaonZL) {
513  if (G4UniformRand() >= 0.5) {
514  currentParticle.SetDefinitionAndUpdateE(aKaonZS);
515  incidentHasChanged = true;
516  }
517  }
518 
519  if (targetParticle.GetDefinition() == aKaonZL) {
520  if (G4UniformRand() >= 0.5) {
521  targetParticle.SetDefinitionAndUpdateE(aKaonZS);
522  targetHasChanged = true;
523  }
524  }
525 
526  SetUpPions(np, nneg, nz, vec, vecLen);
527  return;
528 }
529 
530  /* end of file */
531 
532 
533 
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:278
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
void SetUpChange(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
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
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:102
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
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ParticleDefinition * GetDefinition() const
static G4KaonZeroShort * KaonZeroShort()
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
static G4PionZero * PionZero()
Definition: G4PionZero.cc:108
static G4SigmaMinus * SigmaMinus()
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
void SetEnergyChange(G4double anEnergy)
G4double GetPDGMass() const
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 G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static constexpr double GeV
Definition: G4SIunits.hh:217
#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
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
void GetNormalizationConstant(const G4double availableEnergy, G4double &n, G4double &anpn)
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108
double G4double
Definition: G4Types.hh:76
void SetUpPions(const G4int np, const G4int nm, const G4int nz, G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen)
void SetMomentumChange(const G4ThreeVector &aV)
G4double GetMass() const
G4double GetTotalMomentum() const
G4double GetTotalEnergy() const