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

#include <G4eeToTwoGammaModel.hh>

+ Inheritance diagram for G4eeToTwoGammaModel:

Public Member Functions

 G4eeToTwoGammaModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="eplus2gg")
 
 ~G4eeToTwoGammaModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual G4double ComputeCrossSectionPerElectron (G4double kinEnergy)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0., G4double maxEnergy=DBL_MAX) override
 
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4eeToTwoGammaModeloperator= (const G4eeToTwoGammaModel &right)=delete
 
 G4eeToTwoGammaModel (const G4eeToTwoGammaModel &)=delete
 
- 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 *, 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 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 *)
 
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 = 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
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 59 of file G4eeToTwoGammaModel.hh.

Constructor & Destructor Documentation

◆ G4eeToTwoGammaModel() [1/2]

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

Definition at line 90 of file G4eeToTwoGammaModel.cc.

92 : G4VEmModel(nam),
93 pi_rcl2(pi*classic_electr_radius*classic_electr_radius)
94{
95 theGamma = G4Gamma::Gamma();
96 fParticleChange = nullptr;
97}
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85

◆ ~G4eeToTwoGammaModel()

G4eeToTwoGammaModel::~G4eeToTwoGammaModel ( )
overridedefault

◆ G4eeToTwoGammaModel() [2/2]

G4eeToTwoGammaModel::G4eeToTwoGammaModel ( const G4eeToTwoGammaModel )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 155 of file G4eeToTwoGammaModel.cc.

159{
160 // Calculates the cross section per atom of annihilation into two photons
161 return Z*ComputeCrossSectionPerElectron(kineticEnergy);
162}
const G4int Z[17]
virtual G4double ComputeCrossSectionPerElectron(G4double kinEnergy)

◆ ComputeCrossSectionPerElectron()

G4double G4eeToTwoGammaModel::ComputeCrossSectionPerElectron ( G4double  kinEnergy)
virtual

Reimplemented in G4PolarizedAnnihilationModel.

Definition at line 135 of file G4eeToTwoGammaModel.cc.

136{
137 // Calculates the cross section per electron of annihilation into two photons
138 // from the Heilter formula.
139
140 G4double ekin = std::max(eV,kineticEnergy);
141
142 G4double tau = ekin/electron_mass_c2;
143 G4double gam = tau + 1.0;
144 G4double gamma2= gam*gam;
145 G4double bg2 = tau * (tau+2.0);
146 G4double bg = sqrt(bg2);
147
148 G4double cross = pi_rcl2*((gamma2+4*gam+1.)*G4Log(gam+bg) - (gam+3.)*bg)
149 / (bg2*(gam+1.));
150 return cross;
151}
G4double G4Log(G4double x)
Definition: G4Log.hh:227
double G4double
Definition: G4Types.hh:83

Referenced by ComputeCrossSectionPerAtom(), G4PolarizedAnnihilationModel::ComputeCrossSectionPerElectron(), and CrossSectionPerVolume().

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 166 of file G4eeToTwoGammaModel.cc.

171{
172 // Calculates the cross section per volume of annihilation into two photons
173 return material->GetElectronDensity()*ComputeCrossSectionPerElectron(kineticEnergy);
174}
G4double GetElectronDensity() const
Definition: G4Material.hh:212

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 105 of file G4eeToTwoGammaModel.cc.

107{
108 if(IsMaster()) {
110 // redo initialisation for each new run
111 fSampleAtomicPDF = false;
112 const auto& materialTable = G4Material::GetMaterialTable();
113 for (const auto& material: *materialTable) {
114 const G4double meanEnergyPerIonPair = material->GetIonisation()->GetMeanEnergyPerIonPair();
115 if (meanEnergyPerIonPair > 0.) {
116 fSampleAtomicPDF = true;
117 if(verbose > 0) {
118 G4cout << "### G4eeToTwoGammaModel: for " << material->GetName() << " mean energy per ion pair is "
119 << meanEnergyPerIonPair/CLHEP::eV << " eV" << G4endl;
120 }
121 }
122 }
123 }
124 // If no materials have meanEnergyPerIonPair set. This is probably the usual
125 // case, since most applications are not senstive to the slight
126 // non-collinearity of gammas in eeToTwoGamma. Do not issue any warning.
127
128 if(fParticleChange) { return; }
129 fParticleChange = GetParticleChangeForGamma();
130}
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4EmParameters * Instance()
G4int Verbose() const
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4bool IsMaster() const
Definition: G4VEmModel.hh:725

Referenced by G4PolarizedAnnihilationModel::Initialise().

◆ operator=()

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

◆ SampleSecondaries()

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

!! likely problematic direction to be checked

Implements G4VEmModel.

