Geant4-11
Public Member Functions | Private Member Functions | Private Attributes
G4BetaDecayCorrections Class Reference

#include <G4BetaDecayCorrections.hh>

Public Member Functions

G4double FermiFunction (const G4double &W)
 
 G4BetaDecayCorrections (const G4int Z, const G4int A)
 
G4double ShapeFactor (const G4BetaDecayType &, const G4double &p_e, const G4double &e_nu)
 
 ~G4BetaDecayCorrections ()
 

Private Member Functions

G4double Gamma (const G4double &arg)
 
G4double ModSquared (const G4double &x, const G4double &y)
 

Private Attributes

const G4int A
 
G4double alphaZ
 
G4double gamma0
 
G4double gc [6]
 
G4double Rnuc
 
G4double V0
 
const G4int Z
 

Detailed Description

Definition at line 36 of file G4BetaDecayCorrections.hh.

Constructor & Destructor Documentation

◆ G4BetaDecayCorrections()

G4BetaDecayCorrections::G4BetaDecayCorrections ( const G4int  Z,
const G4int  A 
)

Definition at line 32 of file G4BetaDecayCorrections.cc.

33 : Z(theZ), A(theA)
34{
35 // alphaZ = fine_structure_const*std::abs(Z);
37
38 // Nuclear radius in units of hbar/m_e/c
39 Rnuc = 0.5*fine_structure_const*std::pow(A, 0.33333);
40
41 // Electron screening potential in units of electron mass
43 *std::pow(std::abs(Z), 1.33333);
44
45 gamma0 = std::sqrt(1. - alphaZ*alphaZ);
46
47 // Coefficients for gamma function with real argument
48 gc[0] = -0.1010678;
49 gc[1] = 0.4245549;
50 gc[2] = -0.6998588;
51 gc[3] = 0.9512363;
52 gc[4] = -0.5748646;
53 gc[5] = 1.0;
54}
int fine_structure_const
Definition: hepunit.py:286

References A, alphaZ, source.hepunit::fine_structure_const, gamma0, gc, Rnuc, V0, and Z.

◆ ~G4BetaDecayCorrections()

G4BetaDecayCorrections::~G4BetaDecayCorrections ( )
inline

Definition at line 41 of file G4BetaDecayCorrections.hh.

41{};

Member Function Documentation

◆ FermiFunction()

G4double G4BetaDecayCorrections::FermiFunction ( const G4double W)

Definition at line 57 of file G4BetaDecayCorrections.cc.

58{
59 // Calculate the relativistic Fermi function. Argument W is the
60 // total electron energy in units of electron mass.
61
62 G4double Wprime;
63 if (Z < 0) {
64 Wprime = W + V0;
65 } else {
66 Wprime = W - V0;
67// if (Wprime < 1.) Wprime = W;
68 if (Wprime <= 1.00001) Wprime = 1.00001;
69 }
70
71 G4double p_e = std::sqrt(Wprime*Wprime - 1.);
72 G4double eta = alphaZ*Wprime/p_e;
73 G4double epieta = std::exp(pi*eta);
74 G4double realGamma = Gamma(2.*gamma0+1);
75 G4double mod2Gamma = ModSquared(gamma0, eta);
76
77 // Fermi function
78 G4double factor1 = 2*(1+gamma0)*mod2Gamma/realGamma/realGamma;
79 G4double factor2 = epieta*std::pow(2*p_e*Rnuc, 2*(gamma0-1) );
80
81 // Electron screening factor
82 G4double factor3 = (Wprime/W)*std::sqrt( (Wprime*Wprime - 1.)/(W*W - 1.) );
83
84 return factor1*factor2*factor3;
85}
static constexpr double pi
Definition: G4SIunits.hh:55
double G4double
Definition: G4Types.hh:83
G4double Gamma(const G4double &arg)
G4double ModSquared(const G4double &x, const G4double &y)

References alphaZ, Gamma(), gamma0, ModSquared(), pi, Rnuc, V0, and Z.

Referenced by G4BetaMinusDecay::SetUpBetaSpectrumSampler(), and G4BetaPlusDecay::SetUpBetaSpectrumSampler().

◆ Gamma()

G4double G4BetaDecayCorrections::Gamma ( const G4double arg)
private

Definition at line 106 of file G4BetaDecayCorrections.cc.

