46 G4cout <<
"WARNING: model G4LEPionPlusInelastic is being deprecated and will\n"
47 <<
"disappear in Geant4 version 10.0" <<
G4endl;
53 outFile <<
"G4LEPionPlusInelastic is one of the Low Energy Parameterized\n"
54 <<
"(LEP) models used to implement inelastic pi+ 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 pions with initial energies between 0 and 25\n"
84 G4cout <<
"G4LEPionPlusInelastic::ApplyYourself called" <<
G4endl;
86 G4cout <<
"target material = " << targetMaterial->
GetName() <<
", ";
104 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
116 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
127 G4bool incidentHasChanged =
false;
128 G4bool targetHasChanged =
false;
129 G4bool quasiElastic =
false;
136 Cascade(vec, vecLen, originalIncident, currentParticle, targetParticle,
137 incidentHasChanged, targetHasChanged, quasiElastic);
140 originalIncident, originalTarget, modifiedOriginal,
141 targetNucleus, currentParticle, targetParticle,
142 incidentHasChanged, targetHasChanged, quasiElastic);
144 SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged);
148 delete originalTarget;
153void G4LEPionPlusInelastic::Cascade(
159 G4bool& incidentHasChanged,
178 G4double centerofmassEnergy = std::sqrt(mOriginal*mOriginal +
179 targetMass*targetMass +
180 2.0*targetMass*etOriginal);
181 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];
189 G4int counter, nt=0, npos=0, nneg=0, nzero=0;
191 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 ) {
199 for( nneg=std::max(0,npos-2); nneg<=npos; ++nneg ) {
200 for( nzero=0; nzero<numSec/3; ++nzero ) {
201 if( ++counter < numMul ) {
202 nt = npos+nneg+nzero;
204 protmul[counter] =
Pmltpc(npos,nneg,nzero,nt,b[0],c);
205 protnorm[nt-1] += protmul[counter];
211 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
212 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
214 for( npos=0; npos<numSec/3; ++npos ) {
215 for( nneg=std::max(0,npos-1); nneg<=(npos+1); ++nneg ) {
216 for( nzero=0; nzero<numSec/3; ++nzero ) {
217 if( ++counter < numMul ) {
218 nt = npos+nneg+nzero;
219 if( (nt>0) && (nt<=numSec) ) {
220 neutmul[counter] =
Pmltpc(npos,nneg,nzero,nt,b[1],c);
221 neutnorm[nt-1] += neutmul[counter];
227 for( i=0; i<numSec; ++i ) {
228 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
229 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
238 G4int ieab =
static_cast<G4int>(availableEnergy*5.0/GeV);
239 const G4double supp[] = {0.,0.2,0.45,0.55,0.65,0.75,0.85,0.90,0.94,0.98};
241 if( (availableEnergy < 2.0*GeV) && (
G4UniformRand() >= supp[ieab]) )
247 const G4double cech[] = {1.,0.95,0.79,0.32,0.19,0.16,0.14,0.12,0.10,0.08};
248 G4int iplab =
G4int(std::min( 9.0, pOriginal/GeV*5.0 ));
255 incidentHasChanged =
true;
256 targetHasChanged =
true;
266 nneg = npos = nzero = 0;
268 test = std::exp( std::min( expxu, std::max( expxl, -
sqr(1.0+b[0])/(2.0*c*c) ) ) );
276 test = std::exp( std::min( expxu, std::max( expxl, -
sqr(1.0+b[1])/(2.0*c*c) ) ) );
279 test = std::exp( std::min( expxu, std::max( expxl, -
sqr(-1.0+b[1])/(2.0*c*c) ) ) );
286 else if( ran < wp/wt )
303 for( npos=0; (npos<numSec/3) && (ran>=excs); ++npos ) {
304 for( nneg=std::max(0,npos-2); (nneg<=npos) && (ran>=excs); ++nneg ) {
305 for( nzero=0; (nzero<numSec/3) && (ran>=excs); ++nzero ) {
306 if( ++counter < numMul ) {
307 nt = npos+nneg+nzero;
309 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
310 dum = (
pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
311 if( std::fabs(dum) < 1.0 ) {
312 if( test >= 1.0e-10 )excs += dum*test;
326 npos--; nneg--; nzero--;
329 for( npos=0; (npos<numSec/3) && (ran>=excs); ++npos ) {
330 for( nneg=std::max(0,npos-1); (nneg<=(npos+1)) && (ran>=excs); ++nneg ) {
331 for( nzero=0; (nzero<numSec/3) && (ran>=excs); ++nzero ) {
332 if( ++counter < numMul ) {
333 nt = npos+nneg+nzero;
334 if( (nt>=1) && (nt<=numSec) ) {
335 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
336 dum = (
pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
337 if( std::fabs(dum) < 1.0 ) {
338 if( test >= 1.0e-10 )excs += dum*test;
352 npos--; nneg--; nzero--;
356 switch( npos-nneg ) {
360 incidentHasChanged =
true;
363 targetHasChanged =
true;
369 incidentHasChanged =
true;
370 targetHasChanged =
true;
376 switch( npos-nneg ) {
381 incidentHasChanged =
true;
382 targetHasChanged =
true;
387 incidentHasChanged =
true;
391 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)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
G4LEPionPlusInelastic(const G4String &name="G4LEPionPlusInelastic")
virtual void ModelDescription(std::ostream &outFile) const
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 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