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

#include <G4RDModifiedTsai.hh>

Inheritance diagram for G4RDModifiedTsai:
Collaboration diagram for G4RDModifiedTsai:

Public Member Functions

 G4RDModifiedTsai (const G4String &name)
 
 ~G4RDModifiedTsai ()
 
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 ()
 

Detailed Description

Definition at line 60 of file G4RDModifiedTsai.hh.

Constructor & Destructor Documentation

G4RDModifiedTsai::G4RDModifiedTsai ( const G4String name)

Definition at line 59 of file G4RDModifiedTsai.cc.

60 {;}
G4RDVBremAngularDistribution(const G4String &name)
G4RDModifiedTsai::~G4RDModifiedTsai ( )

Definition at line 64 of file G4RDModifiedTsai.cc.

65 {;}

Member Function Documentation

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

Implements G4RDVBremAngularDistribution.

Definition at line 69 of file G4RDModifiedTsai.cc.

72 {
73 
74  // Sample gamma angle (Z - axis along the parent particle).
75  // Universal distribution suggested by L. Urban (Geant3 manual (1993)
76  // Phys211) derived from Tsai distribution (Rev Mod Phys 49,421(1977))
77 
78  G4double totalEnergy = initial_energy + electron_mass_c2;
79 
80  const G4double a1 = 0.625, a2 = 3.*a1, d = 27.;
81  G4double u, theta = 0;
82 
83  do{
84  u = - std::log(G4UniformRand()*G4UniformRand());
85 
86  if (9./(9.+d) > G4UniformRand()) u /= a1;
87  else u /= a2;
88 
89  theta = u*electron_mass_c2/totalEnergy;
90  }while(u > (totalEnergy*pi/electron_mass_c2));
91 
92  return theta;
93 }
#define G4UniformRand()
Definition: Randomize.hh:97
float electron_mass_c2
Definition: hepunit.py:274
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76
void G4RDModifiedTsai::PrintGeneratorInformation ( ) const
virtual

Implements G4RDVBremAngularDistribution.

Definition at line 96 of file G4RDModifiedTsai.cc.

97 {
98 
99  G4cout << "\n" << G4endl;
100  G4cout << "Bremsstrahlung Angular Generator is Modified Tsai" << G4endl;
101  G4cout << "Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211)" << G4endl;
102  G4cout << "Derived from Tsai distribution (Rev Mod Phys 49,421(1977)) \n" << G4endl;
103 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

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