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