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