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