Geant4 9.6.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
 
const G4StringGetName () const
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmAngularDistribution
G4ThreeVector fLocalDirection
 

Detailed Description

Definition at line 61 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:64

◆ ~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 G4double cost;
86
87 G4double n0 = PP6[Z] - 1.0;
88 G4double n1 = PP7[Z] - 1.0;
89 G4double n2 = PP8[Z] - 1.0;
90 G4double b0 = PP3[Z];
91 G4double b1 = PP4[Z];
92 G4double b2 = PP5[Z];
93 G4double w0 = 0.0;
94 G4double w1 = 0.0;
95 G4double w2 = 0.0;
96
97 const G4double numlim = 0.02;
98 G4double x = 2*xx*b0;
99 if(x < numlim) { w0 = n0*x*(1 - 0.5*(n0 - 1)*x*(1 - (n0 - 2)*x/3.)); }
100 else { w0 = (1 - std::pow(1 + x,-n0)); }
101
102 if(PP1[Z] > 0.0) {
103 x = 2*xx*b1;
104 if(x < numlim) { w1 = n1*x*(1 - 0.5*(n1 - 1)*x*(1 - (n1 - 2)*x/3.)); }
105 else { w1 = (1 - std::pow(1 + x,-n1)); }
106 }
107 if(PP2[Z] > 0.0) {
108 x = 2*xx*b2;
109 if(x < numlim) { w2 = n2*x*(1 - 0.5*(n2 - 1)*x*(1 - (n2 - 2)*x/3.)); }
110 else { w2 = (1 - std::pow(1 + x,-n2)); }
111 }
112
113 G4double x0= w0*PP0[Z]/(b0*n0);
114 G4double x1= w1*PP1[Z]/(b1*n1);
115 G4double x2= w2*PP2[Z]/(b2*n2);
116
117 do {
118
119 G4double w = w0;
120 G4double n = n0;
121 G4double b = b0;
122
123 x = G4UniformRand()*(x0 + x1 + x2);
124 if(x > x0) {
125 x -= x0;
126 if(x <= x1 ) {
127 w = w1;
128 n = n1;
129 b = b1;
130 } else {
131 w = w2;
132 n = n2;
133 b = b2;
134 }
135 }
136 n = 1.0/n;
137
138 // sampling of angle
139 G4double y = w*G4UniformRand();
140 if(y < numlim) { x = y*n*( 1 + 0.5*(n + 1)*y*(1 - (n + 2)*y/3.)); }
141 else { x = 1.0/std::pow(1 - y, n) - 1.0; }
142 cost = 1.0 - x/(b*xx);
143 //G4cout << "cost = " << cost << " w= " << w << " n= " << n
144 // << " b= " << b << " x= " << x << " xx= " << xx << G4endl;
145 } while (2*G4UniformRand() > 1.0 + cost*cost || cost < -1.0);
146
147 G4double phi = twopi*G4UniformRand();
148 G4double sint = sqrt((1. - cost)*(1.0 + cost));
149 fLocalDirection.set(sint*cos(phi),sint*sin(phi),cost);
151
152 return fLocalDirection;
153}
#define G4UniformRand()
Definition: Randomize.hh:53
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const

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