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

#include <G4eplusTo2GammaOKVIModel.hh>

+ Inheritance diagram for G4eplusTo2GammaOKVIModel:

Public Member Functions

 G4eplusTo2GammaOKVIModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="eplus2ggOKVI")
 
virtual ~G4eplusTo2GammaOKVIModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double ComputeCrossSectionPerElectron (G4double kinEnergy)
 
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=0.0, G4double maxEnergy=DBL_MAX) final
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX) final
 
void SetDelta (G4double val)
 
- 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 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 CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
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)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 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) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () 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)
 
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
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 56 of file G4eplusTo2GammaOKVIModel.hh.

Constructor & Destructor Documentation

◆ G4eplusTo2GammaOKVIModel()

G4eplusTo2GammaOKVIModel::G4eplusTo2GammaOKVIModel ( const G4ParticleDefinition p = nullptr,
const G4String nam = "eplus2ggOKVI" 
)
explicit

Definition at line 69 of file G4eplusTo2GammaOKVIModel.cc.

71 : G4VEmModel(nam),
72 fDelta(0.001),
73 fGammaTh(MeV)
74{
75 theGamma = G4Gamma::Gamma();
76 fParticleChange = nullptr;
77 fCuts = nullptr;
78 f3GModel = new G4eplusTo3GammaOKVIModel();
79 SetTripletModel(f3GModel);
80}
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void SetTripletModel(G4VEmModel *)
Definition: G4VEmModel.hh:635

◆ ~G4eplusTo2GammaOKVIModel()

G4eplusTo2GammaOKVIModel::~G4eplusTo2GammaOKVIModel ( )
virtual

Definition at line 84 of file G4eplusTo2GammaOKVIModel.cc.

85{}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 150 of file G4eplusTo2GammaOKVIModel.cc.

154{
155 // Calculates the cross section per atom of annihilation into two photons
156 G4double cross = Z*fCrossSection->Value(kineticEnergy);
157 return cross;
158}
double G4double
Definition: G4Types.hh:83
G4double Value(G4double theEnergy, std::size_t &lastidx) const

◆ ComputeCrossSectionPerElectron()

G4double G4eplusTo2GammaOKVIModel::ComputeCrossSectionPerElectron ( G4double  kinEnergy)

Definition at line 127 of file G4eplusTo2GammaOKVIModel.cc.

128{
129 // Calculates the cross section per electron of annihilation into two
130 // photons from the Heilter formula with the radiation correction to 3 gamma
131 // annihilation channel. (A.A.) rho is changed
132
133 G4double ekin = std::max(eV,kinEnergy);
134 G4double tau = ekin/electron_mass_c2;
135 G4double gam = tau + 1.0;
136 G4double gamma2 = gam*gam;
137 G4double bg2 = tau * (tau+2.0);
138 G4double bg = sqrt(bg2);
139 G4double rho = (gamma2+4.*gam+1.)*G4Log(gam+bg)/(gamma2-1.)
140 - (gam+3.)/(sqrt(gam*gam - 1.));
141
142 static const G4double pir2 = pi*classic_electr_radius*classic_electr_radius;
143 G4double cross = (pir2*rho + alpha_rcl2*2.*G4Log(fDelta)*rho*rho)/(gam+1.);
144
145 return cross;
146}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
const G4double pi

Referenced by Initialise().

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 162 of file G4eplusTo2GammaOKVIModel.cc.

167{
168 // Calculates the cross section per volume of annihilation into two photons
169 G4double eDensity = material->GetElectronDensity();
170 G4double cross = eDensity*fCrossSection->Value(kineticEnergy);
171 return cross;
172}
G4double GetElectronDensity() const
Definition: G4Material.hh:215

◆ Initialise()

void G4eplusTo2GammaOKVIModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
overridevirtual

Implements G4VEmModel.

Definition at line 89 of file G4eplusTo2GammaOKVIModel.cc.

