Geant4-11
Data Structures | Functions
G4InuclSpecialFunctions Namespace Reference

Data Structures

class  paraMaker
 

Functions

G4double bindingEnergy (G4int A, G4int Z)
 
G4double bindingEnergyAsymptotic (G4int A, G4int Z)
 
G4double csNN (G4double e)
 
G4double csPN (G4double e)
 
G4double FermiEnergy (G4int A, G4int Z, G4int ntype)
 
G4double G4cbrt (G4double x)
 
G4double G4cbrt (G4int n)
 
G4LorentzVector generateWithFixedTheta (G4double ct, G4double p, G4double mass=0.)
 
G4LorentzVector generateWithRandomAngles (G4double p, G4double mass=0.)
 
G4double getAL (G4int A)
 
G4double inuclRndm ()
 
G4double nucleiLevelDensity (G4int A)
 
std::pair< G4double, G4doublerandomCOS_SIN ()
 
G4double randomGauss (G4double sigma)
 
G4double randomInuclPowers (G4double ekin, const G4double(&coeff)[4][4])
 
G4double randomPHI ()
 

Function Documentation

◆ bindingEnergy()

G4double G4InuclSpecialFunctions::bindingEnergy ( G4int  A,
G4int  Z 
)

Definition at line 36 of file bindingEnergy.cc.

36 {
37 // NOTE: Test condition copied from G4NucleiProperties.cc; not encapsulated
38 if (A < 1 || Z < 0 || Z > A) return 0.;
39
41}
const G4int Z[17]
const G4double A[17]
static G4double GetBindingEnergy(const G4int A, const G4int Z)

References A, G4NucleiProperties::GetBindingEnergy(), and Z.

Referenced by G4VAtomDeexcitation::AlongStepDeexcitation(), G4eIonisationSpectrum::AverageEnergy(), G4QMDNucleus::CalEnergyAndAngularMomentumInCM(), G4DNAPTBAugerModel::ComputeAugerEffect(), G4BigBanger::deExcite(), G4EquilibriumEvaporator::deExcite(), G4Fissioner::deExcite(), G4NonEquilibriumEvaporator::deExcite(), G4DNAPTBAugerModel::DetermineIonisedAtom(), G4EquilibriumEvaporator::explosion(), G4CascadeDeexciteBase::explosion(), G4NucleiModel::fillBindingEnergies(), G4ProtonField::GetBarrier(), G4AntiProtonField::GetField(), G4KaonMinusField::GetField(), G4KaonPlusField::GetField(), G4KaonZeroField::GetField(), G4PionMinusField::GetField(), G4PionPlusField::GetField(), G4PionZeroField::GetField(), G4SigmaMinusField::GetField(), G4SigmaPlusField::GetField(), G4SigmaZeroField::GetField(), G4CascadeRecoilMaker::goodNucleus(), G4AtomicTransitionManager::Initialise(), G4NucleiModel::initializeCascad(), G4hImpactIonisation::PostStepDoIt(), G4eIonisationSpectrum::Probability(), G4DNAPTBIonisationModel::RandomizeEjectedElectronEnergyFromCumulated(), G4PenelopeOscillatorManager::ReadElementData(), G4DeltaAngle::SampleDirection(), G4eIonisationSpectrum::SampleEnergy(), G4DNABornIonisationModel1::SampleSecondaries(), G4DNABornIonisationModel2::SampleSecondaries(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4LivermoreIonisationModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4MicroElecInelasticModel::SampleSecondaries(), G4MicroElecInelasticModel_new::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeIonisationModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), G4DNAPTBIonisationModel::SampleSecondaries(), G4KM_OpticalEqRhs::SetFactor(), G4MicroElecInelasticModel_new::TransferedEnergy(), and G4MicroElecInelasticModel::TransferedEnergy().

◆ bindingEnergyAsymptotic()

G4double G4InuclSpecialFunctions::bindingEnergyAsymptotic ( G4int  A,
G4int  Z 
)

Definition at line 33 of file bindingEnergyAsymptotic.cc.

33 {
34 // calculates the nuclei binding energy
35 // using smooth liquid high energy formula
36 G4double X = (1.0 - 2.0*Z/A); X *= X;
37 G4double X1 = G4cbrt(A);
38 G4double X2 = X1 * X1;
39 G4double X3 = 1.0 / X1;
40 G4double X4 = 1.0 / X2;
41 G4double X5 = (1.0 - 0.62025 * X4); X5 *= X5;
42 G4double Z13 = G4cbrt(Z);
43
44 G4double DM = 17.035 * (1.0 - 1.846 * X) * A -
45 25.8357 * (1.0 - 1.712 * X) * X2 * X5 -
46 0.779 * Z * (Z - 1) * X3 *
47 (1.0 - 1.5849 * X4 + 1.2273 / A + 1.5772 * X4 * X4) +
48 0.4328 * Z13*Z13*Z13*Z13 * X3 *
49 (1.0 - 0.57811 * X3 - 0.14518 * X4 + 0.496 / A);
50
51 return DM;
52}
double G4double
Definition: G4Types.hh:83

