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

#include <G4hBremsstrahlungModel.hh>

+ Inheritance diagram for G4hBremsstrahlungModel:

Public Member Functions

 G4hBremsstrahlungModel (const G4ParticleDefinition *p=0, const G4String &nam="hBrem")
 
virtual ~G4hBremsstrahlungModel ()
 
- Public Member Functions inherited from G4MuBremsstrahlungModel
 G4MuBremsstrahlungModel (const G4ParticleDefinition *p=0, const G4String &nam="MuBrem")
 
virtual ~G4MuBremsstrahlungModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, 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)
 
void SetLowestKineticEnergy (G4double e)
 
- 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
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
 
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 ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Member Functions

virtual G4double ComputeDMicroscopicCrossSection (G4double tkin, G4double Z, G4double gammaEnergy)
 
- Protected Member Functions inherited from G4MuBremsstrahlungModel
G4double ComputMuBremLoss (G4double Z, G4double tkin, G4double cut)
 
G4double ComputeMicroscopicCrossSection (G4double tkin, G4double Z, G4double cut)
 
virtual G4double ComputeDMicroscopicCrossSection (G4double tkin, G4double Z, G4double gammaEnergy)
 
void SetParticle (const G4ParticleDefinition *)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Additional Inherited Members

- Protected Attributes inherited from G4MuBremsstrahlungModel
const G4ParticleDefinitionparticle
 
G4NistManagernist
 
G4double mass
 
G4double rmass
 
G4double cc
 
G4double coeff
 
G4double sqrte
 
G4double bh
 
G4double bh1
 
G4double btf
 
G4double btf1
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Detailed Description

Definition at line 56 of file G4hBremsstrahlungModel.hh.

Constructor & Destructor Documentation

◆ G4hBremsstrahlungModel()

G4hBremsstrahlungModel::G4hBremsstrahlungModel ( const G4ParticleDefinition p = 0,
const G4String nam = "hBrem" 
)

◆ ~G4hBremsstrahlungModel()

G4hBremsstrahlungModel::~G4hBremsstrahlungModel ( )
virtual

Definition at line 64 of file G4hBremsstrahlungModel.cc.

65{}

Member Function Documentation

◆ ComputeDMicroscopicCrossSection()

G4double G4hBremsstrahlungModel::ComputeDMicroscopicCrossSection ( G4double  tkin,
G4double  Z,
G4double  gammaEnergy 
)
protectedvirtual

Reimplemented from G4MuBremsstrahlungModel.

Definition at line 69 of file G4hBremsstrahlungModel.cc.

74{
75 G4double dxsection = 0.;
76
77 if( gammaEnergy > tkin) return dxsection ;
78 // G4cout << "G4hBremsstrahlungModel m= " << mass
79 // << " " << particle->GetParticleName() << G4endl;
80 G4double E = tkin + mass ;
81 G4double v = gammaEnergy/E ;
82 G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
83 G4double rab0=delta*sqrte ;
84
85 G4int iz = G4int(Z);
86 if(iz < 1) iz = 1;
87
88 G4double z13 = 1.0/nist->GetZ13(iz);
89 G4double dn = mass*nist->GetA27(iz)/(70.*MeV);
90
91 G4double b = btf;
92 if(1 == iz) b = bh;
93
94 // nucleus contribution logarithm
95 G4double rab1=b*z13;
96 G4double fn=log(rab1/(dn*(electron_mass_c2+rab0*rab1))*
97 (mass+delta*(dn*sqrte-2.))) ;
98 if(fn <0.) fn = 0. ;
99
100 G4double x = 1.0 - v;
101 if(particle->GetPDGSpin() != 0) x += 0.75*v*v;
102
103 dxsection = coeff*x*Z*Z*fn/gammaEnergy;
104
105 return dxsection;
106}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
const G4ParticleDefinition * particle
G4double GetA27(G4int Z)
G4double GetZ13(G4double Z)

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