Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LEKaonZeroInelastic.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 //
27 // $Id$
28 //
29  // Hadronic Process: Low Energy KaonZeroShort Inelastic Process
30  // J.L. Chuma, TRIUMF, 11-Feb-1997
31  // Last modified: 27-Mar-1997
32  // Modified by J.L.Chuma 30-Apr-97: added originalTarget for CalculateMomenta
33 
34 #include "G4LEKaonZeroInelastic.hh"
35 #include "G4PhysicalConstants.hh"
36 #include "G4SystemOfUnits.hh"
37 #include "Randomize.hh"
38 
40 {
41  outFile << "G4LEKaonZeroInelastic is one of the Low Energy Parameterized\n"
42  << "(LEP) models used to implement K0 scattering from nuclei. It\n"
43  << "is a re-engineered version of the GHEISHA code of\n"
44  << "H. Fesefeldt. It divides the initial collision products\n"
45  << "into backward- and forward-going clusters which are then\n"
46  << "decayed into final state hadrons. The model does not conserve\n"
47  << "energy on an event-by-event basis. It may be applied to\n"
48  << "K0s with initial energies between 0 and 25 GeV.\n";
49 }
50 
51 
54  G4Nucleus& targetNucleus)
55 {
56  const G4HadProjectile *originalIncident = &aTrack;
57 
58  // create the target particle
59  G4DynamicParticle* originalTarget = targetNucleus.ReturnTargetParticle();
60 
61  if (verboseLevel > 1) {
62  const G4Material *targetMaterial = aTrack.GetMaterial();
63  G4cout << "G4LEKaonZeroInelastic::ApplyYourself called" << G4endl;
64  G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
65  G4cout << "target material = " << targetMaterial->GetName() << ", ";
66  G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
67  << G4endl;
68  }
69 
70  // Fermi motion and evaporation
71  // As of Geant3, the Fermi energy calculation had not been Done
72  G4double ek = originalIncident->GetKineticEnergy()/MeV;
73  G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
74  G4ReactionProduct modifiedOriginal;
75  modifiedOriginal = *originalIncident;
76 
77  G4double tkin = targetNucleus.Cinema(ek);
78  ek += tkin;
79  modifiedOriginal.SetKineticEnergy(ek*MeV);
80  G4double et = ek + amas;
81  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
82  G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
83  if (pp > 0.0) {
84  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
85  modifiedOriginal.SetMomentum( momentum * (p/pp) );
86  }
87 
88  // calculate black track energies
89  tkin = targetNucleus.EvaporationEffects(ek);
90  ek -= tkin;
91  modifiedOriginal.SetKineticEnergy(ek*MeV);
92  et = ek + amas;
93  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
94  pp = modifiedOriginal.GetMomentum().mag()/MeV;
95  if (pp > 0.0) {
96  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
97  modifiedOriginal.SetMomentum( momentum * (p/pp) );
98  }
99  G4ReactionProduct currentParticle = modifiedOriginal;
100  G4ReactionProduct targetParticle;
101  targetParticle = *originalTarget;
102  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
103  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
104  G4bool incidentHasChanged = false;
105  G4bool targetHasChanged = false;
106  G4bool quasiElastic = false;
107  G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
108  G4int vecLen = 0;
109  vec.Initialize(0);
110 
111  const G4double cutOff = 0.1;
112  if (currentParticle.GetKineticEnergy()/MeV > cutOff)
113  Cascade(vec, vecLen, originalIncident, currentParticle, targetParticle,
114  incidentHasChanged, targetHasChanged, quasiElastic);
115 
116  CalculateMomenta(vec, vecLen, originalIncident, originalTarget, modifiedOriginal,
117  targetNucleus, currentParticle, targetParticle,
118  incidentHasChanged, targetHasChanged, quasiElastic);
119 
120  SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged);
121 
122  if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
123 
124  delete originalTarget;
125  return &theParticleChange;
126 }
127 
128 
129 void G4LEKaonZeroInelastic::Cascade(
131  G4int& vecLen,
132  const G4HadProjectile* originalIncident,
133  G4ReactionProduct& currentParticle,
134  G4ReactionProduct& targetParticle,
135  G4bool& incidentHasChanged,
136  G4bool& targetHasChanged,
137  G4bool& quasiElastic)
138 {
139  // derived from original FORTRAN code CASK0 by H. Fesefeldt (13-Sep-1987)
140  //
141  // K0Short undergoes interaction with nucleon within a nucleus. Check if it is
142  // energetically possible to produce pions/kaons. In not, assume nuclear excitation
143  // occurs and input particle is degraded in energy. No other particles are produced.
144  // If reaction is possible, find the correct number of pions/protons/neutrons
145  // produced using an interpolation to multiplicity data. Replace some pions or
146  // protons/neutrons by kaons or strange baryons according to the average
147  // multiplicity per Inelastic reaction.
148 
149  const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
150  const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
151  const G4double targetMass = targetParticle.GetMass()/MeV;
152  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
153  targetMass*targetMass +
154  2.0*targetMass*etOriginal );
155  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
156  if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV )
157  {
158  quasiElastic = true;
159  return;
160  }
161  static G4bool first = true;
162  const G4int numMul = 1200;
163  const G4int numSec = 60;
164  static G4double protmul[numMul], protnorm[numSec]; // proton constants
165  static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
166  // npos = number of pi+, nneg = number of pi-, nzero = number of pi0
167  G4int counter, nt=0, npos=0, nneg=0, nzero=0;
168  const G4double c = 1.25;
169  const G4double b[] = { 0.7, 0.7 };
170  if( first ) // compute normalization constants, this will only be Done once
171  {
172  first = false;
173  G4int i;
174  for( i=0; i<numMul; ++i )protmul[i] = 0.0;
175  for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
176  counter = -1;
177  for( npos=0; npos<(numSec/3); ++npos )
178  {
179  for( nneg=std::max(0,npos-1); nneg<=(npos+1); ++nneg )
180  {
181  for( nzero=0; nzero<numSec/3; ++nzero )
182  {
183  if( ++counter < numMul )
184  {
185  nt = npos+nneg+nzero;
186  if( nt>0 && nt<=numSec )
187  {
188  protmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[0],c);
189  protnorm[nt-1] += protmul[counter];
190  }
191  }
192  }
193  }
194  }
195  for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
196  for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
197  counter = -1;
198  for( npos=0; npos<numSec/3; ++npos )
199  {
200  for( nneg=std::max(0,npos-2); nneg<=npos; ++nneg )
201  {
202  for( nzero=0; nzero<numSec/3; ++nzero )
203  {
204  if( ++counter < numMul )
205  {
206  nt = npos+nneg+nzero;
207  if( nt>0 && nt<=numSec )
208  {
209  neutmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[1],c);
210  neutnorm[nt-1] += neutmul[counter];
211  }
212  }
213  }
214  }
215  }
216  for( i=0; i<numSec; ++i )
217  {
218  if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
219  if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
220  }
221  } // end of initialization
222 
223  const G4double expxu = 82.; // upper bound for arg. of exp
224  const G4double expxl = -expxu; // lower bound for arg. of exp
230  G4int ieab = static_cast<G4int>(5.0*availableEnergy*MeV/GeV);
231  const G4double supp[] = {0.,0.4,0.55,0.65,0.75,0.82,0.86,0.90,0.94,0.98};
232  G4double test, w0, wp, wt, wm;
233  if( (availableEnergy*MeV/GeV < 2.0) && (G4UniformRand() >= supp[ieab]) )
234  {
235  //
236  // suppress high multiplicity events at low momentum
237  // only one pion will be produced
238  //
239  nneg = npos = nzero = 0;
240  if( targetParticle.GetDefinition() == aNeutron )
241  {
242  test = std::exp( std::min( expxu, std::max( expxl, -(1.0+b[0])*(1.0+b[0])/(2.0*c*c) ) ) );
243  w0 = test/2.0;
244  test = std::exp( std::min( expxu, std::max( expxl, -(-1.0+b[0])*(1.0+b[0])/(2.0*c*c) ) ) );
245  wm = test*1.5;
246  if( G4UniformRand() < w0/(w0+wm) )
247  nzero = 1;
248  else
249  nneg = 1;
250  }
251  else // target is a proton
252  {
253  test = std::exp( std::min( expxu, std::max( expxl, -(1.0+b[1])*(1.0+b[1])/(2.0*c*c) ) ) );
254  w0 = test;
255  wp = test;
256  test = std::exp( std::min( expxu, std::max( expxl, -(-1.0+b[1])*(-1.0+b[1])/(2.0*c*c) ) ) );
257  wm = test;
258  wt = w0+wp+wm;
259  wp += w0;
260  G4double ran = G4UniformRand();
261  if( ran < w0/wt )
262  nzero = 1;
263  else if( ran < wp/wt )
264  npos = 1;
265  else
266  nneg = 1;
267  }
268  }
269  else // (availableEnergy*MeV/GeV >= 2.0) || (G4UniformRand() < supp[ieab])
270  {
271  G4double n, anpn;
272  GetNormalizationConstant( availableEnergy, n, anpn );
273  G4double ran = G4UniformRand();
274  G4double dum, excs = 0.0;
275  if( targetParticle.GetDefinition() == aProton )
276  {
277  counter = -1;
278  for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
279  {
280  for( nneg=std::max(0,npos-1); nneg<=(npos+1) && ran>=excs; ++nneg )
281  {
282  for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
283  {
284  if( ++counter < numMul )
285  {
286  nt = npos+nneg+nzero;
287  if( nt>0 && nt<=numSec )
288  {
289  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
290  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
291  if( std::fabs(dum) < 1.0 )
292  {
293  if( test >= 1.0e-10 )excs += dum*test;
294  }
295  else
296  excs += dum*test;
297  }
298  }
299  }
300  }
301  }
302  if( ran >= excs ) // 3 previous loops continued to the end
303  {
304  quasiElastic = true;
305  return;
306  }
307  npos--; nneg--; nzero--;
308  }
309  else // target must be a neutron
310  {
311  counter = -1;
312  for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
313  {
314  for( nneg=std::max(0,npos-2); nneg<=npos && ran>=excs; ++nneg )
315  {
316  for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
317  {
318  if( ++counter < numMul )
319  {
320  nt = npos+nneg+nzero;
321  if( nt>0 && nt<=numSec )
322  {
323  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
324  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
325  if( std::fabs(dum) < 1.0 )
326  {
327  if( test >= 1.0e-10 )excs += dum*test;
328  }
329  else
330  excs += dum*test;
331  }
332  }
333  }
334  }
335  }
336  if( ran >= excs ) // 3 previous loops continued to the end
337  {
338  quasiElastic = true;
339  return;
340  }
341  npos--; nneg--; nzero--;
342  }
343  }
344  if( targetParticle.GetDefinition() == aProton )
345  {
346  switch( npos-nneg )
347  {
348  case 0:
349  if( G4UniformRand() < 0.25 )
350  {
351  currentParticle.SetDefinitionAndUpdateE( aKaonPlus );
352  targetParticle.SetDefinitionAndUpdateE( aNeutron );
353  incidentHasChanged = true;
354  targetHasChanged = true;
355  }
356  break;
357  case 1:
358  targetParticle.SetDefinitionAndUpdateE( aNeutron );
359  targetHasChanged = true;
360  break;
361  default:
362  targetParticle.SetDefinitionAndUpdateE( aNeutron );
363  targetHasChanged = true;
364  break;
365  }
366  }
367  else // targetParticle is a neutron
368  {
369  switch( npos-nneg ) // seems wrong, charge not conserved
370  {
371  case 1:
372  if( G4UniformRand() < 0.5 )
373  {
374  currentParticle.SetDefinitionAndUpdateE( aKaonPlus );
375  incidentHasChanged = true;
376  }
377  else
378  {
379  targetParticle.SetDefinitionAndUpdateE( aProton );
380  targetHasChanged = true;
381  }
382  break;
383  case 2:
384  currentParticle.SetDefinitionAndUpdateE( aKaonPlus );
385  incidentHasChanged = true;
386  targetParticle.SetDefinitionAndUpdateE( aProton );
387  targetHasChanged = true;
388  break;
389  default:
390  break;
391  }
392  }
393  if( currentParticle.GetDefinition() == aKaonZS )
394  {
395  if( G4UniformRand() >= 0.5 )
396  {
397  currentParticle.SetDefinitionAndUpdateE( aKaonZL);
398  incidentHasChanged = true;
399  }
400  }
401  if( targetParticle.GetDefinition() == aKaonZS )
402  {
403  if( G4UniformRand() >= 0.5 )
404  {
405  targetParticle.SetDefinitionAndUpdateE( aKaonZL );
406  targetHasChanged = true;
407  }
408  }
409  SetUpPions( npos, nneg, nzero, vec, vecLen );
410  return;
411 }
412 
413  /* end of file */
414