Definition at line 181 of file G4eeToTwoGammaModel.cc.

186{
187 G4double posiKinEnergy = dp->GetKineticEnergy();
188 G4DynamicParticle *aGamma1, *aGamma2;
189
190 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
191
192 // Case at rest
193 if(posiKinEnergy == 0.0) {
194
195 const G4double eGamma = electron_mass_c2;
196
197 // In rest frame of positronium gammas are back to back
198 const G4ThreeVector& dir1 = G4RandomDirection();
199 const G4ThreeVector& dir2 = -dir1;
200 aGamma1 = new G4DynamicParticle(G4Gamma::Gamma(),dir1,eGamma);
201 aGamma2 = new G4DynamicParticle(G4Gamma::Gamma(),dir2,eGamma);
202
203 // In rest frame the gammas are polarised perpendicular to each other - see
204 // Pryce and Ward, Nature No 4065 (1947) p.435.
205 // Snyder et al, Physical Review 73 (1948) p.440.
206 G4ThreeVector pol1 = (G4RandomDirection().cross(dir1)).unit();
207 G4ThreeVector pol2 = (pol1.cross(dir2)).unit();
208
209 // But the positronium is moving...
210 // A positron in matter slows down and combines with an atomic electron to
211 // make a neutral “atom” called positronium, about half the size of a normal
212 // atom. I expect that when the energy of the positron is small enough,
213 // less than the binding energy of positronium (6.8 eV), it is
214 // energetically favourable for an electron from the outer orbitals of a
215 // nearby atom or molecule to transfer and bind to the positron, as in an
216 // ionic bond, leaving behind a mildly ionised nearby atom/molecule. I
217 // would expect the positronium to come away with a kinetic energy of a
218 // few eV on average. In its para (spin 0) state it annihilates into two
219 // photons, which in the rest frame of the positronium are collinear
220 // (back-to-back) due to momentum conservation. Because of the motion of the
221 // positronium, photons will be not quite back-to-back in the laboratory.
222
223 // The positroniuim acquires an energy of order its binding energy and
224 // doesn't have time to thermalise. Nevertheless, here we approximate its
225 // energy distribution by a Maxwell-Boltzman with mean energy <KE>. In terms
226 // of a more familiar concept of temperature, and the law of equipartition
227 // of energy of translational motion, <KE>=3kT/2. Each component of velocity
228 // has a distribution exp(-mv^2/2kT), which is a Gaussian of mean zero
229 // and variance kT/m=2<KE>/3m, where m is the positronium mass.
230
231 // We take <KE> = material->GetIonisation()->GetMeanEnergyPerIonPair().
232
233 if(fSampleAtomicPDF) {
234 const G4Material* material = pCutsCouple->GetMaterial();
235 const G4double meanEnergyPerIonPair = material->GetIonisation()->GetMeanEnergyPerIonPair();
236 const G4double& meanKE = meanEnergyPerIonPair; // Just an alias
237 if (meanKE > 0.) { // Positronium haas motion
238 // Mass of positronium
239 const G4double mass = 2.*electron_mass_c2;
240 // Mean <KE>=3kT/2, as described above
241 // const G4double T = 2.*meanKE/(3.*k_Boltzmann);
242 // Component velocities: Gaussian, variance kT/m=2<KE>/3m.
243 const G4double sigmav = std::sqrt(2.*meanKE/(3.*mass));
244 // This is in units where c=1
245 const G4double vx = G4RandGauss::shoot(0.,sigmav);
246 const G4double vy = G4RandGauss::shoot(0.,sigmav);
247 const G4double vz = G4RandGauss::shoot(0.,sigmav);
248 const G4ThreeVector v(vx,vy,vz); // In unit where c=1
249 const G4ThreeVector& beta = v; // so beta=v/c=v
250
251 aGamma1->Set4Momentum(aGamma1->Get4Momentum().boost(beta));
252 aGamma2->Set4Momentum(aGamma2->Get4Momentum().boost(beta));
253
254 // Rotate polarisation vectors
255 const G4ThreeVector& newDir1 = aGamma1->GetMomentumDirection();
256 const G4ThreeVector& newDir2 = aGamma2->GetMomentumDirection();
257 const G4ThreeVector& axis1 = dir1.cross(newDir1); // No need to be unit
258 const G4ThreeVector& axis2 = dir2.cross(newDir2); // No need to be unit
259 const G4double& angle1 = std::acos(dir1*newDir1);
260 const G4double& angle2 = std::acos(dir2*newDir2);
261 if (axis1 != G4ThreeVector()) pol1.rotate(axis1,angle1);
262 if (axis2 != G4ThreeVector()) pol2.rotate(axis2,angle2);
263 }
264 }
265 aGamma1->SetPolarization(pol1.x(),pol1.y(),pol1.z());
266 aGamma2->SetPolarization(pol2.x(),pol2.y(),pol2.z());
267
268 } else { // Positron interacts in flight
269
270 G4ThreeVector posiDirection = dp->GetMomentumDirection();
271
272 G4double tau = posiKinEnergy/electron_mass_c2;
273 G4double gam = tau + 1.0;
274 G4double tau2 = tau + 2.0;
275 G4double sqgrate = sqrt(tau/tau2)*0.5;
276 G4double sqg2m1 = sqrt(tau*tau2);
277
278 // limits of the energy sampling
279 G4double epsilmin = 0.5 - sqgrate;
280 G4double epsilmax = 0.5 + sqgrate;
281 G4double epsilqot = epsilmax/epsilmin;
282
283 //
284 // sample the energy rate of the created gammas
285 //
286 G4double epsil, greject;
287
288 do {
289 epsil = epsilmin*G4Exp(G4Log(epsilqot)*rndmEngine->flat());
290 greject = 1. - epsil + (2.*gam*epsil-1.)/(epsil*tau2*tau2);
291 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
292 } while( greject < rndmEngine->flat());
293
294 //
295 // scattered Gamma angles. ( Z - axis along the parent positron)
296 //
297
298 G4double cost = (epsil*tau2-1.)/(epsil*sqg2m1);
299 if(std::abs(cost) > 1.0) {
300 G4cout << "### G4eeToTwoGammaModel WARNING cost= " << cost
301 << " positron Ekin(MeV)= " << posiKinEnergy
302 << " gamma epsil= " << epsil
303 << G4endl;
304 if(cost > 1.0) cost = 1.0;
305 else cost = -1.0;
306 }
307 G4double sint = sqrt((1.+cost)*(1.-cost));
308 G4double phi = twopi * rndmEngine->flat();
309
310 //
311 // kinematic of the created pair
312 //
313
314 G4double totalEnergy = posiKinEnergy + 2.0*electron_mass_c2;
315 G4double phot1Energy = epsil*totalEnergy;
316
317 G4ThreeVector phot1Direction(sint*cos(phi), sint*sin(phi), cost);
318 phot1Direction.rotateUz(posiDirection);
319 aGamma1 = new G4DynamicParticle (theGamma,phot1Direction, phot1Energy);
320 phi = twopi * rndmEngine->flat();
321 G4double cosphi = cos(phi);
322 G4double sinphi = sin(phi);
323 G4ThreeVector pol(cosphi, sinphi, 0.0);
324 pol.rotateUz(phot1Direction);
325 aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
326
327 G4double phot2Energy =(1.-epsil)*totalEnergy;
328 G4double posiP= sqrt(posiKinEnergy*(posiKinEnergy+2.*electron_mass_c2));
329 G4ThreeVector dir = posiDirection*posiP - phot1Direction*phot1Energy;
330 G4ThreeVector phot2Direction = dir.unit();
331
332 // create G4DynamicParticle object for the particle2
333 aGamma2 = new G4DynamicParticle (theGamma, phot2Direction, phot2Energy);
334
335 //!!! likely problematic direction to be checked
336 pol.set(-sinphi, cosphi, 0.0);
337 pol.rotateUz(phot1Direction);
338 cost = pol*phot2Direction;
339 pol -= cost*phot2Direction;
340 pol = pol.unit();
341 aGamma2->SetPolarization(pol.x(),pol.y(),pol.z());
342 /*
343 G4cout << "Annihilation on fly: e0= " << posiKinEnergy
344 << " m= " << electron_mass_c2
345 << " e1= " << phot1Energy
346 << " e2= " << phot2Energy << " dir= " << dir
347 << " -> " << phot1Direction << " "
348 << phot2Direction << G4endl;
349 */
350 }
351
352 vdp->push_back(aGamma1);
353 vdp->push_back(aGamma2);
354
355 // kill primary positron
356 fParticleChange->SetProposedKineticEnergy(0.0);
357 fParticleChange->ProposeTrackStatus(fStopAndKill);
358}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4ThreeVector G4RandomDirection()
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
double z() const
Hep3Vector unit() const
double x() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:24
HepLorentzVector & boost(double, double, double)
virtual double flat()=0
void SetPolarization(const G4ThreeVector &)
const G4ThreeVector & GetMomentumDirection() const
G4LorentzVector Get4Momentum() const
G4double GetKineticEnergy() const
void Set4Momentum(const G4LorentzVector &momentum)
G4double GetMeanEnergyPerIonPair() const
const G4Material * GetMaterial() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:221
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeTrackStatus(G4TrackStatus status)
double flat()
Definition: G4AblaRandom.cc:48

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