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

#include <G4MollerBhabhaModel.hh>

+ Inheritance diagram for G4MollerBhabhaModel:

Public Member Functions

 G4MollerBhabhaModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="MollerBhabha")
 
virtual ~G4MollerBhabhaModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
- 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
 

Protected Member Functions

virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) final
 
void SetParticle (const G4ParticleDefinition *p)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

const G4ParticleDefinitionparticle
 
G4ParticleDefinitiontheElectron
 
G4ParticleChangeForLossfParticleChange
 
G4bool isElectron
 
G4double twoln10
 
G4double lowLimit
 
- 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 62 of file G4MollerBhabhaModel.hh.

Constructor & Destructor Documentation

◆ G4MollerBhabhaModel()

G4MollerBhabhaModel::G4MollerBhabhaModel ( const G4ParticleDefinition p = nullptr,
const G4String nam = "MollerBhabha" 
)
explicit

Definition at line 75 of file G4MollerBhabhaModel.cc.

77 : G4VEmModel(nam),
78 particle(nullptr),
79 isElectron(true),
80 twoln10(2.0*G4Log(10.0)),
81 lowLimit(0.02*keV),
82 isInitialised(false)
83{
85 if(p) { SetParticle(p); }
86 fParticleChange = nullptr;
87}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static G4Electron * Electron()
Definition: G4Electron.cc:93
const G4ParticleDefinition * particle
void SetParticle(const G4ParticleDefinition *p)
G4ParticleChangeForLoss * fParticleChange
G4ParticleDefinition * theElectron

◆ ~G4MollerBhabhaModel()

G4MollerBhabhaModel::~G4MollerBhabhaModel ( )
virtual

Definition at line 91 of file G4MollerBhabhaModel.cc.

92{}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4MollerBhabhaModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  Z,
G4double  A,
G4double  cutEnergy,
G4double  maxEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 177 of file G4MollerBhabhaModel.cc.

183{
184 return Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
185}
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

◆ ComputeCrossSectionPerElectron()

G4double G4MollerBhabhaModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Reimplemented in G4PolarizedMollerBhabhaModel.

Definition at line 123 of file G4MollerBhabhaModel.cc.

127{
128 if(!particle) { SetParticle(p); }
129
130 G4double cross = 0.0;
131 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
132 tmax = std::min(maxEnergy, tmax);
133 //G4cout << "E= " << kineticEnergy << " cut= " << cutEnergy
134 // << " Emax= " << tmax << G4endl;
135 if(cutEnergy < tmax) {
136
137 G4double xmin = cutEnergy/kineticEnergy;
138 G4double xmax = tmax/kineticEnergy;
139 G4double tau = kineticEnergy/electron_mass_c2;
140 G4double gam = tau + 1.0;
141 G4double gamma2= gam*gam;
142 G4double beta2 = tau*(tau + 2)/gamma2;
143
144 //Moller (e-e-) scattering
145 if (isElectron) {
146
147 G4double gg = (2.0*gam - 1.0)/gamma2;
148 cross = ((xmax - xmin)*(1.0 - gg + 1.0/(xmin*xmax)
149 + 1.0/((1.0-xmin)*(1.0 - xmax)))
150 - gg*G4Log( xmax*(1.0 - xmin)/(xmin*(1.0 - xmax)) ) ) / beta2;
151
152 //Bhabha (e+e-) scattering
153 } else {
154
155 G4double y = 1.0/(1.0 + gam);
156 G4double y2 = y*y;
157 G4double y12 = 1.0 - 2.0*y;
158 G4double b1 = 2.0 - y2;
159 G4double b2 = y12*(3.0 + y2);
160 G4double y122= y12*y12;
161 G4double b4 = y122*y12;
162 G4double b3 = b4 + y122;
163
164 cross = (xmax - xmin)*(1.0/(beta2*xmin*xmax) + b2
165 - 0.5*b3*(xmin + xmax)
166 + b4*(xmin*xmin + xmin*xmax + xmax*xmax)/3.0)
167 - b1*G4Log(xmax/xmin);
168 }
169
170 cross *= twopi_mc2_rcl2/kineticEnergy;
171 }
172 return cross;
173}
double G4double
Definition: G4Types.hh:83
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final

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

◆ ComputeDEDXPerVolume()

G4double G4MollerBhabhaModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 202 of file G4MollerBhabhaModel.cc.

