46 G4cout <<
"WARNING: model G4LEKaonMinusInelastic is being deprecated and will\n"
47 <<
"disappear in Geant4 version 10.0" <<
G4endl;
53 outFile <<
"G4LEKaonMinusInelastic is one of the Low Energy Parameterized\n"
54 <<
"(LEP) models used to implement inelastic K- scattering\n"
55 <<
"from nuclei. It is a re-engineered version of the GHEISHA\n"
56 <<
"code of H. Fesefeldt. It divides the initial collision\n"
57 <<
"products into backward- and forward-going clusters which are\n"
58 <<
"then decayed into final state hadrons. The model does not\n"
59 <<
"conserve energy on an event-by-event basis. It may be\n"
60 <<
"applied to kaons with initial energies between 0 and 25\n"
83 G4cout <<
"G4LEKaonMinusInelastic::ApplyYourself called" <<
G4endl;
85 G4cout <<
"target material = " << targetMaterial->
GetName() <<
", ";
103 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
115 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
126 G4bool incidentHasChanged =
false;
127 G4bool targetHasChanged =
false;
128 G4bool quasiElastic =
false;
135 Cascade(vec, vecLen, originalIncident, currentParticle, targetParticle,
136 incidentHasChanged, targetHasChanged, quasiElastic);
139 modifiedOriginal, targetNucleus, currentParticle,
140 targetParticle, incidentHasChanged, targetHasChanged,
143 SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged);
147 delete originalTarget;
152void G4LEKaonMinusInelastic::Cascade(
158 G4bool &incidentHasChanged,
176 G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
177 targetMass*targetMass +
178 2.0*targetMass*etOriginal );
179 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
181 static G4bool first =
true;
182 const G4int numMul = 1200;
183 const G4int numSec = 60;
184 static G4double protmul[numMul], protnorm[numSec];
185 static G4double neutmul[numMul], neutnorm[numSec];
187 G4int nt(0), npos(0), nneg(0), nzero(0);
189 const G4double b[] = { 0.70, 0.70 };
194 for( i=0; i<numMul; ++i )protmul[i] = 0.0;
195 for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
197 for( npos=0; npos<(numSec/3); ++npos )
199 for( nneg=std::max(0,npos-1); nneg<=(npos+1); ++nneg )
201 for( nzero=0; nzero<numSec/3; ++nzero )
203 if( ++counter < numMul )
205 nt = npos+nneg+nzero;
206 if( (nt>0) && (nt<=numSec) )
208 protmul[counter] =
Pmltpc(npos,nneg,nzero,nt,b[0],c);
209 protnorm[nt-1] += protmul[counter];
215 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
216 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
218 for( npos=0; npos<numSec/3; ++npos )
220 for( nneg=npos; nneg<=(npos+2); ++nneg )
222 for( nzero=0; nzero<numSec/3; ++nzero )
224 if( ++counter < numMul )
226 nt = npos+nneg+nzero;
227 if( (nt>0) && (nt<=numSec) )
229 neutmul[counter] =
Pmltpc(npos,nneg,nzero,nt,b[1],c);
230 neutnorm[nt-1] += neutmul[counter];
236 for( i=0; i<numSec; ++i )
238 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
239 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
257 const G4double cech[] = {1.,1.,1.,0.70,0.60,0.55,0.35,0.25,0.18,0.15};
258 G4int iplab =
G4int(std::min( 9.0, pOriginal/GeV*5.0 ));
259 if( (pOriginal <= 2.0*GeV) && (
G4UniformRand() < cech[iplab]) )
261 npos = nneg = nzero = nt = 0;
262 iplab =
G4int(std::min( 19.0, pOriginal/GeV*10.0 ));
263 const G4double cnk0[] = {0.17,0.18,0.17,0.24,0.26,0.20,0.22,0.21,0.34,0.45,
264 0.58,0.55,0.36,0.29,0.29,0.32,0.32,0.33,0.33,0.33};
271 incidentHasChanged =
true;
273 targetHasChanged =
true;
285 incidentHasChanged =
true;
286 targetHasChanged =
true;
289 else if( ran < 0.50 )
296 incidentHasChanged =
true;
297 targetHasChanged =
true;
299 else if( ran < 0.75 )
306 incidentHasChanged =
true;
307 targetHasChanged =
true;
316 incidentHasChanged =
true;
317 targetHasChanged =
true;
323 if( availableEnergy < aPiPlus->GetPDGMass() )
335 for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
337 for( nneg=std::max(0,npos-1); nneg<=(npos+1) && ran>=excs; ++nneg )
339 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
341 if( ++counter < numMul )
343 nt = npos+nneg+nzero;
346 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
347 dum = (
pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
348 if( std::fabs(dum) < 1.0 )
350 if( test >= 1.0e-10 )excs += dum*test;
364 npos--; nneg--; nzero--;
371 incidentHasChanged =
true;
372 targetHasChanged =
true;
375 else if( npos == nneg+1 )
378 targetHasChanged =
true;
383 incidentHasChanged =
true;
389 for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
391 for( nneg=npos; nneg<=(npos+2) && ran>=excs; ++nneg )
393 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
395 if( ++counter < numMul )
397 nt = npos+nneg+nzero;
398 if( (nt>=1) && (nt<=numSec) )
400 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
401 dum = (
pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
402 if( std::fabs(dum) < 1.0 )
404 if( test >= 1.0e-10 )excs += dum*test;
418 npos--; nneg--; nzero--;
424 targetHasChanged =
true;
429 incidentHasChanged =
true;
432 else if( npos != nneg )
435 incidentHasChanged =
true;
452 incidentHasChanged =
true;
453 targetHasChanged =
true;
459 incidentHasChanged =
true;
460 targetHasChanged =
true;
463 else if( ran < 0.84 )
469 incidentHasChanged =
true;
470 targetHasChanged =
true;
476 incidentHasChanged =
true;
477 targetHasChanged =
true;
486 incidentHasChanged =
true;
487 targetHasChanged =
true;
493 incidentHasChanged =
true;
494 targetHasChanged =
true;
506 incidentHasChanged =
true;
507 targetHasChanged =
true;
509 else if( ran < 0.78 )
513 incidentHasChanged =
true;
514 targetHasChanged =
true;
516 else if( ran < 0.89 )
520 incidentHasChanged =
true;
521 targetHasChanged =
true;
527 incidentHasChanged =
true;
528 targetHasChanged =
true;
538 incidentHasChanged =
true;
546 targetHasChanged =
true;
G4DLLIMPORT std::ostream G4cout
G4ParticleDefinition * GetDefinition() const
void Initialize(G4int items)
void SetStatusChange(G4HadFinalStateStatus aS)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4Material * GetMaterial() const
G4double GetTotalMomentum() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetTotalEnergy() const
G4HadFinalState theParticleChange
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 ¤tParticle, 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 ¤tParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
static G4KaonMinus * KaonMinus()
static G4KaonZeroLong * KaonZeroLong()
static G4KaonZeroShort * KaonZeroShort()
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
virtual void ModelDescription(std::ostream &outFile) const
G4LEKaonMinusInelastic(const G4String &name="G4LEKaonMinusInelastic")
static G4Lambda * Lambda()
const G4String & GetName() const
static G4Neutron * Neutron()
G4double EvaporationEffects(G4double kineticEnergy)
G4double Cinema(G4double kineticEnergy)
G4DynamicParticle * ReturnTargetParticle() const
G4double GetPDGMass() const
const G4String & GetParticleName() const
static G4PionMinus * PionMinus()
static G4PionPlus * PionPlus()
static G4PionZero * PionZero()
static G4Proton * Proton()
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
static G4SigmaMinus * SigmaMinus()
static G4SigmaPlus * SigmaPlus()
static G4SigmaZero * SigmaZero()