Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros 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: G4RPGSigmaPlusInelastic.cc 94214 2015-11-09 08:18:05Z gcosmo $
27 //
28 
30 #include "G4Exp.hh"
31 #include "G4PhysicalConstants.hh"
32 #include "G4SystemOfUnits.hh"
33 #include "Randomize.hh"
34 
37  G4Nucleus &targetNucleus )
38 {
39  const G4HadProjectile *originalIncident = &aTrack;
40  if (originalIncident->GetKineticEnergy()<= 0.1*MeV)
41  {
45  return &theParticleChange;
46  }
47 
48  // create the target particle
49 
50  G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
51 
52  if( verboseLevel > 1 )
53  {
54  const G4Material *targetMaterial = aTrack.GetMaterial();
55  G4cout << "G4RPGSigmaPlusInelastic::ApplyYourself called" << G4endl;
56  G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
57  G4cout << "target material = " << targetMaterial->GetName() << ", ";
58  G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
59  << G4endl;
60  }
61 
62  // Fermi motion and evaporation
63  // As of Geant3, the Fermi energy calculation had not been Done
64 
65  G4double ek = originalIncident->GetKineticEnergy()/MeV;
66  G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
67  G4ReactionProduct modifiedOriginal;
68  modifiedOriginal = *originalIncident;
69 
70  G4double tkin = targetNucleus.Cinema( ek );
71  ek += tkin;
72  modifiedOriginal.SetKineticEnergy( ek*MeV );
73  G4double et = ek + amas;
74  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
75  G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
76  if( pp > 0.0 )
77  {
78  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
79  modifiedOriginal.SetMomentum( momentum * (p/pp) );
80  }
81  //
82  // calculate black track energies
83  //
84  tkin = targetNucleus.EvaporationEffects( ek );
85  ek -= tkin;
86  modifiedOriginal.SetKineticEnergy( ek*MeV );
87  et = ek + amas;
88  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
89  pp = modifiedOriginal.GetMomentum().mag()/MeV;
90  if( pp > 0.0 )
91  {
92  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
93  modifiedOriginal.SetMomentum( momentum * (p/pp) );
94  }
95  G4ReactionProduct currentParticle = modifiedOriginal;
96  G4ReactionProduct targetParticle;
97  targetParticle = *originalTarget;
98  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
99  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
100  G4bool incidentHasChanged = false;
101  G4bool targetHasChanged = false;
102  G4bool quasiElastic = false;
103  G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
104  G4int vecLen = 0;
105  vec.Initialize( 0 );
106 
107  const G4double cutOff = 0.1;
108  if( currentParticle.GetKineticEnergy()/MeV > cutOff )
109  Cascade( vec, vecLen,
110  originalIncident, currentParticle, targetParticle,
111  incidentHasChanged, targetHasChanged, quasiElastic );
112 
113  CalculateMomenta( vec, vecLen,
114  originalIncident, originalTarget, modifiedOriginal,
115  targetNucleus, currentParticle, targetParticle,
116  incidentHasChanged, targetHasChanged, quasiElastic );
117 
118  SetUpChange( vec, vecLen,
119  currentParticle, targetParticle,
120  incidentHasChanged );
121 
122  delete originalTarget;
123  return &theParticleChange;
124 }
125 
126 void G4RPGSigmaPlusInelastic::Cascade(
128  G4int& vecLen,
129  const G4HadProjectile *originalIncident,
130  G4ReactionProduct &currentParticle,
131  G4ReactionProduct &targetParticle,
132  G4bool &incidentHasChanged,
133  G4bool &targetHasChanged,
134  G4bool &quasiElastic )
135 {
136  // Derived from H. Fesefeldt's original FORTRAN code CASSP
137  //
138  // SigmaPlus undergoes interaction with nucleon within a nucleus. Check if it is
139  // energetically possible to produce pions/kaons. In not, assume nuclear excitation
140  // occurs and input particle is degraded in energy. No other particles are produced.
141  // If reaction is possible, find the correct number of pions/protons/neutrons
142  // produced using an interpolation to multiplicity data. Replace some pions or
143  // protons/neutrons by kaons or strange baryons according to the average
144  // multiplicity per inelastic reaction.
145 
146  const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
147  const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
148  const G4double targetMass = targetParticle.GetMass()/MeV;
149  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
150  targetMass*targetMass +
151  2.0*targetMass*etOriginal );
152  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
153  if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV )
154  {
155  quasiElastic = true;
156  return;
157  }
158  static G4ThreadLocal G4bool first = true;
159  const G4int numMul = 1200;
160  const G4int numSec = 60;
161  static G4ThreadLocal G4double protmul[numMul], protnorm[numSec]; // proton constants
162  static G4ThreadLocal G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
163  // np = number of pi+, nneg = number of pi-, nz = number of pi0
164  G4int counter, nt=0, np=0, nneg=0, nz=0;
165  G4double test;
166  const G4double c = 1.25;
167  const G4double b[] = { 0.7, 0.7 };
168  if( first ) // compute normalization constants, this will only be Done once
169  {
170  first = false;
171  G4int i;
172  for( i=0; i<numMul; ++i )protmul[i] = 0.0;
173  for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
174  counter = -1;
175  for( np=0; np<(numSec/3); ++np )
176  {
177  for( nneg=np; nneg<=(np+2); ++nneg )
178  {
179  for( nz=0; nz<numSec/3; ++nz )
180  {
181  if( ++counter < numMul )
182  {
183  nt = np+nneg+nz;
184  if( nt>0 && nt<=numSec )
185  {
186  protmul[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
187  protnorm[nt-1] += protmul[counter];
188  }
189  }
190  }
191  }
192  }
193  for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
194  for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
195  counter = -1;
196  for( np=0; np<numSec/3; ++np )
197  {
198  for( nneg=std::max(0,np-1); nneg<=(np+1); ++nneg )
199  {
200  for( nz=0; nz<numSec/3; ++nz )
201  {
202  if( ++counter < numMul )
203  {
204  nt = np+nneg+nz;
205  if( nt>0 && nt<=numSec )
206  {
207  neutmul[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
208  neutnorm[nt-1] += neutmul[counter];
209  }
210  }
211  }
212  }
213  }
214  for( i=0; i<numSec; ++i )
215  {
216  if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
217  if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
218  }
219  } // end of initialization
220 
221  const G4double expxu = 82.; // upper bound for arg. of exp
222  const G4double expxl = -expxu; // lower bound for arg. of exp
227  //
228  // energetically possible to produce pion(s) --> inelastic scattering
229  //
230  G4double n, anpn;
231  GetNormalizationConstant( availableEnergy, n, anpn );
232  G4double ran = G4UniformRand();
233  G4double dum, excs = 0.0;
234  if( targetParticle.GetDefinition() == aProton )
235  {
236  counter = -1;
237  for( np=0; np<numSec/3 && ran>=excs; ++np )
238  {
239  for( nneg=np; nneg<=(np+2) && ran>=excs; ++nneg )
240  {
241  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
242  {
243  if( ++counter < numMul )
244  {
245  nt = np+nneg+nz;
246  if( nt>0 && nt<=numSec )
247  {
248  test = G4Exp( 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  np--; nneg--; nz--;
267  switch( std::min( 3, std::max( 1, np-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( np=0; np<numSec/3 && ran>=excs; ++np )
301  {
302  for( nneg=std::max(0,np-1); nneg<=(np+1) && ran>=excs; ++nneg )
303  {
304  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
305  {
306  if( ++counter < numMul )
307  {
308  nt = np+nneg+nz;
309  if( nt>0 && nt<=numSec )
310  {
311  test = G4Exp( 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  np--; nneg--; nz--;
330  switch( std::min( 3, std::max( 1, np-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 
373  SetUpPions(np, nneg, nz, vec, vecLen);
374  return;
375 }
376 
377  /* end of file */
378 
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:278
void SetUpChange(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
void SetKineticEnergy(const G4double en)
void SetMomentum(const G4double x, const G4double y, const G4double z)
const char * p
Definition: xmltok.h:285
const G4String & GetName() const
Definition: G4Material.hh:178
void SetSide(const G4int sid)
void CalculateMomenta(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, const G4DynamicParticle *originalTarget, G4ReactionProduct &modifiedOriginal, G4Nucleus &targetNucleus, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool quasiElastic)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
G4ParticleDefinition * GetDefinition() const
#define G4ThreadLocal
Definition: tls.hh:89
void Initialize(G4int items)
Definition: G4FastVector.hh:63
int G4int
Definition: G4Types.hh:78
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:241
const G4String & GetParticleName() const
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:102
void SetStatusChange(G4HadFinalStateStatus aS)
G4double Pmltpc(G4int np, G4int nm, G4int nz, G4int n, G4double b, G4double c)
const G4ParticleDefinition * GetDefinition() const
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ParticleDefinition * GetDefinition() const
bool G4bool
Definition: G4Types.hh:79
G4double GetKineticEnergy() const
G4double ek
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
void SetDefinitionAndUpdateE(const G4ParticleDefinition *aParticleDefinition)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
const G4LorentzVector & Get4Momentum() const
G4double GetKineticEnergy() const
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
void SetEnergyChange(G4double anEnergy)
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
Hep3Vector unit() const
G4double Cinema(G4double kineticEnergy)
Definition: G4Nucleus.cc:382
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4ThreeVector GetMomentum() const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
const G4Material * GetMaterial() const
static constexpr double pi
Definition: G4SIunits.hh:75
void GetNormalizationConstant(const G4double availableEnergy, G4double &n, G4double &anpn)
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108
double G4double
Definition: G4Types.hh:76
void SetUpPions(const G4int np, const G4int nm, const G4int nz, G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen)
double mag() const
void SetMomentumChange(const G4ThreeVector &aV)
G4double GetMass() const
G4double GetTotalEnergy() const