Geant4 11.2.2
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 &) final
 
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)
 
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 *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (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 53 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 fParticleChange = nullptr;
66}
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~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 202 of file G4eplusTo3GammaOKVIModel.cc.

206{
207 // Calculates the cross section per atom of annihilation into two photons
208
209
210
211 G4double cross = Z*ComputeCrossSectionPerElectron(kineticEnergy);
212 return cross;
213}
double G4double
Definition G4Types.hh:83
G4double ComputeCrossSectionPerElectron(G4double kinEnergy)

◆ ComputeCrossSectionPerElectron()

G4double G4eplusTo3GammaOKVIModel::ComputeCrossSectionPerElectron ( G4double kinEnergy)

Definition at line 181 of file G4eplusTo3GammaOKVIModel.cc.

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

Referenced by ComputeCrossSectionPerAtom(), CrossSectionPerVolume(), G4eplusTo2GammaOKVIModel::Initialise(), and G4eplusTo2GammaOKVIModel::SampleSecondaries().

◆ ComputeF()

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

Definition at line 85 of file G4eplusTo3GammaOKVIModel.cc.

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

Referenced by ComputeFS().

◆ ComputeF0()

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

Definition at line 129 of file G4eplusTo3GammaOKVIModel.cc.

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

◆ ComputeFS()

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

Definition at line 162 of file G4eplusTo3GammaOKVIModel.cc.

164{
165 G4double ekin = std::max(eV,kinEnergy);
166 G4double tau = ekin/electron_mass_c2;
167 G4double gam = tau + 1.0;
168
169 G4double fsum = fr1*fr1*(ComputeF(fr1,fr2,fr3,ekin) +
170 ComputeF(fr3,fr1,fr2,ekin) +
171 ComputeF(fr2,fr3,fr1,ekin));
172
173 G4double dcross = fsum/((3*fr1*fr1*(gam+1.)));
174
175 return dcross;
176}
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 217 of file G4eplusTo3GammaOKVIModel.cc.

222{
223 // Calculates the cross section per volume of annihilation into two photons
224
225 G4double eDensity = material->GetElectronDensity();
226 G4double cross = eDensity*ComputeCrossSectionPerElectron(kineticEnergy);
227 return cross;
228}
G4double GetElectronDensity() const

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 74 of file G4eplusTo3GammaOKVIModel.cc.

76{
77 // here particle change is set for the triplet model
78 if(fParticleChange) { return; }
79 fParticleChange = GetParticleChangeForGamma();
80}
G4ParticleChangeForGamma * GetParticleChangeForGamma()

Referenced by G4eplusTo2GammaOKVIModel::Initialise().

◆ 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

Implements G4VEmModel.

Definition at line 236 of file G4eplusTo3GammaOKVIModel.cc.

