![]() |
|
|||
File indexing completed on 2025-02-23 09:21:03
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 /* 0027 * \file electromagnetic/TestEm7/src/G4LindhardPartition.cc 0028 * \brief Implementation of the G4LindhardPartition class 0029 * 0030 * Created by Marcus Mendenhall on 1/14/08. 0031 * 2008 Vanderbilt University, Nashville, TN, USA. 0032 * 0033 */ 0034 0035 // 0036 0037 #include "G4LindhardPartition.hh" 0038 0039 #include "G4Element.hh" 0040 #include "G4Material.hh" 0041 #include "G4PhysicalConstants.hh" 0042 #include "G4SystemOfUnits.hh" 0043 0044 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0045 /* 0046 for a first cut, we will compute NIEL from a Lindhard-Robinson partition 0047 based on the most abundant element in the material. 0048 0049 this is from IEEE Trans. Nucl Science Vol. 48 No.1 February 2001 page 162++ 0050 Insoo Jun, "Effects of Secondary Particles on the Total Dose..." 0051 and, by reference, 0052 Lindhard, Nielsen, Scharff & Thompson, 0053 "Integral Equations Governing Radiation Efects...", 0054 Mat. Fys. Medd. Dan. Vid. Selsk. vol 33 #10, pp1-42, 1963 0055 and 0056 Robinson, "The dependence of radiation effects on primary recoil energy", 0057 in Proc. Int. Conf. Radiation-Induced Voids in Metal, 0058 Albany, NY 1972 pp. 397-439 0059 def lindhard_robinson(z1, a1, z2, a2, ke): 0060 el=30.724*z1*z2*math.sqrt(z1**0.6667+z2**0.6667)*(a1+a2)/a2 0061 fl=0.0793*z1**0.6667*math.sqrt(z2)*(a1+a2)**1.5/ 0062 ((z1**0.6667+z2**0.6667)**0.75*a1**1.5*math.sqrt(a2)) 0063 eps=ke*(1.0/el) 0064 return 1.0/(1+fl*(3.4008*eps**0.16667+0.40244*eps**0.75+eps)) 0065 */ 0066 0067 G4LindhardRobinsonPartition::G4LindhardRobinsonPartition() 0068 { 0069 max_z = 120; 0070 for (size_t i = 1; i < max_z; i++) { 0071 z23[i] = std::pow((G4double)i, 2. / 3.); 0072 } 0073 } 0074 0075 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0076 0077 G4double G4LindhardRobinsonPartition::PartitionNIEL(G4int z1, G4double a1, 0078 const G4Material* material, 0079 G4double energy) const 0080 { 0081 size_t nMatElements = material->GetNumberOfElements(); 0082 0083 const G4double* atomDensities = material->GetVecNbOfAtomsPerVolume(); 0084 G4double maxdens = 0.0; 0085 size_t maxindex = 0; 0086 for (size_t k = 0; k < nMatElements; k++) { 0087 if (atomDensities[k] > maxdens) { 0088 maxdens = atomDensities[k]; 0089 maxindex = k; 0090 } 0091 } 0092 const G4Element* element = material->GetElement(maxindex); 0093 0094 G4int z2 = G4int(element->GetZ()); 0095 0096 G4double a2 = element->GetA() / (Avogadro * amu); 0097 0098 G4double zpow = z23[z1] + z23[z2]; 0099 G4double asum = a1 + a2; 0100 0101 G4double el = 30.724 * z1 * z2 * std::sqrt(zpow) * asum / a2; 0102 G4double fl = 0.0793 * z23[z1] * std::sqrt(z2 * asum * asum * asum / (a1 * a1 * a1 * a2)) 0103 / std::pow(zpow, 0.75); 0104 G4double eps = (energy / eV) * (1.0 / el); 0105 0106 return 1.0 / (1 + fl * (3.4008 * std::pow(eps, 0.16667) + 0.40244 * std::pow(eps, 0.75) + eps)); 0107 }
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
![]() ![]() |