Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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: G4ModifiedTsai.cc 74581 2013-10-15 12:03:25Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4ModifiedTsai
34 //
35 // Author: Andreia Trindade (andreia@lip.pt)
36 // Pedro Rodrigues (psilva@lip.pt)
37 // Luis Peralta (luis@lip.pt)
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"
60 #include "G4PhysicalConstants.hh"
61 #include "Randomize.hh"
62 #include "G4Log.hh"
63 
65  : G4VEmAngularDistribution("AngularGenUrban")
66 {}
67 
69 {}
70 
73  G4double, G4int, const G4Material*)
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  } while(u > uMax);
93 
94  G4double cost = 1.0 - 2*u*u/(uMax*uMax);
95  G4double sint = std::sqrt((1 - cost)*(1 + cost));
96  G4double phi = CLHEP::twopi*G4UniformRand();
97 
98  fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi), cost);
100 
101  return fLocalDirection;
102 }
103 
105 {
106  G4cout << "\n" << G4endl;
107  G4cout << "Bremsstrahlung Angular Generator is Modified Tsai" << G4endl;
108  G4cout << "Distribution suggested by L.Urban (Geant3 manual (1993) Phys211)"
109  << G4endl;
110  G4cout << "Derived from Tsai distribution (Rev Mod Phys 49,421(1977)) \n"
111  << G4endl;
112 }
void set(double x, double y, double z)
void PrintGeneratorInformation() const
G4double GetKineticEnergy() const
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
float electron_mass_c2
Definition: hepunit.py:274
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4ModifiedTsai(const G4String &name="")
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=0)
virtual ~G4ModifiedTsai()
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76