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

#include <G4LivermoreNuclearGammaConversionModel.hh>

+ Inheritance diagram for G4LivermoreNuclearGammaConversionModel:

Public Member Functions

 G4LivermoreNuclearGammaConversionModel (const G4ParticleDefinition *p=0, const G4String &nam="LivermoreNuclearConversion")
 
virtual ~G4LivermoreNuclearGammaConversionModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double)
 
- 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 38 of file G4LivermoreNuclearGammaConversionModel.hh.

Constructor & Destructor Documentation

◆ G4LivermoreNuclearGammaConversionModel()

G4LivermoreNuclearGammaConversionModel::G4LivermoreNuclearGammaConversionModel ( const G4ParticleDefinition p = 0,
const G4String nam = "LivermoreNuclearConversion" 
)

Definition at line 46 of file G4LivermoreNuclearGammaConversionModel.cc.

48:G4VEmModel(nam),isInitialised(false),smallEnergy(2.*MeV)
49{
50 fParticleChange = 0;
51
52 lowEnergyLimit = 2.0*electron_mass_c2;
53
54 verboseLevel= 0;
55 // Verbosity scale for debugging purposes:
56 // 0 = nothing
57 // 1 = calculation of cross sections, file openings...
58 // 2 = entering in methods
59
60 if(verboseLevel > 0)
61 {
62 G4cout << "G4LivermoreNuclearGammaConversionModel is constructed " << G4endl;
63 }
64}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ ~G4LivermoreNuclearGammaConversionModel()

G4LivermoreNuclearGammaConversionModel::~G4LivermoreNuclearGammaConversionModel ( )
virtual

Definition at line 68 of file G4LivermoreNuclearGammaConversionModel.cc.

69{
70 if(IsMaster()) {
71 for(G4int i=0; i<maxZ; ++i) {
72 if(data[i]) {
73 delete data[i];
74 data[i] = 0;
75 }
76 }
77 }
78}
int G4int
Definition: G4Types.hh:85
G4bool IsMaster() const
Definition: G4VEmModel.hh:736

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4LivermoreNuclearGammaConversionModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 217 of file G4LivermoreNuclearGammaConversionModel.cc.

221{
222 if (verboseLevel > 1)
223 {
224 G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreNuclearGammaConversionModel"
225 << G4endl;
226 }
227
228 if (GammaEnergy < lowEnergyLimit) { return 0.0; }
229
230 G4double xs = 0.0;
231
232 G4int intZ=G4int(Z);
233
234 if(intZ < 1 || intZ > maxZ) { return xs; }
235
236 G4LPhysicsFreeVector* pv = data[intZ];
237
238 // if element was not initialised
239 // do initialisation safely for MT mode
240 if(!pv)
241 {
242 InitialiseForElement(0, intZ);
243 pv = data[intZ];
244 if(!pv) { return xs; }
245 }
246 // x-section is taken from the table
247 xs = pv->Value(GammaEnergy);
248
249 if(verboseLevel > 0)
250 {
251 G4int n = pv->GetVectorLength() - 1;
252 G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)="
253 << GammaEnergy/MeV << G4endl;
254 G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
255 G4cout << " -> first cs value in EADL data file (iu) =" << (*pv)[0] << G4endl;
256 G4cout << " -> last cs value in EADL data file (iu) =" << (*pv)[n] << G4endl;
257 G4cout << "*********************************************************" << G4endl;
258 }
259
260 return xs;
261
262}
double G4double
Definition: G4Types.hh:83
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
G4double Value(G4double theEnergy, std::size_t &lastidx) const
std::size_t GetVectorLength() const

◆ Initialise()

void G4LivermoreNuclearGammaConversionModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 82 of file G4LivermoreNuclearGammaConversionModel.cc.

