Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DipBustGenerator Class Reference

#include <G4DipBustGenerator.hh>

Inheritance diagram for G4DipBustGenerator:
Collaboration diagram for G4DipBustGenerator:

Public Member Functions

 G4DipBustGenerator (const G4String &name="")
 
virtual ~G4DipBustGenerator ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=nullptr) final
 
G4double PolarAngle (const G4double initial_energy, const G4double final_energy, const G4int Z)
 
virtual void PrintGeneratorInformation () const final
 
- Public Member Functions inherited from G4VEmAngularDistribution
 G4VEmAngularDistribution (const G4String &name)
 
virtual ~G4VEmAngularDistribution ()
 
virtual G4ThreeVectorSampleDirectionForShell (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
 
const G4StringGetName () const
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmAngularDistribution
G4ThreeVector fLocalDirection
 

Detailed Description

Definition at line 55 of file G4DipBustGenerator.hh.

Constructor & Destructor Documentation

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

Definition at line 59 of file G4DipBustGenerator.cc.

60  : G4VEmAngularDistribution("DipBustGen")
61 {}
G4VEmAngularDistribution(const G4String &name)
G4DipBustGenerator::~G4DipBustGenerator ( )
virtual

Definition at line 63 of file G4DipBustGenerator.cc.

64 {}

Member Function Documentation

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

Definition at line 105 of file G4DipBustGenerator.cc.

108 {
109  G4double c, cosTheta, delta, cofA, signc = 1., a;
110  G4double gamma, beta, theta;
111 
112  c = 4. - 8.*G4UniformRand();
113  a = c;
114 
115  if( c < 0. )
116  {
117  signc = -1.;
118  a = -c;
119  }
120  delta = std::sqrt(a*a+4.);
121  delta += a;
122  delta *= 0.5;
123 
124  cofA = -signc*G4Exp(G4Log(delta)/3.0);
125 
126  cosTheta = cofA - 1./cofA;
127 
128  gamma = 1. + eTkin/electron_mass_c2;
129  beta = std::sqrt(1. - 1./gamma/gamma);
130 
131  cosTheta = (cosTheta + beta)/(1 + cosTheta*beta);
132 
133  theta = std::acos(cosTheta);
134 
135  if( theta < 0. ) theta = 0.;
136  if( theta > pi ) theta = pi;
137  // G4cout <<"theta = "<<theta<<"; ";
138 
139  return theta;
140 }
static constexpr double electron_mass_c2
#define G4UniformRand()
Definition: Randomize.hh:97
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4DipBustGenerator::PrintGeneratorInformation ( ) const
finalvirtual

Definition at line 142 of file G4DipBustGenerator.cc.

143 {
144  G4cout << "\n" << G4endl;
145  G4cout << "Angular Generator based on classical formula from" << G4endl;
146  G4cout << "J.D. Jackson, Classical Electrodynamics, Wiley, New York 1975"
147  << G4endl;
148 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4ThreeVector & G4DipBustGenerator::SampleDirection ( const G4DynamicParticle dp,
G4double  out_energy,
G4int  Z,
const G4Material mat = nullptr 
)
finalvirtual

Implements G4VEmAngularDistribution.

Definition at line 67 of file G4DipBustGenerator.cc.

69 {
70  G4double a, c, cosTheta, delta, cofA, signc = 1.;
71 
72  G4double eTkin = dp->GetKineticEnergy();
73 
74  c = 4. - 8.*G4UniformRand();
75  a = c;
76 
77  if( c < 0. )
78  {
79  signc = -1.;
80  a = -c;
81  }
82  delta = std::sqrt(a*a+4.);
83  delta += a;
84  delta *= 0.5;
85 
86  cofA = -signc*G4Exp(G4Log(delta)/3.0);
87 
88  cosTheta = cofA - 1./cofA;
89 
90  G4double tau = eTkin/electron_mass_c2;
91  G4double beta = std::sqrt(tau*(tau + 2.))/(tau + 1.);
92 
93  cosTheta = (cosTheta + beta)/(1 + cosTheta*beta);
94 
95  G4double sinTheta = std::sqrt((1 - cosTheta)*(1 + cosTheta));
96  G4double phi = twopi*G4UniformRand();
97 
98  fLocalDirection.set(sinTheta*std::cos(phi), sinTheta*std::sin(phi),cosTheta);
100 
101  return fLocalDirection;
102 
103 }
void set(double x, double y, double z)
G4double GetKineticEnergy() const
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double electron_mass_c2
#define G4UniformRand()
Definition: Randomize.hh:97
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
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: