Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4RPGKPlusInelastic.cc
이 파일의 문서화 페이지로 가기
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 //
27 // $Id: G4RPGKPlusInelastic.cc 94214 2015-11-09 08:18:05Z gcosmo $
28 //
29 
30 #include "G4RPGKPlusInelastic.hh"
31 #include "G4Exp.hh"
32 #include "G4PhysicalConstants.hh"
33 #include "G4SystemOfUnits.hh"
34 #include "Randomize.hh"
35 
38  G4Nucleus &targetNucleus )
39 {
40  const G4HadProjectile *originalIncident = &aTrack;
41  if (originalIncident->GetKineticEnergy()<= 0.1*MeV)
42  {
46  return &theParticleChange;
47  }
48 
49  // create the target particle
50 
51  G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
52  G4ReactionProduct targetParticle( originalTarget->GetDefinition() );
53 
54  if( verboseLevel > 1 )
55  {
56  const G4Material *targetMaterial = aTrack.GetMaterial();
57  G4cout << "G4RPGKPlusInelastic::ApplyYourself called" << G4endl;
58  G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy() << "MeV, ";
59  G4cout << "target material = " << targetMaterial->GetName() << ", ";
60  G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
61  << G4endl;
62  }
63 
64  G4ReactionProduct currentParticle(originalIncident->GetDefinition() );
65  currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
66  currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
67 
68  // Fermi motion and evaporation
69  // As of Geant3, the Fermi energy calculation had not been Done
70 
71  G4double ek = originalIncident->GetKineticEnergy();
72  G4double amas = originalIncident->GetDefinition()->GetPDGMass();
73 
74  G4double tkin = targetNucleus.Cinema( ek );
75  ek += tkin;
76  currentParticle.SetKineticEnergy( ek );
77  G4double et = ek + amas;
78  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
79  G4double pp = currentParticle.GetMomentum().mag();
80  if( pp > 0.0 )
81  {
82  G4ThreeVector momentum = currentParticle.GetMomentum();
83  currentParticle.SetMomentum( momentum * (p/pp) );
84  }
85 
86  // calculate black track energies
87 
88  tkin = targetNucleus.EvaporationEffects( ek );
89  ek -= tkin;
90  currentParticle.SetKineticEnergy( ek );
91  et = ek + amas;
92  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
93  pp = currentParticle.GetMomentum().mag();
94  if( pp > 0.0 )
95  {
96  G4ThreeVector momentum = currentParticle.GetMomentum();
97  currentParticle.SetMomentum( momentum * (p/pp) );
98  }
99 
100  G4ReactionProduct modifiedOriginal = currentParticle;
101 
102  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
103  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
104  G4bool incidentHasChanged = false;
105  G4bool targetHasChanged = false;
106  G4bool quasiElastic = false;
107  G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
108  G4int vecLen = 0;
109  vec.Initialize( 0 );
110 
111  const G4double cutOff = 0.1*MeV;
112  if( currentParticle.GetKineticEnergy() > cutOff )
113  Cascade( vec, vecLen,
114  originalIncident, currentParticle, targetParticle,
115  incidentHasChanged, targetHasChanged, quasiElastic );
116 
117  CalculateMomenta( vec, vecLen,
118  originalIncident, originalTarget, modifiedOriginal,
119  targetNucleus, currentParticle, targetParticle,
120  incidentHasChanged, targetHasChanged, quasiElastic );
121 
122  SetUpChange( vec, vecLen,
123  currentParticle, targetParticle,
124  incidentHasChanged );
125 
126  delete originalTarget;
127 
128  return &theParticleChange;
129 }
130 
131 
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 CASKP
143  //
144  // K+ undergoes interaction with nucleon within a nucleus. Check if it is
145  // energetically possible to produce pions/kaons. In not, assume nuclear excitation
146  // occurs and input particle is degraded in energy. No other particles are produced.
147  // If reaction is possible, find the correct number of pions/protons/neutrons
148  // produced using an interpolation to multiplicity data. Replace some pions or
149  // protons/neutrons by kaons or strange baryons according to the average
150  // multiplicity per Inelastic reaction.
151  //
152  const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass();
153  const G4double etOriginal = originalIncident->GetTotalEnergy();
154  const G4double targetMass = targetParticle.GetMass();
155  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
156  targetMass*targetMass +
157  2.0*targetMass*etOriginal );
158  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
159  if( availableEnergy < G4PionPlus::PionPlus()->GetPDGMass() )
160  {
161  quasiElastic = true;
162  return;
163  }
164  static G4ThreadLocal G4bool first = true;
165  const G4int numMul = 1200;
166  const G4int numSec = 60;
167  static G4ThreadLocal G4double protmul[numMul], protnorm[numSec]; // proton constants
168  static G4ThreadLocal G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
169 
170  // np = number of pi+, nneg = number of pi-, nz = number of pi0
171 
172  G4int nt=0, np=0, nneg=0, nz=0;
173  const G4double c = 1.25;
174  const G4double b[] = { 0.70, 0.70 };
175  if( first ) // compute normalization constants, this will only be Done once
176  {
177  first = false;
178  G4int i;
179  for( i=0; i<numMul; ++i )protmul[i] = 0.0;
180  for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
181  G4int counter = -1;
182  for( np=0; np<(numSec/3); ++np )
183  {
184  for( nneg=std::max(0,np-2); nneg<=np; ++nneg )
185  {
186  for( nz=0; nz<numSec/3; ++nz )
187  {
188  if( ++counter < numMul )
189  {
190  nt = np+nneg+nz;
191  if( nt > 0 )
192  {
193  protmul[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
194  protnorm[nt-1] += protmul[counter];
195  }
196  }
197  }
198  }
199  }
200  for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
201  for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
202  counter = -1;
203  for( np=0; np<numSec/3; ++np )
204  {
205  for( nneg=std::max(0,np-1); nneg<=(np+1); ++nneg )
206  {
207  for( nz=0; nz<numSec/3; ++nz )
208  {
209  if( ++counter < numMul )
210  {
211  nt = np+nneg+nz;
212  if( (nt>0) && (nt<=numSec) )
213  {
214  neutmul[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
215  neutnorm[nt-1] += neutmul[counter];
216  }
217  }
218  }
219  }
220  }
221  for( i=0; i<numSec; ++i )
222  {
223  if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
224  if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
225  }
226  } // end of initialization
227 
228  const G4double expxu = 82.; // upper bound for arg. of exp
229  const G4double expxl = -expxu; // lower bound for arg. of exp
234  G4int ieab = static_cast<G4int>(availableEnergy*5.0/GeV);
235  const G4double supp[] = {0.,0.4,0.55,0.65,0.75,0.82,0.86,0.90,0.94,0.98};
236  G4double test, w0, wp, wt, wm;
237  if( (availableEnergy < 2.0*GeV) && (G4UniformRand() >= supp[ieab]) )
238  {
239  // suppress high multiplicity events at low momentum
240  // only one pion will be produced
241 
242  nneg = np = nz = 0;
243  if( targetParticle.GetDefinition() == aProton )
244  {
245  test = G4Exp( std::min( expxu, std::max( expxl, -sqr(1.0+b[0])/(2.0*c*c) ) ) );
246  w0 = test;
247  wp = test*2.0;
248  if( G4UniformRand() < w0/(w0+wp) )
249  nz = 1;
250  else
251  np = 1;
252  }
253  else // target is a neutron
254  {
255  test = G4Exp( std::min( expxu, std::max( expxl, -sqr(1.0+b[1])/(2.0*c*c) ) ) );
256  w0 = test;
257  wp = test;
258  test = G4Exp( std::min( expxu, std::max( expxl, -sqr(-1.0+b[1])/(2.0*c*c) ) ) );
259  wm = test;
260  wt = w0+wp+wm;
261  wp += w0;
262  G4double ran = G4UniformRand();
263  if( ran < w0/wt )
264  nz = 1;
265  else if( ran < wp/wt )
266  np = 1;
267  else
268  nneg = 1;
269  }
270  }
271  else
272  {
273  G4double n, anpn;
274  GetNormalizationConstant( availableEnergy, n, anpn );
275  G4double ran = G4UniformRand();
276  G4double dum, excs = 0.0;
277  if( targetParticle.GetDefinition() == aProton )
278  {
279  G4int counter = -1;
280  for( np=0; (np<numSec/3) && (ran>=excs); ++np )
281  {
282  for( nneg=std::max(0,np-2); (nneg<=np) && (ran>=excs); ++nneg )
283  {
284  for( nz=0; (nz<numSec/3) && (ran>=excs); ++nz )
285  {
286  if( ++counter < numMul )
287  {
288  nt = np+nneg+nz;
289  if( nt > 0 )
290  {
291  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
292  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
293  if( std::fabs(dum) < 1.0 )
294  {
295  if( test >= 1.0e-10 )excs += dum*test;
296  }
297  else
298  excs += dum*test;
299  }
300  }
301  }
302  }
303  }
304  if( ran >= excs )return; // 3 previous loops continued to the end
305  np--; nneg--; nz--;
306  }
307  else // target must be a neutron
308  {
309  G4int counter = -1;
310  for( np=0; (np<numSec/3) && (ran>=excs); ++np )
311  {
312  for( nneg=std::max(0,np-1); (nneg<=(np+1)) && (ran>=excs); ++nneg )
313  {
314  for( nz=0; (nz<numSec/3) && (ran>=excs); ++nz )
315  {
316  if( ++counter < numMul )
317  {
318  nt = np+nneg+nz;
319  if( (nt>=1) && (nt<=numSec) )
320  {
321  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
322  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
323  if( std::fabs(dum) < 1.0 )
324  {
325  if( test >= 1.0e-10 )excs += dum*test;
326  }
327  else
328  excs += dum*test;
329  }
330  }
331  }
332  }
333  }
334  if( ran >= excs )return; // 3 previous loops continued to the end
335  np--; nneg--; nz--;
336  }
337  }
338 
339  if( targetParticle.GetDefinition() == aProton )
340  {
341  switch( np-nneg )
342  {
343  case 1:
344  if( G4UniformRand() < 0.5 )
345  {
346  if( G4UniformRand() < 0.5 )
347  currentParticle.SetDefinitionAndUpdateE( aKaonZS );
348  else
349  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
350  incidentHasChanged = true;
351  }
352  else
353  {
354  targetParticle.SetDefinitionAndUpdateE( aNeutron );
355  targetHasChanged = true;
356  }
357  break;
358  case 2:
359  if( G4UniformRand() < 0.5 )
360  currentParticle.SetDefinitionAndUpdateE( aKaonZS );
361  else
362  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
363  incidentHasChanged = true;
364  targetParticle.SetDefinitionAndUpdateE( aNeutron );
365  incidentHasChanged = true;
366  targetHasChanged = true;
367  break;
368  default:
369  break;
370  }
371  }
372  else // target is a neutron
373  {
374  switch( np-nneg )
375  {
376  case 0:
377  if( G4UniformRand() < 0.25 )
378  {
379  if( G4UniformRand() < 0.5 )
380  currentParticle.SetDefinitionAndUpdateE( aKaonZS );
381  else
382  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
383  targetParticle.SetDefinitionAndUpdateE( aProton );
384  incidentHasChanged = true;
385  targetHasChanged = true;
386  }
387  break;
388  case 1:
389  if( G4UniformRand() < 0.5 )
390  currentParticle.SetDefinitionAndUpdateE( aKaonZS );
391  else
392  currentParticle.SetDefinitionAndUpdateE( aKaonZL );
393  incidentHasChanged = true;
394  break;
395  default: // assumes nneg = np+1 so charge is conserved
396  targetParticle.SetDefinitionAndUpdateE( aProton );
397  targetHasChanged = true;
398  break;
399  }
400  }
401 
402  SetUpPions(np, nneg, nz, vec, vecLen);
403  return;
404 }
405 
406  /* end of file */
407 
void Initialize(G4int items)
Definition: G4FastVector.hh:63
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void SetUpChange(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
G4double ek
static constexpr double MeV
Definition: G4SIunits.hh:214
G4double Cinema(G4double kineticEnergy)
Definition: G4Nucleus.cc:382
void SetMomentumChange(const G4ThreeVector &aV)
void SetMomentum(const G4double x, const G4double y, const G4double z)
#define G4endl
Definition: G4ios.hh:61
const char * p
Definition: xmltok.h:285
void SetDefinitionAndUpdateE(const G4ParticleDefinition *aParticleDefinition)
const G4String & GetParticleName() const
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
void GetNormalizationConstant(const G4double availableEnergy, G4double &n, G4double &anpn)
const G4ParticleDefinition * GetDefinition() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double GetTotalEnergy() const
#define G4ThreadLocal
Definition: tls.hh:69
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)
G4double GetPDGMass() const
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:241
static unsigned wp
Definition: csz_inflate.cc:294
const G4String & GetName() const
Definition: G4Material.hh:179
void SetEnergyChange(G4double anEnergy)
G4double Pmltpc(G4int np, G4int nm, G4int nz, G4int n, G4double b, G4double c)
static G4KaonZeroLong * KaonZeroLong()
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4ParticleDefinition * GetDefinition() const
const G4Material * GetMaterial() const
#define G4UniformRand()
Definition: Randomize.hh:53
void SetUpPions(const G4int np, const G4int nm, const G4int nz, G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen)
void SetSide(const G4int sid)
const G4LorentzVector & Get4Momentum() const
Hep3Vector unit() const
G4double GetKineticEnergy() const
static G4KaonZeroShort * KaonZeroShort()
G4double GetMass() const
const G4ParticleDefinition * GetDefinition() const
void Cascade(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool &quasiElastic)
int G4int
Definition: G4Types.hh:78
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
G4GLOB_DLL std::ostream G4cout
Hep3Vector vect() const
Char_t n[5]
T sqr(const T &x)
Definition: templates.hh:145
static constexpr double pi
Definition: G4SIunits.hh:75
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:278
static constexpr double GeV
Definition: G4SIunits.hh:217
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
void SetStatusChange(G4HadFinalStateStatus aS)