Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LESigmaPlusInelastic.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: SigmaPlus 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 << "G4LESigmaPlusInelastic 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 << "G4LESigmaPlusInelastic::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 (currentParticle.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 G4LESigmaPlusInelastic::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 CASSP by H. Fesefeldt (30-Nov-1987)
145  //
146  // SigmaPlus 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  quasiElastic = true;
163  return;
164  }
165  static G4bool first = true;
166  const G4int numMul = 1200;
167  const G4int numSec = 60;
168  static G4double protmul[numMul], protnorm[numSec]; // proton constants
169  static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
170 
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.7, 0.7 };
176  if (first) { // Computation of normalization constants will only be done once
177  first = false;
178  G4int i;
179  for( i=0; i<numMul; ++i )protmul[i] = 0.0;
180  for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
181  counter = -1;
182  for (npos = 0; npos < (numSec/3); ++npos) {
183  for (nneg = npos; nneg <= (npos+2); ++nneg) {
184  for (nzero = 0; nzero < numSec/3; ++nzero) {
185  if (++counter < numMul) {
186  nt = npos+nneg+nzero;
187  if (nt > 0 && nt <= numSec) {
188  protmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[0],c);
189  protnorm[nt-1] += protmul[counter];
190  }
191  }
192  }
193  }
194  }
195 
196  for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
197  for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
198  counter = -1;
199  for( npos=0; npos<numSec/3; ++npos )
200  {
201  for( nneg=std::max(0,npos-1); nneg<=(npos+1); ++nneg )
202  {
203  for( nzero=0; nzero<numSec/3; ++nzero )
204  {
205  if( ++counter < numMul )
206  {
207  nt = npos+nneg+nzero;
208  if( nt>0 && nt<=numSec )
209  {
210  neutmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[1],c);
211  neutnorm[nt-1] += neutmul[counter];
212  }
213  }
214  }
215  }
216  }
217  for (i = 0; i < numSec; ++i) {
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
229 
230  // energetically possible to produce pion(s) --> inelastic scattering
231  G4double n, anpn;
232  GetNormalizationConstant(availableEnergy, n, anpn);
233  G4double ran = G4UniformRand();
234  G4double dum, excs = 0.0;
235  if (targetParticle.GetDefinition() == aProton) {
236  counter = -1;
237  for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
238  {
239  for( nneg=npos; nneg<=(npos+2) && ran>=excs; ++nneg )
240  {
241  for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
242  {
243  if( ++counter < numMul )
244  {
245  nt = npos+nneg+nzero;
246  if( nt>0 && nt<=numSec )
247  {
248  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
249  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
250  if( std::fabs(dum) < 1.0 )
251  {
252  if( test >= 1.0e-10 )excs += dum*test;
253  }
254  else
255  excs += dum*test;
256  }
257  }
258  }
259  }
260  }
261  if( ran >= excs ) // 3 previous loops continued to the end
262  {
263  quasiElastic = true;
264  return;
265  }
266  npos--; nneg--; nzero--;
267  switch( std::min( 3, std::max( 1, npos-nneg+3 ) ) )
268  {
269  case 1:
270  if( G4UniformRand() < 0.5 )
271  currentParticle.SetDefinitionAndUpdateE( aLambda );
272  else
273  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
274  incidentHasChanged = true;
275  targetParticle.SetDefinitionAndUpdateE( aNeutron );
276  targetHasChanged = true;
277  break;
278  case 2:
279  if( G4UniformRand() < 0.5 )
280  {
281  targetParticle.SetDefinitionAndUpdateE( aNeutron );
282  targetHasChanged = true;
283  }
284  else
285  {
286  if( G4UniformRand() < 0.5 )
287  currentParticle.SetDefinitionAndUpdateE( aLambda );
288  else
289  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
290  incidentHasChanged = true;
291  }
292  break;
293  case 3:
294  break;
295  }
296  }
297  else // target must be a neutron
298  {
299  counter = -1;
300  for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
301  {
302  for( nneg=std::max(0,npos-1); nneg<=(npos+1) && ran>=excs; ++nneg )
303  {
304  for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
305  {
306  if( ++counter < numMul )
307  {
308  nt = npos+nneg+nzero;
309  if( nt>0 && nt<=numSec )
310  {
311  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
312  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
313  if( std::fabs(dum) < 1.0 )
314  {
315  if( test >= 1.0e-10 )excs += dum*test;
316  }
317  else
318  excs += dum*test;
319  }
320  }
321  }
322  }
323  }
324  if( ran >= excs ) // 3 previous loops continued to the end
325  {
326  quasiElastic = true;
327  return;
328  }
329  npos--; nneg--; nzero--;
330  switch( std::min( 3, std::max( 1, npos-nneg+2 ) ) )
331  {
332  case 1:
333  targetParticle.SetDefinitionAndUpdateE( aProton );
334  targetHasChanged = true;
335  break;
336  case 2:
337  if( G4UniformRand() < 0.5 )
338  {
339  if( G4UniformRand() < 0.5 )
340  {
341  currentParticle.SetDefinitionAndUpdateE( aLambda );
342  incidentHasChanged = true;
343  targetParticle.SetDefinitionAndUpdateE( aProton );
344  targetHasChanged = true;
345  }
346  else
347  {
348  targetParticle.SetDefinitionAndUpdateE( aNeutron );
349  targetHasChanged = true;
350  }
351  }
352  else
353  {
354  if( G4UniformRand() < 0.5 )
355  {
356  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
357  incidentHasChanged = true;
358  targetParticle.SetDefinitionAndUpdateE( aProton );
359  targetHasChanged = true;
360  }
361  }
362  break;
363  case 3:
364  if( G4UniformRand() < 0.5 )
365  currentParticle.SetDefinitionAndUpdateE( aLambda );
366  else
367  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
368  incidentHasChanged = true;
369  break;
370  }
371  }
372  SetUpPions(npos, nneg, nzero, vec, vecLen);
373  return;
374 }
375 
376  /* end of file */
377