Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DipBustGenerator.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$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4DipBustGenerator
34//
35// Author: Vladimir Grichine
36//
37// Creation date: 17 May 2011
38//
39// Modifications:
40//
41//
42// Class Description:
43//
44// Bremsstrahlung Angular Distribution Generation
45// suggested the dipole approximation in the rest frame of electron
46// busted in the laboratory frame.
47//
48// Class Description: End
49//
50// -------------------------------------------------------------------
51//
52
53#include "G4DipBustGenerator.hh"
55#include "Randomize.hh"
56
58 : G4VEmAngularDistribution("DipBustGen")
59{}
60
62{}
63
66 G4double, G4int, const G4Material*)
67{
68 G4double c, cosTheta, delta, cofA, signc = 1., a, power = 1./3.;
69
70 G4double eTkin = dp->GetKineticEnergy();
71
72 c = 4. - 8.*G4UniformRand();
73 a = c;
74
75 if( c < 0. )
76 {
77 signc = -1.;
78 a = -c;
79 }
80 delta = std::sqrt(a*a+4.);
81 delta += a;
82 delta *= 0.5;
83
84 cofA = -signc*std::pow(delta, power);
85
86 cosTheta = cofA - 1./cofA;
87
88 G4double tau = eTkin/electron_mass_c2;
89 G4double beta = std::sqrt(tau*(tau + 2.))/(tau + 1.);
90
91 cosTheta = (cosTheta + beta)/(1 + cosTheta*beta);
92
93 G4double sinTheta = std::sqrt((1 - cosTheta)*(1 + cosTheta));
94 G4double phi = twopi*G4UniformRand();
95
96 fLocalDirection.set(sinTheta*std::cos(phi), sinTheta*std::sin(phi),cosTheta);
98
99 return fLocalDirection;
100
101}
102
104 const G4double, // final_energy
105 const G4int ) // Z
106{
107 G4double c, cosTheta, delta, cofA, signc = 1., a, power = 1./3.;
108 G4double gamma, beta, theta;
109
110 c = 4. - 8.*G4UniformRand();
111 a = c;
112
113 if( c < 0. )
114 {
115 signc = -1.;
116 a = -c;
117 }
118 delta = std::sqrt(a*a+4.);
119 delta += a;
120 delta *= 0.5;
121
122 cofA = -signc*std::pow(delta, power);
123
124 cosTheta = cofA - 1./cofA;
125
126 gamma = 1. + eTkin/electron_mass_c2;
127 beta = std::sqrt(1. - 1./gamma/gamma);
128
129 cosTheta = (cosTheta + beta)/(1 + cosTheta*beta);
130
131 theta = std::acos(cosTheta);
132
133 if( theta < 0. ) theta = 0.;
134 if( theta > pi ) theta = pi;
135 // G4cout <<"theta = "<<theta<<"; ";
136
137 return theta;
138}
139
141{
142 G4cout << "\n" << G4endl;
143 G4cout << "Angular Generator based on classical formula from" << G4endl;
144 G4cout << "J.D. Jackson, Classical Electrodynamics, Wiley, New York 1975"
145 << G4endl;
146}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
G4double PolarAngle(const G4double initial_energy, const G4double final_energy, const G4int Z)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=0)
void PrintGeneratorInformation() const
G4DipBustGenerator(const G4String &name="")
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const