Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4DipBustGenerator.cc
이 파일의 문서화 페이지로 가기
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 110415 2018-05-23 06:44:31Z 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 "Randomize.hh"
55 #include "G4Log.hh"
56 #include "G4Exp.hh"
58 
59 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
60 
62  : G4VEmAngularDistribution("DipBustGen")
63 {}
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66 
68 {}
69 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
71 
73 {
74  G4double c = 4. - 8.*G4UniformRand();
75  G4double delta = 0.5*(std::sqrt(c*c+4.) + std::abs(c));
76  G4double signc = (c < 0.) ? -1.0 : 1.0;
77  G4double cofA = -signc*G4Exp(G4Log(delta)/3.0);
78  G4double cosTheta = cofA - 1./cofA;
79  G4double tau = kinEnergy/CLHEP::electron_mass_c2;
80  G4double beta = std::sqrt(tau*(tau + 2.))/(tau + 1.);
81 
82  return (cosTheta + beta)/(1 + cosTheta*beta);
83 }
84 
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86 
89  G4double, G4int, const G4Material*)
90 {
91  G4double cosTheta = SampleCosTheta(dp->GetKineticEnergy());
92 
93  G4double sinTheta = std::sqrt((1 - cosTheta)*(1 + cosTheta));
95 
96  fLocalDirection.set(sinTheta*std::cos(phi), sinTheta*std::sin(phi),cosTheta);
98 
99  return fLocalDirection;
100 }
101 
102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
103 
105  G4double, // final_energy
106  G4int ) // Z
107 {
108  G4double cosTheta = SampleCosTheta(eTkin);
109  G4double theta = std::acos(cosTheta);
110 
111  theta = std::min(std::max(theta, 0.), CLHEP::pi);
112  return theta;
113 }
114 
115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116 
118  G4double elecKinEnergy,
119  G4double posiKinEnergy,
120  G4ThreeVector& dirElectron,
121  G4ThreeVector& dirPositron,
122  G4int, const G4Material*)
123 {
125  G4double sinp = std::sin(phi);
126  G4double cosp = std::cos(phi);
127 
128  G4double cost = SampleCosTheta(elecKinEnergy);
129  G4double sint = std::sqrt((1. - cost)*(1. + cost));
130 
131  dirElectron.set(sint*cosp, sint*sinp, cost);
132  dirElectron.rotateUz(dp->GetMomentumDirection());
133 
134  cost = SampleCosTheta(posiKinEnergy);
135  sint = std::sqrt((1. - cost)*(1. + cost));
136 
137  dirPositron.set(-sint*cosp, -sint*sinp, cost);
138  dirPositron.rotateUz(dp->GetMomentumDirection());
139 }
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142 
144 {
145  G4cout << "\n" << G4endl;
146  G4cout << "Angular Generator based on classical formula from" << G4endl;
147  G4cout << "J.D. Jackson, Classical Electrodynamics, Wiley, New York 1975"
148  << G4endl;
149 }
150 
151 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
void set(double x, double y, double z)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4DipBustGenerator(const G4String &name="")
T max(const T t1, const T t2)
brief Return the largest of the two arguments
#define G4endl
Definition: G4ios.hh:61
const G4ThreeVector & GetMomentumDirection() const
G4double SampleCosTheta(G4double kinEnergy)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
Double_t beta
virtual void SamplePairDirections(const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr) final
G4double G4Log(G4double x)
Definition: G4Log.hh:230
double G4double
Definition: G4Types.hh:76
virtual void PrintGeneratorInformation() const final
static constexpr double electron_mass_c2
#define G4UniformRand()
Definition: Randomize.hh:53
G4double PolarAngle(G4double initial_energy, G4double final_energy, G4int Z)
static constexpr double twopi
Definition: SystemOfUnits.h:55
int G4int
Definition: G4Types.hh:78
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=nullptr) final
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
static constexpr double pi
Definition: SystemOfUnits.h:54
T min(const T t1, const T t2)
brief Return the smallest of the two arguments