Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AtomicFormFactor.hh
Go to the documentation of this file.
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"
44#include "globals.hh"
45
46#include <map>
47#include <vector>
48
50{
51 public:
53 {
54 if (s_G4AtomicFormFactorManager == nullptr) {
55 s_G4AtomicFormFactorManager = new G4AtomicFormFactor();
56 }
57 return s_G4AtomicFormFactorManager;
58 }
59
60 G4double Get(G4double kScatteringVector, G4int Z, G4int charge = 0)
61 {
62 if (loadedIndex != GetIndex(Z, charge)) {
63 LoadCoefficiencts(GetIndex(Z, charge));
64 }
65 G4double result = 0.;
66 // Convert k from mm^-1 to Å^-1 and divide by 4π , then square.
67 const G4double kVecOn4Pi = (kScatteringVector * 1.0e-7) * (1.0 / (4.0 * CLHEP::pi));
68 G4double kVecOn4PiSquared = kVecOn4Pi * kVecOn4Pi; // (k/(4π))^2 in Å^-2
69
70 for (unsigned int i0 = 0; i0 < 4; i0++) {
71 result += theCoefficients[i0 * 2] * G4Exp(-theCoefficients[i0 * 2 + 1] * kVecOn4PiSquared);
72 }
73 result += theCoefficients[8];
74 return result;
75 }
76
77 protected:
80
81 private:
82 void InsertCoefficients(G4int index, const std::vector<G4double>& aDoubleVec)
83 {
84 theCoefficientsMap.insert(std::pair<G4int, std::vector<G4double>>(index, aDoubleVec));
85 }
86
87 // LoadCoefficiencts() method allows the evaluation of the atomic form
88 // factor coefficients and the storage in theCoefficients.
89 // If theCoefficients are already correct, no need to get new ones
90 // Reference: International Tables for Crystallography (2006).
91 // Vol. C, ch. 6.1, pp. 554-595
92 // doi: 10.1107/97809553602060000600
93 // Chapter 6.1. Intensity of diffracted intensities
94 // IUCr Eq. 6.1.1.15, Coefficients Table 6.1.1.4
95 void LoadCoefficiencts(G4int index)
96 {
97 loadedIndex = index;
98 for (unsigned int i0 = 0; i0 < 9; i0++) {
99 theCoefficients[i0] = theCoefficientsMap[index][i0];
100 }
101 }
102
103 // Get() function gives back the Atomic Form Factor of the Z material
104 inline G4int GetIndex(G4int Z, G4int charge = 0) { return Z * 100 + charge; }
105
106 private:
107 inline static G4AtomicFormFactor* s_G4AtomicFormFactorManager = nullptr;
108
109 // theCoefficientsMap stores the coefficients for the form factor
110 // calculations. It can be loaded only by LoadCoefficiencts()
111 // and accessed by theCoefficients[].
112 std::map<G4int, std::vector<G4double>> theCoefficientsMap;
113 G4double theCoefficients[9];
114 G4int loadedIndex;
115};
116#endif
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:132
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
~G4AtomicFormFactor()=default
G4double Get(G4double kScatteringVector, G4int Z, G4int charge=0)
static G4AtomicFormFactor * GetManager()