Geant4 11.1.1
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 ()
 
G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=nullptr) override
 
G4RayleighAngularGeneratoroperator= (const G4RayleighAngularGenerator &right)=delete
 
 G4RayleighAngularGenerator (const G4RayleighAngularGenerator &)=delete
 
- 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)
 
virtual void PrintGeneratorInformation () const
 
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() [1/2]

G4RayleighAngularGenerator::G4RayleighAngularGenerator ( )
explicit

Definition at line 63 of file G4RayleighAngularGenerator.cc.

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

◆ ~G4RayleighAngularGenerator()

G4RayleighAngularGenerator::~G4RayleighAngularGenerator ( )
virtual

Definition at line 72 of file G4RayleighAngularGenerator.cc.

73{}

◆ G4RayleighAngularGenerator() [2/2]

G4RayleighAngularGenerator::G4RayleighAngularGenerator ( const G4RayleighAngularGenerator )
delete

Member Function Documentation

◆ operator=()

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

◆ SampleDirection()

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

Implements G4VEmAngularDistribution.

Definition at line 78 of file G4RayleighAngularGenerator.cc.

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