|
||||
File indexing completed on 2025-01-18 09:57:56
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 // 0028 // ClassName: G4AtomicFormFactor 0029 // 0030 // Description: Contains the function for the evaluation of the atomic form 0031 // factor. The tabulated data are available on IUCr website 0032 // 0033 // Class description: 0034 // 0035 // XXX 0036 // 0037 // 21-04-16, created by E.Bagli 0038 0039 #ifndef G4ATOMICFORMFACTOR_HH 0040 #define G4ATOMICFORMFACTOR_HH 1 0041 0042 #include "G4Exp.hh" 0043 #include "globals.hh" 0044 0045 #include <map> 0046 #include <vector> 0047 0048 class G4AtomicFormFactor 0049 { 0050 public: 0051 static G4AtomicFormFactor* GetManager() 0052 { 0053 if (s_G4AtomicFormFactorManager == nullptr) { 0054 s_G4AtomicFormFactorManager = new G4AtomicFormFactor(); 0055 } 0056 return s_G4AtomicFormFactorManager; 0057 } 0058 0059 G4double Get(G4double kScatteringVector, G4int Z, G4int charge = 0) 0060 { 0061 if (loadedIndex != GetIndex(Z, charge)) { 0062 LoadCoefficiencts(GetIndex(Z, charge)); 0063 } 0064 G4double result = 0.; 0065 G4double kVecOn4PiSquared = 0066 (kScatteringVector / 1.e-7 / 3.1415926536) * 0.125; // (k/(4pi))/ angstrom 0067 kVecOn4PiSquared = kVecOn4PiSquared * kVecOn4PiSquared; // (k/(4pi))^2 0068 0069 for (unsigned int i0 = 0; i0 < 4; i0++) { 0070 result += theCoefficients[i0 * 2] * G4Exp(-theCoefficients[i0 * 2 + 1] * kVecOn4PiSquared); 0071 } 0072 result += theCoefficients[8]; 0073 return result; 0074 } 0075 0076 protected: 0077 G4AtomicFormFactor(); 0078 ~G4AtomicFormFactor() = default; 0079 0080 private: 0081 void InsertCoefficients(G4int index, const std::vector<G4double>& aDoubleVec) 0082 { 0083 theCoefficientsMap.insert(std::pair<G4int, std::vector<G4double>>(index, aDoubleVec)); 0084 } 0085 0086 // LoadCoefficiencts() method allows the evaluation of the atomic form 0087 // factor coefficients and the storage in theCoefficients. 0088 // If theCoefficients are already correct, no need to get new ones 0089 // Reference: International Tables for Crystallography (2006). 0090 // Vol. C, ch. 6.1, pp. 554-595 0091 // doi: 10.1107/97809553602060000600 0092 // Chapter 6.1. Intensity of diffracted intensities 0093 // IUCr Eq. 6.1.1.15, Coefficients Table 6.1.1.4 0094 void LoadCoefficiencts(G4int index) 0095 { 0096 loadedIndex = index; 0097 for (unsigned int i0 = 0; i0 < 9; i0++) { 0098 theCoefficients[i0] = theCoefficientsMap[index][i0]; 0099 } 0100 } 0101 0102 // Get() function gives back the Atomic Form Factor of the Z material 0103 inline G4int GetIndex(G4int Z, G4int charge = 0) { return Z * 100 + charge; } 0104 0105 private: 0106 inline static G4AtomicFormFactor* s_G4AtomicFormFactorManager = nullptr; 0107 0108 // theCoefficientsMap stores the coefficients for the form factor 0109 // calculations. It can be loaded only by LoadCoefficiencts() 0110 // and accessed by theCoefficients[]. 0111 std::map<G4int, std::vector<G4double>> theCoefficientsMap; 0112 G4double theCoefficients[9]; 0113 G4int loadedIndex; 0114 }; 0115 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |