Geant4-11
Public Member Functions | Protected Attributes | Private Attributes
G4DNARuddAngle Class Reference

#include <G4DNARuddAngle.hh>

Inheritance diagram for G4DNARuddAngle:
G4VEmAngularDistribution

Public Member Functions

 G4DNARuddAngle (const G4DNARuddAngle &)=delete
 
 G4DNARuddAngle (const G4String &name="")
 
const G4StringGetName () const
 
G4DNARuddAngleoperator= (const G4DNARuddAngle &right)=delete
 
void PrintGeneratorInformation () const override
 
G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, const G4Material *mat=nullptr) override
 
G4ThreeVectorSampleDirectionForShell (const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, G4int shellIdx, const G4Material *mat=nullptr) override
 
virtual void SamplePairDirections (const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr)
 
 ~G4DNARuddAngle () override
 

Protected Attributes

G4ThreeVector fLocalDirection
 
G4bool fPolarisation
 

Private Attributes

const G4ParticleDefinitionfElectron
 
G4String fName
 

Detailed Description

Definition at line 56 of file G4DNARuddAngle.hh.

Constructor & Destructor Documentation

◆ G4DNARuddAngle() [1/2]

G4DNARuddAngle::G4DNARuddAngle ( const G4String name = "")

Definition at line 58 of file G4DNARuddAngle.cc.

59 : G4VEmAngularDistribution("deltaRudd")
60{
62}
const G4ParticleDefinition * fElectron
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4VEmAngularDistribution(const G4String &name)

References G4Electron::Electron(), and fElectron.

◆ ~G4DNARuddAngle()

G4DNARuddAngle::~G4DNARuddAngle ( )
override

Definition at line 64 of file G4DNARuddAngle.cc.

65{}

◆ G4DNARuddAngle() [2/2]

G4DNARuddAngle::G4DNARuddAngle ( const G4DNARuddAngle )
delete

Member Function Documentation

◆ GetName()

const G4String & G4VEmAngularDistribution::GetName ( ) const
inlineinherited

Definition at line 111 of file G4VEmAngularDistribution.hh.

112{
113 return fName;
114}

References G4VEmAngularDistribution::fName.

◆ operator=()

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

◆ PrintGeneratorInformation()

void G4DNARuddAngle::PrintGeneratorInformation ( ) const
overridevirtual

Reimplemented from G4VEmAngularDistribution.

Definition at line 112 of file G4DNARuddAngle.cc.

113{}

◆ SampleDirection()

G4ThreeVector & G4DNARuddAngle::SampleDirection ( const G4DynamicParticle dp,
G4double  kinEnergyFinal,
G4int  Z,
const G4Material mat = nullptr 
)
overridevirtual

Implements G4VEmAngularDistribution.

Definition at line 105 of file G4DNARuddAngle.cc.

108{
109 return SampleDirectionForShell(dp, secEkin, Z, 0, mat);
110}
const G4int Z[17]
G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, G4int shellIdx, const G4Material *mat=nullptr) override

References SampleDirectionForShell(), and Z.

◆ SampleDirectionForShell()

G4ThreeVector & G4DNARuddAngle::SampleDirectionForShell ( const G4DynamicParticle dp,
G4double  kinEnergyFinal,
G4int  Z,
G4int  shellIdx,
const G4Material mat = nullptr 
)
overridevirtual

Reimplemented from G4VEmAngularDistribution.

Definition at line 68 of file G4DNARuddAngle.cc.

72{
73 G4double k = dp->GetKineticEnergy();
74 G4double cosTheta = 1.0;
75
76 const G4ParticleDefinition* particle = dp->GetDefinition();
77 G4double mass = particle->GetPDGMass();
78
79 G4double maximumEnergyTransfer = k;
80 if(particle == fElectron) { maximumEnergyTransfer *= 0.5; }
81 else if(mass > MeV) {
82 G4double ratio = electron_mass_c2/mass;
83 G4double tau = k/mass;
84 maximumEnergyTransfer = 2.0*electron_mass_c2*tau*(tau + 2.) /
85 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
86
87 }
88
89 if (secKinetic>100*eV && secKinetic <= maximumEnergyTransfer) {
90 cosTheta = std::sqrt(secKinetic / maximumEnergyTransfer);
91 } else {
92 cosTheta = (2.*G4UniformRand())-1.;
93 }
94
95 G4double sint = sqrt((1.0 - cosTheta)*(1.0 + cosTheta));
97
98 fLocalDirection.set(sint*cos(phi), sint*sin(phi), cosTheta);
100
101 return fLocalDirection;
102}
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double eV
Definition: G4SIunits.hh:201
static constexpr double MeV
Definition: G4SIunits.hh:200
double G4double
Definition: G4Types.hh:83
#define G4UniformRand()
Definition: Randomize.hh:52
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
float electron_mass_c2
Definition: hepunit.py:273

References source.hepunit::electron_mass_c2, eV, fElectron, G4VEmAngularDistribution::fLocalDirection, G4UniformRand, G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetPDGMass(), MeV, CLHEP::Hep3Vector::rotateUz(), CLHEP::Hep3Vector::set(), and twopi.

Referenced by SampleDirection().

◆ SamplePairDirections()

void G4VEmAngularDistribution::SamplePairDirections ( const G4DynamicParticle dp,
G4double  elecKinEnergy,
G4double  posiKinEnergy,
G4ThreeVector dirElectron,
G4ThreeVector dirPositron,
G4int  Z = 0,
const G4Material mat = nullptr 
)
virtualinherited

Field Documentation

◆ fElectron

const G4ParticleDefinition* G4DNARuddAngle::fElectron
private

Definition at line 83 of file G4DNARuddAngle.hh.

Referenced by G4DNARuddAngle(), and SampleDirectionForShell().

◆ fLocalDirection

G4ThreeVector G4VEmAngularDistribution::fLocalDirection
protectedinherited

◆ fName

G4String G4VEmAngularDistribution::fName
privateinherited

Definition at line 108 of file G4VEmAngularDistribution.hh.

Referenced by G4VEmAngularDistribution::GetName().

◆ fPolarisation

G4bool G4VEmAngularDistribution::fPolarisation
protectedinherited

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