Geant4 Cross Reference |
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: G4DipBustGenerator 33 // 34 // Author: Vladimir Grichine 35 // 36 // Creation date: 17 May 2011 37 // 38 // Modifications: 39 // 40 // 17.07.2018 optimisations, using G4Pow in SampleCosTheta() method M.Novak 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" 57 #include "G4Pow.hh" 58 #include <CLHEP/Units/PhysicalConstants.h> 59 60 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 61 62 G4DipBustGenerator::G4DipBustGenerator(const G4String&) 63 : G4VEmAngularDistribution("DipBustGen") 64 {} 65 66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 67 68 G4DipBustGenerator::~G4DipBustGenerator() = default; 69 70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 71 72 G4double G4DipBustGenerator::SampleCosTheta(const G4double kinEnergy) 73 { 74 const G4double c = 4. - 8.*G4UniformRand(); 75 const G4double delta = 0.5*(std::sqrt(c*c+4.) + std::abs(c)); 76 const G4double signc = (c < 0.) ? -1.0 : 1.0; 77 // const G4double cofA = -signc*G4Exp(G4Log(delta)/3.0); 78 const G4double cofA = -signc*G4Pow::GetInstance()->A13(delta); 79 const G4double cosTheta = std::min(1.,std::max(-1.,cofA - 1./cofA)); 80 const G4double tau = kinEnergy/CLHEP::electron_mass_c2; 81 const G4double beta = std::sqrt(tau*(tau + 2.))/(tau + 1.); 82 83 return (cosTheta + beta)/(1. + cosTheta*beta); 84 } 85 86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 87 88 G4ThreeVector& 89 G4DipBustGenerator::SampleDirection(const G4DynamicParticle* dp, G4double, 90 G4int, const G4Material*) 91 { 92 const G4double cosTheta = SampleCosTheta(dp->GetKineticEnergy()); 93 94 const G4double sinTheta = std::sqrt((1. - cosTheta)*(1. + cosTheta)); 95 const G4double phi = CLHEP::twopi*G4UniformRand(); 96 97 fLocalDirection.set(sinTheta*std::cos(phi), sinTheta*std::sin(phi),cosTheta); 98 fLocalDirection.rotateUz(dp->GetMomentumDirection()); 99 100 return fLocalDirection; 101 } 102 103 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 104 105 G4double G4DipBustGenerator::PolarAngle(G4double eTkin, 106 G4double, // final_energy 107 G4int ) // Z 108 { 109 const G4double cosTheta = SampleCosTheta(eTkin); 110 G4double theta = std::acos(cosTheta); 111 theta = std::min(std::max(theta, 0.), CLHEP::pi); 112 return theta; 113 } 114 115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 116 117 void G4DipBustGenerator::SamplePairDirections(const G4DynamicParticle* dp, 118 G4double elecKinEnergy, 119 G4double posiKinEnergy, 120 G4ThreeVector& dirElectron, 121 G4ThreeVector& dirPositron, 122 G4int, const G4Material*) 123 { 124 const G4double phi = CLHEP::twopi * G4UniformRand(); 125 const G4double sinp = std::sin(phi); 126 const 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 143 void G4DipBustGenerator::PrintGeneratorInformation() const 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.... 152