Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
63
65= default;
66
69 G4double secKinetic, G4int,
70 G4int,
71 const G4Material*)
72{
73 G4double k = dp->GetKineticEnergy();
74 G4double cosTheta = 1.0;
75
76 if (dp->GetDefinition() == fElectron)
77 {
78 if (secKinetic < 50.*eV) cosTheta = (2.*G4UniformRand())-1.;
79 else if (secKinetic <= 200.*eV)
80 {
81 if (G4UniformRand() <= 0.1) cosTheta = (2.*G4UniformRand())-1.;
82 else cosTheta = G4UniformRand()*(std::sqrt(2.)/2);
83 }
84 else
85 {
86 G4double sin2O =
87 (1.-secKinetic/k)/(1.+secKinetic/(2.*electron_mass_c2));
88 cosTheta = std::sqrt(1.-sin2O);
89 }
90 }
91 else
92 {
93 G4double tau = k/dp->GetDefinition()->GetPDGMass();
94 G4double maxSecKinetic = 2.* electron_mass_c2*tau*(tau + 2.0);
95
96 // Restriction below 100 eV from Emfietzoglou (2000)
97
98 if (secKinetic>100*eV) cosTheta = std::min(std::sqrt(secKinetic / maxSecKinetic), 1.0);
99 else cosTheta = (2.*G4UniformRand())-1.;
100 }
101
102 G4double sint = sqrt((1.0 - cosTheta)*(1.0 + cosTheta));
103 G4double phi = twopi*G4UniformRand();
104
105 fLocalDirection.set(sint*cos(phi), sint*sin(phi), cosTheta);
107
108 return fLocalDirection;
109}
110
113 G4double secEkin, G4int Z,
114 const G4Material* mat)
115{
116 return SampleDirectionForShell(dp, secEkin, Z, 0, mat);
117}
118
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
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
void PrintGeneratorInformation() const override
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition G4Electron.cc:91
G4VEmAngularDistribution(const G4String &name)