Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4RPGKPlusInelastic.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 
30 #include "G4RPGKPlusInelastic.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 << "G4RPGKPlusInelastic::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  G4ReactionProduct currentParticle( const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition()));
63  currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
64  currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
65 
66  // Fermi motion and evaporation
67  // As of Geant3, the Fermi energy calculation had not been Done
68 
69  G4double ek = originalIncident->GetKineticEnergy();
70  G4double amas = originalIncident->GetDefinition()->GetPDGMass();
71 
72  G4double tkin = targetNucleus.Cinema( ek );
73  ek += tkin;
74  currentParticle.SetKineticEnergy( ek );
75  G4double et = ek + amas;
76  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
77  G4double pp = currentParticle.GetMomentum().mag();
78  if( pp > 0.0 )
79  {
80  G4ThreeVector momentum = currentParticle.GetMomentum();
81  currentParticle.SetMomentum( momentum * (p/pp) );
82  }
83 
84  // calculate black track energies
85 
86  tkin = targetNucleus.EvaporationEffects( ek );
87  ek -= tkin;
88  currentParticle.SetKineticEnergy( ek );
89  et = ek + amas;
90  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
91  pp = currentParticle.GetMomentum().mag();
92  if( pp > 0.0 )
93  {
94  G4ThreeVector momentum = currentParticle.GetMomentum();
95  currentParticle.SetMomentum( momentum * (p/pp) );
96  }
97 
98  G4ReactionProduct modifiedOriginal = currentParticle;
99 
100  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
101  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
102  G4bool incidentHasChanged = false;
103  G4bool targetHasChanged = false;
104  G4bool quasiElastic = false;
105  G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
106  G4int vecLen = 0;
107  vec.Initialize( 0 );
108 
109  const G4double cutOff = 0.1*MeV;
110  if( currentParticle.GetKineticEnergy() > cutOff )
111  Cascade( vec, vecLen,
112  originalIncident, currentParticle, targetParticle,
113  incidentHasChanged, targetHasChanged, quasiElastic );
114 
115  CalculateMomenta( vec, vecLen,
116  originalIncident, originalTarget, modifiedOriginal,
117  targetNucleus, currentParticle, targetParticle,
118  incidentHasChanged, targetHasChanged, quasiElastic );
119 
120  SetUpChange( vec, vecLen,
121  currentParticle, targetParticle,
122  incidentHasChanged );
123 
124  delete originalTarget;
125 
126  return &theParticleChange;
127 }
128 
129 
130 void G4RPGKPlusInelastic::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 CASKP
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 targetMass = targetParticle.GetMass();
153  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
154  targetMass*targetMass +
155  2.0*targetMass*etOriginal );
156  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
157  if( availableEnergy < G4PionPlus::PionPlus()->GetPDGMass() )
158  {
159  quasiElastic = true;
160  return;
161  }
162  static G4bool first = true;
163  const G4int numMul = 1200;
164  const G4int numSec = 60;
165  static G4double protmul[numMul], protnorm[numSec]; // proton constants
166  static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
167 
168  // np = number of pi+, nneg = number of pi-, nz = number of pi0
169 
170  G4int nt=0, np=0, nneg=0, nz=0;
171  const G4double c = 1.25;
172  const G4double b[] = { 0.70, 0.70 };
173  if( first ) // compute normalization constants, this will only be Done once
174  {
175  first = false;
176  G4int i;
177  for( i=0; i<numMul; ++i )protmul[i] = 0.0;
178  for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
179  G4int counter = -1;
180  for( np=0; np<(numSec/3); ++np )
181  {
182  for( nneg=std::max(0,np-2); nneg<=np; ++nneg )
183  {
184  for( nz=0; nz<numSec/3; ++nz )
185  {
186  if( ++counter < numMul )
187  {
188  nt = np+nneg+nz;
189  if( nt > 0 )
190  {
191  protmul[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
192  protnorm[nt-1] += protmul[counter];
193  }
194  }
195  }
196  }
197  }
198  for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
199  for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
200  counter = -1;
201  for( np=0; np<numSec/3; ++np )
202  {
203  for( nneg=std::max(0,np-1); nneg<=(np+1); ++nneg )
204  {
205  for( nz=0; nz<numSec/3; ++nz )
206  {
207  if( ++counter < numMul )
208  {
209  nt = np+nneg+nz;
210  if( (nt>0) && (nt<=numSec) )
211  {
212  neutmul[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
213  neutnorm[nt-1] += neutmul[counter];
214  }
215  }
216  }
217  }
218  }
219  for( i=0; i<numSec; ++i )
220  {
221  if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
222  if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
223  }
224  } // end of initialization
225 
226  const G4double expxu = 82.; // upper bound for arg. of exp
227  const G4double expxl = -expxu; // lower bound for arg. of exp
232  G4int ieab = static_cast<G4int>(availableEnergy*5.0/GeV);
233  const G4double supp[] = {0.,0.4,0.55,0.65,0.75,0.82,0.86,0.90,0.94,0.98};
234  G4double test, w0, wp, wt, wm;
235  if( (availableEnergy < 2.0*GeV) && (G4UniformRand() >= supp[ieab]) )
236  {
237  // suppress high multiplicity events at low momentum
238  // only one pion will be produced
239 
240  nneg = np = nz = 0;
241  if( targetParticle.GetDefinition() == aProton )
242  {
243  test = std::exp( std::min( expxu, std::max( expxl, -sqr(1.0+b[0])/(2.0*c*c) ) ) );
244  w0 = test;
245  wp = test*2.0;
246  if( G4UniformRand() < w0/(w0+wp) )
247  nz = 1;
248  else
249  np = 1;
250  }
251  else // target is a neutron
252  {
253  test = std::exp( std::min( expxu, std::max( expxl, -sqr(1.0+b[1])/(2.0*c*c) ) ) );
254  w0 = test;
255  wp = test;
256  test = std::exp( std::min( expxu, std::max( expxl, -sqr(-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  nz = 1;
263  else if( ran < wp/wt )
264  np = 1;
265  else
266  nneg = 1;
267  }
268  }
269  else
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  G4int counter = -1;
278  for( np=0; (np<numSec/3) && (ran>=excs); ++np )
279  {
280  for( nneg=std::max(0,np-2); (nneg<=np) && (ran>=excs); ++nneg )
281  {
282  for( nz=0; (nz<numSec/3) && (ran>=excs); ++nz )
283  {
284  if( ++counter < numMul )
285  {
286  nt = np+nneg+nz;
287  if( nt > 0 )
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 )return; // 3 previous loops continued to the end
303  np--; nneg--; nz--;
304  }
305  else // target must be a neutron
306  {
307  G4int counter = -1;
308  for( np=0; (np<numSec/3) && (ran>=excs); ++np )
309  {
310  for( nneg=std::max(0,np-1); (nneg<=(np+1)) && (ran>=excs); ++nneg )
311  {
312  for( nz=0; (nz<numSec/3) && (ran>=excs); ++nz )
313  {
314  if( ++counter < numMul )
315  {
316  nt = np+nneg+nz;
317  if( (nt>=1) && (nt<=numSec) )
318  {
319  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
320  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
321  if( std::fabs(dum) < 1.0 )
322  {
323  if( test >= 1.0e-10 )excs += dum*test;
324  }
325  else
326  excs += dum*test;
327  }
328  }
329  }
330  }
331  }
332  if( ran >= excs )return; // 3 previous loops continued to the end
333  np--; nneg--; nz--;
334  }
335  }
336 
337  if( targetParticle.GetDefinition() == aProton )
338  {
339  switch( np-nneg )
340  {
341  case 1:
342  if( G4UniformRand() < 0.5 )
343  {
344  if( G4UniformRand() < 0.5 )
345  currentParticle.SetDefinitionAndUpdateE( aKaonZS );
346  else
347  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
348  incidentHasChanged = true;
349  }
350  else
351  {
352  targetParticle.SetDefinitionAndUpdateE( aNeutron );
353  targetHasChanged = true;
354  }
355  break;
356  case 2:
357  if( G4UniformRand() < 0.5 )
358  currentParticle.SetDefinitionAndUpdateE( aKaonZS );
359  else
360  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
361  incidentHasChanged = true;
362  targetParticle.SetDefinitionAndUpdateE( aNeutron );
363  incidentHasChanged = true;
364  targetHasChanged = true;
365  break;
366  default:
367  break;
368  }
369  }
370  else // target is a neutron
371  {
372  switch( np-nneg )
373  {
374  case 0:
375  if( G4UniformRand() < 0.25 )
376  {
377  if( G4UniformRand() < 0.5 )
378  currentParticle.SetDefinitionAndUpdateE( aKaonZS );
379  else
380  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
381  targetParticle.SetDefinitionAndUpdateE( aProton );
382  incidentHasChanged = true;
383  targetHasChanged = true;
384  }
385  break;
386  case 1:
387  if( G4UniformRand() < 0.5 )
388  currentParticle.SetDefinitionAndUpdateE( aKaonZS );
389  else
390  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
391  incidentHasChanged = true;
392  break;
393  default: // assumes nneg = np+1 so charge is conserved
394  targetParticle.SetDefinitionAndUpdateE( aProton );
395  targetHasChanged = true;
396  break;
397  }
398  }
399 
400  SetUpPions(np, nneg, nz, vec, vecLen);
401  return;
402 }
403 
404  /* end of file */
405