44 G4cout <<
"WARNING: model G4LELambdaInelastic is being deprecated and will\n"
45 <<
"disappear in Geant4 version 10.0" <<
G4endl;
51 outFile <<
"G4LELambdaInelastic is one of the Low Energy Parameterized\n"
52 <<
"(LEP) models used to implement inelastic Lambda 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 lambdas with initial energies between 0 and 25\n"
75 G4cout <<
"G4LELambdaInelastic::ApplyYourself called" <<
G4endl;
77 G4cout <<
"target material = " << targetMaterial->
GetName() <<
", ";
87 modifiedOriginal = *originalIncident;
93 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
105 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
114 targetParticle = *originalTarget;
117 G4bool incidentHasChanged =
false;
118 G4bool targetHasChanged =
false;
119 G4bool quasiElastic =
false;
126 Cascade(vec, vecLen, originalIncident, currentParticle, targetParticle,
127 incidentHasChanged, targetHasChanged, quasiElastic);
130 modifiedOriginal, targetNucleus, currentParticle,
131 targetParticle, incidentHasChanged, targetHasChanged,
134 SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged);
138 delete originalTarget;
143void G4LELambdaInelastic::Cascade(
149 G4bool& incidentHasChanged,
167 G4double centerofmassEnergy = std::sqrt(mOriginal*mOriginal +
168 targetMass*targetMass +
169 2.0*targetMass*etOriginal );
170 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
175 static G4bool first =
true;
176 const G4int numMul = 1200;
177 const G4int numSec = 60;
178 static G4double protmul[numMul], protnorm[numSec];
179 static G4double neutmul[numMul], neutnorm[numSec];
182 G4int counter, nt=0, npos=0, nneg=0, nzero=0;
185 const G4double b[] = { 0.70, 0.35 };
189 for( i=0; i<numMul; ++i )protmul[i] = 0.0;
190 for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
192 for( npos=0; npos<(numSec/3); ++npos ) {
193 for( nneg=std::max(0,npos-2); nneg<=(npos+1); ++nneg ) {
194 for( nzero=0; nzero<numSec/3; ++nzero ) {
195 if( ++counter < numMul ) {
196 nt = npos+nneg+nzero;
197 if( nt>0 && nt<=numSec ) {
198 protmul[counter] =
Pmltpc(npos,nneg,nzero,nt,b[0],c);
199 protnorm[nt-1] += protmul[counter];
205 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
206 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
208 for( npos=0; npos<numSec/3; ++npos ) {
209 for( nneg=std::max(0,npos-1); nneg<=(npos+2); ++nneg ) {
210 for( nzero=0; nzero<numSec/3; ++nzero ) {
211 if( ++counter < numMul ) {
212 nt = npos+nneg+nzero;
213 if( nt>0 && nt<=numSec ) {
214 neutmul[counter] =
Pmltpc(npos,nneg,nzero,nt,b[1],c);
215 neutnorm[nt-1] += neutmul[counter];
221 for( i=0; i<numSec; ++i ) {
222 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
223 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
244 for( npos=0; npos<numSec/3 && ran>=excs; ++npos ) {
245 for( nneg=std::max(0,npos-2); nneg<=(npos+1) && ran>=excs; ++nneg ) {
246 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero ) {
247 if( ++counter < numMul ) {
248 nt = npos+nneg+nzero;
249 if( nt>0 && nt<=numSec ) {
250 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
251 dum = (
pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
252 if( std::fabs(dum) < 1.0 ) {
253 if( test >= 1.0e-10 )excs += dum*test;
267 npos--; nneg--; nzero--;
268 G4int ncht = std::max( 1, npos-nneg );
272 incidentHasChanged =
true;
278 incidentHasChanged =
true;
282 incidentHasChanged =
true;
284 targetHasChanged =
true;
291 incidentHasChanged =
true;
294 targetHasChanged =
true;
297 incidentHasChanged =
true;
302 incidentHasChanged =
true;
304 targetHasChanged =
true;
311 for( npos=0; npos<numSec/3 && ran>=excs; ++npos ) {
312 for( nneg=std::max(0,npos-1); nneg<=(npos+2) && ran>=excs; ++nneg ) {
313 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero ) {
314 if( ++counter < numMul ) {
315 nt = npos+nneg+nzero;
316 if( nt>0 && nt<=numSec ) {
317 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
318 dum = (
pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
319 if( std::fabs(dum) < 1.0 ) {
320 if( test >= 1.0e-10 )excs += dum*test;
334 npos--; nneg--; nzero--;
335 G4int ncht = std::max( 1, npos-nneg+3 );
339 incidentHasChanged =
true;
341 targetHasChanged =
true;
347 incidentHasChanged =
true;
350 targetHasChanged =
true;
353 incidentHasChanged =
true;
360 incidentHasChanged =
true;
364 incidentHasChanged =
true;
366 targetHasChanged =
true;
371 incidentHasChanged =
true;
G4DLLIMPORT std::ostream G4cout
G4ParticleDefinition * GetDefinition() const
void Initialize(G4int items)
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() 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)
G4LELambdaInelastic(const G4String &name="G4LELambdaInelastic")
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
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 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
static G4SigmaMinus * SigmaMinus()
static G4SigmaPlus * SigmaPlus()
static G4SigmaZero * SigmaZero()