Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4RPGOmegaMinusInelastic.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: G4RPGOmegaMinusInelastic.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 // G4double targetMass = originalTarget->GetDefinition()->GetPDGMass();
52  G4ReactionProduct targetParticle( originalTarget->GetDefinition() );
53 
54  if( verboseLevel > 1 )
55  {
56  const G4Material *targetMaterial = aTrack.GetMaterial();
57  G4cout << "G4RPGOmegaMinusInelastic::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  return &theParticleChange;
128 }
129 
130 
133  G4int& vecLen,
134  const G4HadProjectile *originalIncident,
135  G4ReactionProduct &currentParticle,
136  G4ReactionProduct &targetParticle,
137  G4bool &incidentHasChanged,
138  G4bool &targetHasChanged,
139  G4bool &quasiElastic )
140 {
141  // Derived from H. Fesefeldt's original FORTRAN code CASOM
142  // OmegaMinus undergoes interaction with nucleon within a nucleus. Check if it is
143  // energetically possible to produce pions/kaons. In not, assume nuclear excitation
144  // occurs and input particle is degraded in energy. No other particles are produced.
145  // If reaction is possible, find the correct number of pions/protons/neutrons
146  // produced using an interpolation to multiplicity data. Replace some pions or
147  // protons/neutrons by kaons or strange baryons according to the average
148  // multiplicity per Inelastic reaction.
149 
150  const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass();
151  const G4double etOriginal = originalIncident->GetTotalEnergy();
152 // const G4double pOriginal = originalIncident->GetTotalMomentum();
153  const G4double targetMass = targetParticle.GetMass();
154  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
155  targetMass*targetMass +
156  2.0*targetMass*etOriginal );
157  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
158  if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass() )
159  {
160  quasiElastic = true;
161  return;
162  }
163  static G4ThreadLocal G4bool first = true;
164  const G4int numMul = 1200;
165  const G4int numSec = 60;
166  static G4ThreadLocal G4double protmul[numMul], protnorm[numSec]; // proton constants
167  static G4ThreadLocal G4double neutmul[numMul], neutnorm[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.70, 0.70 };
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 )
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  } // end of initialization
225 
226  const G4double expxu = 82.; // upper bound for arg. of exp
227  const G4double expxl = -expxu; // lower bound for arg. of exp
234 
235  // energetically possible to produce pion(s) --> inelastic scattering
236 
237  G4double n, anpn;
238  GetNormalizationConstant( availableEnergy, n, anpn );
239  G4double ran = G4UniformRand();
240  G4double dum, excs = 0.0;
241  if( targetParticle.GetDefinition() == aProton )
242  {
243  counter = -1;
244  for( np=0; np<numSec/3 && ran>=excs; ++np )
245  {
246  for( nneg=std::max(0,np-1); nneg<=(np+1) && ran>=excs; ++nneg )
247  {
248  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
249  {
250  if( ++counter < numMul )
251  {
252  nt = np+nneg+nz;
253  if( nt > 0 )
254  {
255  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
256  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
257  if( std::fabs(dum) < 1.0 )
258  {
259  if( test >= 1.0e-10 )excs += dum*test;
260  }
261  else
262  excs += dum*test;
263  }
264  }
265  }
266  }
267  }
268  if( ran >= excs ) // 3 previous loops continued to the end
269  {
270  quasiElastic = true;
271  return;
272  }
273  np--; nneg--; nz--;
274  }
275  else // target must be a neutron
276  {
277  counter = -1;
278  for( np=0; np<numSec/3 && ran>=excs; ++np )
279  {
280  for( nneg=np; nneg<=(np+2) && ran>=excs; ++nneg )
281  {
282  for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
283  {
284  if( ++counter < numMul )
285  {
286  nt = np+nneg+nz;
287  if( (nt>=1) && (nt<=numSec) )
288  {
289  test = G4Exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
290  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
291  if( std::fabs(dum) < 1.0 )
292  {
293  if( test >= 1.0e-10 )excs += dum*test;
294  }
295  else
296  excs += dum*test;
297  }
298  }
299  }
300  }
301  }
302  if( ran >= excs ) // 3 previous loops continued to the end
303  {
304  quasiElastic = true;
305  return;
306  }
307  np--; nneg--; nz--;
308  }
309  // number of secondary mesons determined by kno distribution
310  // check for total charge of final state mesons to determine
311  // the kind of baryons to be produced, taking into account
312  // charge and strangeness conservation
313  //
314  G4int nvefix = 0;
315  if( targetParticle.GetDefinition() == aProton )
316  {
317  if( nneg > np )
318  {
319  if( nneg == np+1 )
320  {
321  currentParticle.SetDefinitionAndUpdateE( aXiZero );
322  nvefix = 1;
323  }
324  else
325  {
326  currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
327  nvefix = 2;
328  }
329  incidentHasChanged = true;
330  }
331  else if( nneg < np )
332  {
333  targetParticle.SetDefinitionAndUpdateE( aNeutron );
334  targetHasChanged = true;
335  }
336  }
337  else // target is a neutron
338  {
339  if( np+1 < nneg )
340  {
341  if( nneg == np+2 )
342  {
343  currentParticle.SetDefinitionAndUpdateE( aXiZero );
344  incidentHasChanged = true;
345  nvefix = 1;
346  }
347  else // charge mismatch
348  {
349  currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
350  incidentHasChanged = true;
351  nvefix = 2;
352  }
353  targetParticle.SetDefinitionAndUpdateE( aProton );
354  targetHasChanged = true;
355  }
356  else if( nneg == np+1 )
357  {
358  targetParticle.SetDefinitionAndUpdateE( aProton );
359  targetHasChanged = true;
360  }
361  }
362 
363  SetUpPions(np, nneg, nz, vec, vecLen);
364  for (G4int i = 0; i < vecLen && nvefix > 0; ++i) {
365  if (vec[i]->GetDefinition() == aPiMinus) {
366  if( nvefix >= 1 )vec[i]->SetDefinitionAndUpdateE(aKaonMinus);
367  --nvefix;
368  }
369  }
370 
371  return;
372 }
373 
374  /* end of file */
375 
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
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
static G4XiZero * XiZero()
Definition: G4XiZero.cc:106
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)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
const G4ParticleDefinition * GetDefinition() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double GetTotalEnergy() const
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
#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
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
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
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)
void SetSide(const G4int sid)
const G4LorentzVector & Get4Momentum() const
Hep3Vector unit() const
G4double GetKineticEnergy() const
G4double GetMass() const
const G4ParticleDefinition * GetDefinition() 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
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
void SetStatusChange(G4HadFinalStateStatus aS)