References A, G4cbrt(), and Z.

Referenced by G4Fissioner::deExcite().

◆ csNN()

G4double G4InuclSpecialFunctions::csNN ( G4double  e)

Definition at line 82 of file G4InuclSpecialFunctions.cc.

82 {
83 G4double snn;
84
85 if (e < 40.0) {
86 snn = -1174.8 / (e * e) + 3088.5 / e + 5.3107;
87 } else {
88 snn = 93074.0 / (e * e) - 11.148 / e + 22.429;
89 }
90
91 return snn;
92}

Referenced by G4NonEquilibriumEvaporator::deExcite().

◆ csPN()

G4double G4InuclSpecialFunctions::csPN ( G4double  e)

Definition at line 94 of file G4InuclSpecialFunctions.cc.

94 {
95 G4double spn;
96
97 if (e < 40.0) {
98 spn = -5057.4 / (e * e) + 9069.2 / e + 6.9466;
99 } else {
100 spn = 239380.0 / (e * e) + 1802.0 / e + 27.147;
101 }
102
103 return spn;
104}

Referenced by G4NonEquilibriumEvaporator::deExcite().

◆ FermiEnergy()

G4double G4InuclSpecialFunctions::FermiEnergy ( G4int  A,
G4int  Z,
G4int  ntype 
)

Definition at line 108 of file G4InuclSpecialFunctions.cc.

108 {
109 G4Pow* g4pow = G4Pow::GetInstance();
110 const G4double C = 55.4 / g4pow->Z23(A);
111 G4double arg = (ntype==0) ? g4pow->Z23(A-Z) : g4pow->Z23(Z);
112
113 return C * arg;
114}
G4double C(G4double temp)
Definition: G4Pow.hh:49
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double Z23(G4int Z) const
Definition: G4Pow.hh:125

References A, C(), G4Pow::GetInstance(), Z, and G4Pow::Z23().

Referenced by G4NonEquilibriumEvaporator::deExcite().

◆ G4cbrt() [1/2]

G4double G4InuclSpecialFunctions::G4cbrt ( G4double  x)

◆ G4cbrt() [2/2]

G4double G4InuclSpecialFunctions::G4cbrt ( G4int  n)

Definition at line 120 of file G4InuclSpecialFunctions.cc.

120 {
121 return n==0 ? 0. : (n<0?-1.:1.)*G4Pow::GetInstance()->Z13(std::abs(n));
122}
G4double Z13(G4int Z) const
Definition: G4Pow.hh:123

References G4Pow::GetInstance(), CLHEP::detail::n, and G4Pow::Z13().

◆ generateWithFixedTheta()

G4LorentzVector G4InuclSpecialFunctions::generateWithFixedTheta ( G4double  ct,
G4double  p,
G4double  mass = 0. 
)

Definition at line 150 of file G4InuclSpecialFunctions.cc.

151 {
152 G4double phi = randomPHI();
153 G4double pt = p * std::sqrt(std::fabs(1.0 - ct * ct));
154
155 // Buffers to avoid memory thrashing
156 static G4ThreadLocal G4ThreeVector *pvec_G4MT_TLS_ = 0;
157 if (!pvec_G4MT_TLS_) {
158 pvec_G4MT_TLS_ = new G4ThreeVector;
159 G4AutoDelete::Register(pvec_G4MT_TLS_);
160 }
161 G4ThreeVector &pvec = *pvec_G4MT_TLS_;
162
163 static G4ThreadLocal G4LorentzVector *momr_G4MT_TLS_ = 0;
164 if (!momr_G4MT_TLS_) {
165 momr_G4MT_TLS_ = new G4LorentzVector;
166 G4AutoDelete::Register(momr_G4MT_TLS_);
167 }
168 G4LorentzVector &momr = *momr_G4MT_TLS_;
169
170 pvec.set(pt*std::cos(phi), pt*std::sin(phi), p*ct);
171 momr.setVectM(pvec, mass);
172
173 return momr;
174}
CLHEP::HepLorentzVector G4LorentzVector
CLHEP::Hep3Vector G4ThreeVector
void set(double x, double y, double z)
void setVectM(const Hep3Vector &spatial, double mass)
void Register(T *inst)
Definition: G4AutoDelete.hh:65
#define G4ThreadLocal
Definition: tls.hh:77

References G4ThreadLocal, randomPHI(), G4AutoDelete::Register(), CLHEP::Hep3Vector::set(), and CLHEP::HepLorentzVector::setVectM().

