Geant4  10.02.p01
G4RPGXiMinusInelastic.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: G4RPGXiMinusInelastic.cc 94214 2015-11-09 08:18:05Z gcosmo $
27 //
28 
29 #include "G4RPGXiMinusInelastic.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)
41  {
44  theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
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 << "G4RPGXiMinusInelastic::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 
127 void
129  G4int& vecLen,
130  const G4HadProjectile* originalIncident,
131  G4ReactionProduct& currentParticle,
132  G4ReactionProduct& targetParticle,
133  G4bool& incidentHasChanged,
134  G4bool& targetHasChanged,
135  G4bool& quasiElastic)
136 {
137  // Derived from H. Fesefeldt's original FORTRAN code CASXM
138  //
139  // XiMinus undergoes interaction with nucleon within a nucleus. Check if it is
140  // energetically possible to produce pions/kaons. In not, assume nuclear excitation
141  // occurs and input particle is degraded in energy. No other particles are produced.
142  // If reaction is possible, find the correct number of pions/protons/neutrons
143  // produced using an interpolation to multiplicity data. Replace some pions or
144  // protons/neutrons by kaons or strange baryons according to the average
145  // multiplicity per inelastic reaction.
146 
147  const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
148  const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
149  const G4double targetMass = targetParticle.GetMass()/MeV;
150  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
151  targetMass*targetMass +
152  2.0*targetMass*etOriginal );
153  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
154  if (availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV) {
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 
164  // np = number of pi+, nneg = number of pi-, nz = number of pi0
165  G4int counter, nt = 0, np = 0, nneg = 0, nz = 0;
166  G4double test;
167  const G4double c = 1.25;
168  const G4double b[] = { 0.7, 0.7 };
169  if (first) { // Computation of normalization constants will only be done once
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  for (nneg = std::max(0,np-1); nneg <= (np+1); ++nneg) {
177  for (nz = 0; nz < numSec/3; ++nz) {
178  if (++counter < numMul) {
179  nt = np + nneg + nz;
180  if (nt > 0 && nt <= numSec) {
181  protmul[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
182  protnorm[nt-1] += protmul[counter];
183  }
184  }
185  }
186  }
187  }
188 
189  for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
190  for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
191  counter = -1;
192  for( np=0; np<numSec/3; ++np )
193  {
194  for( nneg=np; nneg<=(np+2); ++nneg )
195  {
196  for( nz=0; nz<numSec/3; ++nz )
197  {
198  if( ++counter < numMul )
199  {
200  nt = np+nneg+nz;
201  if( nt>0 && nt<=numSec )
202  {
203  neutmul[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
204  neutnorm[nt-1] += neutmul[counter];
205  }
206  }
207  }
208  }
209  }
210  for( i=0; i<numSec; ++i )
211  {
212  if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
213  if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
214  }
215  } // end of initialization
216 
217  const G4double expxu = 82.; // upper bound for arg. of exp
218  const G4double expxl = -expxu; // lower bound for arg. of exp
224  //
225  // energetically possible to produce pion(s) --> inelastic scattering
226  //
227  G4double n, anpn;
228  GetNormalizationConstant( availableEnergy, n, anpn );
229  G4double ran = G4UniformRand();
230  G4double dum, excs = 0.0;
231  if( targetParticle.GetDefinition() == aProton )
232  {
233  counter = -1;
234  for( np=0; np<numSec/3 && ran>=excs; ++np )
235  {
236  for( nneg=std::max(0,np-1); nneg<=(np+1) && ran>=excs; ++nneg )
237  {
238  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
239  {
240  if( ++counter < numMul )
241  {
242  nt = np+nneg+nz;
243  if( nt>0 && nt<=numSec )
244  {
245  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
246  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
247  if( std::fabs(dum) < 1.0 )
248  {
249  if( test >= 1.0e-10 )excs += dum*test;
250  }
251  else
252  excs += dum*test;
253  }
254  }
255  }
256  }
257  }
258  if( ran >= excs ) // 3 previous loops continued to the end
259  {
260  quasiElastic = true;
261  return;
262  }
263  np--; nneg--; nz--;
264  //
265  // number of secondary mesons determined by kno distribution
266  // check for total charge of final state mesons to determine
267  // the kind of baryons to be produced, taking into account
268  // charge and strangeness conservation
269  //
270  if( np < nneg )
271  {
272  if( np+1 == nneg )
273  {
274  currentParticle.SetDefinitionAndUpdateE( aXiZero );
275  incidentHasChanged = true;
276  }
277  else // charge mismatch
278  {
279  currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
280  incidentHasChanged = true;
281  //
282  // correct the strangeness by replacing a pi- by a kaon-
283  //
284  vec.Initialize( 1 );
286  p->SetDefinition( aKaonMinus );
287  (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
288  vec.SetElement( vecLen++, p );
289  --nneg;
290  }
291  }
292  else if( np == nneg )
293  {
294  if( G4UniformRand() >= 0.5 )
295  {
296  currentParticle.SetDefinitionAndUpdateE( aXiZero );
297  incidentHasChanged = true;
298  targetParticle.SetDefinitionAndUpdateE( aNeutron );
299  targetHasChanged = true;
300  }
301  }
302  else
303  {
304  targetParticle.SetDefinitionAndUpdateE( aNeutron );
305  targetHasChanged = true;
306  }
307  }
308  else // target must be a neutron
309  {
310  counter = -1;
311  for( np=0; np<numSec/3 && ran>=excs; ++np )
312  {
313  for( nneg=np; nneg<=(np+2) && ran>=excs; ++nneg )
314  {
315  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
316  {
317  if( ++counter < numMul )
318  {
319  nt = np+nneg+nz;
320  if( nt>0 && nt<=numSec )
321  {
322  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
323  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
324  if( std::fabs(dum) < 1.0 )
325  {
326  if( test >= 1.0e-10 )excs += dum*test;
327  }
328  else
329  excs += dum*test;
330  }
331  }
332  }
333  }
334  }
335  if( ran >= excs ) // 3 previous loops continued to the end
336  {
337  quasiElastic = true;
338  return;
339  }
340  np--; nneg--; nz--;
341  if( np+1 < nneg )
342  {
343  if( np+2 == nneg )
344  {
345  currentParticle.SetDefinitionAndUpdateE( aXiZero );
346  incidentHasChanged = true;
347  targetParticle.SetDefinitionAndUpdateE( aProton );
348  targetHasChanged = true;
349  }
350  else // charge mismatch
351  {
352  currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
353  incidentHasChanged = true;
354  targetParticle.SetDefinitionAndUpdateE( aProton );
355  targetHasChanged = true;
356  //
357  // correct the strangeness by replacing a pi- by a kaon-
358  //
359  vec.Initialize( 1 );
361  p->SetDefinition( aKaonMinus );
362  (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
363  vec.SetElement( vecLen++, p );
364  --nneg;
365  }
366  }
367  else if( np+1 == nneg )
368  {
369  if( G4UniformRand() < 0.5 )
370  {
371  currentParticle.SetDefinitionAndUpdateE( aXiZero );
372  incidentHasChanged = true;
373  }
374  else
375  {
376  targetParticle.SetDefinitionAndUpdateE( aProton );
377  targetHasChanged = true;
378  }
379  }
380  }
381  SetUpPions(np, nneg, nz, vec, vecLen);
382  return;
383 }
384 
385  /* end of file */
386 
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 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)
static G4XiZero * XiZero()
Definition: G4XiZero.cc:106
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 G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
void Cascade(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool &quasiElastic)
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
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
G4ThreeVector GetMomentum() const
#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