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

#include <G4RayleighAngularGenerator.hh>

+ Inheritance diagram for G4RayleighAngularGenerator:

Public Member Functions

 G4RayleighAngularGenerator ()
 
virtual ~G4RayleighAngularGenerator ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=0)
 
- 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 59 of file G4RayleighAngularGenerator.hh.

Constructor & Destructor Documentation

◆ G4RayleighAngularGenerator()

G4RayleighAngularGenerator::G4RayleighAngularGenerator ( )

Definition at line 64 of file G4RayleighAngularGenerator.cc.

65 : G4VEmAngularDistribution("CullenGenerator")
66{
67 G4double x = cm/(h_Planck*c_light);
68 fFactor = 0.5*x*x;
69}
double G4double
Definition: G4Types.hh:83

◆ ~G4RayleighAngularGenerator()

G4RayleighAngularGenerator::~G4RayleighAngularGenerator ( )
virtual

Definition at line 73 of file G4RayleighAngularGenerator.cc.

74{}

Member Function Documentation

◆ SampleDirection()

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

Implements G4VEmAngularDistribution.

Definition at line 79 of file G4RayleighAngularGenerator.cc.

82{
83 G4double ekin = dp->GetKineticEnergy();
84 G4double xx = fFactor*ekin*ekin;
85
86 G4double n0 = PP6[Z] - 1.0;
87 G4double n1 = PP7[Z] - 1.0;
88 G4double n2 = PP8[Z] - 1.0;
89 G4double b0 = PP3[Z];
90 G4double b1 = PP4[Z];
91 G4double b2 = PP5[Z];
92
93 static const G4double numlim = 0.02;
94 G4double x = 2.*xx*b0;
95 G4double w0 = (x < numlim) ? n0*x*(1.0 - 0.5*(n0 - 1.0)*x*(1.0 - (n0 - 2.0)*x/3.))
96 : 1.0 - G4Exp(-n0*G4Log(1.0 + x));
97
98 x = 2.*xx*b1;
99 G4double w1 = (x < numlim) ? n1*x*(1.0 - 0.5*(n1 - 1.0)*x*(1.0 - (n1 - 2.0)*x/3.))
100 : 1.0 - G4Exp(-n1*G4Log(1.0 + x));
101
102 x = 2.*xx*b2;
103 G4double w2 = (x < numlim) ? n2*x*(1.0 - 0.5*(n2 - 1.0)*x*(1.0 - (n2 - 2.0)*x/3.))
104 : 1.0 - G4Exp(-n2*G4Log(1.0 + x));
105
106 G4double x0= w0*PP0[Z]/(b0*n0);
107 G4double x1= w1*PP1[Z]/(b1*n1);
108 G4double x2= w2*PP2[Z]/(b2*n2);
109
110 G4double cost;
111 do {
112
113 G4double w = w0;
114 G4double n = n0;
115 G4double b = b0;
116
117 x = G4UniformRand()*(x0 + x1 + x2);
118 if(x > x0) {
119 x -= x0;
120 if(x <= x1 ) {
121 w = w1;
122 n = n1;
123 b = b1;
124 } else {
125 w = w2;
126 n = n2;
127 b = b2;
128 }
129 }
130 n = 1.0/n;
131
132 // sampling of angle
133 G4double y = w*G4UniformRand();
134 if(y < numlim) { x = y*n*( 1. + 0.5*(n + 1.)*y*(1. - (n + 2.)*y/3.)); }
135 //else { x = 1.0/std::pow(1 - y, n) - 1.0; }
136 else { x = G4Exp(-n*G4Log(1. - y)) - 1.0; }
137 cost = 1.0 - x/(b*xx);
138 //G4cout << "cost = " << cost << " w= " << w << " n= " << n
139 // << " b= " << b << " x= " << x << " xx= " << xx << G4endl;
140 } while (2*G4UniformRand() > 1.0 + cost*cost || cost < -1.0);
141
142 G4double phi = twopi*G4UniformRand();
143 G4double sint = sqrt((1. - cost)*(1.0 + cost));
144 fLocalDirection.set(sint*cos(phi),sint*sin(phi),cost);
146
147 return fLocalDirection;
148}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
#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 GetKineticEnergy() const

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