Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4RDGenerator2BS Class Reference

#include <G4RDGenerator2BS.hh>

Inheritance diagram for G4RDGenerator2BS:
Collaboration diagram for G4RDGenerator2BS:

Public Member Functions

 G4RDGenerator2BS (const G4String &name)
 
 ~G4RDGenerator2BS ()
 
G4double PolarAngle (const G4double initial_energy, const G4double final_energy, const G4int Z)
 
void PrintGeneratorInformation () const
 
- Public Member Functions inherited from G4RDVBremAngularDistribution
 G4RDVBremAngularDistribution (const G4String &name)
 
virtual ~G4RDVBremAngularDistribution ()
 

Protected Member Functions

G4double RejectionFunction (G4double value) const
 

Detailed Description

Definition at line 59 of file G4RDGenerator2BS.hh.

Constructor & Destructor Documentation

G4RDGenerator2BS::G4RDGenerator2BS ( const G4String name)

Definition at line 60 of file G4RDGenerator2BS.cc.

61 {;}
G4RDVBremAngularDistribution(const G4String &name)
G4RDGenerator2BS::~G4RDGenerator2BS ( )

Definition at line 65 of file G4RDGenerator2BS.cc.

66 {;}

Member Function Documentation

G4double G4RDGenerator2BS::PolarAngle ( const G4double  initial_energy,
const G4double  final_energy,
const G4int  Z 
)
virtual

Implements G4RDVBremAngularDistribution.

Definition at line 70 of file G4RDGenerator2BS.cc.

73 {
74 
75  // Adapted from "Improved bremsstrahlung photon angular sampling in the EGS4 code system"
76  // by Alex F. Bielajew, Rahde Mohan anc Chen-Shou Chui, PIRS-0203
77  // Ionizing Radiation Standards
78  // Institute for National Measurement Standards
79  // National Research Council of Canada
80  // Departement of Medical Physics, Memorial Sloan-Kettering Cancer Center, New York
81 
82 
83  G4double theta = 0;
84 
85  G4double initialTotalEnergy = (initial_energy+electron_mass_c2)/electron_mass_c2;
86  G4double finalTotalEnergy = (final_energy+electron_mass_c2)/electron_mass_c2;
87  EnergyRatio = finalTotalEnergy/initialTotalEnergy;
88  G4double gMaxEnergy = (pi*initialTotalEnergy)*(pi*initialTotalEnergy);
89 
90  G4double Zeff = std::sqrt(static_cast<G4double>(Z) * (static_cast<G4double>(Z) + 1.0));
91  z = (0.00008116224*(std::pow(Zeff,0.3333333)));
92 
93  // Rejection arguments
94  rejection_argument1 = (1.0+EnergyRatio*EnergyRatio);
95  rejection_argument2 = - 2.0*EnergyRatio + 3.0*rejection_argument1;
96  rejection_argument3 = ((1-EnergyRatio)/(2.0*initialTotalEnergy*EnergyRatio))*
97  ((1-EnergyRatio)/(2.0*initialTotalEnergy*EnergyRatio));
98 
99  // Calculate rejection function at 0, 1 and Emax
100  G4double gfunction0 = RejectionFunction(0);
101  G4double gfunction1 = RejectionFunction(1);
102  G4double gfunctionEmax = RejectionFunction(gMaxEnergy);
103 
104 
105  // Calculate Maximum value
106  G4double gMaximum = std::max(gfunction0,gfunction1);
107  gMaximum = std::max(gMaximum,gfunctionEmax);
108 
109  G4double rand, gfunctionTest, randTest;
110 
111  do{
112  rand = G4UniformRand();
113  rand = rand/(1-rand+1.0/gMaxEnergy);
114  gfunctionTest = RejectionFunction(rand);
115  randTest = G4UniformRand();
116 
117  }while(randTest > (gfunctionTest/gMaximum));
118 
119  theta = std::sqrt(rand)/initialTotalEnergy;
120 
121 
122  return theta;
123 }
#define G4UniformRand()
Definition: Randomize.hh:97
float electron_mass_c2
Definition: hepunit.py:274
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76
G4double RejectionFunction(G4double value) const

Here is the call graph for this function:

void G4RDGenerator2BS::PrintGeneratorInformation ( ) const
virtual

Implements G4RDVBremAngularDistribution.

Definition at line 138 of file G4RDGenerator2BS.cc.

139 {
140  G4cout << "\n" << G4endl;
141  G4cout << "Bremsstrahlung Angular Generator is 2BS Generator from 2BS Koch & Motz distribution (Rev Mod Phys 31(4), 920 (1959))" << G4endl;
142  G4cout << "Sampling algorithm adapted from PIRS-0203" << G4endl;
143  G4cout << "\n" << G4endl;
144 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4double G4RDGenerator2BS::RejectionFunction ( G4double  value) const
protected

Definition at line 126 of file G4RDGenerator2BS.cc.

127 {
128 
129  G4double argument = (1+value)*(1+value);
130 
131  G4double gfunction = (4+std::log(rejection_argument3+(z/argument)))*
132  ((4*EnergyRatio*value/argument)-rejection_argument1)+rejection_argument2;
133 
134  return gfunction;
135 
136 }
const XML_Char int const XML_Char * value
Definition: expat.h:331
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:


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