Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DeltaAngleFreeScat.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: G4DeltaAngleFreeScat
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
51#include "Randomize.hh"
52
53using namespace std;
54
56 : G4VEmAngularDistribution("deltaFree")
57{}
58
60
63 G4double kinEnergyFinal, G4int,
64 const G4Material*)
65{
66 G4double deltaMomentum =
67 sqrt(kinEnergyFinal*(kinEnergyFinal + 2*electron_mass_c2));
68
69 G4double costet = kinEnergyFinal*(dp->GetTotalEnergy() + electron_mass_c2) /
70 (deltaMomentum * dp->GetTotalMomentum());
71
72 G4double phi = G4UniformRand()*twopi;
73 G4double sintet = sqrt((1 - costet)*(1 + costet));
74
75 fLocalDirection.set(sintet*cos(phi), sintet*sin(phi), costet);
77
78 return fLocalDirection;
79}
80
82{}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
G4DeltaAngleFreeScat(const G4String &name="")
~G4DeltaAngleFreeScat() override
G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, const G4Material *mat=nullptr) final
void PrintGeneratorInformation() const final
const G4ThreeVector & GetMomentumDirection() const
G4double GetTotalEnergy() const
G4double GetTotalMomentum() const