42#define MAX_SECONDARIES 100
46G4AntiProtonAnnihilationAtRest::G4AntiProtonAnnihilationAtRest(
const G4String& processName,
49 massPionMinus(
G4PionMinus::PionMinus()->GetPDGMass()/GeV),
50 massProton(
G4Proton::Proton()->GetPDGMass()/GeV),
51 massPionZero(
G4PionZero::PionZero()->GetPDGMass()/GeV),
52 massAntiProton(
G4AntiProton::AntiProton()->GetPDGMass()/GeV),
53 massPionPlus(
G4PionPlus::PionPlus()->GetPDGMass()/GeV),
54 massGamma(
G4Gamma::Gamma()->GetPDGMass()/GeV),
104 return ( &particle == pdefAntiProton );
137 G4cout <<
"G4AntiProtonAnnihilationAtRestProcess::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();
199 G4cout <<
"G4AntiProtonAnnihilationAtRest::AtRestDoIt is invoked " <<
G4endl;
207 GenerateSecondaries();
211 for (
G4int isec = 0; isec < ngkine; isec++ ) {
214 aNewParticle->
SetMomentum( gkin[isec].GetMomentum() * GeV );
216 localtime = globalTime + gkin[isec].
GetTOF();
237void G4AntiProtonAnnihilationAtRest::GenerateSecondaries()
252 result.
SetMass( massAntiProton );
257 AntiProtonAnnihilation(&nopt);
276 for (l = 1; l <= ntot; ++l) {
282 gkin[ngkine] = eve[index];
283 gkin[ngkine].
SetTOF( eve[index].GetTOF() * 5e-11 );
302void G4AntiProtonAnnihilationAtRest::Poisso(
G4float xav,
G4int *iran)
307 static G4float rr, ran, rrr, ran1;
317 ran1 = xav + ran1 * std::sqrt(xav);
331 for (i = 1; i <= fivex; ++i) {
334 rrr = std::pow(xav,
G4float(i)) / NFac(i);
338 rrr = std::exp(i * std::log(xav) -
351 p1 = xav * std::exp(-
G4double(xav));
377G4int G4AntiProtonAnnihilationAtRest::NFac(
G4int n)
392 for (i = 2; i <= j; ++i) {
401void G4AntiProtonAnnihilationAtRest::Normal(
G4float *ran)
409 for (i = 1; i <= 12; ++i) {
416void G4AntiProtonAnnihilationAtRest::AntiProtonAnnihilation(
G4int *nopt)
422 static G4int i, ii, kk;
425 static G4int ika, nbl;
430 static G4float ran1, ran2, ekin, tkin;
433 static G4float ekin1, ekin2, black;
434 static G4float pnrat, rmnve1, rmnve2;
448 pv[1].
SetMass( massAntiProton );
468 rmnve1 = massPionPlus;
469 rmnve2 = massPionMinus;
472 rmnve1 = massPionZero;
473 rmnve2 = massPionZero;
476 rmnve1 = massPionMinus;
477 rmnve2 = massPionZero;
483 ek = massProton + massAntiProton - rmnve1 - rmnve2;
490 en = ek + (rmnve1 + rmnve2) *
G4float(.5);
491 r__1 = en * en - rmnve1 * rmnve2;
492 pcm = r__1 > 0 ? std::sqrt(r__1) : 0;
503 pv[3].
SetEnergy( std::sqrt(pv[3].GetMass()*pv[3].GetMass()+pcm*pcm) );
504 pv[3].
SetMomentumAndUpdate( -pv[2].GetMomentum().x(), -pv[2].GetMomentum().y(), -pv[2].GetMomentum().z() );
527 if (targetAtomicMass >=
G4float(1.5)) {
528 cfa = (targetAtomicMass -
G4float(1.)) /
534 black = (targ *
G4float(1.25) +
535 G4float(1.5)) * evapEnergy1 / (evapEnergy1 + evapEnergy3);
538 nbl =
G4int(targetAtomicMass - targ);
546 for (i = 1; i <= nbl; ++i) {
555 ekin1 = -
G4double(ekin) * std::log(ran1) -
558 ekin1 = std::log(ran1) *
G4float(-.01);
563 ekin1 = tex - (ekin2 - ekin1);
569 pnrat =
G4float(1.) - targetCharge / targetAtomicMass;
591 for (i = 1; i <= nt; ++i) {
593 if (pv[ii].GetParticleDef() != pdefProton) {
613 (evapEnergy1 + evapEnergy3);
621 for (i = 1; i <= nbl; ++i) {
630 ekin1 = -
G4double(ekin) * std::log(ran1) -
633 ekin1 = std::log(ran1) *
G4float(-.01);
638 ekin1 = tex - (ekin2 - ekin1);
665 for (i = 3; i <= nt; ++i) {
679 static G4float cfa, gfa, ran1, ran2, ekin1, atno3;
689 if (targetAtomicMass >=
G4float(1.5)) {
703 cfa =
G4float(.13043478260869565);
704 cfa = cfa * std::log(ekin1) +
G4float(.35);
709 atno3 = targetAtomicMass;
721 gfa = (targetAtomicMass -
G4float(1.)) /
724 evapEnergy1 = ret_val * fpdiv;
725 evapEnergy3 = ret_val - evapEnergy1;
732 evapEnergy1 *= ran1 * gfa +
G4float(1.);
733 if (evapEnergy1 <
G4float(0.)) {
736 evapEnergy3 *= ran2 * gfa +
G4float(1.);
737 if (evapEnergy3 <
G4float(0.)) {
740 while ((ret_val = evapEnergy1 + evapEnergy3) >= ek1) {
std::vector< G4Element * > G4ElementVector
G4double condition(const G4ErrorSymMatrix &m)
#define G4HadronicDeprecate(name)
G4DLLIMPORT std::ostream G4cout
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
void PreparePhysicsTable(const G4ParticleDefinition &)
~G4AntiProtonAnnihilationAtRest()
G4bool IsApplicable(const G4ParticleDefinition &)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4GHEKinematicsVector * GetSecondaryKinematics()
G4double GetMeanLifeTime(const G4Track &, G4ForceCondition *)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4int GetNumberOfSecondaries()
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