44 G4cout <<
"WARNING: model G4LEKaonPlusInelastic is being deprecated and will\n"
45 <<
"disappear in Geant4 version 10.0" <<
G4endl;
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"
81 G4cout <<
"G4LEKaonPlusInelastic::ApplyYourself called" <<
G4endl;
83 G4cout <<
"target material = " << targetMaterial->
GetName() <<
", ";
100 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
112 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
123 G4bool incidentHasChanged =
false;
124 G4bool targetHasChanged =
false;
125 G4bool quasiElastic =
false;
132 Cascade(vec, vecLen, originalIncident, currentParticle,
133 targetParticle, incidentHasChanged, targetHasChanged,
137 modifiedOriginal, targetNucleus, currentParticle,
138 targetParticle, incidentHasChanged, targetHasChanged,
141 SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged);
145 delete originalTarget;
150void G4LEKaonPlusInelastic::Cascade(
156 G4bool &incidentHasChanged,
173 G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
174 targetMass*targetMass +
175 2.0*targetMass*etOriginal );
176 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
182 static G4bool first =
true;
183 const G4int numMul = 1200;
184 const G4int numSec = 60;
185 static G4double protmul[numMul], protnorm[numSec];
186 static G4double neutmul[numMul], neutnorm[numSec];
188 G4int nt=0, npos=0, nneg=0, nzero=0;
190 const G4double b[] = { 0.70, 0.70 };
195 for( i=0; i<numMul; ++i )protmul[i] = 0.0;
196 for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
198 for( npos=0; npos<(numSec/3); ++npos )
200 for( nneg=std::max(0,npos-2); nneg<=npos; ++nneg )
202 for( nzero=0; nzero<numSec/3; ++nzero )
204 if( ++counter < numMul )
206 nt = npos+nneg+nzero;
209 protmul[counter] =
Pmltpc(npos,nneg,nzero,nt,b[0],c);
210 protnorm[nt-1] += protmul[counter];
216 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
217 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
219 for( npos=0; npos<numSec/3; ++npos )
221 for( nneg=std::max(0,npos-1); nneg<=(npos+1); ++nneg )
223 for( nzero=0; nzero<numSec/3; ++nzero )
225 if( ++counter < numMul )
227 nt = npos+nneg+nzero;
228 if( (nt>0) && (nt<=numSec) )
230 neutmul[counter] =
Pmltpc(npos,nneg,nzero,nt,b[1],c);
231 neutnorm[nt-1] += neutmul[counter];
237 for( i=0; i<numSec; ++i )
239 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
240 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
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};
253 if( (availableEnergy < 2.0*GeV) && (
G4UniformRand() >= supp[ieab]) )
258 nneg = npos = nzero = 0;
261 test = std::exp( std::min( expxu, std::max( expxl, -
sqr(1.0+b[0])/(2.0*c*c) ) ) );
271 test = std::exp( std::min( expxu, std::max( expxl, -
sqr(1.0+b[1])/(2.0*c*c) ) ) );
274 test = std::exp( std::min( expxu, std::max( expxl, -
sqr(-1.0+b[1])/(2.0*c*c) ) ) );
281 else if( ran < wp/wt )
296 for( npos=0; (npos<numSec/3) && (ran>=excs); ++npos )
298 for( nneg=std::max(0,npos-2); (nneg<=npos) && (ran>=excs); ++nneg )
300 for( nzero=0; (nzero<numSec/3) && (ran>=excs); ++nzero )
302 if( ++counter < numMul )
304 nt = npos+nneg+nzero;
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 )
311 if( test >= 1.0e-10 )excs += dum*test;
320 if( ran >= excs )
return;
321 npos--; nneg--; nzero--;
326 for( npos=0; (npos<numSec/3) && (ran>=excs); ++npos )
328 for( nneg=std::max(0,npos-1); (nneg<=(npos+1)) && (ran>=excs); ++nneg )
330 for( nzero=0; (nzero<numSec/3) && (ran>=excs); ++nzero )
332 if( ++counter < numMul )
334 nt = npos+nneg+nzero;
335 if( (nt>=1) && (nt<=numSec) )
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 )
341 if( test >= 1.0e-10 )excs += dum*test;
350 if( ran >= excs )
return;
351 npos--; nneg--; nzero--;
365 incidentHasChanged =
true;
370 targetHasChanged =
true;
378 incidentHasChanged =
true;
380 incidentHasChanged =
true;
381 targetHasChanged =
true;
399 incidentHasChanged =
true;
400 targetHasChanged =
true;
408 incidentHasChanged =
true;
412 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
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 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
static G4Neutron * Neutron()
G4double EvaporationEffects(G4double kineticEnergy)
G4double Cinema(G4double kineticEnergy)
G4DynamicParticle * ReturnTargetParticle() const
G4double GetPDGMass() const
const G4String & GetParticleName() const
static G4PionPlus * PionPlus()
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