00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058 #include "G4Generator2BS.hh"
00059 #include "Randomize.hh"
00060 #include "G4PhysicalConstants.hh"
00061 #include "G4SystemOfUnits.hh"
00062 #include "G4Pow.hh"
00063
00064 G4Generator2BS::G4Generator2BS(const G4String&)
00065 : G4VEmAngularDistribution("AngularGen2BS"),fz(1),ratio(1),
00066 ratio1(1),ratio2(1),delta(0)
00067 {
00068 g4pow = G4Pow::GetInstance();
00069 nwarn = 0;
00070 }
00071
00072 G4Generator2BS::~G4Generator2BS()
00073 {}
00074
00075 G4ThreeVector& G4Generator2BS::SampleDirection(const G4DynamicParticle* dp,
00076 G4double final_energy,
00077 G4int Z,
00078 const G4Material*)
00079 {
00080
00081
00082
00083
00084
00085
00086
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
00097
00098
00099
00100 fz = 0.00008116224*g4pow->Z13(Z)*g4pow->Z13(Z+1);
00101
00102
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
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 }
00139
00140 void G4Generator2BS::PrintGeneratorInformation() const
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 }
00148