85{
86
87 if (verboseLevel > 1)
88 {
89 G4cout << "Calling Initialise() of G4LivermoreNuclearGammaConversionModel."
90 << G4endl
91 << "Energy range: "
92 << LowEnergyLimit() / MeV << " MeV - "
93 << HighEnergyLimit() / GeV << " GeV"
94 << G4endl;
95 }
96
97 if(IsMaster())
98 {
99
100 // Initialise element selector
101
102 InitialiseElementSelectors(particle, cuts);
103
104 // Access to elements
105
106 char* path = std::getenv("G4LEDATA");
107
108 G4ProductionCutsTable* theCoupleTable =
110
111 G4int numOfCouples = theCoupleTable->GetTableSize();
112
113 for(G4int i=0; i<numOfCouples; ++i)
114 {
115 const G4Material* material =
116 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
117 const G4ElementVector* theElementVector = material->GetElementVector();
118 G4int nelm = material->GetNumberOfElements();
119
120 for (G4int j=0; j<nelm; ++j)
121 {
122 G4int Z = (G4int)(*theElementVector)[j]->GetZ();
123 if(Z < 1) { Z = 1; }
124 else if(Z > maxZ) { Z = maxZ; }
125 if(!data[Z]) { ReadData(Z, path); }
126 }
127 }
128 }
129 if(isInitialised) { return; }
130 fParticleChange = GetParticleChangeForGamma();
131 isInitialised = true;
132}
std::vector< G4Element * > G4ElementVector
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148

◆ InitialiseForElement()

void G4LivermoreNuclearGammaConversionModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 476 of file G4LivermoreNuclearGammaConversionModel.cc.

479{
480 G4AutoLock l(&LivermoreNuclearGammaConversionModelMutex);
481 // G4cout << "G4LivermoreNuclearGammaConversionModel::InitialiseForElement Z= "
482 // << Z << G4endl;
483 if(!data[Z]) { ReadData(Z); }
484 l.unlock();
485}

Referenced by ComputeCrossSectionPerAtom().

◆ InitialiseLocal()

void G4LivermoreNuclearGammaConversionModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 136 of file G4LivermoreNuclearGammaConversionModel.cc.

138{
140}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:842
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:834

◆ MinPrimaryEnergy()

G4double G4LivermoreNuclearGammaConversionModel::MinPrimaryEnergy ( const G4Material ,
const G4ParticleDefinition ,
G4double   
)
virtual

Reimplemented from G4VEmModel.

Definition at line 145 of file G4LivermoreNuclearGammaConversionModel.cc.

148{
149 return lowEnergyLimit;
150}

◆ SampleSecondaries()

void G4LivermoreNuclearGammaConversionModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 266 of file G4LivermoreNuclearGammaConversionModel.cc.