Referenced by G4CascadeFinalStateAlgorithm::FillDirManyBody(), G4CascadeFinalStateAlgorithm::FillDirThreeBody(), G4BigBanger::generateBangInSCM(), and G4NucleiModel::initializeCascad().

◆ generateWithRandomAngles()

G4LorentzVector G4InuclSpecialFunctions::generateWithRandomAngles ( G4double  p,
G4double  mass = 0. 
)

Definition at line 177 of file G4InuclSpecialFunctions.cc.

177 {
178 std::pair<G4double, G4double> COS_SIN = randomCOS_SIN();
179 G4double phi = randomPHI();
180 G4double pt = p * COS_SIN.second;
181
182 // Buffers to avoid memory thrashing
183 static G4ThreadLocal G4ThreeVector *pvec_G4MT_TLS_ = 0;
184 if (!pvec_G4MT_TLS_) {
185 pvec_G4MT_TLS_ = new G4ThreeVector;
186 G4AutoDelete::Register(pvec_G4MT_TLS_);
187 }
188 G4ThreeVector &pvec = *pvec_G4MT_TLS_;
189
190 static G4ThreadLocal G4LorentzVector *momr_G4MT_TLS_ = 0;
191 if (!momr_G4MT_TLS_) {
192 momr_G4MT_TLS_ = new G4LorentzVector;
193 G4AutoDelete::Register(momr_G4MT_TLS_);
194 }
195 G4LorentzVector &momr = *momr_G4MT_TLS_;
196
197 pvec.set(pt*std::cos(phi), pt*std::sin(phi), p*COS_SIN.first);
198 momr.setVectM(pvec, mass);
199
200 return momr;
201}
std::pair< G4double, G4double > randomCOS_SIN()

References G4ThreadLocal, randomCOS_SIN(), randomPHI(), G4AutoDelete::Register(), CLHEP::Hep3Vector::set(), and CLHEP::HepLorentzVector::setVectM().

Referenced by G4EquilibriumEvaporator::deExcite(), G4Fissioner::deExcite(), G4NonEquilibriumEvaporator::deExcite(), G4BigBanger::generateBangInSCM(), G4NucleiModel::generateNucleonMomentum(), G4ElementaryParticleCollider::generateSCMpionAbsorption(), G4ElementaryParticleCollider::generateSCMpionNAbsorption(), and G4NucleiModel::initializeCascad().

◆ getAL()

G4double G4InuclSpecialFunctions::getAL ( G4int  A)

Definition at line 78 of file G4InuclSpecialFunctions.cc.

78 {
79 return 0.76 + 2.2 / G4cbrt(A);
80}

References A, and G4cbrt().

Referenced by G4EquilibriumEvaporator::deExcite(), and G4NonEquilibriumEvaporator::deExcite().

◆ inuclRndm()

G4double G4InuclSpecialFunctions::inuclRndm ( )

◆ nucleiLevelDensity()

G4double G4InuclSpecialFunctions::nucleiLevelDensity ( G4int  A)

Definition at line 32 of file nucleiLevelDensity.cc.

