Geant4  10.02.p02
G4eBremParametrizedModel Class Reference

#include <G4eBremParametrizedModel.hh>

+ Inheritance diagram for G4eBremParametrizedModel:
+ Collaboration diagram for G4eBremParametrizedModel:

Public Member Functions

 G4eBremParametrizedModel (const G4ParticleDefinition *p=0, const G4String &nam="eBremParam")
 
virtual ~G4eBremParametrizedModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
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 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 Attributes

G4NistManagernist
 
const G4ParticleDefinitionparticle
 
G4ParticleDefinitiontheGamma
 
G4ParticleChangeForLossfParticleChange
 
G4double minThreshold
 
G4double particleMass
 
G4double kinEnergy
 
G4double totalEnergy
 
G4double currentZ
 
G4double z13
 
G4double z23
 
G4double lnZ
 
G4double densityFactor
 
G4double densityCorr
 
G4double Fel
 
G4double Finel
 
G4double facFel
 
G4double facFinel
 
G4double fMax
 
G4double fCoulomb
 
G4bool isElectron
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Static Protected Attributes

static const G4double xgi [8]
 
static const G4double wgi [8]
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Private Member Functions

void InitialiseConstants ()
 
G4double ComputeBremLoss (G4double cutEnergy)
 
G4double ComputeXSectionPerAtom (G4double cutEnergy)
 
G4double ComputeDXSectionPerAtom (G4double gammaEnergy)
 
G4double ComputeParametrizedDXSectionPerAtom (G4double kineticEnergy, G4double gammaEnergy, G4double Z)
 
void SetParticle (const G4ParticleDefinition *p)
 
G4double ScreenFunction1 (G4double ScreenVariable)
 
G4double ScreenFunction2 (G4double ScreenVariable)
 
void SetCurrentElement (const G4double)
 
G4eBremParametrizedModeloperator= (const G4eBremParametrizedModel &right)
 
 G4eBremParametrizedModel (const G4eBremParametrizedModel &)
 

Private Attributes

G4double lowKinEnergy
 
G4double fMigdalConstant
 
G4double bremFactor
 
G4bool isInitialised
 

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 *)
 

Detailed Description

Definition at line 61 of file G4eBremParametrizedModel.hh.

Constructor & Destructor Documentation

G4eBremParametrizedModel::G4eBremParametrizedModel ( const G4ParticleDefinition p = 0,
const G4String nam = "eBremParam" 
)

Definition at line 100 of file G4eBremParametrizedModel.cc.

References currentZ, densityCorr, densityFactor, fCoulomb, Fel, Finel, fMax, G4Gamma::Gamma(), InitialiseConstants(), G4NistManager::Instance(), keV, kinEnergy, lnZ, lowKinEnergy, MeV, minThreshold, nist, particleMass, G4VEmModel::SetAngularDistribution(), G4VEmModel::SetLowEnergyLimit(), SetParticle(), theGamma, totalEnergy, z13, and z23.

+ Here is the call graph for this function:

G4eBremParametrizedModel::~G4eBremParametrizedModel ( )
virtual

Definition at line 136 of file G4eBremParametrizedModel.cc.

G4eBremParametrizedModel::G4eBremParametrizedModel ( const G4eBremParametrizedModel )
private

Member Function Documentation

G4double G4eBremParametrizedModel::ComputeBremLoss ( G4double  cutEnergy)
private

Definition at line 234 of file G4eBremParametrizedModel.cc.

References ComputeDXSectionPerAtom(), densityCorr, n, totalEnergy, wgi, and xgi.

Referenced by ComputeDEDXPerVolume().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4double G4eBremParametrizedModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  tkin,
G4double  Z,
G4double  ,
G4double  cutEnergy,
G4double  maxEnergy = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 267 of file G4eBremParametrizedModel.cc.

References bremFactor, ComputeXSectionPerAtom(), kinEnergy, lowKinEnergy, G4INCL::Math::min(), particle, SetCurrentElement(), and SetParticle().

+ Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 201 of file G4eBremParametrizedModel.cc.

