Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4RPGSigmaMinusInelastic.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: G4RPGSigmaMinusInelastic.cc 94214 2015-11-09 08:18:05Z gcosmo $
27 //
28 
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  {
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 << "G4RPGSigmaMinusInelastic::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( originalIncident->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 
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 CASSM
138  //
139  // SigmaMinus 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  {
156  quasiElastic = true;
157  return;
158  }
159  static G4ThreadLocal G4bool first = true;
160  const G4int numMul = 1200;
161  const G4int numSec = 60;
162  static G4ThreadLocal G4double protmul[numMul], protnorm[numSec]; // proton constants
163  static G4ThreadLocal G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
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.70, 0.70 };
169  if( first ) // compute normalization constants, this will only be Done once
170  {
171  first = false;
172  G4int i;
173  for( i=0; i<numMul; ++i )protmul[i] = 0.0;
174  for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
175  counter = -1;
176  for( np=0; np<(numSec/3); ++np )
177  {
178  for( nneg=std::max(0,np-1); nneg<=(np+1); ++nneg )
179  {
180  for( nz=0; nz<numSec/3; ++nz )
181  {
182  if( ++counter < numMul )
183  {
184  nt = np+nneg+nz;
185  if( nt>0 && nt<=numSec )
186  {
187  protmul[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
188  protnorm[nt-1] += protmul[counter];
189  }
190  }
191  }
192  }
193  }
194  for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
195  for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
196  counter = -1;
197  for( np=0; np<numSec/3; ++np )
198  {
199  for( nneg=np; nneg<=(np+2); ++nneg )
200  {
201  for( nz=0; nz<numSec/3; ++nz )
202  {
203  if( ++counter < numMul )
204  {
205  nt = np+nneg+nz;
206  if( nt>0 && nt<=numSec )
207  {
208  neutmul[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
209  neutnorm[nt-1] += neutmul[counter];
210  }
211  }
212  }
213  }
214  }
215  for( i=0; i<numSec; ++i )
216  {
217  if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
218  if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
219  }
220  } // end of initialization
221 
222  const G4double expxu = 82.; // upper bound for arg. of exp
223  const G4double expxl = -expxu; // lower bound for arg. of exp
228 
229  // energetically possible to produce pion(s) --> inelastic scattering
230 
231  G4double n, anpn;
232  GetNormalizationConstant( availableEnergy, n, anpn );
233  G4double ran = G4UniformRand();
234  G4double dum, excs = 0.0;
235  if( targetParticle.GetDefinition() == aProton )
236  {
237  counter = -1;
238  for( np=0; np<numSec/3 && ran>=excs; ++np )
239  {
240  for( nneg=std::max(0,np-1); nneg<=(np+1) && ran>=excs; ++nneg )
241  {
242  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
243  {
244  if( ++counter < numMul )
245  {
246  nt = np+nneg+nz;
247  if( nt>0 && nt<=numSec )
248  {
249  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
250  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
251  if( std::fabs(dum) < 1.0 )
252  {
253  if( test >= 1.0e-10 )excs += dum*test;
254  }
255  else
256  excs += dum*test;
257  }
258  }
259  }
260  }
261  }
262  if( ran >= excs ) // 3 previous loops continued to the end
263  {
264  quasiElastic = true;
265  return;
266  }
267  np--; nneg--; nz--;
268  G4int ncht = std::max( 1, np-nneg+2 );
269  switch( ncht )
270  {
271  case 1:
272  if( G4UniformRand() < 0.5 )
273  currentParticle.SetDefinitionAndUpdateE( aLambda );
274  else
275  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
276  incidentHasChanged = true;
277  break;
278  case 2:
279  if( G4UniformRand() >= 0.5 )
280  {
281  if( G4UniformRand() < 0.5 )
282  currentParticle.SetDefinitionAndUpdateE( aLambda );
283  else
284  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
285  incidentHasChanged = true;
286  targetParticle.SetDefinitionAndUpdateE( aNeutron );
287  targetHasChanged = true;
288  }
289  break;
290  default:
291  targetParticle.SetDefinitionAndUpdateE( aNeutron );
292  targetHasChanged = true;
293  break;
294  }
295  }
296  else // target must be a neutron
297  {
298  counter = -1;
299  for( np=0; np<numSec/3 && ran>=excs; ++np )
300  {
301  for( nneg=np; nneg<=(np+2) && ran>=excs; ++nneg )
302  {
303  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
304  {
305  if( ++counter < numMul )
306  {
307  nt = np+nneg+nz;
308  if( nt>0 && nt<=numSec )
309  {
310  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
311  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
312  if( std::fabs(dum) < 1.0 )
313  {
314  if( test >= 1.0e-10 )excs += dum*test;
315  }
316  else
317  excs += dum*test;
318  }
319  }
320  }
321  }
322  }
323  if( ran >= excs ) // 3 previous loops continued to the end
324  {
325  quasiElastic = true;
326  return;
327  }
328  np--; nneg--; nz--;
329  G4int ncht = std::max( 1, np-nneg+3 );
330  switch( ncht )
331  {
332  case 1:
333  if( G4UniformRand() < 0.5 )
334  currentParticle.SetDefinitionAndUpdateE( aLambda );
335  else
336  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
337  incidentHasChanged = true;
338  targetParticle.SetDefinitionAndUpdateE( aProton );
339  targetHasChanged = true;
340  break;
341  case 2:
342  if( G4UniformRand() < 0.5 )
343  {
344  if( G4UniformRand() < 0.5 )
345  currentParticle.SetDefinitionAndUpdateE( aLambda );
346  else
347  currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
348  incidentHasChanged = true;
349  }
350  else
351  {
352  targetParticle.SetDefinitionAndUpdateE( aProton );
353  targetHasChanged = true;
354  }
355  break;
356  default:
357  break;
358  }
359  }
360 
361  SetUpPions(np, nneg, nz, vec, vecLen);
362  return;
363 }
364 
365  /* end of file */
366 
void Initialize(G4int items)
Definition: G4FastVector.hh:63
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108
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
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
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 G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:102
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
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)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
void SetKineticEnergy(const G4double en)
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)
G4ThreeVector GetMomentum() const
void SetSide(const G4int sid)
const G4LorentzVector & Get4Momentum() const
Hep3Vector unit() const
G4double GetKineticEnergy() const
G4double GetMass() const
const G4ParticleDefinition * GetDefinition() const
double mag() const
int G4int
Definition: G4Types.hh:78
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4GLOB_DLL std::ostream G4cout
Hep3Vector vect() const
Char_t n[5]
static constexpr double pi
Definition: G4SIunits.hh:75
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:278
void Cascade(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool &quasiElastic)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
void SetStatusChange(G4HadFinalStateStatus aS)