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

#include <G4eeToTwoGammaModel.hh>

+ Inheritance diagram for G4eeToTwoGammaModel:

Public Member Functions

 G4eeToTwoGammaModel (const G4ParticleDefinition *p=0, const G4String &nam="eplus2gg")
 
virtual ~G4eeToTwoGammaModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kinEnergy, G4double cutEnergy=0., G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0., G4double maxEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Detailed Description

Definition at line 60 of file G4eeToTwoGammaModel.hh.

Constructor & Destructor Documentation

◆ G4eeToTwoGammaModel()

G4eeToTwoGammaModel::G4eeToTwoGammaModel ( const G4ParticleDefinition p = 0,
const G4String nam = "eplus2gg" 
)

Definition at line 85 of file G4eeToTwoGammaModel.cc.

87 : G4VEmModel(nam),
88 pi_rcl2(pi*classic_electr_radius*classic_electr_radius),
89 isInitialised(false)
90{
91 theGamma = G4Gamma::Gamma();
92 fParticleChange = 0;
93}
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86

◆ ~G4eeToTwoGammaModel()

G4eeToTwoGammaModel::~G4eeToTwoGammaModel ( )
virtual

Definition at line 97 of file G4eeToTwoGammaModel.cc.

98{}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4eeToTwoGammaModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0.,
G4double  cutEnergy = 0.,
G4double  maxEnergy = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 135 of file G4eeToTwoGammaModel.cc.

139{
140 // Calculates the cross section per atom of annihilation into two photons
141
142 G4double cross = Z*ComputeCrossSectionPerElectron(p,kineticEnergy);
143 return cross;
144}
double G4double
Definition: G4Types.hh:64
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kinEnergy, G4double cutEnergy=0., G4double maxEnergy=DBL_MAX)

◆ ComputeCrossSectionPerElectron()

G4double G4eeToTwoGammaModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  cutEnergy = 0.,
G4double  maxEnergy = DBL_MAX 
)
virtual

Reimplemented in G4PolarizedAnnihilationModel.

Definition at line 112 of file G4eeToTwoGammaModel.cc.

116{
117 // Calculates the cross section per electron of annihilation into two photons
118 // from the Heilter formula.
119
120 G4double ekin = std::max(eV,kineticEnergy);
121
122 G4double tau = ekin/electron_mass_c2;
123 G4double gam = tau + 1.0;
124 G4double gamma2= gam*gam;
125 G4double bg2 = tau * (tau+2.0);
126 G4double bg = sqrt(bg2);
127
128 G4double cross = pi_rcl2*((gamma2+4*gam+1.)*log(gam+bg) - (gam+3.)*bg)
129 / (bg2*(gam+1.));
130 return cross;
131}

Referenced by ComputeCrossSectionPerAtom(), G4PolarizedAnnihilationModel::ComputeCrossSectionPerElectron(), and CrossSectionPerVolume().

◆ CrossSectionPerVolume()

G4double G4eeToTwoGammaModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 148 of file G4eeToTwoGammaModel.cc.

153{
154 // Calculates the cross section per volume of annihilation into two photons
155
156 G4double eDensity = material->GetElectronDensity();
157 G4double cross = eDensity*ComputeCrossSectionPerElectron(p,kineticEnergy);
158 return cross;
159}
G4double GetElectronDensity() const
Definition: G4Material.hh:216

◆ Initialise()

void G4eeToTwoGammaModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Reimplemented in G4PolarizedAnnihilationModel.

Definition at line 102 of file G4eeToTwoGammaModel.cc.

104{
105 if(isInitialised) { return; }
106 fParticleChange = GetParticleChangeForGamma();
107 isInitialised = true;
108}
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109

◆ SampleSecondaries()

void G4eeToTwoGammaModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  vdp,
const G4MaterialCutsCouple ,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
virtual

!! likely problematic direction to be checked

Implements G4VEmModel.

Reimplemented in G4PolarizedAnnihilationModel.

Definition at line 163 of file G4eeToTwoGammaModel.cc.