References bremFactor, ComputeBremLoss(), currentZ, G4Material::GetAtomicNumDensityVector(), G4Material::GetElementVector(), G4Material::GetNumberOfElements(), lowKinEnergy, G4INCL::Math::min(), particle, SetCurrentElement(), G4VEmModel::SetCurrentElement(), SetParticle(), and SetupForMaterial().

+ Here is the call graph for this function:

G4double G4eBremParametrizedModel::ComputeDXSectionPerAtom ( G4double  gammaEnergy)
private

Definition at line 441 of file G4eBremParametrizedModel.cc.

References ComputeParametrizedDXSectionPerAtom(), currentZ, fCoulomb, Fel, Finel, kinEnergy, main(), ScreenFunction1(), ScreenFunction2(), and totalEnergy.

Referenced by ComputeBremLoss(), ComputeXSectionPerAtom(), and SampleSecondaries().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4double G4eBremParametrizedModel::ComputeParametrizedDXSectionPerAtom ( G4double  kineticEnergy,
G4double  gammaEnergy,
G4double  Z 
)
private

Definition at line 363 of file G4eBremParametrizedModel.cc.

References ah10, ah11, ah12, ah20, ah21, ah22, ah30, ah31, ah32, al00, al01, al02, al10, al11, al12, al20, al21, al22, bh10, bh11, bh12, bh20, bh21, bh22, bh30, bh31, bh32, bl00, bl01, bl02, bl10, bl11, bl12, bl20, bl21, bl22, G4Log(), G4lrint(), G4NistManager::GetZ13(), lnZ, G4INCL::Math::max(), nist, ScreenFunction1(), ScreenFunction2(), SetCurrentElement(), tlow, totalEnergy, x, and z13.

Referenced by ComputeDXSectionPerAtom().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4double G4eBremParametrizedModel::ComputeXSectionPerAtom ( G4double  cutEnergy)
private

Definition at line 295 of file G4eBremParametrizedModel.cc.

References ComputeDXSectionPerAtom(), densityCorr, G4Exp(), G4Log(), kinEnergy, n, totalEnergy, wgi, and xgi.

Referenced by ComputeCrossSectionPerAtom().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void G4eBremParametrizedModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 175 of file G4eBremParametrizedModel.cc.

References currentZ, fParticleChange, G4VEmModel::GetParticleChangeForLoss(), G4VEmModel::InitialiseElementSelectors(), isInitialised, G4VEmModel::IsMaster(), G4VEmModel::LowEnergyLimit(), lowKinEnergy, and SetParticle().

+ Here is the call graph for this function:

void G4eBremParametrizedModel::InitialiseConstants ( )
private

Definition at line 128 of file G4eBremParametrizedModel.cc.

References facFel, facFinel, and G4Log().

Referenced by G4eBremParametrizedModel().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void G4eBremParametrizedModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 193 of file G4eBremParametrizedModel.cc.

References G4VEmModel::GetElementSelectors(), and G4VEmModel::SetElementSelectors().

+ Here is the call graph for this function:

G4double G4eBremParametrizedModel::MinEnergyCut ( const G4ParticleDefinition ,
const G4MaterialCutsCouple  
)
virtual

Reimplemented from G4VEmModel.

Definition at line 152 of file G4eBremParametrizedModel.cc.

References minThreshold.

G4eBremParametrizedModel& G4eBremParametrizedModel::operator= ( const G4eBremParametrizedModel right)
private
G4double G4eBremParametrizedModel::ScreenFunction1 ( G4double  ScreenVariable)
private

Definition at line 332 of file G4eBremParametrizedModel.cc.

References G4Log().

Referenced by ComputeDXSectionPerAtom(), and ComputeParametrizedDXSectionPerAtom().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4double G4eBremParametrizedModel::ScreenFunction2 ( G4double  ScreenVariable)
private

Definition at line 348 of file G4eBremParametrizedModel.cc.

References G4Log().

Referenced by ComputeDXSectionPerAtom(), and ComputeParametrizedDXSectionPerAtom().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void G4eBremParametrizedModel::SetCurrentElement ( const G4double  Z)
inlineprivate

Definition at line 164 of file G4eBremParametrizedModel.hh.