91{
92 f3GModel->Initialise(p, cuts);
93 fCuts = &cuts;
95 f3GModel->SetDelta(fDelta);
96
97 if(IsMaster()) {
98 if(!fCrossSection) {
99 G4double emin = 10*eV;
100 G4double emax = 100*TeV;
101 G4int nbins = 20*G4lrint(std::log10(emax/emin));
102 fCrossSection = new G4PhysicsLogVector(emin, emax, nbins);
103 fCrossSection3G = new G4PhysicsLogVector(emin, emax, nbins);
104 f3GProbability = new G4PhysicsLogVector(emin, emax, nbins);
105 fCrossSection->SetSpline(true);
106 fCrossSection3G->SetSpline(true);
107 f3GProbability->SetSpline(true);
108 for(G4int i=0; i<= nbins; ++i) {
109 G4double e = fCrossSection->Energy(i);
111 G4double cs3 = f3GModel->ComputeCrossSectionPerElectron(e);
112 cs2 += cs3;
113 fCrossSection->PutValue(i, cs2);
114 fCrossSection3G->PutValue(i, cs3);
115 f3GProbability->PutValue(i, cs3/cs2);
116 }
117 }
118 }
119 // here particle change is set for the triplet model
120 if(fParticleChange) { return; }
121 fParticleChange = GetParticleChangeForGamma();
122}
int G4int
Definition: G4Types.hh:85
static G4EmParameters * Instance()
G4double LowestTripletEnergy() const
G4double Energy(std::size_t index) const
void PutValue(std::size_t index, G4double theValue)
void SetSpline(G4bool)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
G4double ComputeCrossSectionPerElectron(G4double kinEnergy)
G4double ComputeCrossSectionPerElectron(G4double kinEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) final
int G4lrint(double ad)
Definition: templates.hh:134

◆ SampleSecondaries()

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

!! likely problematic direction to be checked

Implements G4VEmModel.

Definition at line 180 of file G4eplusTo2GammaOKVIModel.cc.

