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