168{
169 G4double PositKinEnergy = dp->GetKineticEnergy();
170 G4DynamicParticle *aGamma1, *aGamma2;
171
172 // Case at rest
173 if(PositKinEnergy == 0.0) {
174 G4double cost = 2.*G4UniformRand()-1.;
175 G4double sint = sqrt((1. - cost)*(1. + cost));
176 G4double phi = twopi * G4UniformRand();
177 G4ThreeVector dir(sint*cos(phi), sint*sin(phi), cost);
178 phi = twopi * G4UniformRand();
179 G4ThreeVector pol(cos(phi), sin(phi), 0.0);
180 pol.rotateUz(dir);
181 aGamma1 = new G4DynamicParticle(theGamma, dir, electron_mass_c2);
182 aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
183 aGamma2 = new G4DynamicParticle(theGamma,-dir, electron_mass_c2);
184 aGamma1->SetPolarization(-pol.x(),-pol.y(),-pol.z());
185
186 } else {
187
188 G4ThreeVector PositDirection = dp->GetMomentumDirection();
189
190 G4double tau = PositKinEnergy/electron_mass_c2;
191 G4double gam = tau + 1.0;
192 G4double tau2 = tau + 2.0;
193 G4double sqgrate = sqrt(tau/tau2)*0.5;
194 G4double sqg2m1 = sqrt(tau*tau2);
195
196 // limits of the energy sampling
197 G4double epsilmin = 0.5 - sqgrate;
198 G4double epsilmax = 0.5 + sqgrate;
199 G4double epsilqot = epsilmax/epsilmin;
200
201 //
202 // sample the energy rate of the created gammas
203 //
204 G4double epsil, greject;
205
206 do {
207 epsil = epsilmin*pow(epsilqot,G4UniformRand());
208 greject = 1. - epsil + (2.*gam*epsil-1.)/(epsil*tau2*tau2);
209 } while( greject < G4UniformRand() );
210
211 //
212 // scattered Gamma angles. ( Z - axis along the parent positron)
213 //
214
215 G4double cost = (epsil*tau2-1.)/(epsil*sqg2m1);
216 if(std::abs(cost) > 1.0) {
217 G4cout << "### G4eeToTwoGammaModel WARNING cost= " << cost
218 << " positron Ekin(MeV)= " << PositKinEnergy
219 << " gamma epsil= " << epsil
220 << G4endl;
221 if(cost > 1.0) cost = 1.0;
222 else cost = -1.0;
223 }
224 G4double sint = sqrt((1.+cost)*(1.-cost));
225 G4double phi = twopi * G4UniformRand();
226
227 //
228 // kinematic of the created pair
229 //
230
231 G4double TotalAvailableEnergy = PositKinEnergy + 2.0*electron_mass_c2;
232 G4double Phot1Energy = epsil*TotalAvailableEnergy;
233
234 G4ThreeVector Phot1Direction(sint*cos(phi), sint*sin(phi), cost);
235 Phot1Direction.rotateUz(PositDirection);
236 aGamma1 = new G4DynamicParticle (theGamma,Phot1Direction, Phot1Energy);
237 phi = twopi * G4UniformRand();
238 G4ThreeVector pol1(cos(phi), sin(phi), 0.0);
239 pol1.rotateUz(Phot1Direction);
240 aGamma1->SetPolarization(pol1.x(),pol1.y(),pol1.z());
241
242 G4double Phot2Energy =(1.-epsil)*TotalAvailableEnergy;
243 G4double PositP= sqrt(PositKinEnergy*(PositKinEnergy+2.*electron_mass_c2));
244 G4ThreeVector dir = PositDirection*PositP - Phot1Direction*Phot1Energy;
245 G4ThreeVector Phot2Direction = dir.unit();
246
247 // create G4DynamicParticle object for the particle2
248 aGamma2 = new G4DynamicParticle (theGamma,Phot2Direction, Phot2Energy);
249
250 //!!! likely problematic direction to be checked
251 aGamma2->SetPolarization(-pol1.x(),-pol1.y(),-pol1.z());
252 }
253 /*
254 G4cout << "Annihilation in fly: e0= " << PositKinEnergy
255 << " m= " << electron_mass_c2
256 << " e1= " << Phot1Energy
257 << " e2= " << Phot2Energy << " dir= " << dir
258 << " -> " << Phot1Direction << " "
259 << Phot2Direction << G4endl;
260 */
261
262 vdp->push_back(aGamma1);
263 vdp->push_back(aGamma2);
264 fParticleChange->SetProposedKineticEnergy(0.);
265 fParticleChange->ProposeTrackStatus(fStopAndKill);
266}
@ fStopAndKill
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
const G4ThreeVector & GetMomentumDirection() const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeTrackStatus(G4TrackStatus status)

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