Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4RPGLambdaInelastic.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: G4RPGLambdaInelastic.cc 94214 2015-11-09 08:18:05Z gcosmo $
27 //
28 
29 #include "G4RPGLambdaInelastic.hh"
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 
41  // create the target particle
42 
43  G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
44 
45  if( verboseLevel > 1 )
46  {
47  const G4Material *targetMaterial = aTrack.GetMaterial();
48  G4cout << "G4RPGLambdaInelastic::ApplyYourself called" << G4endl;
49  G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
50  G4cout << "target material = " << targetMaterial->GetName() << ", ";
51  G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
52  << G4endl;
53  }
54 
55  // Fermi motion and evaporation
56  // As of Geant3, the Fermi energy calculation had not been Done
57 
58  G4double ek = originalIncident->GetKineticEnergy()/MeV;
59  G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
60  G4ReactionProduct modifiedOriginal;
61  modifiedOriginal = *originalIncident;
62 
63  G4double tkin = targetNucleus.Cinema( ek );
64  ek += tkin;
65  modifiedOriginal.SetKineticEnergy( ek*MeV );
66  G4double et = ek + amas;
67  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
68  G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
69  if( pp > 0.0 )
70  {
71  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
72  modifiedOriginal.SetMomentum( momentum * (p/pp) );
73  }
74  //
75  // calculate black track energies
76  //
77  tkin = targetNucleus.EvaporationEffects( ek );
78  ek -= tkin;
79  modifiedOriginal.SetKineticEnergy( ek*MeV );
80  et = ek + amas;
81  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
82  pp = modifiedOriginal.GetMomentum().mag()/MeV;
83  if( pp > 0.0 )
84  {
85  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
86  modifiedOriginal.SetMomentum( momentum * (p/pp) );
87  }
88 
89  G4ReactionProduct currentParticle = modifiedOriginal;
90  G4ReactionProduct targetParticle;
91  targetParticle = *originalTarget;
92  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
93  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
94  G4bool incidentHasChanged = false;
95  G4bool targetHasChanged = false;
96  G4bool quasiElastic = false;
97  G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
98  G4int vecLen = 0;
99  vec.Initialize( 0 );
100 
101  const G4double cutOff = 0.1;
102  if( currentParticle.GetKineticEnergy()/MeV > cutOff )
103  Cascade( vec, vecLen,
104  originalIncident, currentParticle, targetParticle,
105  incidentHasChanged, targetHasChanged, quasiElastic );
106 
107  CalculateMomenta( vec, vecLen,
108  originalIncident, originalTarget, modifiedOriginal,
109  targetNucleus, currentParticle, targetParticle,
110  incidentHasChanged, targetHasChanged, quasiElastic );
111 
112  SetUpChange( vec, vecLen,
113  currentParticle, targetParticle,
114  incidentHasChanged );
115 
116  delete originalTarget;
117  return &theParticleChange;
118 }
119 
120 
121 void G4RPGLambdaInelastic::Cascade(
123  G4int& vecLen,
124  const G4HadProjectile *originalIncident,
125  G4ReactionProduct &currentParticle,
126  G4ReactionProduct &targetParticle,
127  G4bool &incidentHasChanged,
128  G4bool &targetHasChanged,
129  G4bool &quasiElastic )
130 {
131  // Derived from H. Fesefeldt's original FORTRAN code CASL0
132  //
133  // Lambda undergoes interaction with nucleon within a nucleus. Check if it is
134  // energetically possible to produce pions/kaons. In not, assume nuclear excitation
135  // occurs and input particle is degraded in energy. No other particles are produced.
136  // If reaction is possible, find the correct number of pions/protons/neutrons
137  // produced using an interpolation to multiplicity data. Replace some pions or
138  // protons/neutrons by kaons or strange baryons according to the average
139  // multiplicity per Inelastic reaction.
140 
141  const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
142  const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
143  const G4double targetMass = targetParticle.GetMass()/MeV;
144  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
145  targetMass*targetMass +
146  2.0*targetMass*etOriginal );
147  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
148  if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV )
149  {
150  quasiElastic = true;
151  return;
152  }
153  static G4ThreadLocal G4bool first = true;
154  const G4int numMul = 1200;
155  const G4int numSec = 60;
156  static G4ThreadLocal G4double protmul[numMul], protnorm[numSec]; // proton constants
157  static G4ThreadLocal G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
158 
159  // np = number of pi+, nneg = number of pi-, nz = number of pi0
160 
161  G4int counter, nt=0, np=0, nneg=0, nz=0;
162  G4double test;
163  const G4double c = 1.25;
164  const G4double b[] = { 0.70, 0.35 };
165  if( first ) { // compute normalization constants, this will only be Done once
166  first = false;
167  G4int i;
168  for( i=0; i<numMul; ++i )protmul[i] = 0.0;
169  for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
170  counter = -1;
171  for( np=0; np<(numSec/3); ++np ) {
172  for( nneg=std::max(0,np-2); nneg<=(np+1); ++nneg ) {
173  for( nz=0; nz<numSec/3; ++nz ) {
174  if( ++counter < numMul ) {
175  nt = np+nneg+nz;
176  if( nt>0 && nt<=numSec ) {
177  protmul[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
178  protnorm[nt-1] += protmul[counter];
179  }
180  }
181  }
182  }
183  }
184  for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
185  for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
186  counter = -1;
187  for( np=0; np<numSec/3; ++np ) {
188  for( nneg=std::max(0,np-1); nneg<=(np+2); ++nneg ) {
189  for( nz=0; nz<numSec/3; ++nz ) {
190  if( ++counter < numMul ) {
191  nt = np+nneg+nz;
192  if( nt>0 && nt<=numSec ) {
193  neutmul[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
194  neutnorm[nt-1] += neutmul[counter];
195  }
196  }
197  }
198  }
199  }
200  for( i=0; i<numSec; ++i ) {
201  if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
202  if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
203  }
204  } // end of initialization
205 
206  const G4double expxu = 82.; // upper bound for arg. of exp
207  const G4double expxl = -expxu; // lower bound for arg. of exp
213 
214  // energetically possible to produce pion(s) --> inelastic scattering
215  // otherwise quasi-elastic scattering
216 
217  G4double n, anpn;
218  GetNormalizationConstant( availableEnergy, n, anpn );
219  G4double ran = G4UniformRand();
220  G4double dum, excs = 0.0;
221  if( targetParticle.GetDefinition() == aProton ) {
222  counter = -1;
223  for( np=0; np<numSec/3 && ran>=excs; ++np ) {
224  for( nneg=std::max(0,np-2); nneg<=(np+1) && ran>=excs; ++nneg ) {
225  for( nz=0; nz<numSec/3 && ran>=excs; ++nz ) {
226  if( ++counter < numMul ) {
227  nt = np+nneg+nz;
228  if( nt>0 && nt<=numSec ) {
229  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
230  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
231  if( std::fabs(dum) < 1.0 ) {
232  if( test >= 1.0e-10 )excs += dum*test;
233  } else {
234  excs += dum*test;
235  }
236  }
237  }
238  }
239  }
240  }
241  if( ran >= excs ) // 3 previous loops continued to the end
242  {
243  quasiElastic = true;
244  return;
245  }
246  np--; nneg--; nz--;
247  G4int ncht = std::max( 1, np-nneg );
248  switch( ncht ) {
249  case 1:
250  currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
251  incidentHasChanged = true;
252  break;
253  case 2:
254  if( G4UniformRand() < 0.5 ) {
255  if( G4UniformRand() < 0.5 ) {
256  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
257  incidentHasChanged = true;
258  }
259  } else {
260  currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
261  incidentHasChanged = true;
262  targetParticle.SetDefinitionAndUpdateE( aNeutron );
263  targetHasChanged = true;
264  }
265  break;
266  case 3:
267  if( G4UniformRand() < 0.5 ) {
268  if( G4UniformRand() < 0.5 ) {
269  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
270  incidentHasChanged = true;
271  }
272  targetParticle.SetDefinitionAndUpdateE( aNeutron );
273  targetHasChanged = true;
274  } else {
275  currentParticle.SetDefinitionAndUpdateE( aSigmaMinus );
276  incidentHasChanged = true;
277  }
278  break;
279  default:
280  currentParticle.SetDefinitionAndUpdateE( aSigmaMinus );
281  incidentHasChanged = true;
282  targetParticle.SetDefinitionAndUpdateE( aNeutron );
283  targetHasChanged = true;
284  break;
285  }
286  }
287  else // target must be a neutron
288  {
289  counter = -1;
290  for( np=0; np<numSec/3 && ran>=excs; ++np ) {
291  for( nneg=std::max(0,np-1); nneg<=(np+2) && ran>=excs; ++nneg ) {
292  for( nz=0; nz<numSec/3 && ran>=excs; ++nz ) {
293  if( ++counter < numMul ) {
294  nt = np+nneg+nz;
295  if( nt>0 && nt<=numSec ) {
296  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
297  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
298  if( std::fabs(dum) < 1.0 ) {
299  if( test >= 1.0e-10 )excs += dum*test;
300  } else {
301  excs += dum*test;
302  }
303  }
304  }
305  }
306  }
307  }
308  if( ran >= excs ) // 3 previous loops continued to the end
309  {
310  quasiElastic = true;
311  return;
312  }
313  np--; nneg--; nz--;
314  G4int ncht = std::max( 1, np-nneg+3 );
315  switch( ncht ) {
316  case 1:
317  currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
318  incidentHasChanged = true;
319  targetParticle.SetDefinitionAndUpdateE( aProton );
320  targetHasChanged = true;
321  break;
322  case 2:
323  if( G4UniformRand() < 0.5 ) {
324  if( G4UniformRand() < 0.5 ) {
325  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
326  incidentHasChanged = true;
327  }
328  targetParticle.SetDefinitionAndUpdateE( aProton );
329  targetHasChanged = true;
330  } else {
331  currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
332  incidentHasChanged = true;
333  }
334  break;
335  case 3:
336  if( G4UniformRand() < 0.5 ) {
337  if( G4UniformRand() < 0.5 ) {
338  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
339  incidentHasChanged = true;
340  }
341  } else {
342  currentParticle.SetDefinitionAndUpdateE( aSigmaMinus );
343  incidentHasChanged = true;
344  targetParticle.SetDefinitionAndUpdateE( aProton );
345  targetHasChanged = true;
346  }
347  break;
348  default:
349  currentParticle.SetDefinitionAndUpdateE( aSigmaMinus );
350  incidentHasChanged = true;
351  break;
352  }
353  }
354 
355  SetUpPions(np, nneg, nz, vec, vecLen);
356  return;
357 }
358 
359  /* end of file */
360 
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)
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
G4double Pmltpc(G4int np, G4int nm, G4int nz, G4int n, G4double b, G4double c)
const G4ParticleDefinition * GetDefinition() 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
static G4SigmaMinus * SigmaMinus()
G4double GetKineticEnergy() const
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
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
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
void GetNormalizationConstant(const G4double availableEnergy, G4double &n, G4double &anpn)
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
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
G4double GetMass() const
G4double GetTotalEnergy() const