![]() |
|
|||
File indexing completed on 2025-02-23 09:22:16
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 /// \file medical/fanoCavity2/src/MyMollerBhabhaModel.cc 0027 /// \brief Implementation of the MyMollerBhabhaModel class 0028 // 0029 // 0030 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0031 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0032 0033 #include "MyMollerBhabhaModel.hh" 0034 0035 #include "G4PhysicalConstants.hh" 0036 #include "G4SystemOfUnits.hh" 0037 0038 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0039 0040 using namespace std; 0041 0042 MyMollerBhabhaModel::MyMollerBhabhaModel(const G4ParticleDefinition* p, const G4String& nam) 0043 : G4MollerBhabhaModel(p, nam) 0044 {} 0045 0046 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0047 0048 MyMollerBhabhaModel::~MyMollerBhabhaModel() {} 0049 0050 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0051 0052 G4double MyMollerBhabhaModel::ComputeDEDXPerVolume(const G4Material* material, 0053 const G4ParticleDefinition* p, 0054 G4double kineticEnergy, G4double cutEnergy) 0055 { 0056 if (!particle) SetParticle(p); 0057 // calculate the dE/dx due to the ionization by Seltzer-Berger formula 0058 0059 G4double electronDensity = material->GetElectronDensity(); 0060 G4double Zeff = electronDensity / material->GetTotNbOfAtomsPerVolume(); 0061 G4double th = 0.25 * sqrt(Zeff) * keV; 0062 G4double tkin = kineticEnergy; 0063 G4bool lowEnergy = false; 0064 if (kineticEnergy < th) { 0065 tkin = th; 0066 lowEnergy = true; 0067 } 0068 G4double tau = tkin / electron_mass_c2; 0069 G4double gam = tau + 1.0; 0070 G4double gamma2 = gam * gam; 0071 G4double beta2 = 1. - 1. / gamma2; 0072 // G4double bg2 = beta2*gamma2; 0073 0074 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy(); 0075 eexc /= electron_mass_c2; 0076 G4double eexc2 = eexc * eexc; 0077 0078 G4double d = min(cutEnergy, MaxSecondaryEnergy(p, tkin)) / electron_mass_c2; 0079 G4double dedx; 0080 0081 // electron 0082 if (isElectron) { 0083 dedx = log(2.0 * (tau + 2.0) / eexc2) - 1.0 - beta2 + log((tau - d) * d) + tau / (tau - d) 0084 + (0.5 * d * d + (2.0 * tau + 1.) * log(1. - d / tau)) / gamma2; 0085 0086 // positron 0087 } 0088 else { 0089 G4double d2 = d * d * 0.5; 0090 G4double d3 = d2 * d / 1.5; 0091 G4double d4 = d3 * d * 3.75; 0092 G4double y = 1.0 / (1.0 + gam); 0093 dedx = 0094 log(2.0 * (tau + 2.0) / eexc2) + log(tau * d) 0095 - beta2 * (tau + 2.0 * d - y * (3.0 * d2 + y * (d - d3 + y * (d2 - tau * d3 + d4)))) / tau; 0096 } 0097 0098 // do not apply density correction 0099 // G4double cden = material->GetIonisation()->GetCdensity(); 0100 // G4double mden = material->GetIonisation()->GetMdensity(); 0101 // G4double aden = material->GetIonisation()->GetAdensity(); 0102 // G4double x0den = material->GetIonisation()->GetX0density(); 0103 // G4double x1den = material->GetIonisation()->GetX1density(); 0104 // G4double x = log(bg2)/twoln10; 0105 0106 // if (x >= x0den) { 0107 // dedx -= twoln10*x - cden; 0108 // if (x < x1den) dedx -= aden*pow(x1den-x, mden); 0109 // } 0110 0111 // now you can compute the total ionization loss 0112 dedx *= twopi_mc2_rcl2 * electronDensity / beta2; 0113 if (dedx < 0.0) dedx = 0.0; 0114 0115 // lowenergy extrapolation 0116 0117 if (lowEnergy) { 0118 if (kineticEnergy >= lowLimit) 0119 dedx *= sqrt(tkin / kineticEnergy); 0120 else 0121 dedx *= sqrt(tkin * kineticEnergy) / lowLimit; 0122 } 0123 return dedx; 0124 } 0125 0126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
![]() ![]() |