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

#include <G4DNARuddAngle.hh>

+ Inheritance diagram for G4DNARuddAngle:

Public Member Functions

 G4DNARuddAngle (const G4String &name="")
 
virtual ~G4DNARuddAngle ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, const G4Material *mat=0)
 
virtual G4ThreeVectorSampleDirectionForShell (const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, G4int shellIdx, const G4Material *mat=0)
 
void PrintGeneratorInformation () const
 
- Public Member Functions inherited from G4VEmAngularDistribution
 G4VEmAngularDistribution (const G4String &name)
 
virtual ~G4VEmAngularDistribution ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
 
virtual G4ThreeVectorSampleDirectionForShell (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
 
virtual void SamplePairDirections (const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr)
 
const G4StringGetName () const
 
G4VEmAngularDistributionoperator= (const G4VEmAngularDistribution &right)=delete
 
 G4VEmAngularDistribution (const G4VEmAngularDistribution &)=delete
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmAngularDistribution
G4ThreeVector fLocalDirection
 
G4bool fPolarisation
 

Detailed Description

Definition at line 56 of file G4DNARuddAngle.hh.

Constructor & Destructor Documentation

◆ G4DNARuddAngle()

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

Definition at line 58 of file G4DNARuddAngle.cc.

59 : G4VEmAngularDistribution("deltaRudd")
60{
61 fElectron = G4Electron::Electron();
62}
static G4Electron * Electron()
Definition: G4Electron.cc:93

◆ ~G4DNARuddAngle()

G4DNARuddAngle::~G4DNARuddAngle ( )
virtual

Definition at line 64 of file G4DNARuddAngle.cc.

65{}

Member Function Documentation

◆ PrintGeneratorInformation()

void G4DNARuddAngle::PrintGeneratorInformation ( ) const

Definition at line 112 of file G4DNARuddAngle.cc.

113{}

◆ SampleDirection()

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

Implements G4VEmAngularDistribution.

Definition at line 105 of file G4DNARuddAngle.cc.

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

◆ SampleDirectionForShell()

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

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));
96 G4double phi = twopi*G4UniformRand();
97
98 fLocalDirection.set(sint*cos(phi), sint*sin(phi), cosTheta);
100
101 return fLocalDirection;
102}
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

Referenced by SampleDirection().


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