Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MyMollerBhabhaModel Class Reference

#include <MyMollerBhabhaModel.hh>

Inheritance diagram for MyMollerBhabhaModel:
Collaboration diagram for MyMollerBhabhaModel:

Public Member Functions

 MyMollerBhabhaModel (const G4ParticleDefinition *p=0, const G4String &nam="myMollerBhabha")
 
 ~MyMollerBhabhaModel ()
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
 MyMollerBhabhaModel (const G4ParticleDefinition *p=0, const G4String &nam="myMollerBhabha")
 
 ~MyMollerBhabhaModel ()
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
- Public Member Functions inherited from G4MollerBhabhaModel
 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 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 SetFluctuationFlag (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)
 

Additional Inherited Members

- Protected Member Functions inherited from G4MollerBhabhaModel
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 inherited from G4MollerBhabhaModel
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
 
G4bool lossFlucFlag
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 41 of file MyMollerBhabhaModel.hh.

Constructor & Destructor Documentation

MyMollerBhabhaModel::MyMollerBhabhaModel ( const G4ParticleDefinition p = 0,
const G4String nam = "myMollerBhabha" 
)

Definition at line 43 of file MyMollerBhabhaModel.cc.

45  : G4MollerBhabhaModel(p,nam)
46 {}
G4MollerBhabhaModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="MollerBhabha")
MyMollerBhabhaModel::~MyMollerBhabhaModel ( )

Definition at line 50 of file MyMollerBhabhaModel.cc.

51 {}
MyMollerBhabhaModel::MyMollerBhabhaModel ( const G4ParticleDefinition p = 0,
const G4String nam = "myMollerBhabha" 
)
MyMollerBhabhaModel::~MyMollerBhabhaModel ( )

Member Function Documentation

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

Reimplemented from G4MollerBhabhaModel.

Definition at line 55 of file MyMollerBhabhaModel.cc.

60 {
61  if(!particle) SetParticle(p);
62  // calculate the dE/dx due to the ionization by Seltzer-Berger formula
63 
64  G4double electronDensity = material->GetElectronDensity();
65  G4double Zeff = electronDensity/material->GetTotNbOfAtomsPerVolume();
66  G4double th = 0.25*sqrt(Zeff)*keV;
67  G4double tkin = kineticEnergy;
68  G4bool lowEnergy = false;
69  if (kineticEnergy < th) {
70  tkin = th;
71  lowEnergy = true;
72  }
73  G4double tau = tkin/electron_mass_c2;
74  G4double gam = tau + 1.0;
75  G4double gamma2= gam*gam;
76  G4double beta2 = 1. - 1./gamma2;
77  //G4double bg2 = beta2*gamma2;
78 
79  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
80  eexc /= electron_mass_c2;
81  G4double eexc2 = eexc*eexc;
82 
83  G4double d = min(cutEnergy, MaxSecondaryEnergy(p, tkin))/electron_mass_c2;
84  G4double dedx;
85 
86  // electron
87  if (isElectron) {
88 
89  dedx = log(2.0*(tau + 2.0)/eexc2) - 1.0 - beta2
90  + log((tau-d)*d) + tau/(tau-d)
91  + (0.5*d*d + (2.0*tau + 1.)*log(1. - d/tau))/gamma2;
92 
93  //positron
94  } else {
95 
96  G4double d2 = d*d*0.5;
97  G4double d3 = d2*d/1.5;
98  G4double d4 = d3*d*3.75;
99  G4double y = 1.0/(1.0 + gam);
100  dedx = log(2.0*(tau + 2.0)/eexc2) + log(tau*d)
101  - beta2*(tau + 2.0*d - y*(3.0*d2
102  + y*(d - d3 + y*(d2 - tau*d3 + d4))))/tau;
103  }
104 
105  //do not apply density correction
106  //G4double cden = material->GetIonisation()->GetCdensity();
107  //G4double mden = material->GetIonisation()->GetMdensity();
108  //G4double aden = material->GetIonisation()->GetAdensity();
109  //G4double x0den = material->GetIonisation()->GetX0density();
110  //G4double x1den = material->GetIonisation()->GetX1density();
111  //G4double x = log(bg2)/twoln10;
112 
113  //if (x >= x0den) {
114  // dedx -= twoln10*x - cden;
115  // if (x < x1den) dedx -= aden*pow(x1den-x, mden);
116  //}
117 
118  // now you can compute the total ionization loss
119  dedx *= twopi_mc2_rcl2*electronDensity/beta2;
120  if (dedx < 0.0) dedx = 0.0;
121 
122  // lowenergy extrapolation
123 
124  if (lowEnergy) {
125 
126  if (kineticEnergy >= lowLimit) dedx *= sqrt(tkin/kineticEnergy);
127  else dedx *= sqrt(tkin*kineticEnergy)/lowLimit;
128 
129  }
130  return dedx;
131 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
static const G4double d2
void SetParticle(const G4ParticleDefinition *p)
G4double GetElectronDensity() const
Definition: G4Material.hh:217
const G4ParticleDefinition * particle
bool G4bool
Definition: G4Types.hh:79
float electron_mass_c2
Definition: hepunit.py:274
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:209
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:

virtual G4double MyMollerBhabhaModel::ComputeDEDXPerVolume ( const G4Material ,
const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  cutEnergy 
)
virtual

Reimplemented from G4MollerBhabhaModel.


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