Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4MollerBhabhaModel Class Reference

#include <G4MollerBhabhaModel.hh>

Inheritance diagram for G4MollerBhabhaModel:
Collaboration 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 InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
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 *, 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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, 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 *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () 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)
 

Protected Member Functions

virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) final
 
void SetParticle (const G4ParticleDefinition *p)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
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 std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Additional Inherited Members

- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 63 of file G4MollerBhabhaModel.hh.

Constructor & Destructor Documentation

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

Definition at line 76 of file G4MollerBhabhaModel.cc.

78  : G4VEmModel(nam),
79  particle(nullptr),
80  isElectron(true),
81  twoln10(2.0*G4Log(10.0)),
82  lowLimit(0.02*keV),
83  isInitialised(false)
84 {
86  if(p) { SetParticle(p); }
87  fParticleChange = nullptr;
88 }
void SetParticle(const G4ParticleDefinition *p)
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
G4ParticleChangeForLoss * fParticleChange
const G4ParticleDefinition * particle
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4ParticleDefinition * theElectron
static G4Electron * Electron()
Definition: G4Electron.cc:94
static constexpr double keV
Definition: G4SIunits.hh:216

Here is the call graph for this function:

G4MollerBhabhaModel::~G4MollerBhabhaModel ( )
virtual

Definition at line 92 of file G4MollerBhabhaModel.cc.

93 {}

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 178 of file G4MollerBhabhaModel.cc.

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

Here is the call graph for this function:

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

Reimplemented in G4PolarizedMollerBhabhaModel.

Definition at line 124 of file G4MollerBhabhaModel.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 203 of file G4MollerBhabhaModel.cc.

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

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 190 of file G4MollerBhabhaModel.cc.

196 {
197  G4double eDensity = material->GetElectronDensity();
198  return eDensity*ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
199 }
G4double GetElectronDensity() const
Definition: G4Material.hh:217
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 107 of file G4MollerBhabhaModel.cc.

109 {
110  if(!particle) { SetParticle(p); }
111 
112  if(isInitialised) { return; }
113 
114  isInitialised = true;
118  }
119 }
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
void SetParticle(const G4ParticleDefinition *p)
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:616
G4ParticleChangeForLoss * fParticleChange
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:696
const G4ParticleDefinition * particle
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:623

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 97 of file G4MollerBhabhaModel.cc.

99 {
100  G4double tmax = kinEnergy;
101  if(isElectron) { tmax *= 0.5; }
102  return tmax;
103 }
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

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 271 of file G4MollerBhabhaModel.cc.

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

Here is the call graph for this function:

void G4MollerBhabhaModel::SetParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 133 of file G4MollerBhabhaModel.hh.

134 {
135  particle = p;
136  if(p != theElectron) { isElectron = false; }
137 }
const char * p
Definition: xmltok.h:285
const G4ParticleDefinition * particle
G4ParticleDefinition * theElectron

Here is the caller graph for this function:

Member Data Documentation

G4ParticleChangeForLoss* G4MollerBhabhaModel::fParticleChange
protected

Definition at line 115 of file G4MollerBhabhaModel.hh.

G4bool G4MollerBhabhaModel::isElectron
protected

Definition at line 117 of file G4MollerBhabhaModel.hh.

G4double G4MollerBhabhaModel::lowLimit
protected

Definition at line 119 of file G4MollerBhabhaModel.hh.

const G4ParticleDefinition* G4MollerBhabhaModel::particle
protected

Definition at line 113 of file G4MollerBhabhaModel.hh.

G4ParticleDefinition* G4MollerBhabhaModel::theElectron
protected

Definition at line 114 of file G4MollerBhabhaModel.hh.

G4double G4MollerBhabhaModel::twoln10
protected

Definition at line 118 of file G4MollerBhabhaModel.hh.


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