107{
108 // Use recursion relation to get argument < 1
109 G4double fac = 1.0;
110 G4double x = arg - 1.;
111
112 G4int loop = 0;
114 ed << " While count exceeded " << G4endl;
115 while (x > 1.0) { /* Loop checking, 01.09.2015, D.Wright */
116 fac *= x;
117 x -= 1.0;
118 loop++;
119 if (loop > 1000) {
120 G4Exception("G4BetaDecayCorrections::Gamma()", "HAD_RDM_100", JustWarning, ed);
121 break;
122 }
123 }
124
125 // Calculation of Gamma function with real argument
126 // 0 < arg < 1 using polynomial from Abramowitz and Stegun
127 G4double sum = gc[0];
128 for (G4int i = 1; i < 6; i++) sum = sum*x + gc[i];
129
130 return sum*fac;
131}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
static const G4double fac
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57

References fac, G4endl, G4Exception(), gc, and JustWarning.

Referenced by FermiFunction(), and ShapeFactor().

◆ ModSquared()

G4double G4BetaDecayCorrections::ModSquared ( const G4double x,
const G4double y 
)
private

Definition at line 89 of file G4BetaDecayCorrections.cc.

90{
91 // Calculate the squared modulus of the Gamma function
92 // with complex argument (re, im) using approximation B
93 // of Wilkinson, Nucl. Instr. & Meth. 82, 122 (1970).
94 // Here, choose N = 1 in Wilkinson's notation for approximation B
95
96 G4double factor1 = std::pow( (1+re)*(1+re) + im*im, re+0.5);
97 G4double factor2 = std::exp(2*im * std::atan(im/(1+re)));
98 G4double factor3 = std::exp(2*(1+re));
99 G4double factor4 = 2.*pi;
100 G4double factor5 = std::exp( (1+re)/( (1+re)*(1+re) + im*im)/6 );
101 G4double factor6 = re*re + im*im;
102 return factor1*factor4*factor5/factor2/factor3/factor6;
103}

References pi.

Referenced by FermiFunction(), and ShapeFactor().

◆ ShapeFactor()

G4double G4BetaDecayCorrections::ShapeFactor ( const G4BetaDecayType bdt,
const G4double p_e,
const G4double e_nu 
)

Definition at line 135 of file G4BetaDecayCorrections.cc.

