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

#include <G4ModifiedTsai.hh>

Inheritance diagram for G4ModifiedTsai:
Collaboration diagram for G4ModifiedTsai:

Public Member Functions

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

Additional Inherited Members

- Protected Attributes inherited from G4VEmAngularDistribution
G4ThreeVector fLocalDirection
 

Detailed Description

Definition at line 63 of file G4ModifiedTsai.hh.

Constructor & Destructor Documentation

G4ModifiedTsai::G4ModifiedTsai ( const G4String name = "")
explicit

Definition at line 64 of file G4ModifiedTsai.cc.

65  : G4VEmAngularDistribution("AngularGenUrban")
66 {}
G4VEmAngularDistribution(const G4String &name)
G4ModifiedTsai::~G4ModifiedTsai ( )
virtual

Definition at line 68 of file G4ModifiedTsai.cc.

69 {}

Member Function Documentation

void G4ModifiedTsai::PrintGeneratorInformation ( ) const
finalvirtual

Definition at line 105 of file G4ModifiedTsai.cc.

106 {
107  G4cout << "\n" << G4endl;
108  G4cout << "Bremsstrahlung Angular Generator is Modified Tsai" << G4endl;
109  G4cout << "Distribution suggested by L.Urban (Geant3 manual (1993) Phys211)"
110  << G4endl;
111  G4cout << "Derived from Tsai distribution (Rev Mod Phys 49,421(1977)) \n"
112  << G4endl;
113 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4ThreeVector & G4ModifiedTsai::SampleDirection ( const G4DynamicParticle dp,
G4double  out_energy,
G4int  Z,
const G4Material mat = nullptr 
)
finalvirtual

Implements G4VEmAngularDistribution.

Definition at line 72 of file G4ModifiedTsai.cc.

74 {
75  // Sample gamma angle (Z - axis along the parent particle).
76  // Universal distribution suggested by L. Urban (Geant3 manual (1993)
77  // Phys211) derived from Tsai distribution (Rev Mod Phys 49,421(1977))
78 
79  G4double uMax = 2*(1. + dp->GetKineticEnergy()/electron_mass_c2);
80 
81  static const G4double a1 = 0.625;
82  static const G4double a2 = 1.875;
83  static const G4double border = 0.25;
84  G4double u;
85 
86  do {
88 
89  if ( border > G4UniformRand() ) { u /= a1; }
90  else { u /= a2; }
91 
92  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
93  } while(u > uMax);
94 
95  G4double cost = 1.0 - 2*u*u/(uMax*uMax);
96  G4double sint = std::sqrt((1 - cost)*(1 + cost));
98 
99  fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi), cost);
101 
102  return fLocalDirection;
103 }
void set(double x, double y, double z)
static unsigned border[]
Definition: csz_inflate.cc:298
G4double GetKineticEnergy() const
static constexpr double electron_mass_c2
#define G4UniformRand()
Definition: Randomize.hh:97
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
G4double G4Log(G4double x)
Definition: G4Log.hh:230
double G4double
Definition: G4Types.hh:76
static constexpr double twopi
Definition: SystemOfUnits.h:55

Here is the call graph for this function:


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