Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AllisonPositronAtRestModel Class Reference

#include <G4AllisonPositronAtRestModel.hh>

+ Inheritance diagram for G4AllisonPositronAtRestModel:

Public Member Functions

 G4AllisonPositronAtRestModel ()
 
 ~G4AllisonPositronAtRestModel () override=default
 
void SampleSecondaries (std::vector< G4DynamicParticle * > &secParticles, G4double &, const G4Material *) const override
 
void PrintGeneratorInformation () const override
 
G4AllisonPositronAtRestModeloperator= (const G4AllisonPositronAtRestModel &right)=delete
 
 G4AllisonPositronAtRestModel (const G4AllisonPositronAtRestModel &)=delete
 
- Public Member Functions inherited from G4VPositronAtRestModel
 G4VPositronAtRestModel (const G4String &name)
 
virtual ~G4VPositronAtRestModel ()=default
 
const G4StringGetName () const
 
G4VPositronAtRestModeloperator= (const G4VPositronAtRestModel &right)=delete
 
 G4VPositronAtRestModel (const G4VPositronAtRestModel &)=delete
 

Detailed Description

Definition at line 53 of file G4AllisonPositronAtRestModel.hh.

Constructor & Destructor Documentation

◆ G4AllisonPositronAtRestModel() [1/2]

G4AllisonPositronAtRestModel::G4AllisonPositronAtRestModel ( )

Definition at line 51 of file G4AllisonPositronAtRestModel.cc.

52 : G4VPositronAtRestModel("Allison")
53{}
G4VPositronAtRestModel(const G4String &name)

Referenced by G4AllisonPositronAtRestModel(), and operator=().

◆ ~G4AllisonPositronAtRestModel()

G4AllisonPositronAtRestModel::~G4AllisonPositronAtRestModel ( )
overridedefault

◆ G4AllisonPositronAtRestModel() [2/2]

G4AllisonPositronAtRestModel::G4AllisonPositronAtRestModel ( const G4AllisonPositronAtRestModel & )
delete

Member Function Documentation

◆ operator=()

G4AllisonPositronAtRestModel & G4AllisonPositronAtRestModel::operator= ( const G4AllisonPositronAtRestModel & right)
delete

◆ PrintGeneratorInformation()

void G4AllisonPositronAtRestModel::PrintGeneratorInformation ( ) const
overridevirtual

Implements G4VPositronAtRestModel.

Definition at line 136 of file G4AllisonPositronAtRestModel.cc.

137{
138 G4cout << "\n" << G4endl;
139 G4cout << "Allison AtRest positron 2-gamma annihilation model." << G4endl;
140 G4cout << "Takes into account positronium motion in the media." << G4endl;
141}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout

◆ SampleSecondaries()

void G4AllisonPositronAtRestModel::SampleSecondaries ( std::vector< G4DynamicParticle * > & secParticles,
G4double & ,
const G4Material * material ) const
overridevirtual

Implements G4VPositronAtRestModel.

Definition at line 57 of file G4AllisonPositronAtRestModel.cc.

60{
61 const G4double eGamma = CLHEP::electron_mass_c2;
62
63 // In rest frame of positronium gammas are back to back
64 const G4ThreeVector& dir1 = G4RandomDirection();
65 const G4ThreeVector& dir2 = -dir1;
66 auto aGamma1 = new G4DynamicParticle(G4Gamma::Gamma(),dir1,eGamma);
67 auto aGamma2 = new G4DynamicParticle(G4Gamma::Gamma(),dir2,eGamma);
68
69 // In rest frame the gammas are polarised perpendicular to each other - see
70 // Pryce and Ward, Nature No 4065 (1947) p.435.
71 // Snyder et al, Physical Review 73 (1948) p.440.
72 G4ThreeVector pol1 = (G4RandomDirection().cross(dir1)).unit();
73 G4ThreeVector pol2 = (pol1.cross(dir2)).unit();
74
75 // A positron in matter slows down and combines with an atomic electron to
76 // make a neutral atom called positronium, about half the size of a normal
77 // atom. I expect that when the energy of the positron is small enough,
78 // less than the binding energy of positronium (6.8 eV), it is
79 // energetically favourable for an electron from the outer orbitals of a
80 // nearby atom or molecule to transfer and bind to the positron, as in an
81 // ionic bond, leaving behind a mildly ionised nearby atom/molecule. I
82 // would expect the positronium to come away with a kinetic energy of a
83 // few eV on average. In its para (spin 0) state it annihilates into two
84 // photons, which in the rest frame of the positronium are collinear
85 // (back-to-back) due to momentum conservation. Because of the motion of the
86 // positronium, photons will be not quite back-to-back in the laboratory.
87
88 // The positroniuim acquires an energy of order its binding energy and
89 // doesn't have time to thermalise. Nevertheless, here we approximate its
90 // energy distribution by a Maxwell-Boltzman with mean energy <KE>. In terms
91 // of a more familiar concept of temperature, and the law of equipartition
92 // of energy of translational motion, <KE>=3kT/2. Each component of velocity
93 // has a distribution exp(-mv^2/2kT), which is a Gaussian of mean zero
94 // and variance kT/m=2<KE>/3m, where m is the positronium mass.
95
96 const G4double meanEnergyPerIonPair = material->GetIonisation()->GetMeanEnergyPerIonPair();
97 const G4double& meanKE = meanEnergyPerIonPair; // Just an alias
98
99 if (meanKE > 0.) { // Positronium has motion
100 // Mass of positronium
101 const G4double mass = 2.*CLHEP::electron_mass_c2;
102 // Mean <KE>=3kT/2, as described above
103 // const G4double T = 2.*meanKE/(3.*k_Boltzmann);
104 // Component velocities: Gaussian, variance kT/m=2<KE>/3m.
105 const G4double sigmav = std::sqrt(2.*meanKE/(3.*mass));
106 // This is in units where c=1
107 const G4double vx = G4RandGauss::shoot(0.,sigmav);
108 const G4double vy = G4RandGauss::shoot(0.,sigmav);
109 const G4double vz = G4RandGauss::shoot(0.,sigmav);
110 const G4ThreeVector v(vx,vy,vz); // In unit where c=1
111 const G4ThreeVector& beta = v; // so beta=v/c=v
112 aGamma1->Set4Momentum(aGamma1->Get4Momentum().boost(beta));
113 aGamma2->Set4Momentum(aGamma2->Get4Momentum().boost(beta));
114
115 // Rotate polarisation vectors
116 const G4ThreeVector& newDir1 = aGamma1->GetMomentumDirection();
117 const G4ThreeVector& newDir2 = aGamma2->GetMomentumDirection();
118 const G4ThreeVector& axis1 = dir1.cross(newDir1); // No need to be unit
119 const G4ThreeVector& axis2 = dir2.cross(newDir2); // No need to be unit
120 const G4double& angle1 = std::acos(dir1*newDir1);
121 const G4double& angle2 = std::acos(dir2*newDir2);
122 pol1.rotate(axis1, angle1);
123 pol2.rotate(axis2, angle2);
124 }
125
126 // use constructors optimal for massless particle
127 aGamma1->SetPolarization(pol1);
128 aGamma2->SetPolarization(pol2);
129
130 secParticles.push_back(aGamma1);
131 secParticles.push_back(aGamma2);
132}
G4ThreeVector G4RandomDirection()
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
Hep3Vector unit() const
Hep3Vector cross(const Hep3Vector &) const
Hep3Vector & rotate(double, const Hep3Vector &)
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
G4double GetMeanEnergyPerIonPair() const
G4IonisParamMat * GetIonisation() const

The documentation for this class was generated from the following files: