Geant4 11.1.1
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=nullptr, const G4String &nam="hBrem")
 
 ~G4hBremsstrahlungModel () override
 
G4hBremsstrahlungModeloperator= (const G4hBremsstrahlungModel &right)=delete
 
 G4hBremsstrahlungModel (const G4hBremsstrahlungModel &)=delete
 
- Public Member Functions inherited from G4MuBremsstrahlungModel
 G4MuBremsstrahlungModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="MuBrem")
 
 ~G4MuBremsstrahlungModel ()=default
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SetLowestKineticEnergy (G4double e)
 
G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double) override
 
G4MuBremsstrahlungModeloperator= (const G4MuBremsstrahlungModel &right)=delete
 
 G4MuBremsstrahlungModel (const G4MuBremsstrahlungModel &)=delete
 
- 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 void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
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 GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
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 *, const G4double &length, G4double &eloss)
 
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 FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
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)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
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) const
 
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
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

G4double ComputeDMicroscopicCrossSection (G4double tkin, G4double Z, G4double gammaEnergy) override
 
- 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 = nullptr
 
G4ParticleDefinitiontheGamma = nullptr
 
G4ParticleChangeForLossfParticleChange = nullptr
 
G4NistManagernist = nullptr
 
G4double mass = 1.0
 
G4double rmass = 1.0
 
G4double cc = 1.0
 
G4double coeff = 1.0
 
G4double sqrte
 
G4double bh = 202.4
 
G4double bh1 = 446.
 
G4double btf = 183.
 
G4double btf1 = 1429.
 
G4double lowestKinEnergy
 
G4double minThreshold
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 
- Static Protected Attributes inherited from G4MuBremsstrahlungModel
static const G4double xgi [6]
 
static const G4double wgi [6]
 
static G4double fDN [93] = {0.0}
 

Detailed Description

Definition at line 52 of file G4hBremsstrahlungModel.hh.

Constructor & Destructor Documentation

◆ G4hBremsstrahlungModel() [1/2]

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

◆ ~G4hBremsstrahlungModel()

G4hBremsstrahlungModel::~G4hBremsstrahlungModel ( )
override

Definition at line 65 of file G4hBremsstrahlungModel.cc.

66{}

◆ G4hBremsstrahlungModel() [2/2]

G4hBremsstrahlungModel::G4hBremsstrahlungModel ( const G4hBremsstrahlungModel )
delete

Member Function Documentation

◆ ComputeDMicroscopicCrossSection()

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

Reimplemented from G4MuBremsstrahlungModel.

Definition at line 70 of file G4hBremsstrahlungModel.cc.

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

◆ operator=()

G4hBremsstrahlungModel & G4hBremsstrahlungModel::operator= ( const G4hBremsstrahlungModel right)
delete

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