Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4Generator2BS Class Reference

#include <G4Generator2BS.hh>

Inheritance diagram for G4Generator2BS:
Collaboration diagram for G4Generator2BS:

Public Member Functions

 G4Generator2BS (const G4String &name="")
 
virtual ~G4Generator2BS ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=0)
 
void PrintGeneratorInformation () const
 
- 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
 

Protected Member Functions

G4double RejectionFunction (G4double value) const
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmAngularDistribution
G4ThreeVector fLocalDirection
 

Detailed Description

Definition at line 64 of file G4Generator2BS.hh.

Constructor & Destructor Documentation

G4Generator2BS::G4Generator2BS ( const G4String name = "")

Definition at line 64 of file G4Generator2BS.cc.

65  : G4VEmAngularDistribution("AngularGen2BS"),fz(1),ratio(1),
66  ratio1(1),ratio2(1),delta(0)
67 {
68  g4pow = G4Pow::GetInstance();
69  nwarn = 0;
70 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4VEmAngularDistribution(const G4String &name)

Here is the call graph for this function:

G4Generator2BS::~G4Generator2BS ( )
virtual

Definition at line 72 of file G4Generator2BS.cc.

73 {}

Member Function Documentation

void G4Generator2BS::PrintGeneratorInformation ( ) const

Definition at line 140 of file G4Generator2BS.cc.

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

Definition at line 101 of file G4Generator2BS.hh.

102 {
103  G4double y2 = (1 + y)*(1 + y);
104  G4double x = 4*y*ratio/y2;
105  return 4*x - ratio1 - (ratio2 - x)*std::log(delta + fz/y2);
106 }
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

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

Implements G4VEmAngularDistribution.

Definition at line 75 of file G4Generator2BS.cc.

79 {
80 
81  // Adapted from "Improved bremsstrahlung photon angular sampling in the EGS4 code system"
82  // by Alex F. Bielajew, Rahde Mohan anc Chen-Shou Chui, PIRS-0203
83  // Ionizing Radiation Standards
84  // Institute for National Measurement Standards
85  // National Research Council of Canada
86  // Departement of Medical Physics, Memorial Sloan-Kettering Cancer Center, New York
87 
89  ratio = final_energy/energy;
90  ratio1 = (1 + ratio)*(1 + ratio);
91  ratio2 = 1 + ratio*ratio;
92 
93  G4double gamma = energy/electron_mass_c2;
94  G4double beta = std::sqrt((gamma - 1)*(gamma + 1))/gamma;
95 
96  //G4double Zeff = std::sqrt(static_cast<G4double>(Z) * (static_cast<G4double>(Z) + 1.0));
97  //z = (0.00008116224*(std::pow(Zeff,0.3333333)));
98 
99  // VI speadup
100  fz = 0.00008116224*g4pow->Z13(Z)*g4pow->Z13(Z+1);
101 
102  // majoranta
103  G4double ymax = 2*beta*(1 + beta)*gamma*gamma;
104  G4double gMax = RejectionFunction(0.0);
105  gMax = std::max(gMax,RejectionFunction(ymax));
106 
107  G4double y, gfun;
108 
109  do{
110  G4double q = G4UniformRand();
111  y = q*ymax/(1 + ymax*(1 - q));
112  gfun = RejectionFunction(y);
113 
114  // violation point
115  if(gfun > gMax && nwarn >= 20) {
116  ++nwarn;
117  G4cout << "### WARNING in G4Generator2BS: Etot(MeV)= " << energy/MeV
118  << " Egamma(MeV)" << (energy - final_energy)/MeV
119  << " gMax= " << gMax << " < " << gfun
120  << " results are not reliable!"
121  << G4endl;
122  if(20 == nwarn) {
123  G4cout << " WARNING in G4Generator2BS is closed" << G4endl;
124  }
125  }
126 
127  } while(G4UniformRand()*gMax > gfun || y > ymax);
128 
129 
130  G4double cost = 1 - 2*y/ymax;
131  G4double sint = std::sqrt((1 - cost)*(1 + cost));
132  G4double phi = twopi*G4UniformRand();
133 
134  fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi),cost);
136 
137  return fLocalDirection;
138 }
void set(double x, double y, double z)
G4double RejectionFunction(G4double value) const
G4double GetTotalEnergy() const
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double electron_mass_c2
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4double Z13(G4int Z) const
Definition: G4Pow.hh:127
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double energy(const ThreeVector &p, const G4double m)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:


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