Geant4  10.02.p03
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 &)
 
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)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle *> *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
- 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, const G4ParticleDefinition *, G4double)
 
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 *> *)
 
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=0)
 
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)
 
void SetParticle (const G4ParticleDefinition *p)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

const G4ParticleDefinitionparticle
 
G4ParticleDefinitiontheElectron
 
G4ParticleChangeForLoss * fParticleChange
 
G4bool isElectron
 
G4double twoln10
 
G4double lowLimit
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Private Member Functions

G4MollerBhabhaModeloperator= (const G4MollerBhabhaModel &right)
 
 G4MollerBhabhaModel (const G4MollerBhabhaModel &)
 

Private Attributes

G4bool isInitialised
 

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() [1/2]

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

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:69
G4ParticleChangeForLoss * fParticleChange
const G4ParticleDefinition * particle
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4ParticleDefinition * theElectron
static G4Electron * Electron()
Definition: G4Electron.cc:94
static const double keV
Definition: G4SIunits.hh:213
Here is the call graph for this function:

◆ ~G4MollerBhabhaModel()

G4MollerBhabhaModel::~G4MollerBhabhaModel ( )
virtual

Definition at line 92 of file G4MollerBhabhaModel.cc.

93 {}

◆ G4MollerBhabhaModel() [2/2]

G4MollerBhabhaModel::G4MollerBhabhaModel ( const G4MollerBhabhaModel )
private

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 178 of file G4MollerBhabhaModel.cc.

184 {
185  return Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
186 }
Float_t Z
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
Here is the call graph for this function:

◆ ComputeCrossSectionPerElectron()

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 }
Double_t y2[nxs]
int twopi_mc2_rcl2
Definition: hepunit.py:294
void SetParticle(const G4ParticleDefinition *p)
Double_t y
static const G4double b3
const G4ParticleDefinition * particle
static const G4double b2
float electron_mass_c2
Definition: hepunit.py:274
G4double G4Log(G4double x)
Definition: G4Log.hh:230
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
static const G4double b1
double G4double
Definition: G4Types.hh:76
static const G4double b4
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Reimplemented in MyMollerBhabhaModel, and MyMollerBhabhaModel.

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 cut;
217  // if(isElectron) { cut = std::max(th*0.5, cutEnergy); }
218  // else { cut = std::max(th, cutEnergy); }
219  G4double tkin = kineticEnergy;
220  if (kineticEnergy < th) { tkin = th; }
221 
222  G4double tau = tkin/electron_mass_c2;
223  G4double gam = tau + 1.0;
224  G4double gamma2= gam*gam;
225  G4double bg2 = tau*(tau + 2);
226  G4double beta2 = bg2/gamma2;
227 
228  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
229  eexc /= electron_mass_c2;
230  G4double eexc2 = eexc*eexc;
231 
233  G4double dedx;
234 
235  // electron
236  if (isElectron) {
237 
238  dedx = G4Log(2.0*(tau + 2.0)/eexc2) - 1.0 - beta2
239  + G4Log((tau-d)*d) + tau/(tau-d)
240  + (0.5*d*d + (2.0*tau + 1.)*G4Log(1. - d/tau))/gamma2;
241 
242  //positron
243  } else {
244 
245  G4double d2 = d*d*0.5;
246  G4double d3 = d2*d/1.5;
247  G4double d4 = d3*d*0.75;
248  G4double y = 1.0/(1.0 + gam);
249  dedx = G4Log(2.0*(tau + 2.0)/eexc2) + G4Log(tau*d)
250  - beta2*(tau + 2.0*d - y*(3.0*d2
251  + y*(d - d3 + y*(d2 - tau*d3 + d4))))/tau;
252  }
253 
254  //density correction
255  G4double x = G4Log(bg2)/twoln10;
256  dedx -= material->GetIonisation()->DensityCorrection(x);
257 
258  // now you can compute the total ionization loss
259  dedx *= twopi_mc2_rcl2*electronDensity/beta2;
260  if (dedx < 0.0) { dedx = 0.0; }
261 
262  // lowenergy extrapolation
263 
264  if (kineticEnergy < th) {
265  x = kineticEnergy/th;
266  if(x > 0.25) { dedx /= sqrt(x); }
267  else { dedx *= 1.4*sqrt(x)/(0.1 + x); }
268  }
269  return dedx;
270 }
static const G4double d3
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
int twopi_mc2_rcl2
Definition: hepunit.py:294
Float_t d
G4double GetZeffective() const
void SetParticle(const G4ParticleDefinition *p)
Double_t y
G4double GetMeanExcitationEnergy() const
const G4ParticleDefinition * particle
float electron_mass_c2
Definition: hepunit.py:274
G4double G4Log(G4double x)
Definition: G4Log.hh:230
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
G4double DensityCorrection(G4double x)
G4double GetElectronDensity() const
Definition: G4Material.hh:217
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
static const G4double d2
static const G4double d4
Here is the call graph for this function:

◆ CrossSectionPerVolume()

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

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:

◆ Initialise()

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

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:120
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:697
void SetParticle(const G4ParticleDefinition *p)
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:617
G4ParticleChangeForLoss * fParticleChange
const G4ParticleDefinition * particle
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:624
Here is the call graph for this function:
Here is the caller graph for this function:

◆ MaxSecondaryEnergy()

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

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:

◆ operator=()

G4MollerBhabhaModel& G4MollerBhabhaModel::operator= ( const G4MollerBhabhaModel right)
private

◆ SampleSecondaries()

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

Implements G4VEmModel.

Reimplemented in G4PolarizedMollerBhabhaModel.

Definition at line 275 of file G4MollerBhabhaModel.cc.

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

◆ 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 }
const G4ParticleDefinition * particle
G4ParticleDefinition * theElectron
Here is the caller graph for this function:

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForLoss* G4MollerBhabhaModel::fParticleChange
protected

Definition at line 114 of file G4MollerBhabhaModel.hh.

◆ isElectron

G4bool G4MollerBhabhaModel::isElectron
protected

Definition at line 116 of file G4MollerBhabhaModel.hh.

◆ isInitialised

G4bool G4MollerBhabhaModel::isInitialised
private

Definition at line 126 of file G4MollerBhabhaModel.hh.

◆ lowLimit

G4double G4MollerBhabhaModel::lowLimit
protected

Definition at line 118 of file G4MollerBhabhaModel.hh.

◆ particle

const G4ParticleDefinition* G4MollerBhabhaModel::particle
protected

Definition at line 112 of file G4MollerBhabhaModel.hh.

◆ theElectron

G4ParticleDefinition* G4MollerBhabhaModel::theElectron
protected

Definition at line 113 of file G4MollerBhabhaModel.hh.

◆ twoln10

G4double G4MollerBhabhaModel::twoln10
protected

Definition at line 117 of file G4MollerBhabhaModel.hh.


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