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

#include <G4DipBustGenerator.hh>

Inheritance diagram for G4DipBustGenerator:
G4VEmAngularDistribution

Public Member Functions

 G4DipBustGenerator (const G4DipBustGenerator &)=delete
 
 G4DipBustGenerator (const G4String &name="")
 
const G4StringGetName () const
 
G4DipBustGeneratoroperator= (const G4DipBustGenerator &right)=delete
 
G4double PolarAngle (G4double initial_energy, G4double final_energy, G4int Z)
 
void PrintGeneratorInformation () const final
 
G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=nullptr) final
 
virtual G4ThreeVectorSampleDirectionForShell (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
 
void SamplePairDirections (const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr) final
 
 ~G4DipBustGenerator () override
 

Protected Attributes

G4ThreeVector fLocalDirection
 
G4bool fPolarisation
 

Private Member Functions

G4double SampleCosTheta (G4double kinEnergy)
 

Private Attributes

G4String fName
 

Detailed Description

Definition at line 54 of file G4DipBustGenerator.hh.

Constructor & Destructor Documentation

◆ G4DipBustGenerator() [1/2]

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

Definition at line 62 of file G4DipBustGenerator.cc.

63 : G4VEmAngularDistribution("DipBustGen")
64{}
G4VEmAngularDistribution(const G4String &name)

◆ ~G4DipBustGenerator()

G4DipBustGenerator::~G4DipBustGenerator ( )
override

Definition at line 68 of file G4DipBustGenerator.cc.

69{}

◆ G4DipBustGenerator() [2/2]

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

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

◆ PolarAngle()

G4double G4DipBustGenerator::PolarAngle ( G4double  initial_energy,
G4double  final_energy,
G4int  Z 
)

Definition at line 106 of file G4DipBustGenerator.cc.

109{
110 const G4double cosTheta = SampleCosTheta(eTkin);
111 G4double theta = std::acos(cosTheta);
112 theta = std::min(std::max(theta, 0.), CLHEP::pi);
113 return theta;
114}
double G4double
Definition: G4Types.hh:83
G4double SampleCosTheta(G4double kinEnergy)
static constexpr double pi
Definition: SystemOfUnits.h:55
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References G4INCL::Math::max(), G4INCL::Math::min(), CLHEP::pi, and SampleCosTheta().

◆ PrintGeneratorInformation()

void G4DipBustGenerator::PrintGeneratorInformation ( ) const
finalvirtual

Reimplemented from G4VEmAngularDistribution.

Definition at line 144 of file G4DipBustGenerator.cc.

145{
146 G4cout << "\n" << G4endl;
147 G4cout << "Angular Generator based on classical formula from" << G4endl;
148 G4cout << "J.D. Jackson, Classical Electrodynamics, Wiley, New York 1975"
149 << G4endl;
150}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

References G4cout, and G4endl.

◆ SampleCosTheta()

G4double G4DipBustGenerator::SampleCosTheta ( G4double  kinEnergy)
private

Definition at line 73 of file G4DipBustGenerator.cc.

74{
75 const G4double c = 4. - 8.*G4UniformRand();
76 const G4double delta = 0.5*(std::sqrt(c*c+4.) + std::abs(c));
77 const G4double signc = (c < 0.) ? -1.0 : 1.0;
78 // const G4double cofA = -signc*G4Exp(G4Log(delta)/3.0);
79 const G4double cofA = -signc*G4Pow::GetInstance()->A13(delta);
80 const G4double cosTheta = std::min(1.,std::max(-1.,cofA - 1./cofA));
81 const G4double tau = kinEnergy/CLHEP::electron_mass_c2;
82 const G4double beta = std::sqrt(tau*(tau + 2.))/(tau + 1.);
83
84 return (cosTheta + beta)/(1. + cosTheta*beta);
85}
#define G4UniformRand()
Definition: Randomize.hh:52
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double A13(G4double A) const
Definition: G4Pow.cc:120
static constexpr double electron_mass_c2

References G4Pow::A13(), anonymous_namespace{G4PionRadiativeDecayChannel.cc}::beta, CLHEP::electron_mass_c2, G4UniformRand, G4Pow::GetInstance(), G4INCL::Math::max(), and G4INCL::Math::min().

Referenced by PolarAngle(), SampleDirection(), and SamplePairDirections().

◆ SampleDirection()

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

Implements G4VEmAngularDistribution.

Definition at line 90 of file G4DipBustGenerator.cc.

92{
93 const G4double cosTheta = SampleCosTheta(dp->GetKineticEnergy());
94
95 const G4double sinTheta = std::sqrt((1. - cosTheta)*(1. + cosTheta));
97
98 fLocalDirection.set(sinTheta*std::cos(phi), sinTheta*std::sin(phi),cosTheta);
100
101 return fLocalDirection;
102}
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static constexpr double twopi
Definition: SystemOfUnits.h:56

References G4VEmAngularDistribution::fLocalDirection, G4UniformRand, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), CLHEP::Hep3Vector::rotateUz(), SampleCosTheta(), CLHEP::Hep3Vector::set(), and CLHEP::twopi.

◆ SampleDirectionForShell()

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

◆ SamplePairDirections()

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

Reimplemented from G4VEmAngularDistribution.

Definition at line 118 of file G4DipBustGenerator.cc.

124{
125 const G4double phi = CLHEP::twopi * G4UniformRand();
126 const G4double sinp = std::sin(phi);
127 const G4double cosp = std::cos(phi);
128
129 G4double cost = SampleCosTheta(elecKinEnergy);
130 G4double sint = std::sqrt((1. - cost)*(1. + cost));
131
132 dirElectron.set(sint*cosp, sint*sinp, cost);
133 dirElectron.rotateUz(dp->GetMomentumDirection());
134
135 cost = SampleCosTheta(posiKinEnergy);
136 sint = std::sqrt((1. - cost)*(1. + cost));
137
138 dirPositron.set(-sint*cosp, -sint*sinp, cost);
139 dirPositron.rotateUz(dp->GetMomentumDirection());
140}

References G4UniformRand, G4DynamicParticle::GetMomentumDirection(), CLHEP::Hep3Vector::rotateUz(), SampleCosTheta(), CLHEP::Hep3Vector::set(), and CLHEP::twopi.

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: