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

#include <G4DeltaAngleFreeScat.hh>

Inheritance diagram for G4DeltaAngleFreeScat:
G4VEmAngularDistribution

Public Member Functions

 G4DeltaAngleFreeScat (const G4DeltaAngleFreeScat &)=delete
 
 G4DeltaAngleFreeScat (const G4String &name="")
 
const G4StringGetName () const
 
G4DeltaAngleFreeScatoperator= (const G4DeltaAngleFreeScat &right)=delete
 
void PrintGeneratorInformation () const final
 
G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, const G4Material *mat=nullptr) final
 
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)
 
 ~G4DeltaAngleFreeScat () override
 

Protected Attributes

G4ThreeVector fLocalDirection
 
G4bool fPolarisation
 

Private Attributes

G4String fName
 

Detailed Description

Definition at line 53 of file G4DeltaAngleFreeScat.hh.

Constructor & Destructor Documentation

◆ G4DeltaAngleFreeScat() [1/2]

G4DeltaAngleFreeScat::G4DeltaAngleFreeScat ( const G4String name = "")
explicit

Definition at line 55 of file G4DeltaAngleFreeScat.cc.

56 : G4VEmAngularDistribution("deltaFree")
57{}
G4VEmAngularDistribution(const G4String &name)

◆ ~G4DeltaAngleFreeScat()

G4DeltaAngleFreeScat::~G4DeltaAngleFreeScat ( )
override

Definition at line 59 of file G4DeltaAngleFreeScat.cc.

60{}

◆ G4DeltaAngleFreeScat() [2/2]

G4DeltaAngleFreeScat::G4DeltaAngleFreeScat ( const G4DeltaAngleFreeScat )
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=()

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

◆ PrintGeneratorInformation()

void G4DeltaAngleFreeScat::PrintGeneratorInformation ( ) const
finalvirtual

Reimplemented from G4VEmAngularDistribution.

Definition at line 82 of file G4DeltaAngleFreeScat.cc.

83{}

◆ SampleDirection()

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

Implements G4VEmAngularDistribution.

Definition at line 63 of file G4DeltaAngleFreeScat.cc.

66{
67 G4double deltaMomentum =
68 sqrt(kinEnergyFinal*(kinEnergyFinal + 2*electron_mass_c2));
69
70 G4double costet = kinEnergyFinal*(dp->GetTotalEnergy() + electron_mass_c2) /
71 (deltaMomentum * dp->GetTotalMomentum());
72
74 G4double sintet = sqrt((1 - costet)*(1 + costet));
75
76 fLocalDirection.set(sintet*cos(phi), sintet*sin(phi), costet);
78
79 return fLocalDirection;
80}
static constexpr double twopi
Definition: G4SIunits.hh:56
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
G4double GetTotalEnergy() const
G4double GetTotalMomentum() const
float electron_mass_c2
Definition: hepunit.py:273

References source.hepunit::electron_mass_c2, G4VEmAngularDistribution::fLocalDirection, G4UniformRand, G4DynamicParticle::GetMomentumDirection(), G4DynamicParticle::GetTotalEnergy(), G4DynamicParticle::GetTotalMomentum(), CLHEP::Hep3Vector::rotateUz(), CLHEP::Hep3Vector::set(), and twopi.

◆ SampleDirectionForShell()

G4ThreeVector & G4VEmAngularDistribution::SampleDirectionForShell ( const G4DynamicParticle dp,
G4double  finalTotalEnergy,
G4int  Z,
G4int  shellID,
const G4Material mat 
)
virtualinherited

◆ 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

◆ 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: