Geant4_10
G4RPGAntiProtonInelastic.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$
27 //
28 
30 #include "G4PhysicalConstants.hh"
31 #include "G4SystemOfUnits.hh"
32 #include "Randomize.hh"
33 
36  G4Nucleus &targetNucleus )
37 {
38  const G4HadProjectile* originalIncident = &aTrack;
39 
40  if (originalIncident->GetKineticEnergy() <= 0.1*MeV) {
44  return &theParticleChange;
45  }
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 << "G4RPGAntiProtonInelastic::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  const G4double anni = std::min( 1.3*originalIncident->GetTotalMomentum()/GeV, 0.4 );
109 
110  if( (currentParticle.GetKineticEnergy()/MeV > cutOff) ||
111  (G4UniformRand() > anni) )
112  Cascade( vec, vecLen,
113  originalIncident, currentParticle, targetParticle,
114  incidentHasChanged, targetHasChanged, quasiElastic );
115  else
116  quasiElastic = true;
117 
118  CalculateMomenta( vec, vecLen,
119  originalIncident, originalTarget, modifiedOriginal,
120  targetNucleus, currentParticle, targetParticle,
121  incidentHasChanged, targetHasChanged, quasiElastic );
122 
123  SetUpChange( vec, vecLen,
124  currentParticle, targetParticle,
125  incidentHasChanged );
126 
127  delete originalTarget;
128  return &theParticleChange;
129 }
130 
131 
132 void G4RPGAntiProtonInelastic::Cascade(
134  G4int& vecLen,
135  const G4HadProjectile *originalIncident,
136  G4ReactionProduct &currentParticle,
137  G4ReactionProduct &targetParticle,
138  G4bool &incidentHasChanged,
139  G4bool &targetHasChanged,
140  G4bool &quasiElastic )
141 {
142  // Derived from H. Fesefeldt's original FORTRAN code CASPB
143  // AntiProton undergoes interaction with nucleon within a nucleus. Check if it is
144  // energetically possible to produce pions/kaons. In not, assume nuclear excitation
145  // occurs and input particle is degraded in energy. No other particles are produced.
146  // If reaction is possible, find the correct number of pions/protons/neutrons
147  // produced using an interpolation to multiplicity data. Replace some pions or
148  // protons/neutrons by kaons or strange baryons according to the average
149  // multiplicity per Inelastic reaction.
150 
151  const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
152  const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
153  const G4double pOriginal = originalIncident->GetTotalMomentum()/MeV;
154  const G4double targetMass = targetParticle.GetMass()/MeV;
155  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
156  targetMass*targetMass +
157  2.0*targetMass*etOriginal );
158  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
159 
160  static G4ThreadLocal G4bool first = true;
161  const G4int numMul = 1200;
162  const G4int numMulA = 400;
163  const G4int numSec = 60;
164  static G4ThreadLocal G4double protmul[numMul], protnorm[numSec]; // proton constants
165  static G4ThreadLocal G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
166  static G4ThreadLocal G4double protmulA[numMulA], protnormA[numSec]; // proton constants
167  static G4ThreadLocal G4double neutmulA[numMulA], neutnormA[numSec]; // neutron constants
168  // np = number of pi+, nneg = number of pi-, nz = number of pi0
169  G4int counter, nt=0, np=0, nneg=0, nz=0;
170  G4double test;
171  const G4double c = 1.25;
172  const G4double b[] = { 0.7, 0.7 };
173  if( first ) // compute normalization constants, this will only be Done once
174  {
175  first = false;
176  G4int i;
177  for( i=0; i<numMul; ++i )protmul[i] = 0.0;
178  for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
179  counter = -1;
180  for( np=0; np<(numSec/3); ++np )
181  {
182  for( nneg=std::max(0,np-1); nneg<=(np+1); ++nneg )
183  {
184  for( nz=0; nz<numSec/3; ++nz )
185  {
186  if( ++counter < numMul )
187  {
188  nt = np+nneg+nz;
189  if( nt>0 && nt<=numSec )
190  {
191  protmul[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
192  protnorm[nt-1] += protmul[counter];
193  }
194  }
195  }
196  }
197  }
198  for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
199  for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
200  counter = -1;
201  for( np=0; np<numSec/3; ++np )
202  {
203  for( nneg=np; nneg<=(np+2); ++nneg )
204  {
205  for( nz=0; nz<numSec/3; ++nz )
206  {
207  if( ++counter < numMul )
208  {
209  nt = np+nneg+nz;
210  if( (nt>0) && (nt<=numSec) )
211  {
212  neutmul[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
213  neutnorm[nt-1] += neutmul[counter];
214  }
215  }
216  }
217  }
218  }
219  for( i=0; i<numSec; ++i )
220  {
221  if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
222  if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
223  }
224  //
225  // do the same for annihilation channels
226  //
227  for( i=0; i<numMulA; ++i )protmulA[i] = 0.0;
228  for( i=0; i<numSec; ++i )protnormA[i] = 0.0;
229  counter = -1;
230  for( np=0; np<(numSec/3); ++np )
231  {
232  nneg = np;
233  for( nz=0; nz<numSec/3; ++nz )
234  {
235  if( ++counter < numMulA )
236  {
237  nt = np+nneg+nz;
238  if( nt>1 && nt<=numSec )
239  {
240  protmulA[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
241  protnormA[nt-1] += protmulA[counter];
242  }
243  }
244  }
245  }
246  for( i=0; i<numMulA; ++i )neutmulA[i] = 0.0;
247  for( i=0; i<numSec; ++i )neutnormA[i] = 0.0;
248  counter = -1;
249  for( np=0; np<numSec/3; ++np )
250  {
251  nneg = np+1;
252  for( nz=0; nz<numSec/3; ++nz )
253  {
254  if( ++counter < numMulA )
255  {
256  nt = np+nneg+nz;
257  if( (nt>1) && (nt<=numSec) )
258  {
259  neutmulA[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
260  neutnormA[nt-1] += neutmulA[counter];
261  }
262  }
263  }
264  }
265  for( i=0; i<numSec; ++i )
266  {
267  if( protnormA[i] > 0.0 )protnormA[i] = 1.0/protnormA[i];
268  if( neutnormA[i] > 0.0 )neutnormA[i] = 1.0/neutnormA[i];
269  }
270  } // end of initialization
271  //
272  // initialization is OK
273  //
274  const G4double expxu = 82.; // upper bound for arg. of exp
275  const G4double expxl = -expxu; // lower bound for arg. of exp
276 
281 
282  const G4double anhl[] = {1.0,1.0,1.0,1.0,1.0,1.0,1.0,1.0,1.0,0.90,
283  0.6,0.52,0.47,0.44,0.41,0.39,0.37,0.35,0.34,0.24,
284  0.19,0.15,0.12,0.10,0.09,0.07,0.06,0.05,0.0};
285 
286  G4int iplab = G4int( pOriginal/GeV*10.0 );
287  if( iplab > 9 )iplab = G4int( pOriginal/GeV ) + 9;
288  if( iplab > 18 )iplab = G4int( pOriginal/GeV/10.0 ) + 18;
289  if( iplab > 27 )iplab = 28;
290  if( G4UniformRand() > anhl[iplab] )
291  {
292  if( availableEnergy <= aPiPlus->GetPDGMass()/MeV )
293  {
294  quasiElastic = true;
295  return;
296  }
297  G4int ieab = static_cast<G4int>(availableEnergy*5.0/GeV);
298  const G4double supp[] = {0.,0.4,0.55,0.65,0.75,0.82,0.86,0.90,0.94,0.98};
299  G4double w0, wp, wt, wm;
300  if( (availableEnergy < 2.0*GeV) && (G4UniformRand() >= supp[ieab]) )
301  {
302  //
303  // suppress high multiplicity events at low momentum
304  // only one pion will be produced
305  //
306  np = nneg = nz = 0;
307  if( targetParticle.GetDefinition() == aProton )
308  {
309  test = std::exp( std::min( expxu, std::max( expxl, -(1.0+b[1])*(1.0+b[1])/(2.0*c*c) ) ) );
310  w0 = test;
311  wp = test;
312  test = std::exp( std::min( expxu, std::max( expxl, -(-1.0+b[1])*(-1.0+b[1])/(2.0*c*c) ) ) );
313  wm = test;
314  wt = w0+wp+wm;
315  wp += w0;
316  G4double ran = G4UniformRand();
317  if( ran < w0/wt )
318  nz = 1;
319  else if( ran < wp/wt )
320  np = 1;
321  else
322  nneg = 1;
323  }
324  else // target is a neutron
325  {
326  test = std::exp( std::min( expxu, std::max( expxl, -(1.0+b[0])*(1.0+b[0])/(2.0*c*c) ) ) );
327  w0 = test;
328  test = std::exp( std::min( expxu, std::max( expxl, -(-1.0+b[0])*(-1.0+b[0])/(2.0*c*c) ) ) );
329  wm = test;
330  G4double ran = G4UniformRand();
331  if( ran < w0/(w0+wm) )
332  nz = 1;
333  else
334  nneg = 1;
335  }
336  }
337  else // (availableEnergy >= 2.0*GeV) || (random number < supp[ieab])
338  {
339  G4double n, anpn;
340  GetNormalizationConstant( availableEnergy, n, anpn );
341  G4double ran = G4UniformRand();
342  G4double dum, excs = 0.0;
343  if( targetParticle.GetDefinition() == aProton )
344  {
345  counter = -1;
346  for( np=0; np<numSec/3 && ran>=excs; ++np )
347  {
348  for( nneg=std::max(0,np-1); nneg<=(np+1) && ran>=excs; ++nneg )
349  {
350  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
351  {
352  if( ++counter < numMul )
353  {
354  nt = np+nneg+nz;
355  if( (nt>0) && (nt<=numSec) )
356  {
357  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
358  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
359  if( std::fabs(dum) < 1.0 )
360  {
361  if( test >= 1.0e-10 )excs += dum*test;
362  }
363  else
364  excs += dum*test;
365  }
366  }
367  }
368  }
369  }
370  if( ran >= excs ) // 3 previous loops continued to the end
371  {
372  quasiElastic = true;
373  return;
374  }
375  }
376  else // target must be a neutron
377  {
378  counter = -1;
379  for( np=0; np<numSec/3 && ran>=excs; ++np )
380  {
381  for( nneg=np; nneg<=(np+2) && ran>=excs; ++nneg )
382  {
383  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
384  {
385  if( ++counter < numMul )
386  {
387  nt = np+nneg+nz;
388  if( (nt>0) && (nt<=numSec) )
389  {
390  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
391  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
392  if( std::fabs(dum) < 1.0 )
393  {
394  if( test >= 1.0e-10 )excs += dum*test;
395  }
396  else
397  excs += dum*test;
398  }
399  }
400  }
401  }
402  }
403  if( ran >= excs ) // 3 previous loops continued to the end
404  {
405  quasiElastic = true;
406  return;
407  }
408  }
409  np--; nneg--; nz--;
410  }
411  if( targetParticle.GetDefinition() == aProton )
412  {
413  switch( np-nneg )
414  {
415  case 0:
416  if( G4UniformRand() < 0.33 )
417  {
418  currentParticle.SetDefinitionAndUpdateE( anAntiNeutron );
419  targetParticle.SetDefinitionAndUpdateE( aNeutron );
420  incidentHasChanged = true;
421  targetHasChanged = true;
422  }
423  break;
424  case 1:
425  targetParticle.SetDefinitionAndUpdateE( aNeutron );
426  targetHasChanged = true;
427  break;
428  default:
429  currentParticle.SetDefinitionAndUpdateE( anAntiNeutron );
430  incidentHasChanged = true;
431  break;
432  }
433  }
434  else // target must be a neutron
435  {
436  switch( np-nneg )
437  {
438  case -1:
439  if( G4UniformRand() < 0.5 )
440  {
441  targetParticle.SetDefinitionAndUpdateE( aProton );
442  targetHasChanged = true;
443  }
444  else
445  {
446  currentParticle.SetDefinitionAndUpdateE( anAntiNeutron );
447  incidentHasChanged = true;
448  }
449  break;
450  case 0:
451  break;
452  default:
453  currentParticle.SetDefinitionAndUpdateE( anAntiNeutron );
454  targetParticle.SetDefinitionAndUpdateE( aProton );
455  incidentHasChanged = true;
456  targetHasChanged = true;
457  break;
458  }
459  }
460  }
461  else // random number <= anhl[iplab]
462  {
463  if( centerofmassEnergy <= 2*aPiPlus->GetPDGMass()/MeV )
464  {
465  quasiElastic = true;
466  return;
467  }
468  //
469  // annihilation channels
470  //
471  G4double n, anpn;
472  GetNormalizationConstant( -centerofmassEnergy, n, anpn );
473  G4double ran = G4UniformRand();
474  G4double dum, excs = 0.0;
475  if( targetParticle.GetDefinition() == aProton )
476  {
477  counter = -1;
478  for( np=0; (np<numSec/3) && (ran>=excs); ++np )
479  {
480  nneg = np;
481  for( nz=0; (nz<numSec/3) && (ran>=excs); ++nz )
482  {
483  if( ++counter < numMulA )
484  {
485  nt = np+nneg+nz;
486  if( (nt>1) && (nt<=numSec) )
487  {
488  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
489  dum = (pi/anpn)*nt*protmulA[counter]*protnormA[nt-1]/(2.0*n*n);
490  if( std::abs(dum) < 1.0 )
491  {
492  if( test >= 1.0e-10 )excs += dum*test;
493  }
494  else
495  excs += dum*test;
496  }
497  }
498  }
499  }
500  }
501  else // target must be a neutron
502  {
503  counter = -1;
504  for( np=0; (np<numSec/3) && (ran>=excs); ++np )
505  {
506  nneg = np+1;
507  for( nz=0; (nz<numSec/3) && (ran>=excs); ++nz )
508  {
509  if( ++counter < numMulA )
510  {
511  nt = np+nneg+nz;
512  if( (nt>1) && (nt<=numSec) )
513  {
514  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
515  dum = (pi/anpn)*nt*neutmulA[counter]*neutnormA[nt-1]/(2.0*n*n);
516  if( std::fabs(dum) < 1.0 )
517  {
518  if( test >= 1.0e-10 )excs += dum*test;
519  }
520  else
521  excs += dum*test;
522  }
523  }
524  }
525  }
526  }
527  if( ran >= excs ) // 3 previous loops continued to the end
528  {
529  quasiElastic = true;
530  return;
531  }
532  np--; nz--;
533  currentParticle.SetMass( 0.0 );
534  targetParticle.SetMass( 0.0 );
535  }
536 
537  while(np + nneg + nz < 3) nz++;
538  SetUpPions( np, nneg, nz, vec, vecLen );
539  return;
540 }
541 
542  /* end of file */
543 
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:264
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:176
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:52
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
void Initialize(G4int items)
Definition: G4FastVector.hh:63
int G4int
Definition: G4Types.hh:78
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:227
void SetDefinitionAndUpdateE(G4ParticleDefinition *aParticleDefinition)
const G4String & GetParticleName() const
G4ParticleDefinition * GetDefinition() const
void SetStatusChange(G4HadFinalStateStatus aS)
G4double Pmltpc(G4int np, G4int nm, G4int nz, G4int n, G4double b, G4double c)
tuple b
Definition: test.py:12
Char_t n[5]
void SetMass(const G4double mas)
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
const G4ParticleDefinition * GetDefinition() const
G4double ek
bool G4bool
Definition: G4Types.hh:79
G4double GetKineticEnergy() const
TTree * nt
Definition: plotHisto.C:21
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
const G4LorentzVector & Get4Momentum() const
G4double GetKineticEnergy() const
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:368
G4int first
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
const G4Material * GetMaterial() const
def test
Definition: mcscore.py:117
void GetNormalizationConstant(const G4double availableEnergy, G4double &n, G4double &anpn)
double G4double
Definition: G4Types.hh:76
tuple c
Definition: test.py:13
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
static G4AntiNeutron * AntiNeutron()
G4double GetTotalMomentum() const
G4double GetTotalEnergy() const