184{
185 G4double posiKinEnergy = dp->GetKineticEnergy();
186 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
187
188 if(rndmEngine->flat() < f3GProbability->Value(posiKinEnergy)) {
189 G4double cutd = std::max(fGammaTh,(*fCuts)[mcc->GetIndex()])
190 /(posiKinEnergy + electron_mass_c2);
191 // check cut to avoid production of 3d gamma below
192 if(cutd > fDelta) {
193 G4double cs30 = fCrossSection3G->Value(posiKinEnergy);
194 f3GModel->SetDelta(cutd);
195 G4double cs3 = f3GModel->ComputeCrossSectionPerElectron(posiKinEnergy);
196 if(rndmEngine->flat()*cs30 < cs3) {
197 f3GModel->SampleSecondaries(vdp, mcc, dp);
198 return;
199 }
200 } else {
201 f3GModel->SampleSecondaries(vdp, mcc, dp);
202 return;
203 }
204 }
205
206 G4DynamicParticle *aGamma1, *aGamma2;
207
208 // Case at rest
209 if(posiKinEnergy == 0.0) {
210 G4double cost = 2.*rndmEngine->flat()-1.;
211 G4double sint = sqrt((1. - cost)*(1. + cost));
212 G4double phi = twopi * rndmEngine->flat();
213 G4ThreeVector dir(sint*cos(phi), sint*sin(phi), cost);
214 phi = twopi * rndmEngine->flat();
215 G4double cosphi = cos(phi);
216 G4double sinphi = sin(phi);
217 G4ThreeVector pol(cosphi, sinphi, 0.0);
218 pol.rotateUz(dir);
219 aGamma1 = new G4DynamicParticle(theGamma, dir, electron_mass_c2);
220 aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
221 aGamma2 = new G4DynamicParticle(theGamma,-dir, electron_mass_c2);
222 pol.set(-sinphi, cosphi, 0.0);
223 pol.rotateUz(dir);
224 aGamma2->SetPolarization(pol.x(),pol.y(),pol.z());
225
226 } else {
227
228 G4ThreeVector posiDirection = dp->GetMomentumDirection();
229
230 G4double tau = posiKinEnergy/electron_mass_c2;
231 G4double gam = tau + 1.0;
232 G4double tau2 = tau + 2.0;
233 G4double sqgrate = sqrt(tau/tau2)*0.5;
234 G4double sqg2m1 = sqrt(tau*tau2);
235
236 // limits of the energy sampling
237 G4double epsilmin = 0.5 - sqgrate;
238 G4double epsilmax = 0.5 + sqgrate;
239 G4double epsilqot = epsilmax/epsilmin;
240
241 //
242 // sample the energy rate of the created gammas
243 //
244 G4double epsil, greject;
245
246 do {
247 epsil = epsilmin*G4Exp(G4Log(epsilqot)*rndmEngine->flat());
248 greject = 1. - epsil + (2.*gam*epsil-1.)/(epsil*tau2*tau2);
249 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
250 } while( greject < rndmEngine->flat());
251
252 //
253 // scattered Gamma angles. ( Z - axis along the parent positron)
254 //
255
256 G4double cost = (epsil*tau2-1.)/(epsil*sqg2m1);
257 if(std::abs(cost) > 1.0) {
258 G4cout << "### G4eplusTo2GammaOKVIModel WARNING cost= " << cost
259 << " positron Ekin(MeV)= " << posiKinEnergy
260 << " gamma epsil= " << epsil
261 << G4endl;
262 if(cost > 1.0) cost = 1.0;
263 else cost = -1.0;
264 }
265 G4double sint = sqrt((1.+cost)*(1.-cost));
266 G4double phi = twopi * rndmEngine->flat();
267
268 //
269 // kinematic of the created pair
270 //
271
272 G4double TotalAvailableEnergy = posiKinEnergy + 2.0*electron_mass_c2;
273 G4double phot1Energy = epsil*TotalAvailableEnergy;
274
275 G4ThreeVector phot1Direction(sint*cos(phi), sint*sin(phi), cost);
276 phot1Direction.rotateUz(posiDirection);
277 aGamma1 = new G4DynamicParticle (theGamma,phot1Direction, phot1Energy);
278 phi = twopi * rndmEngine->flat();
279 G4double cosphi = cos(phi);
280 G4double sinphi = sin(phi);
281 G4ThreeVector pol(cosphi, sinphi, 0.0);
282 pol.rotateUz(phot1Direction);
283 aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
284
285 G4double phot2Energy =(1.-epsil)*TotalAvailableEnergy;
286 G4double posiP= sqrt(posiKinEnergy*(posiKinEnergy+2.*electron_mass_c2));
287 G4ThreeVector dir = posiDirection*posiP - phot1Direction*phot1Energy;
288 G4ThreeVector phot2Direction = dir.unit();
289
290 // create G4DynamicParticle object for the particle2
291 aGamma2 = new G4DynamicParticle (theGamma,phot2Direction, phot2Energy);
292
293 //!!! likely problematic direction to be checked
294 pol.set(-sinphi, cosphi, 0.0);
295 pol.rotateUz(phot1Direction);
296 cost = pol*phot2Direction;
297 pol -= cost*phot2Direction;
298 pol = pol.unit();
299 aGamma2->SetPolarization(pol.x(),pol.y(),pol.z());
300 }
301 /*
302 G4cout << "Annihilation in fly: e0= " << posiKinEnergy
303 << " m= " << electron_mass_c2
304 << " e1= " << phot1Energy
305 << " e2= " << phot2Energy << " dir= " << dir
306 << " -> " << phot1Direction << " "
307 << phot2Direction << G4endl;
308 */
309
310 vdp->push_back(aGamma1);
311 vdp->push_back(aGamma2);
312
313 // kill primary positron
314 fParticleChange->SetProposedKineticEnergy(0.0);
315 fParticleChange->ProposeTrackStatus(fStopAndKill);
316}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
@ fStopAndKill
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
virtual double flat()=0
void SetPolarization(const G4ThreeVector &)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeTrackStatus(G4TrackStatus status)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX) final
double flat()
Definition: G4AblaRandom.cc:48

◆ SetDelta()

void G4eplusTo2GammaOKVIModel::SetDelta ( G4double  val)
inline

Definition at line 91 of file G4eplusTo2GammaOKVIModel.hh.

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

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