207{
208 if(nullptr == particle) { SetParticle(p); }
209 // calculate the dE/dx due to the ionization by Seltzer-Berger formula
210 // checl low-energy limit
211 G4double electronDensity = material->GetElectronDensity();
212
213 G4double Zeff = material->GetIonisation()->GetZeffective();
214 G4double th = 0.25*sqrt(Zeff)*keV;
215 G4double tkin = std::max(kineticEnergy, th);
216
217 G4double tau = tkin/electron_mass_c2;
218 G4double gam = tau + 1.0;
219 G4double gamma2= gam*gam;
220 G4double bg2 = tau*(tau + 2);
221 G4double beta2 = bg2/gamma2;
222
223 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
224 eexc /= electron_mass_c2;
225 G4double eexc2 = eexc*eexc;
226
227 G4double d = std::min(cut, MaxSecondaryEnergy(p, tkin))/electron_mass_c2;
228 G4double dedx;
229
230 // electron
231 if (isElectron) {
232
233 dedx = G4Log(2.0*(tau + 2.0)/eexc2) - 1.0 - beta2
234 + G4Log((tau-d)*d) + tau/(tau-d)
235 + (0.5*d*d + (2.0*tau + 1.)*G4Log(1. - d/tau))/gamma2;
236
237 //positron
238 } else {
239
240 G4double d2 = d*d*0.5;
241 G4double d3 = d2*d/1.5;
242 G4double d4 = d3*d*0.75;
243 G4double y = 1.0/(1.0 + gam);
244 dedx = G4Log(2.0*(tau + 2.0)/eexc2) + G4Log(tau*d)
245 - beta2*(tau + 2.0*d - y*(3.0*d2
246 + y*(d - d3 + y*(d2 - tau*d3 + d4))))/tau;
247 }
248
249 //density correction
250 G4double x = G4Log(bg2)/twoln10;
251 dedx -= material->GetIonisation()->DensityCorrection(x);
252
253 // now you can compute the total ionization loss
254 dedx *= twopi_mc2_rcl2*electronDensity/beta2;
255 if (dedx < 0.0) { dedx = 0.0; }
256
257 // lowenergy extrapolation
258
259 if (kineticEnergy < th) {
260 x = kineticEnergy/th;
261 if(x > 0.25) { dedx /= sqrt(x); }
262 else { dedx *= 1.4*sqrt(x)/(0.1 + x); }
263 }
264 return dedx;
265}
G4double GetMeanExcitationEnergy() const
G4double DensityCorrection(G4double x)
G4double GetZeffective() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
G4double GetElectronDensity() const
Definition: G4Material.hh:215

◆ CrossSectionPerVolume()

G4double G4MollerBhabhaModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 189 of file G4MollerBhabhaModel.cc.

195{
196 G4double eDensity = material->GetElectronDensity();
197 return eDensity*ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
198}

◆ Initialise()

void G4MollerBhabhaModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector  
)
overridevirtual

Implements G4VEmModel.

Definition at line 106 of file G4MollerBhabhaModel.cc.

108{
109 if(!particle) { SetParticle(p); }
110
111 if(isInitialised) { return; }
112
113 isInitialised = true;
117 }
118}
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:611
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:708
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118

◆ MaxSecondaryEnergy()

G4double G4MollerBhabhaModel::MaxSecondaryEnergy ( const G4ParticleDefinition ,
G4double  kinEnergy 
)
finalprotectedvirtual

Reimplemented from G4VEmModel.

Definition at line 96 of file G4MollerBhabhaModel.cc.

98{
99 G4double tmax = kinEnergy;
100 if(isElectron) { tmax *= 0.5; }
101 return tmax;
102}

Referenced by G4PolarizedMollerBhabhaModel::ComputeCrossSectionPerElectron(), ComputeCrossSectionPerElectron(), and ComputeDEDXPerVolume().

◆ SampleSecondaries()

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

Implements G4VEmModel.

Reimplemented in G4PolarizedMollerBhabhaModel.

Definition at line 270 of file G4MollerBhabhaModel.cc.

