43 outFile <<
"G4LEAntiXiZeroInelastic is one of the Low Energy Parameterized\n"
44 <<
"(LEP) models used to implement inelastic antiXi0 scattering\n"
45 <<
"from nuclei. It is a re-engineered version of the GHEISHA\n"
46 <<
"code of H. Fesefeldt. It divides the initial collision\n"
47 <<
"products into backward- and forward-going clusters which are\n"
48 <<
"then decayed into final state hadrons. The model does not\n"
49 <<
"conserve energy on an event-by-event basis. It may be\n"
50 <<
"applied to antiXi0 with initial energies between 0 and 25\n"
65 G4cout <<
"G4LEAntiXiZeroInelastic::ApplyYourself called" <<
G4endl;
67 G4cout <<
"target material = " << targetMaterial->
GetName() <<
", ";
77 modifiedOriginal = *originalIncident;
83 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
95 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
103 targetParticle = *originalTarget;
106 G4bool incidentHasChanged =
false;
107 G4bool targetHasChanged =
false;
108 G4bool quasiElastic =
false;
117 Cascade(vec, vecLen, originalIncident, currentParticle, targetParticle,
118 incidentHasChanged, targetHasChanged, quasiElastic);
121 modifiedOriginal, targetNucleus, currentParticle,
122 targetParticle, incidentHasChanged, targetHasChanged,
125 SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged);
129 delete originalTarget;
134void G4LEAntiXiZeroInelastic::Cascade(
140 G4bool &incidentHasChanged,
158 G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
159 targetMass*targetMass +
160 2.0*targetMass*etOriginal );
161 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
167 static G4bool first =
true;
168 const G4int numMul = 1200;
169 const G4int numSec = 60;
170 static G4double protmul[numMul], protnorm[numSec];
171 static G4double neutmul[numMul], neutnorm[numSec];
173 G4int counter, nt=0, npos=0, nneg=0, nzero=0;
181 for( i=0; i<numMul; ++i )protmul[i] = 0.0;
182 for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
184 for( npos=0; npos<(numSec/3); ++npos )
186 for( nneg=std::max(0,npos-2); nneg<=(npos+1); ++nneg )
188 for( nzero=0; nzero<numSec/3; ++nzero )
190 if( ++counter < numMul )
192 nt = npos+nneg+nzero;
193 if( nt>0 && nt<=numSec )
195 protmul[counter] =
Pmltpc(npos,nneg,nzero,nt,b[0],c);
196 protnorm[nt-1] += protmul[counter];
202 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
203 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
205 for( npos=0; npos<numSec/3; ++npos )
207 for( nneg=std::max(0,npos-1); nneg<=(npos+2); ++nneg )
209 for( nzero=0; nzero<numSec/3; ++nzero )
211 if( ++counter < numMul )
213 nt = npos+nneg+nzero;
214 if( nt>0 && nt<=numSec )
216 neutmul[counter] =
Pmltpc(npos,nneg,nzero,nt,b[1],c);
217 neutnorm[nt-1] += neutmul[counter];
223 for( i=0; i<numSec; ++i )
225 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
226 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
247 for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
249 for( nneg=std::max(0,npos-2); nneg<=(npos+1) && ran>=excs; ++nneg )
251 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
253 if( ++counter < numMul )
255 nt = npos+nneg+nzero;
256 if( nt>0 && nt<=numSec )
258 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
259 dum = (
pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
260 if( std::fabs(dum) < 1.0 )
262 if( test >= 1.0e-10 )excs += dum*test;
276 npos--; nneg--; nzero--;
288 incidentHasChanged =
true;
300 else if( npos == nneg+1 )
305 targetHasChanged =
true;
310 incidentHasChanged =
true;
316 incidentHasChanged =
true;
318 targetHasChanged =
true;
324 for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
326 for( nneg=std::max(0,npos-1); nneg<=(npos+2) && ran>=excs; ++nneg )
328 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
330 if( ++counter < numMul )
332 nt = npos+nneg+nzero;
333 if( nt>0 && 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 )
339 if( test >= 1.0e-10 )excs += dum*test;
353 npos--; nneg--; nzero--;
359 targetHasChanged =
true;
364 incidentHasChanged =
true;
366 targetHasChanged =
true;
378 else if( npos == nneg )
383 incidentHasChanged =
true;
385 targetHasChanged =
true;
391 incidentHasChanged =
true;
G4DLLIMPORT std::ostream G4cout
G4ParticleDefinition * GetDefinition() const
void SetElement(G4int anIndex, Type *anElement)
void Initialize(G4int items)
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4HadFinalState theParticleChange
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()
virtual void ModelDescription(std::ostream &outFile) const
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 GetTotalMomentum() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void SetSide(const G4int sid)
void SetDefinitionAndUpdateE(G4ParticleDefinition *aParticleDefinition)
void SetKineticEnergy(const G4double en)
G4ParticleDefinition * GetDefinition() const
void SetDefinition(G4ParticleDefinition *aParticleDefinition)
static G4SigmaPlus * SigmaPlus()
static G4XiMinus * XiMinus()