Geant4-11
G4DNABornAngle.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: G4DNABornAngle
33//
34// Author: Vladimir Ivantcheko
35//
36// Creation date: 23 August 2013
37//
38// Modifications:
39//
40// Class Description:
41//
42// Delta-electron Angular Distribution Generation
43//
44// Class Description: End
45//
46// -------------------------------------------------------------------
47//
48
49#include "G4DNABornAngle.hh"
51#include "Randomize.hh"
53#include "G4Electron.hh"
54#include "G4SystemOfUnits.hh"
55
56using namespace std;
57
59 : G4VEmAngularDistribution("deltaBorn")
60{
62}
63
65{}
66
69 G4double secKinetic, G4int,
70 G4int,
71 const G4Material*)
72{
73 G4double k = dp->GetKineticEnergy();
74 G4double cosTheta = 1.0;
75 if (dp->GetDefinition() == fElectron)
76 {
77 if (secKinetic < 50.*eV) cosTheta = (2.*G4UniformRand())-1.;
78 else if (secKinetic <= 200.*eV)
79 {
80 if (G4UniformRand() <= 0.1) cosTheta = (2.*G4UniformRand())-1.;
81 else cosTheta = G4UniformRand()*(std::sqrt(2.)/2);
82 }
83 else
84 {
85 G4double sin2O =
86 (1.-secKinetic/k)/(1.+secKinetic/(2.*electron_mass_c2));
87 cosTheta = std::sqrt(1.-sin2O);
88 }
89 }
90 else
91 {
92 G4double mass = dp->GetDefinition()->GetPDGMass();
93 G4double maxSecKinetic = 4.* (electron_mass_c2 / mass) * k;
94
95 // cosTheta = std::sqrt(secKinetic / maxSecKinetic);
96
97 // Restriction below 100 eV from Emfietzoglou (2000)
98
99 if (secKinetic>100*eV) cosTheta = std::sqrt(secKinetic / maxSecKinetic);
100 else cosTheta = (2.*G4UniformRand())-1.;
101 }
102
103 G4double sint = sqrt((1.0 - cosTheta)*(1.0 + cosTheta));
105
106 fLocalDirection.set(sint*cos(phi), sint*sin(phi), cosTheta);
108
109 return fLocalDirection;
110}
111
114 G4double secEkin, G4int Z,
115 const G4Material* mat)
116{
117 return SampleDirectionForShell(dp, secEkin, Z, 0, mat);
118}
119
121{}
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double eV
Definition: G4SIunits.hh:201
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4UniformRand()
Definition: Randomize.hh:52
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, const G4Material *mat=nullptr) override
G4DNABornAngle(const G4String &name="")
~G4DNABornAngle() override
G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, G4int shellIdx, const G4Material *mat=nullptr) override
const G4ParticleDefinition * fElectron
void PrintGeneratorInformation() const override
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
float electron_mass_c2
Definition: hepunit.py:273