47 G4cout <<
"G4RPGAntiLambdaInelastic::ApplyYourself called" <<
G4endl;
49 G4cout <<
"target material = " << targetMaterial->
GetName() <<
", ";
60 modifiedOriginal = *originalIncident;
66 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
80 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
90 targetParticle = *originalTarget;
93 G4bool incidentHasChanged =
false;
94 G4bool targetHasChanged =
false;
95 G4bool quasiElastic =
false;
103 Cascade( vec, vecLen,
104 originalIncident, currentParticle, targetParticle,
105 incidentHasChanged, targetHasChanged, quasiElastic );
108 originalIncident, originalTarget, modifiedOriginal,
109 targetNucleus, currentParticle, targetParticle,
110 incidentHasChanged, targetHasChanged, quasiElastic );
113 currentParticle, targetParticle,
114 incidentHasChanged );
116 delete originalTarget;
121void G4RPGAntiLambdaInelastic::Cascade(
127 G4bool &incidentHasChanged,
144 G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
145 targetMass*targetMass +
146 2.0*targetMass*etOriginal );
147 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
149 static G4bool first =
true;
150 const G4int numMul = 1200;
151 const G4int numMulA = 400;
152 const G4int numSec = 60;
153 static G4double protmul[numMul], protnorm[numSec];
154 static G4double neutmul[numMul], neutnorm[numSec];
155 static G4double protmulA[numMulA], protnormA[numSec];
156 static G4double neutmulA[numMulA], neutnormA[numSec];
158 G4int nt=0, np=0, nneg=0, nz=0;
166 for( i=0; i<numMul; ++i )protmul[i] = 0.0;
167 for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
169 for( np=0; np<(numSec/3); ++np )
171 for( nneg=std::max(0,np-2); nneg<=(np+1); ++nneg )
173 for( nz=0; nz<numSec/3; ++nz )
175 if( ++counter < numMul )
178 if( nt>0 && nt<=numSec )
180 protmul[counter] =
Pmltpc(np,nneg,nz,nt,b[0],c);
181 protnorm[nt-1] += protmul[counter];
187 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
188 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
190 for( np=0; np<numSec/3; ++np )
192 for( nneg=std::max(0,np-1); nneg<=(np+2); ++nneg )
194 for( nz=0; nz<numSec/3; ++nz )
196 if( ++counter < numMul )
199 if( nt>0 && nt<=numSec )
201 neutmul[counter] =
Pmltpc(np,nneg,nz,nt,b[1],c);
202 neutnorm[nt-1] += neutmul[counter];
208 for( i=0; i<numSec; ++i )
210 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
211 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
216 for( i=0; i<numMulA; ++i )protmulA[i] = 0.0;
217 for( i=0; i<numSec; ++i )protnormA[i] = 0.0;
219 for( np=1; np<(numSec/3); ++np )
222 for( nz=0; nz<numSec/3; ++nz )
224 if( ++counter < numMulA )
227 if( nt>1 && nt<=numSec )
229 protmulA[counter] =
Pmltpc(np,nneg,nz,nt,b[0],c);
230 protnormA[nt-1] += protmulA[counter];
235 for( i=0; i<numMulA; ++i )neutmulA[i] = 0.0;
236 for( i=0; i<numSec; ++i )neutnormA[i] = 0.0;
238 for( np=0; np<numSec/3; ++np )
241 for( nz=0; nz<numSec/3; ++nz )
243 if( ++counter < numMulA )
246 if( nt>1 && nt<=numSec )
248 neutmulA[counter] =
Pmltpc(np,nneg,nz,nt,b[1],c);
249 neutnormA[nt-1] += neutmulA[counter];
254 for( i=0; i<numSec; ++i )
256 if( protnormA[i] > 0.0 )protnormA[i] = 1.0/protnormA[i];
257 if( neutnormA[i] > 0.0 )neutnormA[i] = 1.0/neutnormA[i];
274 const G4double anhl[] = {1.00,1.00,1.00,1.00,1.00,1.00,1.00,1.00,0.97,0.88,
275 0.85,0.81,0.75,0.64,0.64,0.55,0.55,0.45,0.47,0.40,
276 0.39,0.36,0.33,0.10,0.01};
278 if( iplab > 9 )iplab =
G4int( (pOriginal- 1.0)*5.0 ) + 10;
279 if( iplab > 14 )iplab =
G4int( pOriginal- 2.0 ) + 15;
280 if( iplab > 22 )iplab =
G4int( (pOriginal-10.0)/10.0 ) + 23;
281 if( iplab > 24 )iplab = 24;
284 if( availableEnergy <= aPiPlus->GetPDGMass()/MeV )
296 for( np=0; np<numSec/3 && ran>=excs; ++np )
298 for( nneg=std::max(0,np-2); nneg<=(np+1) && ran>=excs; ++nneg )
300 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
302 if( ++counter < numMul )
305 if( nt>0 && nt<=numSec )
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;
326 G4int ncht = std::min( 4, std::max( 1, np-nneg+2 ) );
331 incidentHasChanged =
true;
339 incidentHasChanged =
true;
344 incidentHasChanged =
true;
346 targetHasChanged =
true;
354 incidentHasChanged =
true;
356 targetHasChanged =
true;
366 incidentHasChanged =
true;
368 targetHasChanged =
true;
373 incidentHasChanged =
true;
381 targetHasChanged =
true;
386 incidentHasChanged =
true;
392 incidentHasChanged =
true;
394 targetHasChanged =
true;
401 for( np=0; np<numSec/3 && ran>=excs; ++np )
403 for( nneg=std::max(0,np-1); nneg<=(np+2) && ran>=excs; ++nneg )
405 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
407 if( ++counter < numMul )
410 if( nt>0 && nt<=numSec )
412 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
413 dum = (
pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
414 if( std::fabs(dum) < 1.0 )
416 if( test >= 1.0e-10 )excs += dum*test;
431 G4int ncht = std::min( 4, std::max( 1, np-nneg+3 ) );
436 incidentHasChanged =
true;
438 targetHasChanged =
true;
446 incidentHasChanged =
true;
448 targetHasChanged =
true;
453 incidentHasChanged =
true;
461 targetHasChanged =
true;
466 incidentHasChanged =
true;
476 incidentHasChanged =
true;
481 incidentHasChanged =
true;
483 targetHasChanged =
true;
491 incidentHasChanged =
true;
493 targetHasChanged =
true;
499 incidentHasChanged =
true;
506 if( centerofmassEnergy <= aPiPlus->GetPDGMass()/MeV+aKaonPlus->
GetPDGMass()/MeV )
521 for( np=1; np<numSec/3 && ran>=excs; ++np )
524 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
526 if( ++counter < numMulA )
529 if( nt>1 && nt<=numSec )
531 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
532 dum = (
pi/anpn)*nt*protmulA[counter]*protnormA[nt-1]/(2.0*n*n);
533 if( std::fabs(dum) < 1.0 )
535 if( test >= 1.0e-10 )excs += dum*test;
547 for( np=0; np<numSec/3 && ran>=excs; ++np )
550 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
552 if( ++counter < numMulA )
555 if( nt>1 && nt<=numSec )
557 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
558 dum = (
pi/anpn)*nt*neutmulA[counter]*neutnormA[nt-1]/(2.0*n*n);
559 if( std::fabs(dum) < 1.0 )
561 if( test >= 1.0e-10 )excs += dum*test;
576 currentParticle.
SetMass( 0.0 );
580 if( currentParticle.
GetMass() == 0.0 )
586 for(
G4int i=0; i<vecLen; ++i )
588 if( vec[i]->GetDefinition() == aPiMinus )
590 vec[i]->SetDefinitionAndUpdateE( aKaonMinus );
600 for(
G4int i=0; i<vecLen; ++i )
602 if( vec[i]->GetDefinition() == aPiZero )
604 vec[i]->SetDefinitionAndUpdateE( aKaonZL );
615 for(
G4int i=0; i<vecLen; ++i )
617 if( vec[i]->GetDefinition() == aPiMinus )
619 vec[i]->SetDefinitionAndUpdateE( aKaonMinus );
627 for(
G4int i=0; i<vecLen; ++i )
629 if( vec[i]->GetDefinition() == aPiZero )
631 vec[i]->SetDefinitionAndUpdateE( aKaonZL );
G4DLLIMPORT std::ostream G4cout
static G4AntiSigmaMinus * AntiSigmaMinus()
static G4AntiSigmaPlus * AntiSigmaPlus()
static G4AntiSigmaZero * AntiSigmaZero()
G4ParticleDefinition * GetDefinition() const
void Initialize(G4int items)
const G4Material * GetMaterial() const
G4double GetTotalMomentum() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4HadFinalState theParticleChange
static G4KaonMinus * KaonMinus()
static G4KaonPlus * KaonPlus()
static G4KaonZeroLong * KaonZeroLong()
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()
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
void SetUpPions(const G4int np, const G4int nm, const G4int nz, G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen)
void GetNormalizationConstant(const G4double availableEnergy, G4double &n, G4double &anpn)
void CalculateMomenta(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, const G4DynamicParticle *originalTarget, G4ReactionProduct &modifiedOriginal, G4Nucleus &targetNucleus, G4ReactionProduct ¤tParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool quasiElastic)
void SetUpChange(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, G4ReactionProduct ¤tParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
G4double Pmltpc(G4int np, G4int nm, G4int nz, G4int n, G4double b, G4double c)
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4double GetTotalMomentum() const
G4ThreeVector GetMomentum() const
void SetSide(const G4int sid)
void SetDefinitionAndUpdateE(G4ParticleDefinition *aParticleDefinition)
void SetKineticEnergy(const G4double en)
G4ParticleDefinition * GetDefinition() const
void SetMass(const G4double mas)