Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LEKaonPlusInelastic.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// Hadronic Process: Low Energy KaonPlus Inelastic Process
29// J.L. Chuma, TRIUMF, 05-Feb-1997
30// Modified by J.L.Chuma 30-Apr-97: added originalTarget for CalculateMomenta
31
32#include <iostream>
33
36#include "G4SystemOfUnits.hh"
37#include "Randomize.hh"
38
41{
42 SetMinEnergy(0.0);
43 SetMaxEnergy(25.*GeV);
44 G4cout << "WARNING: model G4LEKaonPlusInelastic is being deprecated and will\n"
45 << "disappear in Geant4 version 10.0" << G4endl;
46}
47
48
49void G4LEKaonPlusInelastic::ModelDescription(std::ostream& outFile) const
50{
51 outFile << "G4LEKaonPlusInelastic is one of the Low Energy Parameterized\n"
52 << "(LEP) models used to implement inelastic K+ scattering\n"
53 << "from nuclei. It is a re-engineered version of the GHEISHA\n"
54 << "code of H. Fesefeldt. It divides the initial collision\n"
55 << "products into backward- and forward-going clusters which are\n"
56 << "then decayed into final state hadrons. The model does not\n"
57 << "conserve energy on an event-by-event basis. It may be\n"
58 << "applied to kaons with initial energies between 0 and 25\n"
59 << "GeV.\n";
60}
61
62
65 G4Nucleus& targetNucleus)
66{
67 const G4HadProjectile *originalIncident = &aTrack;
68 if (originalIncident->GetKineticEnergy()<= 0.1*MeV) {
72 return &theParticleChange;
73 }
74
75 // create the target particle
76 G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
77 G4ReactionProduct targetParticle( originalTarget->GetDefinition() );
78
79 if (verboseLevel > 1) {
80 const G4Material *targetMaterial = aTrack.GetMaterial();
81 G4cout << "G4LEKaonPlusInelastic::ApplyYourself called" << G4endl;
82 G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy() << "MeV, ";
83 G4cout << "target material = " << targetMaterial->GetName() << ", ";
84 G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
85 << G4endl;
86 }
87 G4ReactionProduct currentParticle( const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition()));
88 currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
89 currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
90
91 // Fermi motion and evaporation
92 // As of Geant3, the Fermi energy calculation had not been done
93 G4double ek = originalIncident->GetKineticEnergy();
94 G4double amas = originalIncident->GetDefinition()->GetPDGMass();
95
96 G4double tkin = targetNucleus.Cinema(ek);
97 ek += tkin;
98 currentParticle.SetKineticEnergy( ek );
99 G4double et = ek + amas;
100 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
101 G4double pp = currentParticle.GetMomentum().mag();
102 if (pp > 0.0) {
103 G4ThreeVector momentum = currentParticle.GetMomentum();
104 currentParticle.SetMomentum( momentum * (p/pp) );
105 }
106
107 // calculate black track energies
108 tkin = targetNucleus.EvaporationEffects( ek );
109 ek -= tkin;
110 currentParticle.SetKineticEnergy( ek );
111 et = ek + amas;
112 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
113 pp = currentParticle.GetMomentum().mag();
114 if (pp > 0.0) {
115 G4ThreeVector momentum = currentParticle.GetMomentum();
116 currentParticle.SetMomentum( momentum * (p/pp) );
117 }
118
119 G4ReactionProduct modifiedOriginal = currentParticle;
120
121 currentParticle.SetSide(1); // incident always goes in forward hemisphere
122 targetParticle.SetSide(-1); // target always goes in backward hemisphere
123 G4bool incidentHasChanged = false;
124 G4bool targetHasChanged = false;
125 G4bool quasiElastic = false;
126 G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
127 G4int vecLen = 0;
128 vec.Initialize( 0 );
129
130 const G4double cutOff = 0.1*MeV;
131 if (currentParticle.GetKineticEnergy() > cutOff)
132 Cascade(vec, vecLen, originalIncident, currentParticle,
133 targetParticle, incidentHasChanged, targetHasChanged,
134 quasiElastic);
135
136 CalculateMomenta(vec, vecLen, originalIncident, originalTarget,
137 modifiedOriginal, targetNucleus, currentParticle,
138 targetParticle, incidentHasChanged, targetHasChanged,
139 quasiElastic);
140
141 SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged);
142
143 if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
144
145 delete originalTarget;
146 return &theParticleChange;
147}
148
149
150void G4LEKaonPlusInelastic::Cascade(
152 G4int &vecLen,
153 const G4HadProjectile *originalIncident,
154 G4ReactionProduct &currentParticle,
155 G4ReactionProduct &targetParticle,
156 G4bool &incidentHasChanged,
157 G4bool &targetHasChanged,
158 G4bool &quasiElastic)
159{
160 // derived from original FORTRAN code CASKP by H. Fesefeldt (13-Sep-1987)
161 //
162 // K+ undergoes interaction with nucleon within a nucleus. Check if it is
163 // energetically possible to produce pions/kaons. In not, assume nuclear excitation
164 // occurs and input particle is degraded in energy. No other particles are produced.
165 // If reaction is possible, find the correct number of pions/protons/neutrons
166 // produced using an interpolation to multiplicity data. Replace some pions or
167 // protons/neutrons by kaons or strange baryons according to the average
168 // multiplicity per Inelastic reaction.
169 //
170 const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass();
171 const G4double etOriginal = originalIncident->GetTotalEnergy();
172 const G4double targetMass = targetParticle.GetMass();
173 G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
174 targetMass*targetMass +
175 2.0*targetMass*etOriginal );
176 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
177 if( availableEnergy < G4PionPlus::PionPlus()->GetPDGMass() )
178 {
179 quasiElastic = true;
180 return;
181 }
182 static G4bool first = true;
183 const G4int numMul = 1200;
184 const G4int numSec = 60;
185 static G4double protmul[numMul], protnorm[numSec]; // proton constants
186 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
187 // npos = number of pi+, nneg = number of pi-, nzero = number of pi0
188 G4int nt=0, npos=0, nneg=0, nzero=0;
189 const G4double c = 1.25;
190 const G4double b[] = { 0.70, 0.70 };
191 if( first ) // compute normalization constants, this will only be Done once
192 {
193 first = false;
194 G4int i;
195 for( i=0; i<numMul; ++i )protmul[i] = 0.0;
196 for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
197 G4int counter = -1;
198 for( npos=0; npos<(numSec/3); ++npos )
199 {
200 for( nneg=std::max(0,npos-2); nneg<=npos; ++nneg )
201 {
202 for( nzero=0; nzero<numSec/3; ++nzero )
203 {
204 if( ++counter < numMul )
205 {
206 nt = npos+nneg+nzero;
207 if( nt > 0 )
208 {
209 protmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[0],c);
210 protnorm[nt-1] += protmul[counter];
211 }
212 }
213 }
214 }
215 }
216 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
217 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
218 counter = -1;
219 for( npos=0; npos<numSec/3; ++npos )
220 {
221 for( nneg=std::max(0,npos-1); nneg<=(npos+1); ++nneg )
222 {
223 for( nzero=0; nzero<numSec/3; ++nzero )
224 {
225 if( ++counter < numMul )
226 {
227 nt = npos+nneg+nzero;
228 if( (nt>0) && (nt<=numSec) )
229 {
230 neutmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[1],c);
231 neutnorm[nt-1] += neutmul[counter];
232 }
233 }
234 }
235 }
236 }
237 for( i=0; i<numSec; ++i )
238 {
239 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
240 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
241 }
242 } // end of initialization
243
244 const G4double expxu = 82.; // upper bound for arg. of exp
245 const G4double expxl = -expxu; // lower bound for arg. of exp
250 G4int ieab = static_cast<G4int>(availableEnergy*5.0/GeV);
251 const G4double supp[] = {0.,0.4,0.55,0.65,0.75,0.82,0.86,0.90,0.94,0.98};
252 G4double test, w0, wp, wt, wm;
253 if( (availableEnergy < 2.0*GeV) && (G4UniformRand() >= supp[ieab]) )
254 {
255 // suppress high multiplicity events at low momentum
256 // only one pion will be produced
257
258 nneg = npos = nzero = 0;
259 if( targetParticle.GetDefinition() == aProton )
260 {
261 test = std::exp( std::min( expxu, std::max( expxl, -sqr(1.0+b[0])/(2.0*c*c) ) ) );
262 w0 = test;
263 wp = test*2.0;
264 if( G4UniformRand() < w0/(w0+wp) )
265 nzero = 1;
266 else
267 npos = 1;
268 }
269 else // target is a neutron
270 {
271 test = std::exp( std::min( expxu, std::max( expxl, -sqr(1.0+b[1])/(2.0*c*c) ) ) );
272 w0 = test;
273 wp = test;
274 test = std::exp( std::min( expxu, std::max( expxl, -sqr(-1.0+b[1])/(2.0*c*c) ) ) );
275 wm = test;
276 wt = w0+wp+wm;
277 wp += w0;
278 G4double ran = G4UniformRand();
279 if( ran < w0/wt )
280 nzero = 1;
281 else if( ran < wp/wt )
282 npos = 1;
283 else
284 nneg = 1;
285 }
286 }
287 else
288 {
289 G4double n, anpn;
290 GetNormalizationConstant( availableEnergy, n, anpn );
291 G4double ran = G4UniformRand();
292 G4double dum, excs = 0.0;
293 if( targetParticle.GetDefinition() == aProton )
294 {
295 G4int counter = -1;
296 for( npos=0; (npos<numSec/3) && (ran>=excs); ++npos )
297 {
298 for( nneg=std::max(0,npos-2); (nneg<=npos) && (ran>=excs); ++nneg )
299 {
300 for( nzero=0; (nzero<numSec/3) && (ran>=excs); ++nzero )
301 {
302 if( ++counter < numMul )
303 {
304 nt = npos+nneg+nzero;
305 if( nt > 0 )
306 {
307 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
308 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
309 if( std::fabs(dum) < 1.0 )
310 {
311 if( test >= 1.0e-10 )excs += dum*test;
312 }
313 else
314 excs += dum*test;
315 }
316 }
317 }
318 }
319 }
320 if( ran >= excs )return; // 3 previous loops continued to the end
321 npos--; nneg--; nzero--;
322 }
323 else // target must be a neutron
324 {
325 G4int counter = -1;
326 for( npos=0; (npos<numSec/3) && (ran>=excs); ++npos )
327 {
328 for( nneg=std::max(0,npos-1); (nneg<=(npos+1)) && (ran>=excs); ++nneg )
329 {
330 for( nzero=0; (nzero<numSec/3) && (ran>=excs); ++nzero )
331 {
332 if( ++counter < numMul )
333 {
334 nt = npos+nneg+nzero;
335 if( (nt>=1) && (nt<=numSec) )
336 {
337 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
338 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
339 if( std::fabs(dum) < 1.0 )
340 {
341 if( test >= 1.0e-10 )excs += dum*test;
342 }
343 else
344 excs += dum*test;
345 }
346 }
347 }
348 }
349 }
350 if( ran >= excs )return; // 3 previous loops continued to the end
351 npos--; nneg--; nzero--;
352 }
353 }
354 if( targetParticle.GetDefinition() == aProton )
355 {
356 switch( npos-nneg )
357 {
358 case 1:
359 if( G4UniformRand() < 0.5 )
360 {
361 if( G4UniformRand() < 0.5 )
362 currentParticle.SetDefinitionAndUpdateE( aKaonZS );
363 else
364 currentParticle.SetDefinitionAndUpdateE( aKaonZL );
365 incidentHasChanged = true;
366 }
367 else
368 {
369 targetParticle.SetDefinitionAndUpdateE( aNeutron );
370 targetHasChanged = true;
371 }
372 break;
373 case 2:
374 if( G4UniformRand() < 0.5 )
375 currentParticle.SetDefinitionAndUpdateE( aKaonZS );
376 else
377 currentParticle.SetDefinitionAndUpdateE( aKaonZL );
378 incidentHasChanged = true;
379 targetParticle.SetDefinitionAndUpdateE( aNeutron );
380 incidentHasChanged = true;
381 targetHasChanged = true;
382 break;
383 default:
384 break;
385 }
386 }
387 else // target is a neutron
388 {
389 switch( npos-nneg )
390 {
391 case 0:
392 if( G4UniformRand() < 0.25 )
393 {
394 if( G4UniformRand() < 0.5 )
395 currentParticle.SetDefinitionAndUpdateE( aKaonZS );
396 else
397 currentParticle.SetDefinitionAndUpdateE( aKaonZL );
398 targetParticle.SetDefinitionAndUpdateE( aProton );
399 incidentHasChanged = true;
400 targetHasChanged = true;
401 }
402 break;
403 case 1:
404 if( G4UniformRand() < 0.5 )
405 currentParticle.SetDefinitionAndUpdateE( aKaonZS );
406 else
407 currentParticle.SetDefinitionAndUpdateE( aKaonZL );
408 incidentHasChanged = true;
409 break;
410 default: // assumes nneg = npos+1 so charge is conserved
411 targetParticle.SetDefinitionAndUpdateE( aProton );
412 targetHasChanged = true;
413 break;
414 }
415 }
416 SetUpPions( npos, nneg, nzero, vec, vecLen );
417 return;
418 }
419
420 /* end of file */
421
@ isAlive
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
double mag() const
Hep3Vector vect() const
G4ParticleDefinition * GetDefinition() const
void Initialize(G4int items)
Definition: G4FastVector.hh:63
void SetStatusChange(G4HadFinalStateStatus aS)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetTotalEnergy() const
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
G4double Pmltpc(G4int np, G4int nm, G4int nz, G4int n, G4double b, G4double c)
void GetNormalizationConstant(const G4double availableEnergy, G4double &n, G4double &anpn)
void SetUpPions(const G4int np, const G4int nm, const G4int nz, G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen)
void CalculateMomenta(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, const G4DynamicParticle *originalTarget, G4ReactionProduct &modifiedOriginal, G4Nucleus &targetNucleus, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool quasiElastic)
void DoIsotopeCounting(const G4HadProjectile *theProjectile, const G4Nucleus &aNucleus)
void SetUpChange(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
static G4KaonZeroLong * KaonZeroLong()
static G4KaonZeroShort * KaonZeroShort()
virtual void ModelDescription(std::ostream &outFile) const
G4LEKaonPlusInelastic(const G4String &name="G4LEKaonPlusInelastic")
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
const G4String & GetName() const
Definition: G4Material.hh:177
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:264
G4double Cinema(G4double kineticEnergy)
Definition: G4Nucleus.cc:368
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:227
const G4String & GetParticleName() const
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void SetSide(const G4int sid)
void SetDefinitionAndUpdateE(G4ParticleDefinition *aParticleDefinition)
void SetKineticEnergy(const G4double en)
G4ParticleDefinition * GetDefinition() const
G4double GetMass() const
const G4double pi
T sqr(const T &x)
Definition: templates.hh:145