Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4RayleighAngularGenerator Class Reference

#include <G4RayleighAngularGenerator.hh>

Inheritance diagram for G4RayleighAngularGenerator:
G4VEmAngularDistribution

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 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 61 of file G4RayleighAngularGenerator.hh.

Constructor & Destructor Documentation

G4RayleighAngularGenerator::G4RayleighAngularGenerator ( )

Definition at line 66 of file G4RayleighAngularGenerator.cc.

References python.hepunit::c_light, python.hepunit::cm, python.hepunit::h_Planck, and test::x.

67  : G4VEmAngularDistribution("CullenGenerator")
68 {
70  fFactor = 0.5*x*x;
71 }
float h_Planck
Definition: hepunit.py:263
G4VEmAngularDistribution(const G4String &name)
double G4double
Definition: G4Types.hh:76
float c_light
Definition: hepunit.py:257
G4RayleighAngularGenerator::~G4RayleighAngularGenerator ( )
virtual

Definition at line 75 of file G4RayleighAngularGenerator.cc.

76 {}

Member Function Documentation

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

Implements G4VEmAngularDistribution.

Definition at line 81 of file G4RayleighAngularGenerator.cc.

References test::b, G4VEmAngularDistribution::fLocalDirection, G4Exp(), G4Log(), G4UniformRand, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), n, CLHEP::Hep3Vector::rotateUz(), CLHEP::Hep3Vector::set(), python.hepunit::twopi, and test::x.

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

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