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

#include <G4hBremsstrahlungModel.hh>

Inheritance diagram for G4hBremsstrahlungModel:
Collaboration diagram for G4hBremsstrahlungModel:

Public Member Functions

 G4hBremsstrahlungModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="hBrem")
 
virtual ~G4hBremsstrahlungModel ()
 
- Public Member Functions inherited from G4MuBremsstrahlungModel
 G4MuBremsstrahlungModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="MuBrem")
 
virtual ~G4MuBremsstrahlungModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SetLowestKineticEnergy (G4double e)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double) override
 
- 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 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 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)
 

Protected Member Functions

virtual 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)
 
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
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 56 of file G4hBremsstrahlungModel.hh.

Constructor & Destructor Documentation

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

Definition at line 58 of file G4hBremsstrahlungModel.cc.

60  : G4MuBremsstrahlungModel(p, nam)
61 {}
G4MuBremsstrahlungModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="MuBrem")
G4hBremsstrahlungModel::~G4hBremsstrahlungModel ( )
virtual

Definition at line 65 of file G4hBremsstrahlungModel.cc.

66 {}

Member Function Documentation

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 = G4int(Z);
87  if(iz < 1) { iz = 1; }
88 
89  G4double z13 = 1.0/nist->GetZ13(iz);
90  G4double dn = mass*nist->GetA27(iz)/(70.*MeV);
91 
92  G4double b = btf;
93  if(1 == iz) b = bh;
94 
95  // nucleus contribution logarithm
96  G4double rab1=b*z13;
97  G4double fn=G4Log(rab1/(dn*(electron_mass_c2+rab0*rab1))*
98  (mass+delta*(dn*sqrte-2.))) ;
99  if(fn <0.) fn = 0. ;
100 
101  G4double x = 1.0 - v;
102  if(particle->GetPDGSpin() != 0) { x += 0.75*v*v; }
103 
104  dxsection = coeff*x*Z*Z*fn/gammaEnergy;
105 
106  return dxsection;
107 }
G4double GetA27(G4int Z) const
tuple x
Definition: test.py:50
int G4int
Definition: G4Types.hh:78
G4double GetZ13(G4double Z) const
tuple b
Definition: test.py:12
const G4ParticleDefinition * particle
float electron_mass_c2
Definition: hepunit.py:274
tuple v
Definition: test.py:18
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double GetPDGSpin() const
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:


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