Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ModifiedTsai.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4ModifiedTsai
34//
35// Author: Andreia Trindade ([email protected])
36// Pedro Rodrigues ([email protected])
37// Luis Peralta ([email protected])
38//
39// Creation date: 21 March 2003
40//
41// Modifications:
42// 21 Mar 2003 A.Trindade First implementation acording with new design
43// 24 Mar 2003 A.Trindade Fix in Tsai generator in order to prevent theta
44// generation above pi
45// 13 Oct 2010 V.Ivanchenko Moved to standard and improved comment
46// 26.04.2011 V.Grichine Clean-up of PolarAngle method
47//
48// Class Description:
49//
50// Bremsstrahlung Angular Distribution Generation
51// suggested by L.Urban (Geant3 manual (1993) Phys211)
52// Derived from Tsai distribution (Rev Mod Phys 49,421(1977))
53//
54// Class Description: End
55//
56// -------------------------------------------------------------------
57//
58
59#include "G4ModifiedTsai.hh"
61#include "Randomize.hh"
62
64 : G4VEmAngularDistribution("AngularGenUrban")
65{}
66
68{}
69
72 G4double, G4int, const G4Material*)
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 uMax = 2*(1. + dp->GetKineticEnergy()/electron_mass_c2);
79
80 const G4double a1 = 0.625;
81 const G4double a2 = 1.875;
82 const G4double border = 0.25;
83 G4double u;
84
85 do {
86 u = - std::log(G4UniformRand()*G4UniformRand());
87
88 if ( border > G4UniformRand() ) { u /= a1; }
89 else { u /= a2; }
90
91 } while(u > uMax);
92
93 G4double cost = 1.0 - 2*u*u/(uMax*uMax);
94 G4double sint = std::sqrt((1 - cost)*(1 + cost));
95 G4double phi = CLHEP::twopi*G4UniformRand();
96
97 fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi), cost);
99
100 return fLocalDirection;
101}
102
104{
105 G4cout << "\n" << G4endl;
106 G4cout << "Bremsstrahlung Angular Generator is Modified Tsai" << G4endl;
107 G4cout << "Distribution suggested by L.Urban (Geant3 manual (1993) Phys211)"
108 << G4endl;
109 G4cout << "Derived from Tsai distribution (Rev Mod Phys 49,421(1977)) \n"
110 << G4endl;
111}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
virtual ~G4ModifiedTsai()
G4ModifiedTsai(const G4String &name="")
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=0)
void PrintGeneratorInformation() const