Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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: G4DipBustGenerator.cc 74581 2013-10-15 12:03:25Z gcosmo $
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"
54 #include "G4PhysicalConstants.hh"
55 #include "Randomize.hh"
56 #include "G4Log.hh"
57 #include "G4Exp.hh"
58 
60  : G4VEmAngularDistribution("DipBustGen")
61 {}
62 
64 {}
65 
68  G4double, G4int, const G4Material*)
69 {
70  G4double a, c, cosTheta, delta, cofA, signc = 1.;
71 
72  G4double eTkin = dp->GetKineticEnergy();
73 
74  c = 4. - 8.*G4UniformRand();
75  a = c;
76 
77  if( c < 0. )
78  {
79  signc = -1.;
80  a = -c;
81  }
82  delta = std::sqrt(a*a+4.);
83  delta += a;
84  delta *= 0.5;
85 
86  cofA = -signc*G4Exp(G4Log(delta)/3.0);
87 
88  cosTheta = cofA - 1./cofA;
89 
90  G4double tau = eTkin/electron_mass_c2;
91  G4double beta = std::sqrt(tau*(tau + 2.))/(tau + 1.);
92 
93  cosTheta = (cosTheta + beta)/(1 + cosTheta*beta);
94 
95  G4double sinTheta = std::sqrt((1 - cosTheta)*(1 + cosTheta));
96  G4double phi = twopi*G4UniformRand();
97 
98  fLocalDirection.set(sinTheta*std::cos(phi), sinTheta*std::sin(phi),cosTheta);
100 
101  return fLocalDirection;
102 
103 }
104 
106  const G4double, // final_energy
107  const G4int ) // Z
108 {
109  G4double c, cosTheta, delta, cofA, signc = 1., a;
110  G4double gamma, beta, theta;
111 
112  c = 4. - 8.*G4UniformRand();
113  a = c;
114 
115  if( c < 0. )
116  {
117  signc = -1.;
118  a = -c;
119  }
120  delta = std::sqrt(a*a+4.);
121  delta += a;
122  delta *= 0.5;
123 
124  cofA = -signc*G4Exp(G4Log(delta)/3.0);
125 
126  cosTheta = cofA - 1./cofA;
127 
128  gamma = 1. + eTkin/electron_mass_c2;
129  beta = std::sqrt(1. - 1./gamma/gamma);
130 
131  cosTheta = (cosTheta + beta)/(1 + cosTheta*beta);
132 
133  theta = std::acos(cosTheta);
134 
135  if( theta < 0. ) theta = 0.;
136  if( theta > pi ) theta = pi;
137  // G4cout <<"theta = "<<theta<<"; ";
138 
139  return theta;
140 }
141 
143 {
144  G4cout << "\n" << G4endl;
145  G4cout << "Angular Generator based on classical formula from" << G4endl;
146  G4cout << "J.D. Jackson, Classical Electrodynamics, Wiley, New York 1975"
147  << G4endl;
148 }
void set(double x, double y, double z)
G4double GetKineticEnergy() const
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
int G4int
Definition: G4Types.hh:78
static constexpr double twopi
Definition: G4SIunits.hh:76
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
virtual void PrintGeneratorInformation() const final
float electron_mass_c2
Definition: hepunit.py:274
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4DipBustGenerator(const G4String &name="")
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=nullptr) final
#define G4endl
Definition: G4ios.hh:61
static constexpr double pi
Definition: G4SIunits.hh:75
G4double PolarAngle(const G4double initial_energy, const G4double final_energy, const G4int Z)
double G4double
Definition: G4Types.hh:76
tuple c
Definition: test.py:13