275{
276 G4double kineticEnergy = dp->GetKineticEnergy();
277 //G4cout << "G4MollerBhabhaModel::SampleSecondaries: E= " << kineticEnergy
278 // << " in " << couple->GetMaterial()->GetName() << G4endl;
279 G4double tmax;
280 G4double tmin = cutEnergy;
281 if(isElectron) {
282 tmax = 0.5*kineticEnergy;
283 } else {
284 tmax = kineticEnergy;
285 }
286 if(maxEnergy < tmax) { tmax = maxEnergy; }
287 if(tmin >= tmax) { return; }
288
289 G4double energy = kineticEnergy + electron_mass_c2;
290 G4double xmin = tmin/kineticEnergy;
291 G4double xmax = tmax/kineticEnergy;
292 G4double gam = energy/electron_mass_c2;
293 G4double gamma2 = gam*gam;
294 G4double beta2 = 1.0 - 1.0/gamma2;
295 G4double x, z, grej;
296 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
297 G4double rndm[2];
298
299 //Moller (e-e-) scattering
300 if (isElectron) {
301
302 G4double gg = (2.0*gam - 1.0)/gamma2;
303 G4double y = 1.0 - xmax;
304 grej = 1.0 - gg*xmax + xmax*xmax*(1.0 - gg + (1.0 - gg*y)/(y*y));
305
306 do {
307 rndmEngine->flatArray(2, rndm);
308 x = xmin*xmax/(xmin*(1.0 - rndm[0]) + xmax*rndm[0]);
309 y = 1.0 - x;
310 z = 1.0 - gg*x + x*x*(1.0 - gg + (1.0 - gg*y)/(y*y));
311 /*
312 if(z > grej) {
313 G4cout << "G4MollerBhabhaModel::SampleSecondary Warning! "
314 << "Majorant " << grej << " < "
315 << z << " for x= " << x
316 << " e-e- scattering"
317 << G4endl;
318 }
319 */
320 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
321 } while(grej * rndm[1] > z);
322
323 //Bhabha (e+e-) scattering
324 } else {
325
326 G4double y = 1.0/(1.0 + gam);
327 G4double y2 = y*y;
328 G4double y12 = 1.0 - 2.0*y;
329 G4double b1 = 2.0 - y2;
330 G4double b2 = y12*(3.0 + y2);
331 G4double y122= y12*y12;
332 G4double b4 = y122*y12;
333 G4double b3 = b4 + y122;
334
335 y = xmax*xmax;
336 grej = 1.0 + (y*y*b4 - xmin*xmin*xmin*b3 + y*b2 - xmin*b1)*beta2;
337 do {
338 rndmEngine->flatArray(2, rndm);
339 x = xmin*xmax/(xmin*(1.0 - rndm[0]) + xmax*rndm[0]);
340 y = x*x;
341 z = 1.0 + (y*y*b4 - x*y*b3 + y*b2 - x*b1)*beta2;
342 /*
343 if(z > grej) {
344 G4cout << "G4MollerBhabhaModel::SampleSecondary Warning! "
345 << "Majorant " << grej << " < "
346 << z << " for x= " << x
347 << " e+e- scattering"
348 << G4endl;
349 }
350 */
351 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
352 } while(grej * rndm[1] > z);
353 }
354
355 G4double deltaKinEnergy = x * kineticEnergy;
356
357 G4ThreeVector deltaDirection;
358
360 const G4Material* mat = couple->GetMaterial();
362
363 deltaDirection =
364 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
365
366 } else {
367
368 G4double deltaMomentum =
369 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
370 G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
371 (deltaMomentum * dp->GetTotalMomentum());
372 if(cost > 1.0) { cost = 1.0; }
373 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
374
375 G4double phi = twopi * rndmEngine->flat() ;
376
377 deltaDirection.set(sint*cos(phi),sint*sin(phi), cost) ;
378 deltaDirection.rotateUz(dp->GetMomentumDirection());
379 }
380
381 // create G4DynamicParticle object for delta ray
382 G4DynamicParticle* delta =
383 new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
384 vdp->push_back(delta);
385
386 // primary change
387 kineticEnergy -= deltaKinEnergy;
388 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
389 finalP = finalP.unit();
390
393}
int G4int
Definition: G4Types.hh:85
Hep3Vector unit() const
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetTotalMomentum() const
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.cc:315
G4double energy(const ThreeVector &p, const G4double m)

◆ SetParticle()

void G4MollerBhabhaModel::SetParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 132 of file G4MollerBhabhaModel.hh.

133{
134 particle = p;
135 if(p != theElectron) { isElectron = false; }
136}

Referenced by ComputeCrossSectionPerElectron(), ComputeDEDXPerVolume(), G4MollerBhabhaModel(), and Initialise().

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForLoss* G4MollerBhabhaModel::fParticleChange
protected

◆ isElectron

◆ lowLimit

G4double G4MollerBhabhaModel::lowLimit
protected

Definition at line 118 of file G4MollerBhabhaModel.hh.

◆ particle

const G4ParticleDefinition* G4MollerBhabhaModel::particle
protected

◆ theElectron

◆ twoln10

G4double G4MollerBhabhaModel::twoln10
protected

Definition at line 117 of file G4MollerBhabhaModel.hh.

Referenced by ComputeDEDXPerVolume().


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