42#define MAX_SECONDARIES 100
46G4AntiNeutronAnnihilationAtRest::G4AntiNeutronAnnihilationAtRest(
const G4String& processName,
49 massPionMinus(
G4PionMinus::PionMinus()->GetPDGMass()/GeV),
50 massPionZero(
G4PionZero::PionZero()->GetPDGMass()/GeV),
51 massPionPlus(
G4PionPlus::PionPlus()->GetPDGMass()/GeV),
52 massGamma(
G4Gamma::Gamma()->GetPDGMass()/GeV),
53 massAntiNeutron(
G4AntiNeutron::AntiNeutron()->GetPDGMass()/GeV),
54 massNeutron(
G4Neutron::Neutron()->GetPDGMass()/GeV),
104 return ( &particle == pdefAntiNeutron );
137 G4cout <<
"G4AntiNeutronAnnihilationAtRestProcess::AtRestGetPhysicalInteractionLength ";
174 for (
G4int i1=0; i1 < numberOfElements; i1++ )
176 normalization += theAtomicNumberDensity[i1] ;
181 for (
G4int i2=0; i2 < numberOfElements; i2++ )
183 runningSum += theAtomicNumberDensity[i2];
185 if (random<=runningSum)
187 targetCharge =
G4double( ((*theElementVector)[i2])->GetZ());
188 targetAtomicMass = (*theElementVector)[i2]->GetN();
191 if (random>runningSum)
193 targetCharge =
G4double( ((*theElementVector)[numberOfElements-1])->GetZ());
194 targetAtomicMass = (*theElementVector)[numberOfElements-1]->GetN();
198 G4cout <<
"G4AntiNeutronAnnihilationAtRest::AtRestDoIt is invoked " <<
G4endl;
206 GenerateSecondaries();
210 for (
G4int isec = 0; isec < ngkine; isec++ ) {
213 aNewParticle->
SetMomentum( gkin[isec].GetMomentum() * GeV );
215 localtime = globalTime + gkin[isec].
GetTOF();
236void G4AntiNeutronAnnihilationAtRest::GenerateSecondaries()
252 result.
SetMass( massAntiNeutron );
259 AntiNeutronAnnihilation(&nopt);
278 for (l = 1; l <= ntot; ++l) {
284 gkin[ngkine] = eve[index];
285 gkin[ngkine].
SetTOF( eve[index].GetTOF() * 5e-11 );
304void G4AntiNeutronAnnihilationAtRest::Poisso(
G4float xav,
G4int *iran)
309 static G4float rr, ran, rrr, ran1;
319 ran1 = xav + ran1 * std::sqrt(xav);
333 for (i = 1; i <= fivex; ++i) {
336 rrr = std::pow(xav,
G4float(i)) / NFac(i);
340 rrr = std::exp(i * std::log(xav) -
353 p1 = xav * std::exp(-
G4double(xav));
379G4int G4AntiNeutronAnnihilationAtRest::NFac(
G4int n)
394 for (i = 2; i <= j; ++i) {
403void G4AntiNeutronAnnihilationAtRest::Normal(
G4float *ran)
411 for (i = 1; i <= 12; ++i) {
418void G4AntiNeutronAnnihilationAtRest::AntiNeutronAnnihilation(
G4int *nopt)
424 static G4int i, ii, kk;
427 static G4int ika, nbl;
432 static G4float ran1, ran2, ekin, tkin;
435 static G4float ekin1, ekin2, black;
436 static G4float pnrat, rmnve1, rmnve2;
450 pv[1].
SetMass( massAntiNeutron );
469 rmnve1 = massPionPlus;
470 rmnve2 = massPionMinus;
472 rmnve1 = massPionZero;
475 rmnve2 = massPionZero;
478 rmnve2 = massPionZero;
486 ek = massNeutron + massAntiNeutron - rmnve1 - rmnve2;
493 en = ek + (rmnve1 + rmnve2) /
G4float(2.);
494 r__1 = en * en - rmnve1 * rmnve2;
495 pcm = r__1 > 0 ? std::sqrt(r__1) : 0;
506 pv[3].
SetEnergy( std::sqrt(pv[3].GetMass()*pv[3].GetMass()+pcm*pcm) );
507 pv[3].
SetMomentumAndUpdate( -pv[2].GetMomentum().x(), -pv[2].GetMomentum().y(), -pv[2].GetMomentum().z() );
530 if (targetAtomicMass >=
G4float(1.5)) {
537 black = (targ *
G4float(1.25) +
538 G4float(1.5)) * evapEnergy1 / (evapEnergy1 + evapEnergy3);
541 nbl =
G4int(targetAtomicMass - targ);
549 for (i = 1; i <= nbl; ++i) {
558 ekin1 = -
G4double(ekin) * std::log(ran1) -
561 ekin1 = std::log(ran1) *
G4float(-.01);
566 ekin1 = tex - (ekin2 - ekin1);
572 pnrat =
G4float(1.) - targetCharge / targetAtomicMass;
594 for (i = 1; i <= nt; ++i) {
596 if (pv[ii].GetParticleDef() != pdefProton) {
616 (evapEnergy1 + evapEnergy3);
624 for (i = 1; i <= nbl; ++i) {
633 ekin1 = -
G4double(ekin) * std::log(ran1) -
636 ekin1 = std::log(ran1) *
G4float(-.01);
641 ekin1 = tex - (ekin2 - ekin1);
668 for (i = 3; i <= nt; ++i) {
682 static G4float cfa, gfa, ran1, ran2, ekin1, atno3;
692 if (targetAtomicMass >=
G4float(1.5)) {
706 cfa =
G4float(.13043478260869565);
707 cfa = cfa * std::log(ekin1) +
G4float(.35);
712 atno3 = targetAtomicMass;
724 gfa = (targetAtomicMass -
G4float(1.)) /
727 evapEnergy1 = ret_val * fpdiv;
728 evapEnergy3 = ret_val - evapEnergy1;
735 evapEnergy1 *= ran1 * gfa +
G4float(1.);
736 if (evapEnergy1 <
G4float(0.)) {
739 evapEnergy3 *= ran2 * gfa +
G4float(1.);
740 if (evapEnergy3 <
G4float(0.)) {
743 while ((ret_val = evapEnergy1 + evapEnergy3) >= ek1) {
std::vector< G4Element * > G4ElementVector
G4double condition(const G4ErrorSymMatrix &m)
#define G4HadronicDeprecate(name)
G4DLLIMPORT std::ostream G4cout
void PreparePhysicsTable(const G4ParticleDefinition &)
G4GHEKinematicsVector * GetSecondaryKinematics()
G4double GetMeanLifeTime(const G4Track &, G4ForceCondition *)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4int GetNumberOfSecondaries()
void BuildPhysicsTable(const G4ParticleDefinition &)
~G4AntiNeutronAnnihilationAtRest()
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
G4bool IsApplicable(const G4ParticleDefinition &)
void DumpInfo(G4int mode=0) const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMomentum(const G4ThreeVector &momentum)
void SetEnergy(G4double e)
void SetEnergyAndUpdate(G4double e)
G4ParticleDefinition * GetParticleDef()
void SetParticleDef(G4ParticleDefinition *c)
void SetMomentumAndUpdate(G4ParticleMomentum mom)
void SetMass(G4double mas)
void SetKineticEnergyAndUpdate(G4double ekin)
void DeRegisterExtraProcess(G4VProcess *)
void RegisterExtraProcess(G4VProcess *)
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
static G4HadronicProcessStore * Instance()
void PrintInfo(const G4ParticleDefinition *)
const G4ElementVector * GetElementVector() const
size_t GetNumberOfElements() const
const G4double * GetAtomicNumDensityVector() const
const G4String & GetName() const
void AddSecondary(G4Track *aSecondary)
virtual void Initialize(const G4Track &)
G4double GetPDGMass() const
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4double currentInteractionLength
virtual void ResetNumberOfInteractionLengthLeft()
G4ParticleChange aParticleChange
G4double theNumberOfInteractionLengthLeft
void SetProcessSubType(G4int)
const G4String & GetProcessName() const