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 /// \file medical/fanoCavity/src/MyKleinNishinaCompton.cc 27 /// \brief Implementation of the MyKleinNishinaCompton class 28 // 29 // 30 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 32 33 #include "MyKleinNishinaCompton.hh" 34 35 #include "DetectorConstruction.hh" 36 #include "MyKleinNishinaMessenger.hh" 37 38 #include "G4DataVector.hh" 39 #include "G4Electron.hh" 40 #include "G4Gamma.hh" 41 #include "G4ParticleChangeForGamma.hh" 42 #include "G4PhysicalConstants.hh" 43 #include "Randomize.hh" 44 45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 46 47 using namespace std; 48 49 MyKleinNishinaCompton::MyKleinNishinaCompton(DetectorConstruction* det, const G4ParticleDefinition*, 50 const G4String& nam) 51 : G4KleinNishinaCompton(0, nam), fDetector(det), fMessenger(0) 52 { 53 fCrossSectionFactor = 1.; 54 fMessenger = new MyKleinNishinaMessenger(this); 55 } 56 57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 58 59 MyKleinNishinaCompton::~MyKleinNishinaCompton() 60 { 61 delete fMessenger; 62 } 63 64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 65 66 G4double MyKleinNishinaCompton::CrossSectionPerVolume(const G4Material* mat, 67 const G4ParticleDefinition* part, 68 G4double GammaEnergy, G4double, G4double) 69 { 70 G4double xsection = G4VEmModel::CrossSectionPerVolume(mat, part, GammaEnergy); 71 72 return xsection * fCrossSectionFactor; 73 } 74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 75 76 void MyKleinNishinaCompton::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect, 77 const G4MaterialCutsCouple*, 78 const G4DynamicParticle* aDynamicGamma, G4double, 79 G4double) 80 { 81 // The scattered gamma energy is sampled according to Klein - Nishina formula. 82 // The random number techniques of Butcher & Messel are used 83 // (Nuc Phys 20(1960),15). 84 // Note : Effects due to binding of atomic electrons are negliged. 85 86 G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy(); 87 G4double E0_m = gamEnergy0 / electron_mass_c2; 88 89 G4ThreeVector gamDirection0 = aDynamicGamma->GetMomentumDirection(); 90 91 // 92 // sample the energy rate of the scattered gamma 93 // 94 95 G4double epsilon, epsilonsq, onecost, sint2, greject; 96 97 G4double eps0 = 1. / (1. + 2. * E0_m); 98 G4double eps0sq = eps0 * eps0; 99 G4double alpha1 = -log(eps0); 100 G4double alpha2 = 0.5 * (1. - eps0sq); 101 102 do { 103 if (alpha1 / (alpha1 + alpha2) > G4UniformRand()) { 104 epsilon = exp(-alpha1 * G4UniformRand()); // eps0**r 105 epsilonsq = epsilon * epsilon; 106 } 107 else { 108 epsilonsq = eps0sq + (1. - eps0sq) * G4UniformRand(); 109 epsilon = sqrt(epsilonsq); 110 }; 111 112 onecost = (1. - epsilon) / (epsilon * E0_m); 113 sint2 = onecost * (2. - onecost); 114 greject = 1. - epsilon * sint2 / (1. + epsilonsq); 115 116 } while (greject < G4UniformRand()); 117 118 // 119 // scattered gamma angles. ( Z - axis along the parent gamma) 120 // 121 122 G4double cosTeta = 1. - onecost; 123 G4double sinTeta = sqrt(sint2); 124 G4double Phi = twopi * G4UniformRand(); 125 G4double dirx = sinTeta * cos(Phi), diry = sinTeta * sin(Phi), dirz = cosTeta; 126 127 // 128 // update G4VParticleChange for the scattered gamma 129 // 130 // beam regeneration trick : restore incident beam 131 132 G4ThreeVector gamDirection1(dirx, diry, dirz); 133 gamDirection1.rotateUz(gamDirection0); 134 G4double gamEnergy1 = epsilon * gamEnergy0; 135 fParticleChange->SetProposedKineticEnergy(gamEnergy0); 136 fParticleChange->ProposeMomentumDirection(gamDirection0); 137 138 // 139 // kinematic of the scattered electron 140 // 141 142 G4double eKinEnergy = gamEnergy0 - gamEnergy1; 143 144 if (eKinEnergy > DBL_MIN) { 145 G4ThreeVector eDirection = gamEnergy0 * gamDirection0 - gamEnergy1 * gamDirection1; 146 eDirection = eDirection.unit(); 147 148 // create G4DynamicParticle object for the electron. 149 G4DynamicParticle* dp = new G4DynamicParticle(theElectron, eDirection, eKinEnergy); 150 fvect->push_back(dp); 151 } 152 } 153 154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 155