240{
241
242 G4double posiKinEnergy = dp->GetKineticEnergy();
243 G4DynamicParticle *aGamma1, *aGamma2;
244 G4DynamicParticle* aGamma3 = nullptr;
245 G4double border;
246
247 if(posiKinEnergy < 500*MeV) {
248 border = 1. - (electron_mass_c2)/(2*(posiKinEnergy + electron_mass_c2));
249 } else {
250 border = 1. - (100*electron_mass_c2)/(2*(posiKinEnergy + electron_mass_c2));
251 }
252 border = std::min(border, 0.9999);
253
254 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
255
256 // Case at rest
257 if(posiKinEnergy == 0.0) {
258 G4double cost = 2.*rndmEngine->flat()-1.;
259 G4double sint = sqrt((1. - cost)*(1. + cost));
260 G4double phi = twopi * rndmEngine->flat();
261 G4ThreeVector dir(sint*cos(phi), sint*sin(phi), cost);
262 phi = twopi * rndmEngine->flat();
263 G4double cosphi = cos(phi);
264 G4double sinphi = sin(phi);
265 G4ThreeVector pol(cosphi, sinphi, 0.0);
266 pol.rotateUz(dir);
267 aGamma1 = new G4DynamicParticle(theGamma, dir, electron_mass_c2);
268 aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
269 aGamma2 = new G4DynamicParticle(theGamma,-dir, electron_mass_c2);
270 pol.set(-sinphi, cosphi, 0.0);
271 pol.rotateUz(dir);
272 aGamma2->SetPolarization(pol.x(),pol.y(),pol.z());
273
274 } else {
275
276 G4ThreeVector posiDirection = dp->GetMomentumDirection();
277
278 // (A.A.) LIMITS FOR 1st GAMMA
279 G4double xmin = 0.01;
280 G4double xmax = 0.667; // CHANGE to 3/2
281
282 G4double d1, d0, x1, x2, dmax, x2min;
283
284 // (A.A.) sampling of x1 x2 x3 (whole cycle of rejection)
285 do {
286 x1 = 1/((1/xmin) - ((1/xmin)-(1/xmax))*rndmEngine->flat());
287 dmax = ComputeFS(posiKinEnergy, x1,1.-x1,border);
288 x2min = 1.-x1;
289 x2 = 1 - rndmEngine->flat()*(1-x2min);
290 d1 = dmax*rndmEngine->flat();
291 d0 = ComputeFS(posiKinEnergy,x1,x2,2-x1-x2);
292 }
293 while(d0 < d1);
294
295 G4double x3 = 2 - x1 - x2;
296
297 //
298 // angles between Gammas
299 //
300
301 G4double psi13 = 2*asin(sqrt(std::abs((x1+x3-1)/(x1*x3))));
302 G4double psi12 = 2*asin(sqrt(std::abs((x1+x2-1)/(x1*x2))));
303
304 // sin^t
305
306 //G4double phi = twopi * rndmEngine->flat();
307 //G4double psi = acos(x3); // Angle of the plane
308
309 //
310 // kinematic of the created pair
311 //
312
313 G4double TotalAvailableEnergy = posiKinEnergy + 2.0*electron_mass_c2;
314
315 G4double phot1Energy = 0.5*x1*TotalAvailableEnergy;
316 G4double phot2Energy = 0.5*x2*TotalAvailableEnergy;
317 G4double phot3Energy = 0.5*x3*TotalAvailableEnergy;
318
319
320 // DIRECTIONS
321
322 // The azimuthal angles of ql and q3 with respect to some plane
323 // through the beam axis are generated at random.
324
325 G4ThreeVector phot1Direction(0, 0, 1);
326 G4ThreeVector phot2Direction(0, sin(psi12), cos(psi12));
327 G4ThreeVector phot3Direction(0, sin(psi13), cos(psi13));
328
329 phot1Direction.rotateUz(posiDirection);
330 phot2Direction.rotateUz(posiDirection);
331 phot3Direction.rotateUz(posiDirection);
332
333 aGamma1 = new G4DynamicParticle (theGamma,phot1Direction, phot1Energy);
334 aGamma2 = new G4DynamicParticle (theGamma,phot2Direction, phot2Energy);
335 aGamma3 = new G4DynamicParticle (theGamma,phot3Direction, phot3Energy);
336
337
338 //POLARIZATION - ???
339 /*
340
341
342 phi = twopi * rndmEngine->flat();
343 G4double cosphi = cos(phi);
344 G4double sinphi = sin(phi);
345 G4ThreeVector pol(cosphi, sinphi, 0.0);
346 pol.rotateUz(phot1Direction);
347 aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
348
349 G4double phot2Energy =(1.-epsil)*TotalAvailableEnergy;
350 G4double posiP= sqrt(posiKinEnergy*(posiKinEnergy+2.*electron_mass_c2));
351 G4ThreeVector dir = posiDirection*posiP - phot1Direction*phot1Energy;
352 G4ThreeVector phot2Direction = dir.unit();
353
354 // create G4DynamicParticle object for the particle2
355 aGamma2 = new G4DynamicParticle (theGamma,phot2Direction, phot2Energy);
356
357 //!!! likely problematic direction to be checked
358 pol.set(-sinphi, cosphi, 0.0);
359 pol.rotateUz(phot1Direction);
360 cost = pol*phot2Direction;
361 pol -= cost*phot2Direction;
362 pol = pol.unit();
363 aGamma2->SetPolarization(pol.x(),pol.y(),pol.z());
364
365 */
366
367 }
368 /*
369 G4cout << "Annihilation in fly: e0= " << posiKinEnergy
370 << " m= " << electron_mass_c2
371 << " e1= " << phot1Energy
372 << " e2= " << phot2Energy << " dir= " << dir
373 << " -> " << phot1Direction << " "
374 << phot2Direction << G4endl;
375 */
376
377 vdp->push_back(aGamma1);
378 vdp->push_back(aGamma2);
379 if(aGamma3 != nullptr) { vdp->push_back(aGamma3); }
380
381 // kill primary positron
382 fParticleChange->SetProposedKineticEnergy(0.0);
383 fParticleChange->ProposeTrackStatus(fStopAndKill);
384}
@ fStopAndKill
virtual double flat()=0
void SetPolarization(const G4ThreeVector &)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeTrackStatus(G4TrackStatus status)
G4double ComputeFS(G4double fr1, G4double fr2, G4double fr3, G4double kinEnergy)

Referenced by G4eplusTo2GammaOKVIModel::SampleSecondaries().

◆ SetDelta()

void G4eplusTo3GammaOKVIModel::SetDelta ( G4double val)
inline

Definition at line 93 of file G4eplusTo3GammaOKVIModel.hh.

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

Referenced by G4eplusTo2GammaOKVIModel::Initialise(), and G4eplusTo2GammaOKVIModel::SampleSecondaries().


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