Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LEKaonMinusInelastic.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 // Hadronic Process: Low Energy KaonMinus Inelastic Process
29 // J.L. Chuma, TRIUMF, 12-Feb-1997
30 // J.P.Wellisch 23-Apr-97: bug-hunting (missing initialization of npos,nneg,nzero
31 // fixed)
32 // Modified by J.L.Chuma 30-Apr-97: added originalTarget for CalculateMomenta
33 
34 #include <iostream>
35 
37 #include "G4PhysicalConstants.hh"
38 #include "G4SystemOfUnits.hh"
39 #include "Randomize.hh"
40 
43 {
44  SetMinEnergy(0.0);
45  SetMaxEnergy(25.*GeV);
46  G4cout << "WARNING: model G4LEKaonMinusInelastic is being deprecated and will\n"
47  << "disappear in Geant4 version 10.0" << G4endl;
48 }
49 
50 
52 {
53  outFile << "G4LEKaonMinusInelastic is one of the Low Energy Parameterized\n"
54  << "(LEP) models used to implement inelastic K- scattering\n"
55  << "from nuclei. It is a re-engineered version of the GHEISHA\n"
56  << "code of H. Fesefeldt. It divides the initial collision\n"
57  << "products into backward- and forward-going clusters which are\n"
58  << "then decayed into final state hadrons. The model does not\n"
59  << "conserve energy on an event-by-event basis. It may be\n"
60  << "applied to kaons with initial energies between 0 and 25\n"
61  << "GeV.\n";
62 }
63 
64 
67  G4Nucleus& targetNucleus)
68 {
69  const G4HadProjectile *originalIncident = &aTrack;
70  if (originalIncident->GetKineticEnergy()<= 0.1*MeV) {
74  return &theParticleChange;
75  }
76 
77  // create the target particle
78  G4DynamicParticle* originalTarget = targetNucleus.ReturnTargetParticle();
79  G4ReactionProduct targetParticle( originalTarget->GetDefinition() );
80 
81  if (verboseLevel > 1) {
82  const G4Material *targetMaterial = aTrack.GetMaterial();
83  G4cout << "G4LEKaonMinusInelastic::ApplyYourself called" << G4endl;
84  G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy() << "MeV, ";
85  G4cout << "target material = " << targetMaterial->GetName() << ", ";
86  G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
87  << G4endl;
88  }
89 
90  G4ReactionProduct currentParticle(const_cast<G4ParticleDefinition*>(originalIncident->GetDefinition()) );
91  currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
92  currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
93 
94  // Fermi motion and evaporation
95  // As of Geant3, the Fermi energy calculation had not been done
96  G4double ek = originalIncident->GetKineticEnergy();
97  G4double amas = originalIncident->GetDefinition()->GetPDGMass();
98 
99  G4double tkin = targetNucleus.Cinema( ek );
100  ek += tkin;
101  currentParticle.SetKineticEnergy( ek );
102  G4double et = ek + amas;
103  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
104  G4double pp = currentParticle.GetMomentum().mag();
105  if (pp > 0.0) {
106  G4ThreeVector momentum = currentParticle.GetMomentum();
107  currentParticle.SetMomentum( momentum * (p/pp) );
108  }
109 
110  // calculate black track energies
111  tkin = targetNucleus.EvaporationEffects( ek );
112  ek -= tkin;
113  currentParticle.SetKineticEnergy( ek );
114  et = ek + amas;
115  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
116  pp = currentParticle.GetMomentum().mag();
117  if (pp > 0.0) {
118  G4ThreeVector momentum = currentParticle.GetMomentum();
119  currentParticle.SetMomentum( momentum * (p/pp) );
120  }
121 
122  G4ReactionProduct modifiedOriginal = currentParticle;
123 
124  currentParticle.SetSide(1); // incident always goes in forward hemisphere
125  targetParticle.SetSide(-1); // target always goes in backward hemisphere
126  G4bool incidentHasChanged = false;
127  G4bool targetHasChanged = false;
128  G4bool quasiElastic = false;
129  G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
130  G4int vecLen = 0;
131  vec.Initialize(0);
132 
133  const G4double cutOff = 0.1*MeV;
134  if (currentParticle.GetKineticEnergy() > cutOff)
135  Cascade(vec, vecLen, originalIncident, currentParticle, targetParticle,
136  incidentHasChanged, targetHasChanged, quasiElastic);
137 
138  CalculateMomenta(vec, vecLen, originalIncident, originalTarget,
139  modifiedOriginal, targetNucleus, currentParticle,
140  targetParticle, incidentHasChanged, targetHasChanged,
141  quasiElastic);
142 
143  SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged);
144 
145  if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
146 
147  delete originalTarget;
148  return &theParticleChange;
149 }
150 
151 
152 void G4LEKaonMinusInelastic::Cascade(
154  G4int& vecLen,
155  const G4HadProjectile *originalIncident,
156  G4ReactionProduct &currentParticle,
157  G4ReactionProduct &targetParticle,
158  G4bool &incidentHasChanged,
159  G4bool &targetHasChanged,
160  G4bool &quasiElastic )
161 {
162  // derived from original FORTRAN code CASKM by H. Fesefeldt (13-Sep-1987)
163  //
164  // K- undergoes interaction with nucleon within a nucleus. Check if it is
165  // energetically possible to produce pions/kaons. In not, assume nuclear excitation
166  // occurs and input particle is degraded in energy. No other particles are produced.
167  // If reaction is possible, find the correct number of pions/protons/neutrons
168  // produced using an interpolation to multiplicity data. Replace some pions or
169  // protons/neutrons by kaons or strange baryons according to the average
170  // multiplicity per Inelastic reaction.
171  //
172  const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass();
173  const G4double etOriginal = originalIncident->GetTotalEnergy();
174  const G4double pOriginal = originalIncident->GetTotalMomentum();
175  const G4double targetMass = targetParticle.GetMass();
176  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
177  targetMass*targetMass +
178  2.0*targetMass*etOriginal );
179  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
180 
181  static G4bool first = true;
182  const G4int numMul = 1200;
183  const G4int numSec = 60;
184  static G4double protmul[numMul], protnorm[numSec]; // proton constants
185  static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
186  // npos = number of pi+, nneg = number of pi-, nzero = number of pi0
187  G4int nt(0), npos(0), nneg(0), nzero(0);
188  const G4double c = 1.25;
189  const G4double b[] = { 0.70, 0.70 };
190  if( first ) // compute normalization constants, this will only be Done once
191  {
192  first = false;
193  G4int i;
194  for( i=0; i<numMul; ++i )protmul[i] = 0.0;
195  for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
196  G4int counter = -1;
197  for( npos=0; npos<(numSec/3); ++npos )
198  {
199  for( nneg=std::max(0,npos-1); nneg<=(npos+1); ++nneg )
200  {
201  for( nzero=0; nzero<numSec/3; ++nzero )
202  {
203  if( ++counter < numMul )
204  {
205  nt = npos+nneg+nzero;
206  if( (nt>0) && (nt<=numSec) )
207  {
208  protmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[0],c);
209  protnorm[nt-1] += protmul[counter];
210  }
211  }
212  }
213  }
214  }
215  for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
216  for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
217  counter = -1;
218  for( npos=0; npos<numSec/3; ++npos )
219  {
220  for( nneg=npos; nneg<=(npos+2); ++nneg )
221  {
222  for( nzero=0; nzero<numSec/3; ++nzero )
223  {
224  if( ++counter < numMul )
225  {
226  nt = npos+nneg+nzero;
227  if( (nt>0) && (nt<=numSec) )
228  {
229  neutmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[1],c);
230  neutnorm[nt-1] += neutmul[counter];
231  }
232  }
233  }
234  }
235  }
236  for( i=0; i<numSec; ++i )
237  {
238  if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
239  if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
240  }
241  } // end of initialization
242 
243  const G4double expxu = 82.; // upper bound for arg. of exp
244  const G4double expxl = -expxu; // lower bound for arg. of exp
257  const G4double cech[] = {1.,1.,1.,0.70,0.60,0.55,0.35,0.25,0.18,0.15};
258  G4int iplab = G4int(std::min( 9.0, pOriginal/GeV*5.0 ));
259  if( (pOriginal <= 2.0*GeV) && (G4UniformRand() < cech[iplab]) )
260  {
261  npos = nneg = nzero = nt = 0;
262  iplab = G4int(std::min( 19.0, pOriginal/GeV*10.0 ));
263  const G4double cnk0[] = {0.17,0.18,0.17,0.24,0.26,0.20,0.22,0.21,0.34,0.45,
264  0.58,0.55,0.36,0.29,0.29,0.32,0.32,0.33,0.33,0.33};
265  if( G4UniformRand() <= cnk0[iplab] )
266  {
267  quasiElastic = true;
268  if( targetParticle.GetDefinition() == aProton )
269  {
270  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
271  incidentHasChanged = true;
272  targetParticle.SetDefinitionAndUpdateE( aNeutron );
273  targetHasChanged = true;
274  }
275  }
276  else // random number > cnk0[iplab]
277  {
278  G4double ran = G4UniformRand();
279  if( ran < 0.25 ) // k- p --> pi- s+
280  {
281  if( targetParticle.GetDefinition() == aProton )
282  {
283  currentParticle.SetDefinitionAndUpdateE( aPiMinus );
284  targetParticle.SetDefinitionAndUpdateE( aSigmaPlus );
285  incidentHasChanged = true;
286  targetHasChanged = true;
287  }
288  }
289  else if( ran < 0.50 ) // k- p --> pi0 s0 or k- n --> pi- s0
290  {
291  if( targetParticle.GetDefinition() == aNeutron )
292  currentParticle.SetDefinitionAndUpdateE( aPiMinus );
293  else
294  currentParticle.SetDefinitionAndUpdateE( aPiZero );
295  targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
296  incidentHasChanged = true;
297  targetHasChanged = true;
298  }
299  else if( ran < 0.75 ) // k- p --> pi+ s- or k- n --> pi0 s-
300  {
301  if( targetParticle.GetDefinition() == aNeutron )
302  currentParticle.SetDefinitionAndUpdateE( aPiZero );
303  else
304  currentParticle.SetDefinitionAndUpdateE( aPiPlus );
305  targetParticle.SetDefinitionAndUpdateE( aSigmaMinus );
306  incidentHasChanged = true;
307  targetHasChanged = true;
308  }
309  else // k- p --> pi0 L or k- n --> pi- L
310  {
311  if( targetParticle.GetDefinition() == aNeutron )
312  currentParticle.SetDefinitionAndUpdateE( aPiMinus );
313  else
314  currentParticle.SetDefinitionAndUpdateE( aPiZero );
315  targetParticle.SetDefinitionAndUpdateE( aLambda );
316  incidentHasChanged = true;
317  targetHasChanged = true;
318  }
319  }
320  }
321  else // (pOriginal > 2.0*GeV) || (random number >= cech[iplab])
322  {
323  if( availableEnergy < aPiPlus->GetPDGMass() )
324  { // not energetically possible to produce pion(s)
325  quasiElastic = true;
326  return;
327  }
328  G4double n, anpn;
329  GetNormalizationConstant( availableEnergy, n, anpn );
330  G4double ran = G4UniformRand();
331  G4double dum, test, excs = 0.0;
332  if( targetParticle.GetDefinition() == aProton )
333  {
334  G4int counter = -1;
335  for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
336  {
337  for( nneg=std::max(0,npos-1); nneg<=(npos+1) && ran>=excs; ++nneg )
338  {
339  for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
340  {
341  if( ++counter < numMul )
342  {
343  nt = npos+nneg+nzero;
344  if( nt > 0 )
345  {
346  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
347  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
348  if( std::fabs(dum) < 1.0 )
349  {
350  if( test >= 1.0e-10 )excs += dum*test;
351  }
352  else
353  excs += dum*test;
354  }
355  }
356  }
357  }
358  }
359  if( ran >= excs ) // 3 previous loops continued to the end
360  {
361  quasiElastic = true;
362  return;
363  }
364  npos--; nneg--; nzero--;
365  if( npos == nneg )
366  {
367  if( G4UniformRand() >= 0.75 )
368  {
369  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
370  targetParticle.SetDefinitionAndUpdateE( aNeutron );
371  incidentHasChanged = true;
372  targetHasChanged = true;
373  }
374  }
375  else if( npos == nneg+1 )
376  {
377  targetParticle.SetDefinitionAndUpdateE( aNeutron );
378  targetHasChanged = true;
379  }
380  else
381  {
382  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
383  incidentHasChanged = true;
384  }
385  }
386  else // target must be a neutron
387  {
388  G4int counter = -1;
389  for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
390  {
391  for( nneg=npos; nneg<=(npos+2) && ran>=excs; ++nneg )
392  {
393  for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
394  {
395  if( ++counter < numMul )
396  {
397  nt = npos+nneg+nzero;
398  if( (nt>=1) && (nt<=numSec) )
399  {
400  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
401  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
402  if( std::fabs(dum) < 1.0 )
403  {
404  if( test >= 1.0e-10 )excs += dum*test;
405  }
406  else
407  excs += dum*test;
408  }
409  }
410  }
411  }
412  }
413  if( ran >= excs ) // 3 previous loops continued to the end
414  {
415  quasiElastic = true;
416  return;
417  }
418  npos--; nneg--; nzero--;
419  if( npos == nneg-1 )
420  {
421  if( G4UniformRand() < 0.5 )
422  {
423  targetParticle.SetDefinitionAndUpdateE( aProton );
424  targetHasChanged = true;
425  }
426  else
427  {
428  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
429  incidentHasChanged = true;
430  }
431  }
432  else if( npos != nneg )
433  {
434  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
435  incidentHasChanged = true;
436  }
437  }
438  if( G4UniformRand() >= 0.5 )
439  {
440  if( (currentParticle.GetDefinition() == aKaonMinus &&
441  targetParticle.GetDefinition() == aNeutron ) ||
442  (currentParticle.GetDefinition() == aKaonZL &&
443  targetParticle.GetDefinition() == aProton ) )
444  {
445  ran = G4UniformRand();
446  if( ran < 0.68 )
447  {
448  if( targetParticle.GetDefinition() == aProton )
449  {
450  currentParticle.SetDefinitionAndUpdateE( aPiPlus );
451  targetParticle.SetDefinitionAndUpdateE( aLambda );
452  incidentHasChanged = true;
453  targetHasChanged = true;
454  }
455  else
456  {
457  currentParticle.SetDefinitionAndUpdateE( aPiMinus );
458  targetParticle.SetDefinitionAndUpdateE( aLambda );
459  incidentHasChanged = true;
460  targetHasChanged = true;
461  }
462  }
463  else if( ran < 0.84 )
464  {
465  if( targetParticle.GetDefinition() == aProton )
466  {
467  currentParticle.SetDefinitionAndUpdateE( aPiZero );
468  targetParticle.SetDefinitionAndUpdateE( aSigmaPlus );
469  incidentHasChanged = true;
470  targetHasChanged = true;
471  }
472  else
473  {
474  currentParticle.SetDefinitionAndUpdateE( aPiMinus );
475  targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
476  incidentHasChanged = true;
477  targetHasChanged = true;
478  }
479  }
480  else
481  {
482  if( targetParticle.GetDefinition() == aProton )
483  {
484  currentParticle.SetDefinitionAndUpdateE( aPiPlus );
485  targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
486  incidentHasChanged = true;
487  targetHasChanged = true;
488  }
489  else
490  {
491  currentParticle.SetDefinitionAndUpdateE( aPiZero );
492  targetParticle.SetDefinitionAndUpdateE( aSigmaMinus );
493  incidentHasChanged = true;
494  targetHasChanged = true;
495  }
496  }
497  }
498  else // ( current != aKaonMinus || target != aNeutron ) &&
499  // ( current != aKaonZL || target != aProton )
500  {
501  ran = G4UniformRand();
502  if( ran < 0.67 )
503  {
504  currentParticle.SetDefinitionAndUpdateE( aPiZero );
505  targetParticle.SetDefinitionAndUpdateE( aLambda );
506  incidentHasChanged = true;
507  targetHasChanged = true;
508  }
509  else if( ran < 0.78 )
510  {
511  currentParticle.SetDefinitionAndUpdateE( aPiMinus );
512  targetParticle.SetDefinitionAndUpdateE( aSigmaPlus );
513  incidentHasChanged = true;
514  targetHasChanged = true;
515  }
516  else if( ran < 0.89 )
517  {
518  currentParticle.SetDefinitionAndUpdateE( aPiZero );
519  targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
520  incidentHasChanged = true;
521  targetHasChanged = true;
522  }
523  else
524  {
525  currentParticle.SetDefinitionAndUpdateE( aPiPlus );
526  targetParticle.SetDefinitionAndUpdateE( aSigmaMinus );
527  incidentHasChanged = true;
528  targetHasChanged = true;
529  }
530  }
531  }
532  }
533  if( currentParticle.GetDefinition() == aKaonZL )
534  {
535  if( G4UniformRand() >= 0.5 )
536  {
537  currentParticle.SetDefinitionAndUpdateE( aKaonZS );
538  incidentHasChanged = true;
539  }
540  }
541  if( targetParticle.GetDefinition() == aKaonZL )
542  {
543  if( G4UniformRand() >= 0.5 )
544  {
545  targetParticle.SetDefinitionAndUpdateE( aKaonZS );
546  targetHasChanged = true;
547  }
548  }
549 
550  SetUpPions( npos, nneg, nzero, vec, vecLen );
551  return;
552 }