Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4RPGSigmaPlusInelastic.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 
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 
51  if( verboseLevel > 1 )
52  {
53  const G4Material *targetMaterial = aTrack.GetMaterial();
54  G4cout << "G4RPGSigmaPlusInelastic::ApplyYourself called" << G4endl;
55  G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
56  G4cout << "target material = " << targetMaterial->GetName() << ", ";
57  G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
58  << G4endl;
59  }
60 
61  // Fermi motion and evaporation
62  // As of Geant3, the Fermi energy calculation had not been Done
63 
64  G4double ek = originalIncident->GetKineticEnergy()/MeV;
65  G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
66  G4ReactionProduct modifiedOriginal;
67  modifiedOriginal = *originalIncident;
68 
69  G4double tkin = targetNucleus.Cinema( ek );
70  ek += tkin;
71  modifiedOriginal.SetKineticEnergy( ek*MeV );
72  G4double et = ek + amas;
73  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
74  G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
75  if( pp > 0.0 )
76  {
77  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
78  modifiedOriginal.SetMomentum( momentum * (p/pp) );
79  }
80  //
81  // calculate black track energies
82  //
83  tkin = targetNucleus.EvaporationEffects( ek );
84  ek -= tkin;
85  modifiedOriginal.SetKineticEnergy( ek*MeV );
86  et = ek + amas;
87  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
88  pp = modifiedOriginal.GetMomentum().mag()/MeV;
89  if( pp > 0.0 )
90  {
91  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
92  modifiedOriginal.SetMomentum( momentum * (p/pp) );
93  }
94  G4ReactionProduct currentParticle = modifiedOriginal;
95  G4ReactionProduct targetParticle;
96  targetParticle = *originalTarget;
97  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
98  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
99  G4bool incidentHasChanged = false;
100  G4bool targetHasChanged = false;
101  G4bool quasiElastic = false;
102  G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
103  G4int vecLen = 0;
104  vec.Initialize( 0 );
105 
106  const G4double cutOff = 0.1;
107  if( currentParticle.GetKineticEnergy()/MeV > cutOff )
108  Cascade( vec, vecLen,
109  originalIncident, currentParticle, targetParticle,
110  incidentHasChanged, targetHasChanged, quasiElastic );
111 
112  CalculateMomenta( vec, vecLen,
113  originalIncident, originalTarget, modifiedOriginal,
114  targetNucleus, currentParticle, targetParticle,
115  incidentHasChanged, targetHasChanged, quasiElastic );
116 
117  SetUpChange( vec, vecLen,
118  currentParticle, targetParticle,
119  incidentHasChanged );
120 
121  delete originalTarget;
122  return &theParticleChange;
123 }
124 
125 void G4RPGSigmaPlusInelastic::Cascade(
127  G4int& vecLen,
128  const G4HadProjectile *originalIncident,
129  G4ReactionProduct &currentParticle,
130  G4ReactionProduct &targetParticle,
131  G4bool &incidentHasChanged,
132  G4bool &targetHasChanged,
133  G4bool &quasiElastic )
134 {
135  // Derived from H. Fesefeldt's original FORTRAN code CASSP
136  //
137  // SigmaPlus undergoes interaction with nucleon within a nucleus. Check if it is
138  // energetically possible to produce pions/kaons. In not, assume nuclear excitation
139  // occurs and input particle is degraded in energy. No other particles are produced.
140  // If reaction is possible, find the correct number of pions/protons/neutrons
141  // produced using an interpolation to multiplicity data. Replace some pions or
142  // protons/neutrons by kaons or strange baryons according to the average
143  // multiplicity per inelastic reaction.
144 
145  const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
146  const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
147  const G4double targetMass = targetParticle.GetMass()/MeV;
148  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
149  targetMass*targetMass +
150  2.0*targetMass*etOriginal );
151  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
152  if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV )
153  {
154  quasiElastic = true;
155  return;
156  }
157  static G4bool first = true;
158  const G4int numMul = 1200;
159  const G4int numSec = 60;
160  static G4double protmul[numMul], protnorm[numSec]; // proton constants
161  static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
162  // np = number of pi+, nneg = number of pi-, nz = number of pi0
163  G4int counter, nt=0, np=0, nneg=0, nz=0;
164  G4double test;
165  const G4double c = 1.25;
166  const G4double b[] = { 0.7, 0.7 };
167  if( first ) // compute normalization constants, this will only be Done once
168  {
169  first = false;
170  G4int i;
171  for( i=0; i<numMul; ++i )protmul[i] = 0.0;
172  for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
173  counter = -1;
174  for( np=0; np<(numSec/3); ++np )
175  {
176  for( nneg=np; nneg<=(np+2); ++nneg )
177  {
178  for( nz=0; nz<numSec/3; ++nz )
179  {
180  if( ++counter < numMul )
181  {
182  nt = np+nneg+nz;
183  if( nt>0 && nt<=numSec )
184  {
185  protmul[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
186  protnorm[nt-1] += protmul[counter];
187  }
188  }
189  }
190  }
191  }
192  for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
193  for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
194  counter = -1;
195  for( np=0; np<numSec/3; ++np )
196  {
197  for( nneg=std::max(0,np-1); nneg<=(np+1); ++nneg )
198  {
199  for( nz=0; nz<numSec/3; ++nz )
200  {
201  if( ++counter < numMul )
202  {
203  nt = np+nneg+nz;
204  if( nt>0 && nt<=numSec )
205  {
206  neutmul[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
207  neutnorm[nt-1] += neutmul[counter];
208  }
209  }
210  }
211  }
212  }
213  for( i=0; i<numSec; ++i )
214  {
215  if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
216  if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
217  }
218  } // end of initialization
219 
220  const G4double expxu = 82.; // upper bound for arg. of exp
221  const G4double expxl = -expxu; // lower bound for arg. of exp
226  //
227  // energetically possible to produce pion(s) --> inelastic scattering
228  //
229  G4double n, anpn;
230  GetNormalizationConstant( availableEnergy, n, anpn );
231  G4double ran = G4UniformRand();
232  G4double dum, excs = 0.0;
233  if( targetParticle.GetDefinition() == aProton )
234  {
235  counter = -1;
236  for( np=0; np<numSec/3 && ran>=excs; ++np )
237  {
238  for( nneg=np; nneg<=(np+2) && ran>=excs; ++nneg )
239  {
240  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
241  {
242  if( ++counter < numMul )
243  {
244  nt = np+nneg+nz;
245  if( nt>0 && nt<=numSec )
246  {
247  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
248  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
249  if( std::fabs(dum) < 1.0 )
250  {
251  if( test >= 1.0e-10 )excs += dum*test;
252  }
253  else
254  excs += dum*test;
255  }
256  }
257  }
258  }
259  }
260  if( ran >= excs ) // 3 previous loops continued to the end
261  {
262  quasiElastic = true;
263  return;
264  }
265  np--; nneg--; nz--;
266  switch( std::min( 3, std::max( 1, np-nneg+3 ) ) )
267  {
268  case 1:
269  if( G4UniformRand() < 0.5 )
270  currentParticle.SetDefinitionAndUpdateE( aLambda );
271  else
272  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
273  incidentHasChanged = true;
274  targetParticle.SetDefinitionAndUpdateE( aNeutron );
275  targetHasChanged = true;
276  break;
277  case 2:
278  if( G4UniformRand() < 0.5 )
279  {
280  targetParticle.SetDefinitionAndUpdateE( aNeutron );
281  targetHasChanged = true;
282  }
283  else
284  {
285  if( G4UniformRand() < 0.5 )
286  currentParticle.SetDefinitionAndUpdateE( aLambda );
287  else
288  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
289  incidentHasChanged = true;
290  }
291  break;
292  case 3:
293  break;
294  }
295  }
296  else // target must be a neutron
297  {
298  counter = -1;
299  for( np=0; np<numSec/3 && ran>=excs; ++np )
300  {
301  for( nneg=std::max(0,np-1); nneg<=(np+1) && ran>=excs; ++nneg )
302  {
303  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
304  {
305  if( ++counter < numMul )
306  {
307  nt = np+nneg+nz;
308  if( nt>0 && nt<=numSec )
309  {
310  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
311  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
312  if( std::fabs(dum) < 1.0 )
313  {
314  if( test >= 1.0e-10 )excs += dum*test;
315  }
316  else
317  excs += dum*test;
318  }
319  }
320  }
321  }
322  }
323  if( ran >= excs ) // 3 previous loops continued to the end
324  {
325  quasiElastic = true;
326  return;
327  }
328  np--; nneg--; nz--;
329  switch( std::min( 3, std::max( 1, np-nneg+2 ) ) )
330  {
331  case 1:
332  targetParticle.SetDefinitionAndUpdateE( aProton );
333  targetHasChanged = true;
334  break;
335  case 2:
336  if( G4UniformRand() < 0.5 )
337  {
338  if( G4UniformRand() < 0.5 )
339  {
340  currentParticle.SetDefinitionAndUpdateE( aLambda );
341  incidentHasChanged = true;
342  targetParticle.SetDefinitionAndUpdateE( aProton );
343  targetHasChanged = true;
344  }
345  else
346  {
347  targetParticle.SetDefinitionAndUpdateE( aNeutron );
348  targetHasChanged = true;
349  }
350  }
351  else
352  {
353  if( G4UniformRand() < 0.5 )
354  {
355  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
356  incidentHasChanged = true;
357  targetParticle.SetDefinitionAndUpdateE( aProton );
358  targetHasChanged = true;
359  }
360  }
361  break;
362  case 3:
363  if( G4UniformRand() < 0.5 )
364  currentParticle.SetDefinitionAndUpdateE( aLambda );
365  else
366  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
367  incidentHasChanged = true;
368  break;
369  }
370  }
371 
372  SetUpPions(np, nneg, nz, vec, vecLen);
373  return;
374 }
375 
376  /* end of file */
377