G4Generator2BS Class Reference

#include <G4Generator2BS.hh>

Inheritance diagram for G4Generator2BS:

G4VEmAngularDistribution

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

Protected Member Functions

G4double RejectionFunction (G4double value) const

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.

References G4Pow::GetInstance().

00065   : G4VEmAngularDistribution("AngularGen2BS"),fz(1),ratio(1),
00066     ratio1(1),ratio2(1),delta(0)
00067 {
00068   g4pow = G4Pow::GetInstance();
00069   nwarn = 0;
00070 }

G4Generator2BS::~G4Generator2BS (  )  [virtual]

Definition at line 72 of file G4Generator2BS.cc.

00073 {}


Member Function Documentation

void G4Generator2BS::PrintGeneratorInformation (  )  const

Definition at line 140 of file G4Generator2BS.cc.

References G4cout, and G4endl.

00141 {
00142   G4cout << "\n" << G4endl;
00143   G4cout << "Bremsstrahlung Angular Generator is 2BS Generator "
00144          << "from 2BS Koch & Motz distribution (Rev Mod Phys 31(4), 920 (1959))" << G4endl;
00145   G4cout << "Sampling algorithm adapted from PIRS-0203" << G4endl;
00146   G4cout << "\n" << G4endl;
00147 } 

G4double G4Generator2BS::RejectionFunction ( G4double  value  )  const [inline, protected]

Definition at line 101 of file G4Generator2BS.hh.

Referenced by SampleDirection().

00102 {
00103   G4double y2 = (1 + y)*(1 + y);
00104   G4double x  = 4*y*ratio/y2;
00105   return 4*x - ratio1 - (ratio2 - x)*std::log(delta + fz/y2); 
00106 }

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.

References G4VEmAngularDistribution::fLocalDirection, G4cout, G4endl, G4UniformRand, G4DynamicParticle::GetMomentumDirection(), G4DynamicParticle::GetTotalEnergy(), RejectionFunction(), and G4Pow::Z13().

Referenced by G4Generator2BN::SampleDirection().

00079 {
00080 
00081   // Adapted from "Improved bremsstrahlung photon angular sampling in the EGS4 code system"
00082   // by Alex F. Bielajew, Rahde Mohan anc Chen-Shou Chui, PIRS-0203
00083   // Ionizing Radiation Standards
00084   // Institute for National Measurement Standards 
00085   // National Research Council of Canada
00086   // Departement of Medical Physics, Memorial Sloan-Kettering Cancer Center, New York
00087 
00088   G4double energy = dp->GetTotalEnergy();
00089   ratio = final_energy/energy;
00090   ratio1 = (1 + ratio)*(1 + ratio);
00091   ratio2 = 1 + ratio*ratio;
00092 
00093   G4double gamma = energy/electron_mass_c2;
00094   G4double beta  = std::sqrt((gamma - 1)*(gamma + 1))/gamma;
00095 
00096   //G4double Zeff = std::sqrt(static_cast<G4double>(Z) * (static_cast<G4double>(Z) + 1.0));
00097   //z = (0.00008116224*(std::pow(Zeff,0.3333333)));
00098 
00099   // VI speadup
00100   fz = 0.00008116224*g4pow->Z13(Z)*g4pow->Z13(Z+1);
00101 
00102   // majoranta
00103   G4double ymax = 2*beta*(1 + beta)*gamma*gamma;
00104   G4double gMax = RejectionFunction(0.0);
00105   gMax = std::max(gMax,RejectionFunction(ymax));
00106 
00107   G4double y, gfun;
00108 
00109   do{
00110     G4double q = G4UniformRand();
00111     y = q*ymax/(1 + ymax*(1 - q));
00112     gfun = RejectionFunction(y);
00113 
00114     // violation point
00115     if(gfun > gMax && nwarn >= 20) {
00116       ++nwarn;
00117       G4cout << "### WARNING in G4Generator2BS: Etot(MeV)= " << energy/MeV 
00118              << "  Egamma(MeV)" << (energy - final_energy)/MeV
00119              << " gMax= " << gMax << "  < " << gfun
00120              << "  results are not reliable!" 
00121              << G4endl;
00122       if(20 == nwarn) { 
00123         G4cout << "   WARNING in G4Generator2BS is closed" << G4endl; 
00124       }
00125     }
00126 
00127   } while(G4UniformRand()*gMax > gfun || y > ymax);
00128 
00129   
00130   G4double cost = 1 - 2*y/ymax;
00131   G4double sint = std::sqrt((1 - cost)*(1 + cost));
00132   G4double phi  = twopi*G4UniformRand(); 
00133 
00134   fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi),cost);
00135   fLocalDirection.rotateUz(dp->GetMomentumDirection());
00136 
00137   return fLocalDirection;
00138 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:52:04 2013 for Geant4 by  doxygen 1.4.7