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

#include <G4eplusTo3GammaOKVIModel.hh>

+ Inheritance diagram for G4eplusTo3GammaOKVIModel:

Public Member Functions

 G4eplusTo3GammaOKVIModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="eplus3ggOKVI")
 
 ~G4eplusTo3GammaOKVIModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double ComputeCrossSectionPerElectron (G4double kinEnergy)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0., G4double maxEnergy=DBL_MAX) final
 
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) final
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX) final
 
G4double ComputeF (G4double fr1, G4double fr2, G4double fr3, G4double kinEnergy)
 
G4double ComputeF0 (G4double fr1, G4double fr2, G4double fr3)
 
G4double ComputeFS (G4double fr1, G4double fr2, G4double fr3, G4double kinEnergy)
 
void SetDelta (G4double val)
 
G4eplusTo3GammaOKVIModeloperator= (const G4eplusTo3GammaOKVIModel &right)=delete
 
 G4eplusTo3GammaOKVIModel (const G4eplusTo3GammaOKVIModel &)=delete
 
- 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 *, const G4double &length, G4double &eloss)
 
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 FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
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)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () 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) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (G4bool)
 
void SetMasterThread (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

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 = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 56 of file G4eplusTo3GammaOKVIModel.hh.

Constructor & Destructor Documentation

◆ G4eplusTo3GammaOKVIModel() [1/2]

G4eplusTo3GammaOKVIModel::G4eplusTo3GammaOKVIModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "eplus3ggOKVI" )
explicit

Definition at line 60 of file G4eplusTo3GammaOKVIModel.cc.

62 : G4VEmModel(nam), fDelta(0.001)
63{
64 theGamma = G4Gamma::Gamma();
65}
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

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

◆ ~G4eplusTo3GammaOKVIModel()

G4eplusTo3GammaOKVIModel::~G4eplusTo3GammaOKVIModel ( )
overridedefault

◆ G4eplusTo3GammaOKVIModel() [2/2]

G4eplusTo3GammaOKVIModel::G4eplusTo3GammaOKVIModel ( const G4eplusTo3GammaOKVIModel & )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 197 of file G4eplusTo3GammaOKVIModel.cc.

201{
202 G4double cross = Z*ComputeCrossSectionPerElectron(kineticEnergy);
203 return cross;
204}
double G4double
Definition G4Types.hh:83
G4double ComputeCrossSectionPerElectron(G4double kinEnergy)

◆ ComputeCrossSectionPerElectron()

G4double G4eplusTo3GammaOKVIModel::ComputeCrossSectionPerElectron ( G4double kinEnergy)

Definition at line 176 of file G4eplusTo3GammaOKVIModel.cc.

177{
178 // Calculates the cross section per electron of annihilation into 3 photons
179 // from the Heilter formula.
180
181 G4double ekin = std::max(CLHEP::eV, kinEnergy);
182 G4double tau = ekin/CLHEP::electron_mass_c2;
183 G4double gam = tau + 1.0;
184 G4double gamma2 = gam*gam;
185 G4double bg2 = tau * (tau+2.0);
186 G4double bg = sqrt(bg2);
187
188 G4double rho = (gamma2+4*gam+1.)*G4Log(gam+bg)/(gamma2-1.)
189 - (gam+3.)/(sqrt(gam*gam - 1.));
190
191 G4double cross = alpha_rcl2*(4.2 - (2.*G4Log(fDelta)+1.)*rho*rho)/(gam+1.);
192 return cross;
193}
G4double G4Log(G4double x)
Definition G4Log.hh:227

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

◆ ComputeF()

G4double G4eplusTo3GammaOKVIModel::ComputeF ( G4double fr1,
G4double fr2,
G4double fr3,
G4double kinEnergy )

Definition at line 80 of file G4eplusTo3GammaOKVIModel.cc.

