Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Generator2BS Class Reference

#include <G4Generator2BS.hh>

Inheritance diagram for G4Generator2BS:

Public Member Functions

 G4Generator2BS (const G4String &name="")
virtual ~G4Generator2BS ()
G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=nullptr) override
void PrintGeneratorInformation () const override
G4Generator2BSoperator= (const G4Generator2BS &right)=delete
 G4Generator2BS (const G4Generator2BS &)=delete
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 *)
virtual void SamplePairDirections (const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr)
const G4StringGetName () const
G4VEmAngularDistributionoperator= (const G4VEmAngularDistribution &right)=delete
 G4VEmAngularDistribution (const G4VEmAngularDistribution &)=delete

Protected Member Functions

G4double RejectionFunction (G4double value) const

Additional Inherited Members

Protected Attributes inherited from G4VEmAngularDistribution
G4ThreeVector fLocalDirection
G4bool fPolarisation

Detailed Description

Definition at line 64 of file G4Generator2BS.hh.

Constructor & Destructor Documentation

◆ G4Generator2BS() [1/2]

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

Definition at line 63 of file G4Generator2BS.cc.

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

Referenced by G4Generator2BS(), and operator=().

◆ ~G4Generator2BS()

G4Generator2BS::~G4Generator2BS ( )
virtual

Definition at line 71 of file G4Generator2BS.cc.

72{}

◆ G4Generator2BS() [2/2]

G4Generator2BS::G4Generator2BS ( const G4Generator2BS & )
delete

Member Function Documentation

◆ operator=()

G4Generator2BS & G4Generator2BS::operator= ( const G4Generator2BS & right)
delete

◆ PrintGeneratorInformation()

void G4Generator2BS::PrintGeneratorInformation ( ) const
overridevirtual

Reimplemented from G4VEmAngularDistribution.

Definition at line 134 of file G4Generator2BS.cc.

135{
136 G4cout << "\n" << G4endl;
137 G4cout << "Bremsstrahlung Angular Generator is 2BS Generator "
138 << "from 2BS Koch & Motz distribution (Rev Mod Phys 31(4), 920 (1959))" << G4endl;
139 G4cout << "Sampling algorithm adapted from PIRS-0203" << G4endl;
140 G4cout << "\n" << G4endl;
141}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout

◆ RejectionFunction()

G4double G4Generator2BS::RejectionFunction ( G4double value) const
inlineprotected

Definition at line 99 of file G4Generator2BS.hh.

100{
101 G4double y2 = (1 + y)*(1 + y);
102 G4double x = 4*y*ratio/y2;
103 return 4*x - ratio1 - (ratio2 - x)*G4Log(delta + fz/y2);
104}
G4double G4Log(G4double x)
Definition G4Log.hh:169
double G4double
Definition G4Types.hh:83

Referenced by SampleDirection().

◆ SampleDirection()

G4ThreeVector & G4Generator2BS::SampleDirection ( const G4DynamicParticle * dp,
G4double out_energy,
G4int Z,
const G4Material * mat = nullptr )
overridevirtual

Implements G4VEmAngularDistribution.

Definition at line 74 of file G4Generator2BS.cc.

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

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