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