271{
272
273// The energies of the e+ e- secondaries are sampled using the Bethe - Heitler
274// cross sections with Coulomb correction. A modified version of the random
275// number techniques of Butcher & Messel is used (Nuc Phys 20(1960),15).
276
277// Note 1 : Effects due to the breakdown of the Born approximation at low
278// energy are ignored.
279// Note 2 : The differential cross section implicitly takes account of
280// pair creation in both nuclear and atomic electron fields. However triplet
281// prodution is not generated.
282
283 if (verboseLevel > 1) {
284 G4cout << "Calling SampleSecondaries() of G4LivermoreNuclearGammaConversionModel"
285 << G4endl;
286 }
287
288 G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
289 G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
290
292 G4double epsilon0Local = electron_mass_c2 / photonEnergy ;
293
294 // Do it fast if photon energy < 2. MeV
295 if (photonEnergy < smallEnergy )
296 {
297 epsilon = epsilon0Local + (0.5 - epsilon0Local) * G4UniformRand();
298 }
299 else
300 {
301 // Select randomly one element in the current material
302
303 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
304 const G4Element* element = SelectRandomAtom(couple,particle,photonEnergy);
305
306 if (element == 0)
307 {
308 G4cout << "G4LivermoreNuclearGammaConversionModel::SampleSecondaries - element = 0"
309 << G4endl;
310 return;
311 }
312 G4IonisParamElm* ionisation = element->GetIonisation();
313 if (ionisation == 0)
314 {
315 G4cout << "G4LivermoreNuclearGammaConversionModel::SampleSecondaries - ionisation = 0"
316 << G4endl;
317 return;
318 }
319
320 // Extract Coulomb factor for this Elements
321 G4double fZ = 8. * (ionisation->GetlogZ3());
322 if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
323
324 // Limits of the screening variable
325 G4double screenFactor = 136. * epsilon0Local / (element->GetIonisation()->GetZ3()) ;
326 G4double screenMax = G4Exp ((42.24 - fZ)/8.368) - 0.952 ;
327 G4double screenMin = std::min(4.*screenFactor,screenMax) ;
328
329 // Limits of the energy sampling
330 G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
331 G4double epsilonMin = std::max(epsilon0Local,epsilon1);
332 G4double epsilonRange = 0.5 - epsilonMin ;
333
334 // Sample the energy rate of the created electron (or positron)
335 G4double screen;
336 G4double gReject ;
337
338 G4double f10 = ScreenFunction1(screenMin) - fZ;
339 G4double f20 = ScreenFunction2(screenMin) - fZ;
340 G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
341 G4double normF2 = std::max(1.5 * f20,0.);
342
343 do
344 {
345 if (normF1 / (normF1 + normF2) > G4UniformRand() )
346 {
347 epsilon = 0.5 - epsilonRange * std::pow(G4UniformRand(), 0.333333) ;
348 screen = screenFactor / (epsilon * (1. - epsilon));
349 gReject = (ScreenFunction1(screen) - fZ) / f10 ;
350 }
351 else
352 {
353 epsilon = epsilonMin + epsilonRange * G4UniformRand();
354 screen = screenFactor / (epsilon * (1 - epsilon));
355 gReject = (ScreenFunction2(screen) - fZ) / f20 ;
356 }
357 } while ( gReject < G4UniformRand() );
358
359 } // End of epsilon sampling
360
361 // Fix charges randomly
362
363 G4double electronTotEnergy;
364 G4double positronTotEnergy;
365
366 if (G4UniformRand() > 0.5)
367 {
368 electronTotEnergy = (1. - epsilon) * photonEnergy;
369 positronTotEnergy = epsilon * photonEnergy;
370 }
371 else
372 {
373 positronTotEnergy = (1. - epsilon) * photonEnergy;
374 electronTotEnergy = epsilon * photonEnergy;
375 }
376
377 // Scattered electron (positron) angles. ( Z - axis along the parent photon)
378 // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
379 // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
380
381 G4double u;
382 const G4double a1 = 0.625;
383 G4double a2 = 3. * a1;
384 // G4double d = 27. ;
385
386 // if (9. / (9. + d) > G4UniformRand())
387 if (0.25 > G4UniformRand())
388 {
389 u = - G4Log(G4UniformRand() * G4UniformRand()) / a1 ;
390 }
391 else
392 {
393 u = - G4Log(G4UniformRand() * G4UniformRand()) / a2 ;
394 }
395
396 G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
397 G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
398 G4double phi = twopi * G4UniformRand();
399
400 G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
401 G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
402
403
404 // Kinematics of the created pair:
405 // the electron and positron are assumed to have a symetric angular
406 // distribution with respect to the Z axis along the parent photon
407
408 G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
409
410 G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
411 electronDirection.rotateUz(photonDirection);
412
414 electronDirection,
415 electronKineEnergy);
416
417 // The e+ is always created
418 G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
419
420 G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
421 positronDirection.rotateUz(photonDirection);
422
423 // Create G4DynamicParticle object for the particle2
425 positronDirection,
426 positronKineEnergy);
427 // Fill output vector
428 fvect->push_back(particle1);
429 fvect->push_back(particle2);
430
431 // kill incident photon
432 fParticleChange->SetProposedKineticEnergy(0.);
433 fParticleChange->ProposeTrackStatus(fStopAndKill);
434
435}
double epsilon(double density, double temperature)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4double GetfCoulomb() const
Definition: G4Element.hh:190
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:198
G4double GetlogZ3() const
G4double GetZ3() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static G4Positron * Positron()
Definition: G4Positron.cc:93
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:570
void ProposeTrackStatus(G4TrackStatus status)

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