Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PhotoElectricAngularGeneratorSauterGavrila.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4PhotoElectricAngularGeneratorSauterGavrila
33//
34// Creation date: 10 May 2004
35//
36// Modifications:
37// 10 May 2003 P. Rodrigues First implementation acording with new design
38//
39// Class Description:
40//
41// Concrete class for PhotoElectric Electron Angular Distribution Generation
42// This model is a re-implementation of the Photolectric angular distribution
43// developed my M. Maire for the Standard EM Physics G4PhotoElectricEffect
44//
45// Class Description: End
46//
47// -------------------------------------------------------------------
48//
49
52#include "Randomize.hh"
53
55 G4VEmAngularDistribution("AngularGenSauterGavrilaLowE")
56{}
57
59{}
60
63 const G4DynamicParticle* dp,
64 G4double, G4int, const G4Material*)
65{
66
67 // Compute Theta distribution of the emitted electron, with respect to the
68 // incident Gamma.
69 // The Sauter-Gavrila distribution for the K-shell is used.
70
71 G4double costeta = 1.;
72 G4double Phi = twopi * G4UniformRand();
73 G4double cosphi = std::cos(Phi);
74 G4double sinphi = std::sin(Phi);
75 G4double sinteta = 0;
76 G4double gamma = 1. + dp->GetKineticEnergy()/electron_mass_c2;
77
78 if (gamma > 5.) {
80 return fLocalDirection;
81 // Bugzilla 1120
82 // SI on 05/09/2010 as suggested by JG 04/09/10
83 }
84
85 G4double beta = std::sqrt((gamma - 1)*(gamma + 1))/gamma;
86 G4double b = 0.5*gamma*(gamma - 1)*(gamma - 2);
87
88 G4double rndm,term,greject,grejsup;
89 if (gamma < 2.) grejsup = gamma*gamma*(1.+b-beta*b);
90 else grejsup = gamma*gamma*(1.+b+beta*b);
91
92 do { rndm = 1.-2*G4UniformRand();
93 costeta = (rndm+beta)/(rndm*beta+1.);
94 term = 1.-beta*costeta;
95 greject = (1.-costeta*costeta)*(1.+b*term)/(term*term);
96 } while(greject < G4UniformRand()*grejsup);
97
98 sinteta = std::sqrt((1 - costeta)*(1 + costeta));
99 fLocalDirection.set(sinteta*cosphi, sinteta*sinphi, costeta);
101 return fLocalDirection;
102}
103
105{
106 G4cout << "\n" << G4endl;
107 G4cout << "" << G4endl;
108 G4cout << "Re-implementation of the photolectric angular distribution" << G4endl;
109 G4cout << "developed my M. Maire for the Standard EM Physics G4PhotoElectricEffect" << G4endl;
110 G4cout << "It computes the theta distribution of the emitted electron, with respect to the" << G4endl;
111 G4cout << "incident Gamma, using the Sauter-Gavrila distribution for the K-shell\n" << G4endl;
112}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double e=0.0, G4int shellId=0, const G4Material *mat=0)