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