32 {
33 const G4double NLD[226] = {
34 // 20 - 29
35 3.94, 3.84, 3.74, 3.64, 3.55, 4.35, 4.26, 4.09, 3.96, 4.18,
36 // 30 - 39
37 4.39, 4.61, 4.82, 4.44, 4.44, 4.43, 4.42, 5.04, 5.66, 5.8,
38 // 40 - 49
39 5.95, 5.49, 6.18, 7.11, 6.96, 7.2, 7.73, 6.41, 6.85, 6.77,
40 // 50 - 59
41 6.91, 7.3, 7.2, 6.86, 8.06, 7.8, 7.82, 8.41, 8.13, 7.19,
42 // 60 - 69
43 8.35, 8.13, 8.02, 8.93, 8.9, 9.7, 9.65, 10.55, 9.38, 9.72,
44 // 70 - 79
45 10.66, 11.98, 12.76, 12.1, 12.86, 13.0, 12.81, 12.8, 12.65, 12.0,
46 // 80 - 89
47 12.69, 14.05, 13.33, 13.28, 13.22, 13.17, 8.66, 11.03, 10.4, 13.47,
48 // 90 - 99
49 10.17, 12.22, 11.62, 12.95, 13.15, 13.57, 12.87, 16.2, 14.71, 15.69,
50 // 100 - 109
51 14.1, 18.56, 16.22, 16.7, 17.13, 17.0, 16.86, 16.2, 15.61, 16.8,
52 // 110 - 119
53 17.93, 17.5, 16.97, 17.3, 17.6, 15.78, 16.8, 17.49, 16.03, 15.08,
54 // 120 - 129
55 16.74, 17.74, 17.43, 18.1, 17.1, 19.01, 17.02, 17.0, 17.02, 18.51,
56 // 130 - 139
57 17.2, 16.8, 16.97, 16.14, 16.91, 17.69, 15.5, 14.56, 14.35, 16.5,
58 // 140 - 149
59 18.29, 17.8, 17.05, 21.31, 19.15, 19.5, 19.78, 20.3, 20.9, 21.9,
60 // 150 - 159
61 22.89, 25.68, 24.6, 24.91, 23.24, 22.9, 22.46, 21.98, 21.64, 21.8,
62 // 160 - 169
63 21.85, 21.7, 21.69, 23.7, 21.35, 23.03, 20.66, 21.81, 20.77, 22.2,
64 // 170 - 179
65 22.58, 22.55, 21.45, 21.16, 21.02, 20.87, 22.09, 22.0, 21.28, 23.05,
66 // 180 - 189
67 21.7, 21.18, 22.28, 23.0, 22.11, 23.56, 22.83, 24.88, 22.6, 23.5,
68 // 190 - 199
69 23.89, 23.9, 23.94, 21.16, 22.3, 21.7, 21.19, 20.7, 20.29, 21.32,
70 // 200 - 209
71 19.0, 17.93, 17.85, 15.7, 13.54, 11.9, 10.02, 10.48, 10.28, 11.72,
72 // 210 - 219
73 13.81, 14.7, 15.5, 16.3, 17.2, 18.0, 18.9, 19.7, 20.6, 21.4,
74 // 220 - 229
75 22.3, 23.1, 24.0, 24.8, 25.6, 26.5, 27.3, 28.2, 29.0, 29.9,
76 // 230 - 239
77 30.71, 30.53, 31.45, 29.6, 30.2, 30.65, 30.27, 29.52, 30.08, 29.8,
78 // 240 - 245
79 29.87, 30.25, 30.5, 29.8, 29.17, 28.67};
80
81 return (a>=20) ? NLD[a-20] : 0.1*a;
82}

Referenced by G4Fissioner::deExcite().

◆ randomCOS_SIN()

std::pair< G4double, G4double > G4InuclSpecialFunctions::randomCOS_SIN ( )

Definition at line 143 of file G4InuclSpecialFunctions.cc.

143 {
144 G4double CT = 1.0 - 2.0 * inuclRndm();
145
146 return std::pair<G4double, G4double>(CT, std::sqrt(1.0 - CT*CT));
147}

References inuclRndm().

Referenced by generateWithRandomAngles().

◆ randomGauss()

G4double G4InuclSpecialFunctions::randomGauss ( G4double  sigma)

Definition at line 128 of file G4InuclSpecialFunctions.cc.

128 {
129 const G4double eps = 1.0e-6;
130 G4double r1 = inuclRndm();
131 r1 = r1 > eps ? r1 : eps;
132 G4double r2 = inuclRndm();
133 r2 = r2 > eps ? r2 : eps;
134 r2 = r2 < 1.0 - eps ? r2 : 1.0 - eps;
135
136 return sigma * std::sin(twopi * r1) * std::sqrt(-2.0 * G4Log(r2));
137}
static const G4double eps
static constexpr double twopi
Definition: G4SIunits.hh:56

References eps, G4Log(), inuclRndm(), and twopi.

Referenced by G4Fissioner::deExcite().

◆ randomInuclPowers()

G4double G4InuclSpecialFunctions::randomInuclPowers ( G4double  ekin,
const G4double(&)  coeff[4][4] 
)

Definition at line 54 of file G4InuclSpecialFunctions.cc.

55 {
56 G4Pow* theG4Pow = G4Pow::GetInstance();
57
58 G4double S = G4UniformRand(); // Random fraction for expansion
59
60 G4double C, V;
61 G4double PQ=0., PR=0.;
62 for (G4int i=0; i<4; i++) {
63 V = 0.0;
64 for (G4int k=0; k<4; k++) {
65 C = coeff[i][k];
66 V += C * theG4Pow->powN(ekin, k);
67 }
68
69 PQ += V;
70 PR += V * theG4Pow->powN(S, i);
71 }
72
73 return std::sqrt(S) * (PR + (1-PQ)*(S*S*S*S));
74}
G4double S(G4double temp)
int G4int
Definition: G4Types.hh:85
G4double powN(G4double x, G4int n) const
Definition: G4Pow.cc:166

References C(), G4UniformRand, G4Pow::GetInstance(), G4Pow::powN(), and S().

Referenced by G4InuclParamAngDst::GetCosTheta(), and G4InuclParamMomDst::GetMomentum().

◆ randomPHI()

G4double G4InuclSpecialFunctions::randomPHI ( )

Definition at line 139 of file G4InuclSpecialFunctions.cc.

139 {
140 return twopi * inuclRndm();
141}

References inuclRndm(), and twopi.

Referenced by generateWithFixedTheta(), generateWithRandomAngles(), and G4NucleiModel::initializeCascad().