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 //--------------------------------------------------------------------------- 27 // 28 // ClassName: G4AtomicFormFactor 29 // 30 // Description: Contains the function for the evaluation of the atomic form 31 // factor. The tabulated data are available on IUCr website 32 // 33 // Class description: 34 // 35 // XXX 36 // 37 // 21-04-16, created by E.Bagli 38 39 #ifndef G4ATOMICFORMFACTOR_HH 40 #define G4ATOMICFORMFACTOR_HH 1 41 42 #include "G4Exp.hh" 43 #include "globals.hh" 44 45 #include <map> 46 #include <vector> 47 48 class G4AtomicFormFactor 49 { 50 public: 51 static G4AtomicFormFactor* GetManager() 52 { 53 if (s_G4AtomicFormFactorManager == nullptr) { 54 s_G4AtomicFormFactorManager = new G4AtomicFormFactor(); 55 } 56 return s_G4AtomicFormFactorManager; 57 } 58 59 G4double Get(G4double kScatteringVector, G4int Z, G4int charge = 0) 60 { 61 if (loadedIndex != GetIndex(Z, charge)) { 62 LoadCoefficiencts(GetIndex(Z, charge)); 63 } 64 G4double result = 0.; 65 G4double kVecOn4PiSquared = 66 (kScatteringVector / 1.e-7 / 3.1415926536) * 0.125; // (k/(4pi))/ angstrom 67 kVecOn4PiSquared = kVecOn4PiSquared * kVecOn4PiSquared; // (k/(4pi))^2 68 69 for (unsigned int i0 = 0; i0 < 4; i0++) { 70 result += theCoefficients[i0 * 2] * G4Exp(-theCoefficients[i0 * 2 + 1] * kVecOn4PiSquared); 71 } 72 result += theCoefficients[8]; 73 return result; 74 } 75 76 protected: 77 G4AtomicFormFactor(); 78 ~G4AtomicFormFactor() = default; 79 80 private: 81 void InsertCoefficients(G4int index, const std::vector<G4double>& aDoubleVec) 82 { 83 theCoefficientsMap.insert(std::pair<G4int, std::vector<G4double>>(index, aDoubleVec)); 84 } 85 86 // LoadCoefficiencts() method allows the evaluation of the atomic form 87 // factor coefficients and the storage in theCoefficients. 88 // If theCoefficients are already correct, no need to get new ones 89 // Reference: International Tables for Crystallography (2006). 90 // Vol. C, ch. 6.1, pp. 554-595 91 // doi: 10.1107/97809553602060000600 92 // Chapter 6.1. Intensity of diffracted intensities 93 // IUCr Eq. 6.1.1.15, Coefficients Table 6.1.1.4 94 void LoadCoefficiencts(G4int index) 95 { 96 loadedIndex = index; 97 for (unsigned int i0 = 0; i0 < 9; i0++) { 98 theCoefficients[i0] = theCoefficientsMap[index][i0]; 99 } 100 } 101 102 // Get() function gives back the Atomic Form Factor of the Z material 103 inline G4int GetIndex(G4int Z, G4int charge = 0) { return Z * 100 + charge; } 104 105 private: 106 inline static G4AtomicFormFactor* s_G4AtomicFormFactorManager = nullptr; 107 108 // theCoefficientsMap stores the coefficients for the form factor 109 // calculations. It can be loaded only by LoadCoefficiencts() 110 // and accessed by theCoefficients[]. 111 std::map<G4int, std::vector<G4double>> theCoefficientsMap; 112 G4double theCoefficients[9]; 113 G4int loadedIndex; 114 }; 115 #endif 116