82{
83 G4double ekin = std::max(eV,kinEnergy);
84 G4double tau = ekin/electron_mass_c2;
85 G4double gam = tau + 1.0;
86 G4double gamma2 = gam*gam;
87 G4double bg2 = tau * (tau+2.0);
88 G4double bg = sqrt(bg2);
89
90 G4double rho = (gamma2+4.*gam+1.)*G4Log(gam+bg)/(gamma2-1.)
91 - (gam+3.)/(sqrt(gam*gam - 1.)) + 1.;
92 G4double border;
93
94 if(ekin < 500*MeV) {
95 border = 1. - (electron_mass_c2)/(2*(ekin + electron_mass_c2));
96 } else {
97 border = 1. - (100*electron_mass_c2)/(2*(ekin + electron_mass_c2));
98 }
99
100 border = std::min(border, 0.9999);
101
102 if (fr1>border) { fr1 = border; }
103 if (fr2>border) { fr2 = border; }
104 if (fr3>border) { fr3 = border; }
105
106 G4double fr1s = fr1*fr1; // "s" for "squared"
107 G4double fr2s = fr2*fr2;
108 G4double fr3s = fr3*fr3;
109
110 G4double aa = (1.-fr1)*(1.-fr2);
111 G4double ab = fr3s + (fr1-fr2)*(fr1-fr2);
112 G4double add= ((1.-fr1)*(1.-fr1) + (1.-fr2)*(1.-fr2))/(fr3s*aa);
113
114 G4double fres = -rho*(1./fr1s + 1./fr2s)
115 + (ab/(2.*(fr1*fr2*aa)))*(G4Log(2.*gam*aa/(fr1*fr2)))
116 + (ab/(2.*fr1*fr2*(1-fr3)))*G4Log(2.*gam*(1.-fr3)/(fr1*fr2)) - add;
117
118 return fres;
119}

Referenced by ComputeFS().

◆ ComputeF0()

G4double G4eplusTo3GammaOKVIModel::ComputeF0 ( G4double fr1,
G4double fr2,
G4double fr3 )

Definition at line 124 of file G4eplusTo3GammaOKVIModel.cc.

126{
127 G4double tau = 0.0;
128 G4double gam = tau + 1.0;
129 G4double gamma2 = gam*gam;
130 G4double bg2 = tau * (tau+2.0);
131 G4double bg = sqrt(bg2);
132
133 G4double rho = (gamma2+4.*gam+1.)*G4Log(gam+bg)/(gamma2-1.)
134 - (gam+3.)/(sqrt(gam*gam - 1.)) + 1.;
135 G4double border = 0.5;
136
137 if (fr1>border) { fr1 = border; }
138 if (fr2>border) { fr2 = border; }
139 if (fr3>border) { fr3 = border; }
140
141 G4double fr1s = fr1*fr1; // "s" for "squared"
142 G4double fr2s = fr2*fr2;
143 G4double fr3s = fr3*fr3;
144
145 G4double aa = (1.-fr1)*(1.-fr2);
146 G4double ab = fr3s + (fr1-fr2)*(fr1-fr2);
147 G4double add= ((1.-fr1)*(1.-fr1) + (1.-fr2)*(1.-fr2))/(fr3s*aa);
148
149 G4double fres = -rho*(1./fr1s + 1./fr2s)
150 + (ab/(2.*(fr1*fr2*aa)))*(G4Log(2.*gam*aa/(fr1*fr2)))
151 + (ab/(2.*fr1*fr2*(1-fr3)))*G4Log(2.*gam*(1.-fr3)/(fr1*fr2)) - add;
152
153 return fres;
154}

◆ ComputeFS()

G4double G4eplusTo3GammaOKVIModel::ComputeFS ( G4double fr1,
G4double fr2,
G4double fr3,
G4double kinEnergy )

Definition at line 157 of file G4eplusTo3GammaOKVIModel.cc.

159{
160 G4double ekin = std::max(eV,kinEnergy);
161 G4double tau = ekin/electron_mass_c2;
162 G4double gam = tau + 1.0;
163
164 G4double fsum = fr1*fr1*(ComputeF(fr1,fr2,fr3,ekin) +
165 ComputeF(fr3,fr1,fr2,ekin) +
166 ComputeF(fr2,fr3,fr1,ekin));
167
168 G4double dcross = fsum/((3*fr1*fr1*(gam+1.)));
169
170 return dcross;
171}
G4double ComputeF(G4double fr1, G4double fr2, G4double fr3, G4double kinEnergy)

Referenced by SampleSecondaries().

◆ CrossSectionPerVolume()

G4double G4eplusTo3GammaOKVIModel::CrossSectionPerVolume ( const G4Material * material,
const G4ParticleDefinition * ,
G4double kineticEnergy,
G4double cutEnergy = 0.0,
G4double maxEnergy = DBL_MAX )
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 208 of file G4eplusTo3GammaOKVIModel.cc.

213{
214 // Calculates the cross section per volume of annihilation into two photons
215
216 G4double eDensity = material->GetElectronDensity();
217 G4double cross = eDensity*ComputeCrossSectionPerElectron(kineticEnergy);
218 return cross;
219}
G4double GetElectronDensity() const

◆ Initialise()

void G4eplusTo3GammaOKVIModel::Initialise ( const G4ParticleDefinition * ,
const G4DataVector &  )
overridevirtual

Implements G4VEmModel.

Definition at line 73 of file G4eplusTo3GammaOKVIModel.cc.

75{}

◆ operator=()

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

◆ SampleSecondaries()

void G4eplusTo3GammaOKVIModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * vdp,
const G4MaterialCutsCouple * ,
const G4DynamicParticle * dp,
G4double tmin = 0.0,
G4double maxEnergy = DBL_MAX )
finalvirtual

!! POLARIZATION - not yet implemented

Implements G4VEmModel.

Definition at line 227 of file G4eplusTo3GammaOKVIModel.cc.

231{
232 // let us perform sampling in C.M.S. reference frame of e- at rest and e+ on fly
233 G4double posiKinEnergy = dp->GetKineticEnergy();
235 posiKinEnergy + 2*CLHEP::electron_mass_c2);
236 G4double eGammaCMS = 0.5 * lv.mag();
237
238 // the limit value fDelta is defined by a class, which call this method
239 // thickness of border defined by C.M.S. energy
240 G4double border =
241 1.0 - std::min(std::max(CLHEP::electron_mass_c2/eGammaCMS, fDelta), 0.1);
242
243 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
244
245 G4ThreeVector posiDirection = dp->GetMomentumDirection();
246
247 // (A.A.) LIMITS FOR 1st GAMMA
248 G4double xmin = 0.01;
249 G4double xmax = 0.667; // CHANGE to 3/2
250
251 G4double d1, d0, x1, x2, dmax, x2min;
252
253 // (A.A.) sampling of x1 x2 x3 (whole cycle of rejection)
254 do {
255 x1 = 1./((1./xmin) - ((1./xmin)-(1./xmax))*rndmEngine->flat());
256 dmax = ComputeFS(eGammaCMS, x1, 1.-x1, border);
257 x2min = 1. - x1;
258 x2 = 1 - rndmEngine->flat()*(1. - x2min);
259 d1 = dmax*rndmEngine->flat();
260 d0 = ComputeFS(eGammaCMS, x1, x2, 2.-x1-x2);
261 }
262 while(d0 < d1);
263
264 G4double x3 = 2 - x1 - x2;
265 //
266 // angles between Gammas
267 //
268 G4double psi13 = 2*std::asin(std::sqrt(std::abs((x1+x3-1.)/(x1*x3))));
269 G4double psi12 = 2*std::asin(std::sqrt(std::abs((x1+x2-1.)/(x1*x2))));
270 //
271 // kinematic of the created pair
272 //
273 G4double phot1Energy = x1*eGammaCMS;
274 G4double phot2Energy = x2*eGammaCMS;
275 G4double phot3Energy = x3*eGammaCMS;
276
277 // DIRECTIONS
278
279 // The azimuthal angles of q1 and q3 with respect to some plane
280 // through the beam axis are generated at random.
281
282 G4ThreeVector phot1Direction(0, 0, 1);
283 G4ThreeVector phot2Direction(0, std::sin(psi12), std::cos(psi12));
284 G4ThreeVector phot3Direction(0, std::sin(psi13), std::cos(psi13));
285
286 G4LorentzVector lv1(phot1Energy*phot1Direction, phot1Energy);
287 G4LorentzVector lv2(phot2Energy*phot2Direction, phot2Energy);
288 G4LorentzVector lv3(phot3Energy*phot3Direction, phot3Energy);
289
290 auto boostV = lv.boostVector();
291 lv1.boost(boostV);
292 lv2.boost(boostV);
293 lv3.boost(boostV);
294
295 auto aGamma1 = new G4DynamicParticle (theGamma, lv1.vect());
296 auto aGamma2 = new G4DynamicParticle (theGamma, lv2.vect());
297 auto aGamma3 = new G4DynamicParticle (theGamma, lv3.vect());
298
299 //!!! POLARIZATION - not yet implemented
300
301 vdp->push_back(aGamma1);
302 vdp->push_back(aGamma2);
303 vdp->push_back(aGamma3);
304}
CLHEP::HepLorentzVector G4LorentzVector
CLHEP::Hep3Vector G4ThreeVector
virtual double flat()=0
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double ComputeFS(G4double fr1, G4double fr2, G4double fr3, G4double kinEnergy)

◆ SetDelta()

void G4eplusTo3GammaOKVIModel::SetDelta ( G4double val)
inline

Definition at line 96 of file G4eplusTo3GammaOKVIModel.hh.

96{ if(val > 0.0) { fDelta = val; } };

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