Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNARuddAngle.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: G4DNARuddAngle.cc 68380 2013-03-22 18:39:29Z vnivanch $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4DNARuddAngle
34 //
35 // Author: Vladimir Ivantcheko
36 //
37 // Creation date: 23 August 2013
38 //
39 // Modifications:
40 //
41 // Class Description:
42 //
43 // Delta-electron Angular Distribution Generation
44 //
45 // Class Description: End
46 //
47 // -------------------------------------------------------------------
48 //
49 
50 #include "G4DNARuddAngle.hh"
51 #include "G4PhysicalConstants.hh"
52 #include "Randomize.hh"
53 #include "G4ParticleDefinition.hh"
54 #include "G4Electron.hh"
55 #include "G4SystemOfUnits.hh"
56 
57 using namespace std;
58 
60  : G4VEmAngularDistribution("deltaRudd")
61 {
62  fElectron = G4Electron::Electron();
63 }
64 
66 {}
67 
70  G4double secKinetic, G4int,
71  G4int,
72  const G4Material*)
73 {
74  G4double k = dp->GetKineticEnergy();
75  G4double cosTheta = 1.0;
76 
77  const G4ParticleDefinition* particle = dp->GetDefinition();
78  G4double mass = particle->GetPDGMass();
79 
80  G4double maximumEnergyTransfer = k;
81  if(particle == fElectron) { maximumEnergyTransfer *= 0.5; }
82  else if(mass > MeV) {
83  G4double ratio = electron_mass_c2/mass;
84  G4double tau = k/mass;
85  maximumEnergyTransfer = 2.0*electron_mass_c2*tau*(tau + 2.) /
86  (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
87 
88  }
89 
90  if (secKinetic>100*eV && secKinetic <= maximumEnergyTransfer) {
91  cosTheta = std::sqrt(secKinetic / maximumEnergyTransfer);
92  } else {
93  cosTheta = (2.*G4UniformRand())-1.;
94  }
95 
96  G4double sint = sqrt((1.0 - cosTheta)*(1.0 + cosTheta));
97  G4double phi = twopi*G4UniformRand();
98 
99  fLocalDirection.set(sint*cos(phi), sint*sin(phi), cosTheta);
101 
102  return fLocalDirection;
103 }
104 
107  G4double secEkin, G4int Z,
108  const G4Material* mat)
109 {
110  return SampleDirectionForShell(dp, secEkin, Z, 0, mat);
111 }
112 
114 {}
void set(double x, double y, double z)
G4double GetKineticEnergy() const
G4DNARuddAngle(const G4String &name="")
G4ParticleDefinition * GetDefinition() const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, const G4Material *mat=0)
int G4int
Definition: G4Types.hh:78
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double electron_mass_c2
#define G4UniformRand()
Definition: Randomize.hh:97
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
static constexpr double eV
Definition: G4SIunits.hh:215
virtual ~G4DNARuddAngle()
G4double GetPDGMass() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
virtual G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, G4int shellIdx, const G4Material *mat=0)
void PrintGeneratorInformation() const