References currentZ, facFel, facFinel, fCoulomb, Fel, Finel, fMax, G4VEmModel::GetCurrentElement(), G4Element::GetfCoulomb(), G4NistManager::GetLOGZ(), G4NistManager::GetZ13(), iz, lnZ, nist, z13, and z23.

Referenced by ComputeCrossSectionPerAtom(), ComputeDEDXPerVolume(), ComputeParametrizedDXSectionPerAtom(), and SampleSecondaries().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void G4eBremParametrizedModel::SetParticle ( const G4ParticleDefinition p)
private

Definition at line 142 of file G4eBremParametrizedModel.cc.

References G4Electron::Electron(), G4ParticleDefinition::GetPDGMass(), isElectron, particle, and particleMass.

Referenced by ComputeCrossSectionPerAtom(), ComputeDEDXPerVolume(), G4eBremParametrizedModel(), and Initialise().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void G4eBremParametrizedModel::SetupForMaterial ( const G4ParticleDefinition ,
const G4Material mat,
G4double  kineticEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 160 of file G4eBremParametrizedModel.cc.

References densityCorr, densityFactor, fMigdalConstant, G4Material::GetElectronDensity(), kinEnergy, particleMass, and totalEnergy.

Referenced by ComputeDEDXPerVolume(), and SampleSecondaries().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

Member Data Documentation

G4double G4eBremParametrizedModel::bremFactor
private
G4double G4eBremParametrizedModel::densityCorr
protected
G4double G4eBremParametrizedModel::densityFactor
protected
G4double G4eBremParametrizedModel::facFel
protected

Definition at line 147 of file G4eBremParametrizedModel.hh.

Referenced by InitialiseConstants(), and SetCurrentElement().

G4double G4eBremParametrizedModel::facFinel
protected

Definition at line 147 of file G4eBremParametrizedModel.hh.

Referenced by InitialiseConstants(), and SetCurrentElement().

G4double G4eBremParametrizedModel::fCoulomb
protected
G4double G4eBremParametrizedModel::Fel
protected
G4double G4eBremParametrizedModel::Finel
protected
G4double G4eBremParametrizedModel::fMax
protected
G4double G4eBremParametrizedModel::fMigdalConstant
private

Definition at line 156 of file G4eBremParametrizedModel.hh.

Referenced by SetupForMaterial().

G4ParticleChangeForLoss* G4eBremParametrizedModel::fParticleChange
protected

Definition at line 132 of file G4eBremParametrizedModel.hh.

Referenced by Initialise(), and SampleSecondaries().

G4bool G4eBremParametrizedModel::isElectron
protected

Definition at line 150 of file G4eBremParametrizedModel.hh.

Referenced by SetParticle().

G4bool G4eBremParametrizedModel::isInitialised
private

Definition at line 159 of file G4eBremParametrizedModel.hh.

Referenced by Initialise().

G4double G4eBremParametrizedModel::lnZ
protected
G4double G4eBremParametrizedModel::lowKinEnergy
private
G4double G4eBremParametrizedModel::minThreshold
protected

Definition at line 136 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and MinEnergyCut().

G4NistManager* G4eBremParametrizedModel::nist
protected
const G4ParticleDefinition* G4eBremParametrizedModel::particle
protected
G4double G4eBremParametrizedModel::particleMass
protected
G4ParticleDefinition* G4eBremParametrizedModel::theGamma
protected

Definition at line 131 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and SampleSecondaries().

const G4double G4eBremParametrizedModel::wgi
staticprotected
Initial value:
={ 0.0506, 0.1112, 0.1569, 0.1813,
0.1813, 0.1569, 0.1112, 0.0506 }

Definition at line 134 of file G4eBremParametrizedModel.hh.

Referenced by ComputeBremLoss(), and ComputeXSectionPerAtom().

const G4double G4eBremParametrizedModel::xgi
staticprotected
Initial value:
={ 0.0199, 0.1017, 0.2372, 0.4083,
0.5917, 0.7628, 0.8983, 0.9801 }

Definition at line 134 of file G4eBremParametrizedModel.hh.

Referenced by ComputeBremLoss(), and ComputeXSectionPerAtom().

G4double G4eBremParametrizedModel::z13
protected
G4double G4eBremParametrizedModel::z23
protected

Definition at line 143 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and SetCurrentElement().


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