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/MyMollerBhabhaModel.cc 27 /// \brief Implementation of the MyMollerBhabhaModel class 28 // 29 // 30 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 32 33 #include "MyMollerBhabhaModel.hh" 34 35 #include "G4PhysicalConstants.hh" 36 #include "G4SystemOfUnits.hh" 37 38 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 39 40 using namespace std; 41 42 MyMollerBhabhaModel::MyMollerBhabhaModel(const G4ParticleDefinition* p, const G4String& nam) 43 : G4MollerBhabhaModel(p, nam) 44 {} 45 46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 47 48 MyMollerBhabhaModel::~MyMollerBhabhaModel() {} 49 50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 51 52 G4double MyMollerBhabhaModel::ComputeDEDXPerVolume(const G4Material* material, 53 const G4ParticleDefinition* p, 54 G4double kineticEnergy, G4double cutEnergy) 55 { 56 if (!particle) SetParticle(p); 57 // calculate the dE/dx due to the ionization by Seltzer-Berger formula 58 59 G4double electronDensity = material->GetElectronDensity(); 60 G4double Zeff = electronDensity / material->GetTotNbOfAtomsPerVolume(); 61 G4double th = 0.25 * sqrt(Zeff) * keV; 62 G4double tkin = kineticEnergy; 63 G4bool lowEnergy = false; 64 if (kineticEnergy < th) { 65 tkin = th; 66 lowEnergy = true; 67 } 68 G4double tau = tkin / electron_mass_c2; 69 G4double gam = tau + 1.0; 70 G4double gamma2 = gam * gam; 71 G4double beta2 = 1. - 1. / gamma2; 72 // G4double bg2 = beta2*gamma2; 73 74 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy(); 75 eexc /= electron_mass_c2; 76 G4double eexc2 = eexc * eexc; 77 78 G4double d = min(cutEnergy, MaxSecondaryEnergy(p, tkin)) / electron_mass_c2; 79 G4double dedx; 80 81 // electron 82 if (isElectron) { 83 dedx = log(2.0 * (tau + 2.0) / eexc2) - 1.0 - beta2 + log((tau - d) * d) + tau / (tau - d) 84 + (0.5 * d * d + (2.0 * tau + 1.) * log(1. - d / tau)) / gamma2; 85 86 // positron 87 } 88 else { 89 G4double d2 = d * d * 0.5; 90 G4double d3 = d2 * d / 1.5; 91 G4double d4 = d3 * d * 3.75; 92 G4double y = 1.0 / (1.0 + gam); 93 dedx = 94 log(2.0 * (tau + 2.0) / eexc2) + log(tau * d) 95 - beta2 * (tau + 2.0 * d - y * (3.0 * d2 + y * (d - d3 + y * (d2 - tau * d3 + d4)))) / tau; 96 } 97 98 // do not apply density correction 99 // G4double cden = material->GetIonisation()->GetCdensity(); 100 // G4double mden = material->GetIonisation()->GetMdensity(); 101 // G4double aden = material->GetIonisation()->GetAdensity(); 102 // G4double x0den = material->GetIonisation()->GetX0density(); 103 // G4double x1den = material->GetIonisation()->GetX1density(); 104 // G4double x = log(bg2)/twoln10; 105 106 // if (x >= x0den) { 107 // dedx -= twoln10*x - cden; 108 // if (x < x1den) dedx -= aden*pow(x1den-x, mden); 109 // } 110 111 // now you can compute the total ionization loss 112 dedx *= twopi_mc2_rcl2 * electronDensity / beta2; 113 if (dedx < 0.0) dedx = 0.0; 114 115 // lowenergy extrapolation 116 117 if (lowEnergy) { 118 if (kineticEnergy >= lowLimit) 119 dedx *= sqrt(tkin / kineticEnergy); 120 else 121 dedx *= sqrt(tkin * kineticEnergy) / lowLimit; 122 } 123 return dedx; 124 } 125 126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 127