Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4eeToTwoGammaModel Class Reference

#include <G4eeToTwoGammaModel.hh>

Inheritance diagram for G4eeToTwoGammaModel:
Collaboration diagram for G4eeToTwoGammaModel:

Public Member Functions

 G4eeToTwoGammaModel (const G4ParticleDefinition *p=0, const G4String &nam="eplus2gg")
 
virtual ~G4eeToTwoGammaModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
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) final
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) final
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, 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 *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
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 SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

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
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 60 of file G4eeToTwoGammaModel.hh.

Constructor & Destructor Documentation

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

Definition at line 87 of file G4eeToTwoGammaModel.cc.

89  : G4VEmModel(nam),
91  isInitialised(false)
92 {
93  theGamma = G4Gamma::Gamma();
94  fParticleChange = nullptr;
95 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
static constexpr double classic_electr_radius
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static constexpr double pi
Definition: G4SIunits.hh:75

Here is the call graph for this function:

G4eeToTwoGammaModel::~G4eeToTwoGammaModel ( )
virtual

Definition at line 99 of file G4eeToTwoGammaModel.cc.

100 {}

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 137 of file G4eeToTwoGammaModel.cc.

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

Here is the call graph for this function:

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

Reimplemented in G4PolarizedAnnihilationModel.

Definition at line 114 of file G4eeToTwoGammaModel.cc.

118 {
119  // Calculates the cross section per electron of annihilation into two photons
120  // from the Heilter formula.
121 
122  G4double ekin = std::max(eV,kineticEnergy);
123 
124  G4double tau = ekin/electron_mass_c2;
125  G4double gam = tau + 1.0;
126  G4double gamma2= gam*gam;
127  G4double bg2 = tau * (tau+2.0);
128  G4double bg = sqrt(bg2);
129 
130  G4double cross = pi_rcl2*((gamma2+4*gam+1.)*G4Log(gam+bg) - (gam+3.)*bg)
131  / (bg2*(gam+1.));
132  return cross;
133 }
static constexpr double electron_mass_c2
static constexpr double eV
Definition: G4SIunits.hh:215
G4double G4Log(G4double x)
Definition: G4Log.hh:230
T max(const T t1, const T t2)
brief Return the largest of the two arguments
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 150 of file G4eeToTwoGammaModel.cc.

155 {
156  // Calculates the cross section per volume of annihilation into two photons
157 
158  G4double eDensity = material->GetElectronDensity();
159  G4double cross = eDensity*ComputeCrossSectionPerElectron(p,kineticEnergy);
160  return cross;
161 }
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kinEnergy, G4double cutEnergy=0., G4double maxEnergy=DBL_MAX)
G4double GetElectronDensity() const
Definition: G4Material.hh:217
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Implements G4VEmModel.

Reimplemented in G4PolarizedAnnihilationModel.

Definition at line 104 of file G4eeToTwoGammaModel.cc.

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

Here is the call graph for this function:

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

!! likely problematic direction to be checked

Implements G4VEmModel.

Reimplemented in G4PolarizedAnnihilationModel.

Definition at line 168 of file G4eeToTwoGammaModel.cc.

173 {
174  G4double PositKinEnergy = dp->GetKineticEnergy();
175  G4DynamicParticle *aGamma1, *aGamma2;
176 
177  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
178 
179  // Case at rest
180  if(PositKinEnergy == 0.0) {
181  G4double cost = 2.*rndmEngine->flat()-1.;
182  G4double sint = sqrt((1. - cost)*(1. + cost));
183  G4double phi = twopi * rndmEngine->flat();
184  G4ThreeVector dir(sint*cos(phi), sint*sin(phi), cost);
185  phi = twopi * rndmEngine->flat();
186  G4double cosphi = cos(phi);
187  G4double sinphi = sin(phi);
188  G4ThreeVector pol(cosphi, sinphi, 0.0);
189  pol.rotateUz(dir);
190  aGamma1 = new G4DynamicParticle(theGamma, dir, electron_mass_c2);
191  aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
192  aGamma2 = new G4DynamicParticle(theGamma,-dir, electron_mass_c2);
193  pol.set(-sinphi, cosphi, 0.0);
194  pol.rotateUz(dir);
195  aGamma2->SetPolarization(pol.x(),pol.y(),pol.z());
196 
197  } else {
198 
199  G4ThreeVector PositDirection = dp->GetMomentumDirection();
200 
201  G4double tau = PositKinEnergy/electron_mass_c2;
202  G4double gam = tau + 1.0;
203  G4double tau2 = tau + 2.0;
204  G4double sqgrate = sqrt(tau/tau2)*0.5;
205  G4double sqg2m1 = sqrt(tau*tau2);
206 
207  // limits of the energy sampling
208  G4double epsilmin = 0.5 - sqgrate;
209  G4double epsilmax = 0.5 + sqgrate;
210  G4double epsilqot = epsilmax/epsilmin;
211 
212  //
213  // sample the energy rate of the created gammas
214  //
215  G4double epsil, greject;
216 
217  do {
218  epsil = epsilmin*G4Exp(G4Log(epsilqot)*rndmEngine->flat());
219  greject = 1. - epsil + (2.*gam*epsil-1.)/(epsil*tau2*tau2);
220  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
221  } while( greject < rndmEngine->flat());
222 
223  //
224  // scattered Gamma angles. ( Z - axis along the parent positron)
225  //
226 
227  G4double cost = (epsil*tau2-1.)/(epsil*sqg2m1);
228  if(std::abs(cost) > 1.0) {
229  G4cout << "### G4eeToTwoGammaModel WARNING cost= " << cost
230  << " positron Ekin(MeV)= " << PositKinEnergy
231  << " gamma epsil= " << epsil
232  << G4endl;
233  if(cost > 1.0) cost = 1.0;
234  else cost = -1.0;
235  }
236  G4double sint = sqrt((1.+cost)*(1.-cost));
237  G4double phi = twopi * rndmEngine->flat();
238 
239  //
240  // kinematic of the created pair
241  //
242 
243  G4double TotalAvailableEnergy = PositKinEnergy + 2.0*electron_mass_c2;
244  G4double Phot1Energy = epsil*TotalAvailableEnergy;
245 
246  G4ThreeVector Phot1Direction(sint*cos(phi), sint*sin(phi), cost);
247  Phot1Direction.rotateUz(PositDirection);
248  aGamma1 = new G4DynamicParticle (theGamma,Phot1Direction, Phot1Energy);
249  phi = twopi * rndmEngine->flat();
250  G4double cosphi = cos(phi);
251  G4double sinphi = sin(phi);
252  G4ThreeVector pol(cosphi, sinphi, 0.0);
253  pol.rotateUz(Phot1Direction);
254  aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
255 
256  G4double Phot2Energy =(1.-epsil)*TotalAvailableEnergy;
257  G4double PositP= sqrt(PositKinEnergy*(PositKinEnergy+2.*electron_mass_c2));
258  G4ThreeVector dir = PositDirection*PositP - Phot1Direction*Phot1Energy;
259  G4ThreeVector Phot2Direction = dir.unit();
260 
261  // create G4DynamicParticle object for the particle2
262  aGamma2 = new G4DynamicParticle (theGamma,Phot2Direction, Phot2Energy);
263 
265  pol.set(-sinphi, cosphi, 0.0);
266  pol.rotateUz(Phot1Direction);
267  cost = pol*Phot2Direction;
268  pol -= cost*Phot2Direction;
269  pol = pol.unit();
270  aGamma2->SetPolarization(pol.x(),pol.y(),pol.z());
271  }
272  /*
273  G4cout << "Annihilation in fly: e0= " << PositKinEnergy
274  << " m= " << electron_mass_c2
275  << " e1= " << Phot1Energy
276  << " e2= " << Phot2Energy << " dir= " << dir
277  << " -> " << Phot1Direction << " "
278  << Phot2Direction << G4endl;
279  */
280 
281  vdp->push_back(aGamma1);
282  vdp->push_back(aGamma2);
283  fParticleChange->SetProposedKineticEnergy(0.);
284  fParticleChange->ProposeTrackStatus(fStopAndKill);
285 }
G4double GetKineticEnergy() const
virtual double flat()=0
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double electron_mass_c2
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
double flat()
Definition: G4AblaRandom.cc:47
void SetPolarization(G4double polX, G4double polY, G4double polZ)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
Hep3Vector unit() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)

Here is the call graph for this function:


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