Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
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
46#include "Randomize.hh"
47
49 : G4VEmAngularDistribution("AngularGenSauterGavrila")
50{}
51
53{}
54
57 const G4DynamicParticle* dp, G4double, G4int, const G4Material*)
58{
59 G4double tau = dp->GetKineticEnergy()/electron_mass_c2;
60 const G4double taulimit = 30.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
68 G4double invgamma = 1.0/(tau + 1.0);
69 G4double beta = std::sqrt(tau*(tau + 2.0))*invgamma;
70 G4double b = 0.5*tau*(tau*tau - 1.0);
71 G4double invgamma2 = invgamma*invgamma;
72
73 G4double rndm,term,greject,grejsup,costeta,sint2;
74 if (tau < 1.) { grejsup = (1.+b-beta*b)/invgamma2; }
75 else { grejsup = (1.+b+beta*b)/invgamma2; }
76
77 do {
78 rndm = 1 - 2*G4UniformRand();
79 costeta = (rndm + beta)/(rndm*beta + 1);
80 term = invgamma2/(1 + beta*rndm);
81 sint2 = (1 - costeta)*(1 + costeta);
82 greject = sint2*(1 + b*term)/(term*term);
83
84 } while(greject < G4UniformRand()*grejsup);
85
86 G4double sinteta = std::sqrt(sint2);
87 G4double phi = CLHEP::twopi*G4UniformRand();
88
89 fLocalDirection.set(sinteta*std::cos(phi), sinteta*std::sin(phi), costeta);
91 }
92 return fLocalDirection;
93}
94
96{
97 G4cout << "\n" << G4endl;
98 G4cout << "Non-polarized photoelectric effect angular generator." << G4endl;
99 G4cout << "The Sauter-Gavrila distribution for the K-shell is used."<<G4endl;
100 G4cout << "Originally developed by M.Maire for Geant3"
101 << G4endl;
102}
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 G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double e=0.0, G4int shellId=0, const G4Material *mat=0)