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