Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4SauterGavrilaAngularDistribution.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: G4SauterGavrilaAngularDistribution.cc 79188 2014-02-20 09:22:48Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4SauterGavrilaAngularDistribution
34 //
35 // Author: Vladimir Ivanchenko using Michel Maire algorithm
36 // developed for Geant3
37 //
38 // Creation date: 23 July 2012
39 //
40 //
41 // -------------------------------------------------------------------
42 //
43 
45 #include "G4PhysicalConstants.hh"
46 #include "Randomize.hh"
47 
49  : G4VEmAngularDistribution("AngularGenSauterGavrila")
50 {}
51 
53 {}
54 
57  const G4DynamicParticle* dp, G4double, G4int, const G4Material*)
58 {
60  static const G4double taulimit = 50.0;
61 
62  if (tau > taulimit) {
64  // Bugzilla 1120
65  // SI on 05/09/2010 as suggested by JG 04/09/10
66  } else {
67  // algorithm according Penelope 2008 manual and
68  // F.Sauter Ann. Physik 9, 217(1931); 11, 454(1931).
69 
70  G4double gamma = tau + 1;
71  G4double beta = std::sqrt(tau*(tau + 2))/gamma;
72  G4double A = (1 - beta)/beta;
73  G4double Ap2 = A + 2;
74  G4double B = 0.5*beta*gamma*(gamma - 1)*(gamma - 2);
75  G4double grej = 2*(1 + A*B)/A;
76  G4double z, g;
77  do {
78  G4double q = G4UniformRand();
79  z = 2*A*(2*q + Ap2*std::sqrt(q))/(Ap2*Ap2 - 4*q);
80  g = (2 - z)*(1.0/(A + z) + B);
81 
82  } while(g < G4UniformRand()*grej);
83 
84  G4double cost = 1 - z;
85  G4double sint = std::sqrt(z*(2 - z));
86  G4double phi = CLHEP::twopi*G4UniformRand();
87 
88  fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi), cost);
90  }
91  return fLocalDirection;
92 }
93 
95 {
96  G4cout << "\n" << G4endl;
97  G4cout << "Non-polarized photoelectric effect angular generator." << G4endl;
98  G4cout << "The Sauter-Gavrila distribution for the K-shell is used."<<G4endl;
99  G4cout << "Originally developed by M.Maire for Geant3"
100  << G4endl;
101 }
void set(double x, double y, double z)
G4double GetKineticEnergy() const
G4double z
Definition: TRTMaterials.hh:39
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double e=0.0, G4int shellId=0, const G4Material *mat=0)
int G4int
Definition: G4Types.hh:78
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
#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
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76