137{
138 G4double twoPR = 2.*p_e*Rnuc;
139 G4double factor(1.);
140
141 switch (bdt)
142 {
143 case (allowed) :
144 break;
145
146 case (firstForbidden) :
147 {
148 // Parameters for 1st forbidden shape determined from 210Bi data
149 // Not valid for other 1st forbidden nuclei
150 G4double c1 = 0.578;
151 G4double c2 = 28.466;
152 G4double c3 = -0.658;
153
154 G4double w = std::sqrt(1. + p_e*p_e);
155 factor = 1. + c1*w + c2/w + c3*w*w;
156 }
157 break;
158
159 case (uniqueFirstForbidden) :
160 {
161 G4double eta = alphaZ*std::sqrt(1. + p_e*p_e)/p_e;
162 G4double gamma1 = std::sqrt(4. - alphaZ*alphaZ);
163 G4double gamterm1 = Gamma(2.*gamma0+1.)/Gamma(2.*gamma1+1.);
164 G4double term1 = e_nu*e_nu*(1. + gamma0)/6.;
165 G4double term2 = 12.*(2. + gamma1)*p_e*p_e
166 *std::pow(twoPR, 2.*(gamma1-gamma0-1) )
167 *gamterm1*gamterm1
168 *ModSquared(gamma1, eta)/ModSquared(gamma0, eta);
169 factor = term1 + term2;
170 }
171 break;
172
173 case (secondForbidden) :
174 break;
175
176 case (uniqueSecondForbidden) :
177 {
178 G4double eta = alphaZ*std::sqrt(1. + p_e*p_e)/p_e;
179 G4double gamma1 = std::sqrt(4. - alphaZ*alphaZ);
180 G4double gamma2 = std::sqrt(9. - alphaZ*alphaZ);
181 G4double gamterm0 = Gamma(2.*gamma0+1.);
182 G4double gamterm1 = gamterm0/Gamma(2.*gamma1+1.);
183 G4double gamterm2 = gamterm0/Gamma(2.*gamma2+1.);
184 G4double term1 = e_nu*e_nu*e_nu*e_nu*(1. + gamma0)/60.;
185
186 G4double term2 = 4.*(2. + gamma1)*e_nu*e_nu*p_e*p_e
187 *std::pow(twoPR, 2.*(gamma1-gamma0-1.) )
188 *gamterm1*gamterm1
189 *ModSquared(gamma1, eta)/ModSquared(gamma0, eta);
190
191 G4double term3 = 180.*(3.+gamma2)*p_e*p_e*p_e*p_e
192 *std::pow(twoPR, 2.*(gamma2-gamma0-2) )
193 *gamterm2*gamterm2
194 *ModSquared(gamma2, eta)/ModSquared(gamma0, eta);
195
196 factor = term1 + term2 + term3;
197 }
198 break;
199
200 case (thirdForbidden) :
201 break;
202
203 case (uniqueThirdForbidden) :
204 {
205 G4double eta = alphaZ*std::sqrt(1. + p_e*p_e)/p_e;
206 G4double gamma1 = std::sqrt(4. - alphaZ*alphaZ);
207 G4double gamma2 = std::sqrt(9. - alphaZ*alphaZ);
208 G4double gamma3 = std::sqrt(16. - alphaZ*alphaZ);
209 G4double gamterm0 = Gamma(2.*gamma0+1.);
210 G4double gamterm1 = gamterm0/Gamma(2.*gamma1+1.);
211 G4double gamterm2 = gamterm0/Gamma(2.*gamma2+1.);
212 G4double gamterm3 = gamterm0/Gamma(2.*gamma3+1.);
213
214 G4double term1 = e_nu*e_nu*e_nu*e_nu*e_nu*e_nu*(1. + gamma0)/1260.;
215
216 G4double term2 = 2.*(2. + gamma1)*e_nu*e_nu*e_nu*e_nu*p_e*p_e
217 *std::pow(twoPR, 2.*(gamma1-gamma0-1.) )
218 *gamterm1*gamterm1
219 *ModSquared(gamma1, eta)/ModSquared(gamma0, eta)/5.;
220
221 G4double term3 = 60.*(3.+gamma2)*p_e*p_e*p_e*p_e*e_nu*e_nu
222 *std::pow(twoPR, 2.*(gamma2-gamma0-2.) )
223 *gamterm2*gamterm2
224 *ModSquared(gamma2, eta)/ModSquared(gamma0, eta);
225
226 G4double term4 = 2240.*p_e*p_e*p_e*p_e*p_e*p_e*(4. + gamma3)
227 *std::pow(twoPR, 2.*(gamma3-gamma0-3.) )
228 *gamterm3*gamterm3
229 *ModSquared(gamma3, eta)/ModSquared(gamma0, eta);
230
231 factor = term1 + term2 + term3 + term4;
232 }
233 break;
234
235 default:
236 G4Exception("G4BetaDecayCorrections::ShapeFactor()","HAD_RDM_010",
238 "Transition not yet implemented - using allowed shape");
239 break;
240 }
241 return factor;
242}
@ uniqueFirstForbidden
@ uniqueThirdForbidden
@ allowed
@ uniqueSecondForbidden
@ thirdForbidden
@ secondForbidden
@ firstForbidden

References allowed, alphaZ, firstForbidden, G4Exception(), Gamma(), gamma0, JustWarning, ModSquared(), Rnuc, secondForbidden, thirdForbidden, uniqueFirstForbidden, uniqueSecondForbidden, and uniqueThirdForbidden.

Referenced by G4BetaMinusDecay::SetUpBetaSpectrumSampler(), and G4BetaPlusDecay::SetUpBetaSpectrumSampler().

Field Documentation

◆ A

const G4int G4BetaDecayCorrections::A
private

◆ alphaZ

G4double G4BetaDecayCorrections::alphaZ
private

Definition at line 54 of file G4BetaDecayCorrections.hh.

Referenced by FermiFunction(), G4BetaDecayCorrections(), and ShapeFactor().

◆ gamma0

G4double G4BetaDecayCorrections::gamma0
private

Definition at line 57 of file G4BetaDecayCorrections.hh.

Referenced by FermiFunction(), G4BetaDecayCorrections(), and ShapeFactor().

◆ gc

G4double G4BetaDecayCorrections::gc[6]
private

Definition at line 59 of file G4BetaDecayCorrections.hh.

Referenced by G4BetaDecayCorrections(), and Gamma().

◆ Rnuc

G4double G4BetaDecayCorrections::Rnuc
private

Definition at line 55 of file G4BetaDecayCorrections.hh.

Referenced by FermiFunction(), G4BetaDecayCorrections(), and ShapeFactor().

◆ V0

G4double G4BetaDecayCorrections::V0
private

Definition at line 56 of file G4BetaDecayCorrections.hh.

Referenced by FermiFunction(), and G4BetaDecayCorrections().

◆ Z

const G4int G4BetaDecayCorrections::Z
private

The documentation for this class was generated from the following files: