Geant4-11
Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
G4Pow Class Reference

#include <G4Pow.hh>

Public Member Functions

G4double A13 (G4double A) const
 
G4double A23 (G4double A) const
 
G4double expA (G4double A) const
 
G4double factorial (G4int Z) const
 
G4double log10A (G4double A) const
 
G4double log10Z (G4int Z) const
 
G4double logA (G4double A) const
 
G4double logfactorial (G4int Z) const
 
G4double logX (G4double x) const
 
G4double logZ (G4int Z) const
 
G4double powA (G4double A, G4double y) const
 
G4double powN (G4double x, G4int n) const
 
G4double powZ (G4int Z, G4double y) const
 
G4double Z13 (G4int Z) const
 
G4double Z23 (G4int Z) const
 
 ~G4Pow ()
 

Static Public Member Functions

static G4PowGetInstance ()
 

Private Member Functions

G4double A13High (const G4double, const G4bool) const
 
G4double A13Low (const G4double, const G4bool) const
 
 G4Pow ()
 
G4double logBase (G4double x) const
 

Private Attributes

G4DataVector ener
 
G4DataVector fact
 
G4DataVector fexp
 
G4DataVector logen
 
G4DataVector logfact
 
G4DataVector lowa13
 
G4DataVector lz
 
G4DataVector lz2
 
const G4int max2 = 5
 
G4double maxA
 
G4double maxA2
 
G4double maxAexp
 
G4double maxLowA
 
const G4double onethird = 1.0 / 3.0
 
G4DataVector pz13
 

Static Private Attributes

static G4PowfpInstance = nullptr
 

Detailed Description

Definition at line 48 of file G4Pow.hh.

Constructor & Destructor Documentation

◆ ~G4Pow()

G4Pow::~G4Pow ( )

Definition at line 116 of file G4Pow.cc.

116{}

◆ G4Pow()

G4Pow::G4Pow ( )
private

Definition at line 53 of file G4Pow.cc.

54{
55#ifdef G4MULTITHREADED
57 {
58 G4Exception("G4Pow::G4Pow()", "InvalidSetup", FatalException,
59 "Attempt to instantiate G4Pow in worker thread!");
60 }
61#endif
62 const G4int maxZ = 512;
63 const G4int maxZfact = 170;
64 const G4int numLowA = 17;
65
66 maxA = -0.6 + maxZ;
67 maxLowA = 4.0;
68 maxA2 = 1.25 + max2 * 0.2;
69 maxAexp = -0.76 + maxZfact * 0.5;
70
71 ener.resize(max2 + 1, 1.0);
72 logen.resize(max2 + 1, 0.0);
73 lz2.resize(max2 + 1, 0.0);
74 pz13.resize(maxZ, 0.0);
75 lowa13.resize(numLowA, 0.0);
76 lz.resize(maxZ, 0.0);
77 fexp.resize(maxZfact, 0.0);
78 fact.resize(maxZfact, 0.0);
79 logfact.resize(maxZ, 0.0);
80
81 G4double f = 1.0;
82 G4double logf = 0.0;
83 fact[0] = 1.0;
84 fexp[0] = 1.0;
85
86 for(G4int i = 1; i <= max2; ++i)
87 {
88 ener[i] = powN(500., i);
89 logen[i] = G4Log(ener[i]);
90 lz2[i] = G4Log(1.0 + i * 0.2);
91 }
92
93 for(G4int i = 1; i < maxZ; ++i)
94 {
95 G4double x = G4double(i);
96 pz13[i] = std::pow(x, onethird);
97 lz[i] = G4Log(x);
98 if(i < maxZfact)
99 {
100 f *= x;
101 fact[i] = f;
102 fexp[i] = G4Exp(0.5 * i);
103 }
104 logf += lz[i];
105 logfact[i] = logf;
106 }
107
108 for(G4int i = 4; i < numLowA; ++i)
109 {
110 lowa13[i] = std::pow(0.25 * i, onethird);
111 }
112}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
static const G4int maxZ
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4DataVector logfact
Definition: G4Pow.hh:116
G4double maxLowA
Definition: G4Pow.hh:104
const G4int max2
Definition: G4Pow.hh:101
G4DataVector ener
Definition: G4Pow.hh:108
G4double maxAexp
Definition: G4Pow.hh:106
G4DataVector logen
Definition: G4Pow.hh:109
G4DataVector fexp
Definition: G4Pow.hh:114
G4DataVector lz2
Definition: G4Pow.hh:113
G4double powN(G4double x, G4int n) const
Definition: G4Pow.cc:166
G4DataVector pz13
Definition: G4Pow.hh:110
G4DataVector fact
Definition: G4Pow.hh:115
G4double maxA
Definition: G4Pow.hh:103
G4DataVector lowa13
Definition: G4Pow.hh:111
const G4double onethird
Definition: G4Pow.hh:100
G4DataVector lz
Definition: G4Pow.hh:112
G4double maxA2
Definition: G4Pow.hh:105
G4bool IsWorkerThread()
Definition: G4Threading.cc:123

References ener, fact, FatalException, fexp, G4Exception(), G4Exp(), G4Log(), G4Threading::IsWorkerThread(), logen, logfact, lowa13, lz, lz2, max2, maxA, maxA2, maxAexp, maxLowA, maxZ, onethird, powN(), and pz13.

Member Function Documentation

◆ A13()

G4double G4Pow::A13 ( G4double  A) const

Definition at line 120 of file G4Pow.cc.

121{
122 G4double res = 0.;
123 if(A > 0.)
124 {
125 const bool invert = (A < 1.);
126 const G4double a = invert ? 1. / A : A;
127 res = (a < maxLowA) ? A13Low(a, invert) : A13High(a, invert);
128 }
129 return res;
130}
const G4double A[17]
G4double A13Low(const G4double, const G4bool) const
Definition: G4Pow.cc:153
G4double A13High(const G4double, const G4bool) const
Definition: G4Pow.cc:134

References A, A13High(), A13Low(), and maxLowA.

Referenced by A23(), G4LFission::Atomas(), G4AntiNuclElastic::CalculateAm(), G4DiffuseElastic::CalculateAm(), G4DiffuseElasticV2::CalculateAm(), G4NuclNuclDiffuseElastic::CalculateAm(), G4DiffuseElastic::CalculateNuclearRad(), G4DiffuseElasticV2::CalculateNuclearRad(), G4NuclNuclDiffuseElastic::CalculateNuclearRad(), G4StatMFMacroCanonical::CalculateTemperature(), G4FermiMomentum::cbrt(), G4QMDMeanField::DoClusterJudgment(), G4KM_NucleonEqRhs::EvaluateRhsGivenB(), G4IonFluctuations::Factor(), G4StatMFMacroTemperature::FragsExcitEnergy(), G4QMDGroundStateNucleus::G4QMDGroundStateNucleus(), G4QMDParameters::G4QMDParameters(), G4fissionEvent::G4SmpNugDist(), G4WilsonAbrasionModel::GetAbradedNucleons(), G4SigmaMinusField::GetBarrier(), G4SigmaPlusField::GetBarrier(), G4StatMFMicroPartition::GetCoulombEnergy(), G4ComponentAntiNuclNuclearXS::GetInelasticElementCrossSection(), G4IonsShenCrossSection::GetIsoCrossSection(), G4QMDNucleus::GetNuclearMass(), G4StatMFMicroPartition::GetPartitionEnergy(), G4ComponentAntiNuclNuclearXS::GetTotalElementCrossSection(), G4WilsonRadius::GetWilsonRMSRadius(), G4NistManager::GetZ13(), G4mplIonisationModel::Initialise(), G4mplIonisationWithDeltaModel::Initialise(), MCGIDI_KalbachMann_S_a_or_b(), G4StatMFChannel::PlaceFragments(), G4DipBustGenerator::SampleCosTheta(), G4BetheHeitlerModel::SampleSecondaries(), G4PairProductionRelModel::SampleSecondaries(), G4WentzelOKandVIxSection::SampleSingleScattering(), G4QMDGroundStateNucleus::samplingMomentum(), G4ParticleHPKallbachMannSyst::SeparationEnergy(), and G4StatMFChannel::SolveEqOfMotion().

◆ A13High()

G4double G4Pow::A13High ( const G4double  a,
const  G4bool 
) const
private

Definition at line 134 of file G4Pow.cc.

135{
136 G4double res;
137 if(a < maxA)
138 {
139 const G4int i = static_cast<G4int>(a + 0.5);
140 const G4double x = (a / i - 1.) * onethird;
141 res = pz13[i] * (1. + x - x * x * (1. - 1.666667 * x));
142 }
143 else
144 {
145 res = G4Exp(G4Log(a) * onethird);
146 }
147 res = invert ? 1. / res : res;
148 return res;
149}

References G4Exp(), G4Log(), maxA, onethird, and pz13.

Referenced by A13().

◆ A13Low()

G4double G4Pow::A13Low ( const G4double  a,
const G4bool  invert 
) const
private

Definition at line 153 of file G4Pow.cc.

154{
155 G4double res;
156 const G4int i = static_cast<G4int>(4. * (a + 0.125));
157 const G4double y = 0.25 * i;
158 const G4double x = (a / y - 1.) * onethird;
159 res = lowa13[i] * (1. + x - x * x * (1. - 1.666667 * x));
160 res = invert ? 1. / res : res;
161 return res;
162}

References lowa13, and onethird.

Referenced by A13().

◆ A23()

G4double G4Pow::A23 ( G4double  A) const
inline

◆ expA()

G4double G4Pow::expA ( G4double  A) const
inline

Definition at line 203 of file G4Pow.hh.

204{
205 G4double res;
206 G4double a = (0.0 <= A) ? A : -A;
207
208 if(a <= maxAexp)
209 {
210 G4int i = G4int(2 * a + 0.5);
211 G4double x = a - i * 0.5;
212 res = fexp[i] * (1.0 + x * (1.0 + 0.5 * (1.0 + onethird * x) * x));
213 }
214 else
215 {
216 res = G4Exp(a);
217 }
218 if(0.0 > A)
219 {
220 res = 1.0 / res;
221 }
222 return res;
223}

References A, fexp, G4Exp(), maxAexp, and onethird.

Referenced by G4IonFluctuations::Factor(), powA(), and powZ().

◆ factorial()

G4double G4Pow::factorial ( G4int  Z) const
inline

Definition at line 235 of file G4Pow.hh.

235{ return fact[Z]; }
const G4int Z[17]

References fact, and Z.

Referenced by G4LegendrePolynomial::EvalAssocLegendrePoly().

◆ GetInstance()

G4Pow * G4Pow::GetInstance ( )
static

Definition at line 41 of file G4Pow.cc.

42{
43 if(fpInstance == nullptr)
44 {
45 static G4Pow geant4pow;
46 fpInstance = &geant4pow;
47 }
48 return fpInstance;
49}
Definition: G4Pow.hh:49
static G4Pow * fpInstance
Definition: G4Pow.hh:98

References fpInstance.

Referenced by G4ParticleHPKallbachMannSyst::A(), G4FTFAnnihilation::Annihilate(), G4WilsonAbrasionModel::ApplyYourself(), G4LENDElastic::ApplyYourself(), G4QMDReaction::ApplyYourself(), G4ParticleHPElasticFS::ApplyYourself(), G4LFission::Atomas(), G4FissionBarrier::BarashenkovFissionBarrier(), G4CascadeFinalStateAlgorithm::BetaKopylov(), G4HadPhaseSpaceKopylov::BetaKopylov(), G4WilsonAblationModel::BreakItUp(), G4UPiNuclearCrossSection::BuildPhysicsTable(), G4ParticleHPNBodyPhaseSpace::C(), G4QMDMeanField::Cal2BodyQuantities(), G4StatMFMacroMultiplicity::CalcChemicalPotentialMu(), G4StatMFMacroChemicalPotential::CalcChemicalPotentialNu(), G4QuasiElRatios::CalcElTot(), G4StatMFMacroBiNucleon::CalcEnergy(), G4StatMFMacroMultiNucleon::CalcEnergy(), G4StatMFMacroTetraNucleon::CalcEnergy(), G4StatMFMacroTriNucleon::CalcEnergy(), G4StatMFMacroMultiNucleon::CalcEntropy(), G4StatMFMicroCanonical::CalcEntropyOfCompoundNucleus(), G4StatMFFragment::CalcExcitationEnergy(), G4StatMFMicroCanonical::CalcFreeInternalEnergy(), G4StatMFMacroBiNucleon::CalcMeanMultiplicity(), G4StatMFMacroMultiNucleon::CalcMeanMultiplicity(), G4StatMFMacroTetraNucleon::CalcMeanMultiplicity(), G4StatMFMacroTriNucleon::CalcMeanMultiplicity(), G4StatMFMicroPartition::CalcPartitionProbability(), G4QuasiElRatios::CalcQF2IN_Ratio(), G4AntiNuclElastic::CalculateAm(), G4DiffuseElastic::CalculateAm(), G4DiffuseElasticV2::CalculateAm(), G4NuclNuclDiffuseElastic::CalculateAm(), G4MuonVDNuclearModel::CalculateEMVertex(), G4DiffuseElastic::CalculateNuclearRad(), G4DiffuseElasticV2::CalculateNuclearRad(), G4NuclNuclDiffuseElastic::CalculateNuclearRad(), G4StatMFMacroCanonical::CalculateTemperature(), G4RKFieldIntegrator::CalculateTotalEnergy(), G4StatMFMacroMultiNucleon::CalcZARatio(), G4QMDNucleus::CalEnergyAndAngularMomentumInCM(), G4QMDCollision::CalFinalStateOfTheBinaryCollisionJQMD(), G4QMDMeanField::CalGraduate(), G4FermiMomentum::cbrt(), G4Reggeons::Chi_pomeron(), G4Reggeons::Chi_reggeon(), G4DiffractiveExcitation::ChooseP(), G4QGSDiffractiveExcitation::ChooseP(), G4SingleDiffractiveExcitation::ChooseX(), G4StatMFMacroCanonical::ChooseZ(), G4Clebsch::ClebschGordanCoeff(), G4UrbanAdjointMscModel::ComputeCrossSectionPerAtom(), G4UrbanMscModel::ComputeCrossSectionPerAtom(), G4GoudsmitSaundersonMscModel::ComputeGeomPathLength(), G4IonisParamMat::ComputeIonParameters(), G4KalbachCrossSection::ComputePowerParameter(), G4GoudsmitSaundersonMscModel::ComputeTrueStepLength(), G4RToEConvForElectron::ComputeValue(), G4RToEConvForGamma::ComputeValue(), G4RToEConvForPositron::ComputeValue(), G4Fancy3DNucleus::CoulombBarrier(), G4StatMFMicroPartition::CoulombFreeEnergy(), G4DiffractiveExcitation::CreateStrings(), G4XAqmElastic::CrossSection(), G4XAqmTotal::CrossSection(), G4XPDGElastic::CrossSection(), G4XPDGTotal::CrossSection(), G4ChipsAntiBaryonInelasticXS::CrossSectionFormula(), G4QMDMeanField::DoClusterJudgment(), G4LENDCrossSection::DumpPhysicsTable(), G4ParticleHPCaptureData::DumpPhysicsTable(), G4ParticleHPElasticData::DumpPhysicsTable(), G4ParticleHPFissionData::DumpPhysicsTable(), G4ParticleHPInelasticData::DumpPhysicsTable(), G4LegendrePolynomial::EvalAssocLegendrePoly(), G4KM_NucleonEqRhs::EvaluateRhsGivenB(), G4FPYSamplingOps::EvaluateWattConstants(), G4QuarkExchange::ExciteParticipants(), G4NuclearAbrasionGeometry::F(), G4VCrossSectionSource::FcrossX(), G4InuclSpecialFunctions::FermiEnergy(), G4JTPolynomialSolver::FindRoots(), G4CompetitiveFission::FissionKineticEnergy(), G4StatMFMacroTemperature::FragsExcitEnergy(), G4AtimaEnergyLossModel::G4AtimaEnergyLossModel(), G4AtimaFluctuations::G4AtimaFluctuations(), G4BGGPionElasticXS::G4BGGPionElasticXS(), G4BGGPionInelasticXS::G4BGGPionInelasticXS(), G4InuclSpecialFunctions::G4cbrt(), G4EmCorrections::G4EmCorrections(), G4ExcitationHandler::G4ExcitationHandler(), G4FermiPhaseSpaceDecay::G4FermiPhaseSpaceDecay(), G4GEMChannel::G4GEMChannel(), G4GEMProbability::G4GEMProbability(), G4Generator2BS::G4Generator2BS(), G4FissionProductYieldDist::G4GetFission(), G4HadronNucleonXsc::G4HadronNucleonXsc(), G4ICRU49NuclearStoppingModel::G4ICRU49NuclearStoppingModel(), G4IonCoulombCrossSection::G4IonCoulombCrossSection(), G4ionEffectiveCharge::G4ionEffectiveCharge(), G4IonFluctuations::G4IonFluctuations(), G4IonisParamElm::G4IonisParamElm(), G4IonisParamMat::G4IonisParamMat(), G4KM_NucleonEqRhs::G4KM_NucleonEqRhs(), G4LevelManager::G4LevelManager(), G4LevelReader::G4LevelReader(), G4LindhardSorensenData::G4LindhardSorensenData(), G4NistManager::G4NistManager(), G4NuclearFermiDensity::G4NuclearFermiDensity(), G4NuclearLevelData::G4NuclearLevelData(), G4NuclearShellModelDensity::G4NuclearShellModelDensity(), G4PreCompoundEmission::G4PreCompoundEmission(), G4QMDGroundStateNucleus::G4QMDGroundStateNucleus(), G4QMDParameters::G4QMDParameters(), G4FPYSamplingOps::G4SampleWatt(), G4ScreeningMottCrossSection::G4ScreeningMottCrossSection(), G4fissionEvent::G4SmpGEng(), G4fissionEvent::G4SmpNEngCf252(), G4fissionEvent::G4SmpNuDistDataPu239(), G4fissionEvent::G4SmpNuDistDataPu239_241(), G4fissionEvent::G4SmpNuDistDataU232_234_236_238(), G4fissionEvent::G4SmpNuDistDataU233_235(), G4fissionEvent::G4SmpNuDistDataU235(), G4fissionEvent::G4SmpNuDistDataU238(), G4fissionEvent::G4SmpNugDist(), G4VEmissionProbability::G4VEmissionProbability(), G4VPreCompoundFragment::G4VPreCompoundFragment(), G4WentzelOKandVIxSection::G4WentzelOKandVIxSection(), G4ParticleHPMadlandNixSpectrum::Gamma25(), G4QMDParticipant::Get4Momentum(), G4WilsonAbrasionModel::GetAbradedNucleons(), G4AntiProtonField::GetBarrier(), G4KaonMinusField::GetBarrier(), G4KaonPlusField::GetBarrier(), G4PionMinusField::GetBarrier(), G4PionPlusField::GetBarrier(), G4ProtonField::GetBarrier(), G4SigmaMinusField::GetBarrier(), G4SigmaPlusField::GetBarrier(), G4ParticleHPInterpolator::GetBinIntegral(), G4EMDissociationSpectrum::GetClosestApproach(), G4AntiNuclElastic::GetcosTeta1(), G4StatMFMicroPartition::GetCoulombEnergy(), G4EMDissociationCrossSection::GetCrossSectionForProjectile(), G4ChargeExchangeProcess::GetElementCrossSection(), G4StatMFFragment::GetEnergy(), G4ChipsAntiBaryonElasticXS::GetExchangeT(), G4ChipsHyperonElasticXS::GetExchangeT(), G4ChipsPionMinusElasticXS::GetExchangeT(), G4ChipsPionPlusElasticXS::GetExchangeT(), G4NuclearAbrasionGeometry::GetExcitationEnergyOfProjectile(), G4NuclearAbrasionGeometry::GetExcitationEnergyOfTarget(), G4EMDissociationSpectrum::GetGeneralE2Spectrum(), G4NuclNuclDiffuseElastic::GetHadronNucleonXscNS(), G4IonsShenCrossSection::GetIsoCrossSection(), G4LundStringFragmentation::GetLightConeZ(), G4QGSMFragmentation::GetLightConeZ(), G4InuclParamMomDst::GetMomentum(), G4QMDNucleus::GetNuclearMass(), G4StatMFMicroPartition::GetPartitionEnergy(), G4QMDMeanField::GetPotential(), G4ChipsAntiBaryonElasticXS::GetPTables(), G4ChipsHyperonElasticXS::GetPTables(), G4ChipsPionMinusElasticXS::GetPTables(), G4ChipsPionPlusElasticXS::GetPTables(), G4ChipsAntiBaryonElasticXS::GetTabValues(), G4ChipsHyperonElasticXS::GetTabValues(), G4ChipsPionMinusElasticXS::GetTabValues(), G4ChipsPionPlusElasticXS::GetTabValues(), G4QMDMeanField::GetTotalPotential(), G4ParticleHPInterpolator::GetWeightedBinIntegral(), G4WilsonRadius::GetWilsonRMSRadius(), G4ParticleHPJENDLHEData::getXSfromThisIsotope(), G4ParticleHPMadlandNixSpectrum::GIntegral(), G4Bessel::I0(), G4Bessel::I1(), G4FTFParameters::InitForInteraction(), G4mplIonisationModel::Initialise(), G4mplIonisationWithDeltaModel::Initialise(), G4PhotonEvaporation::InitialiseGRData(), G4StatMFMacroCanonical::Initialize(), G4StatMFMicroCanonical::Initialize(), G4ParticleHPVector::Integrate(), G4PolarizationTransition::LnFactorial(), G4ComponentBarNucleonNucleusXsc::LoadData(), G4ParticleHPMadlandNixSpectrum::Madland(), MCGIDI_energy_NBodyPhaseSpacePDF_callback(), MCGIDI_KalbachMann_S_a_or_b(), MCGIDI_sampling_interpolationValues(), MCGIDI_sampling_sampleX_from_pdfsOfXGivenW(), G4NuclearAbrasionGeometry::P(), G4QMDGroundStateNucleus::packNucleons(), G4StatMFChannel::PlaceFragments(), G4ParticleHPNBodyPhaseSpace::Prob(), ptwXY_f_integrate(), ptwXY_LogLogToLinLin(), ptwXY_pow_callback(), G4LEHadronProtonElastic::RandCosThetaDipPen(), G4InuclSpecialFunctions::randomInuclPowers(), G4DipBustGenerator::SampleCosTheta(), G4AntiNuclElastic::SampleInvariantT(), G4HadronElastic::SampleInvariantT(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4ChargeExchange::SampleT(), G4QGSMSplitableHadron::SampleX(), G4QGSParticipants::SampleX(), G4QMDGroundStateNucleus::samplingMomentum(), G4ParticleHPKallbachMannSyst::SeparationEnergy(), G4FissionProductYieldDist::SetNubar(), G4StatMFChannel::SolveEqOfMotion(), G4LundStringFragmentation::Splitup(), stirf(), G4ChipsPionPlusInelasticXS::ThresholdMomentum(), G4ChipsProtonInelasticXS::ThresholdMomentum(), G4Clebsch::TriangleCoeff(), G4Clebsch::Wigner6J(), G4Clebsch::WignerLittleD(), and G4BigBanger::xProbability().

◆ log10A()

G4double G4Pow::log10A ( G4double  A) const
inline

Definition at line 201 of file G4Pow.hh.

201{ return logX(A) / lz[10]; }
G4double logX(G4double x) const
Definition: G4Pow.hh:170

References A, logX(), and lz.

◆ log10Z()

G4double G4Pow::log10Z ( G4int  Z) const
inline

Definition at line 199 of file G4Pow.hh.

199{ return lz[Z] / lz[10]; }

References lz, and Z.

◆ logA()

G4double G4Pow::logA ( G4double  A) const
inline

Definition at line 165 of file G4Pow.hh.

166{
167 return (1.0 <= A ? logBase(A) : -logBase(1. / A));
168}
G4double logBase(G4double x) const
Definition: G4Pow.hh:139

References A, and logBase().

◆ logBase()

G4double G4Pow::logBase ( G4double  x) const
inlineprivate

Definition at line 139 of file G4Pow.hh.

140{
141 G4double res;
142 if(a <= maxA2)
143 {
144 G4int i = G4int(max2 * (a - 1) + 0.5);
145 if(i > max2)
146 {
147 i = max2;
148 }
149 G4double x = a / (G4double(i) / max2 + 1) - 1;
150 res = lz2[i] + x * (1.0 - (0.5 - onethird * x) * x);
151 }
152 else if(a <= maxA)
153 {
154 G4int i = G4int(a + 0.5);
155 G4double x = a / G4double(i) - 1;
156 res = lz[i] + x * (1.0 - (0.5 - onethird * x) * x);
157 }
158 else
159 {
160 res = G4Log(a);
161 }
162 return res;
163}

References G4Log(), lz, lz2, max2, maxA, maxA2, and onethird.

Referenced by logA(), and logX().

◆ logfactorial()

G4double G4Pow::logfactorial ( G4int  Z) const
inline

◆ logX()

G4double G4Pow::logX ( G4double  x) const
inline

Definition at line 170 of file G4Pow.hh.

171{
172 G4double res = 0.0;
173 G4double a = (1.0 <= x) ? x : 1.0 / x;
174
175 if(a <= maxA)
176 {
177 res = logBase(a);
178 }
179 else if(a <= ener[2])
180 {
181 res = logen[1] + logBase(a / ener[1]);
182 }
183 else if(a <= ener[3])
184 {
185 res = logen[2] + logBase(a / ener[2]);
186 }
187 else
188 {
189 res = G4Log(a);
190 }
191
192 if(1.0 > x)
193 {
194 res = -res;
195 }
196 return res;
197}

References ener, G4Log(), logBase(), logen, and maxA.

Referenced by log10A(), and powA().

◆ logZ()

G4double G4Pow::logZ ( G4int  Z) const
inline

◆ powA()

G4double G4Pow::powA ( G4double  A,
G4double  y 
) const
inline

Definition at line 230 of file G4Pow.hh.

231{
232 return (0.0 == A ? 0.0 : expA(y * logX(A)));
233}
G4double expA(G4double A) const
Definition: G4Pow.hh:203

References A, expA(), and logX().

Referenced by G4FTFAnnihilation::Annihilate(), G4WilsonAbrasionModel::ApplyYourself(), G4LENDElastic::ApplyYourself(), G4ParticleHPElasticFS::ApplyYourself(), G4LFission::Atomas(), G4AtimaEnergyLossModel::Bethek_dedx_e(), G4WilsonAblationModel::BreakItUp(), G4UPiNuclearCrossSection::BuildPhysicsTable(), G4ParticleHPNBodyPhaseSpace::C(), G4QuasiElRatios::CalcElTot(), G4QuasiElRatios::CalcQF2IN_Ratio(), G4DiffuseElastic::CalculateNuclearRad(), G4DiffuseElasticV2::CalculateNuclearRad(), G4RKFieldIntegrator::CalculateTotalEnergy(), G4QMDMeanField::CalGraduate(), G4Reggeons::Chi_pomeron(), G4Reggeons::Chi_reggeon(), G4DiffractiveExcitation::ChooseP(), G4QGSDiffractiveExcitation::ChooseP(), G4SingleDiffractiveExcitation::ChooseX(), G4GoudsmitSaundersonMscModel::ComputeGeomPathLength(), G4GoudsmitSaundersonMscModel::ComputeTrueStepLength(), G4XAqmElastic::CrossSection(), G4XPDGElastic::CrossSection(), G4XPDGTotal::CrossSection(), G4ChipsAntiBaryonInelasticXS::CrossSectionFormula(), G4AtimaEnergyLossModel::dedx_n(), G4NonEquilibriumEvaporator::deExcite(), G4AtimaFluctuations::Dispersion(), G4LENDCrossSection::DumpPhysicsTable(), G4ParticleHPCaptureData::DumpPhysicsTable(), G4ParticleHPElasticData::DumpPhysicsTable(), G4ParticleHPFissionData::DumpPhysicsTable(), G4ParticleHPInelasticData::DumpPhysicsTable(), G4FPYSamplingOps::EvaluateWattConstants(), G4QuarkExchange::ExciteParticipants(), G4NuclearAbrasionGeometry::F(), G4IonFluctuations::Factor(), G4VCrossSectionSource::FcrossX(), G4QMDParameters::G4QMDParameters(), G4fissionEvent::G4SmpGEng(), G4fissionEvent::G4SmpNEngCf252(), G4fissionEvent::G4SmpNuDistDataU233_235(), G4fissionEvent::G4SmpNugDist(), G4ParticleHPMadlandNixSpectrum::Gamma25(), G4ParticleHPInterpolator::GetBinIntegral(), G4EMDissociationSpectrum::GetClosestApproach(), G4EMDissociationCrossSection::GetCrossSectionForProjectile(), G4ChargeExchangeProcess::GetElementCrossSection(), G4ChipsAntiBaryonElasticXS::GetExchangeT(), G4ChipsHyperonElasticXS::GetExchangeT(), G4ChipsPionMinusElasticXS::GetExchangeT(), G4ChipsPionPlusElasticXS::GetExchangeT(), G4EMDissociationSpectrum::GetGeneralE2Spectrum(), G4NuclNuclDiffuseElastic::GetHadronNucleonXscNS(), G4ComponentAntiNuclNuclearXS::GetInelasticElementCrossSection(), G4LundStringFragmentation::GetLightConeZ(), G4QGSMFragmentation::GetLightConeZ(), G4QMDMeanField::GetPotential(), G4ChipsAntiBaryonElasticXS::GetPTables(), G4ChipsHyperonElasticXS::GetPTables(), G4ChipsPionMinusElasticXS::GetPTables(), G4ChipsPionPlusElasticXS::GetPTables(), G4ChipsAntiBaryonElasticXS::GetTabValues(), G4ChipsHyperonElasticXS::GetTabValues(), G4ChipsPionMinusElasticXS::GetTabValues(), G4ChipsPionPlusElasticXS::GetTabValues(), G4ComponentAntiNuclNuclearXS::GetTotalElementCrossSection(), G4QMDMeanField::GetTotalPotential(), G4ParticleHPInterpolator::GetWeightedBinIntegral(), G4ParticleHPMadlandNixSpectrum::GIntegral(), G4FTFParameters::InitForInteraction(), G4ParticleHPVector::Integrate(), G4ParticleHPMadlandNixSpectrum::Madland(), MCGIDI_energy_NBodyPhaseSpacePDF_callback(), MCGIDI_sampling_interpolationValues(), MCGIDI_sampling_sampleX_from_pdfsOfXGivenW(), G4QMDGroundStateNucleus::packNucleons(), G4ParticleHPNBodyPhaseSpace::Prob(), ptwXY_f_integrate(), ptwXY_LogLogToLinLin(), ptwXY_pow_callback(), G4LEHadronProtonElastic::RandCosThetaDipPen(), G4QGSMSplitableHadron::SampleX(), G4QGSParticipants::SampleX(), G4ParticleHPKallbachMannSyst::SeparationEnergy(), G4FissionProductYieldDist::SetNubar(), G4AtimaEnergyLossModel::sezi_dedx_e(), G4AtimaEnergyLossModel::sezi_p_se(), G4LundStringFragmentation::Splitup(), stirf(), G4ChipsPionPlusInelasticXS::ThresholdMomentum(), and G4ChipsProtonInelasticXS::ThresholdMomentum().

◆ powN()

G4double G4Pow::powN ( G4double  x,
G4int  n 
) const

Definition at line 166 of file G4Pow.cc.

167{
168 if(0.0 == x)
169 {
170 return 0.0;
171 }
172 if(std::abs(n) > 8)
173 {
174 return std::pow(x, G4double(n));
175 }
176 G4double res = 1.0;
177 if(n >= 0)
178 {
179 for(G4int i = 0; i < n; ++i)
180 {
181 res *= x;
182 }
183 }
184 else if(n < 0)
185 {
186 G4double y = 1.0 / x;
187 G4int nn = -n;
188 for(G4int i = 0; i < nn; ++i)
189 {
190 res *= y;
191 }
192 }
193 return res;
194}

References CLHEP::detail::n, and G4InuclParticleNames::nn.

Referenced by G4ParticleHPKallbachMannSyst::A(), G4QMDReaction::ApplyYourself(), G4CascadeFinalStateAlgorithm::BetaKopylov(), G4HadPhaseSpaceKopylov::BetaKopylov(), G4FermiPhaseSpaceDecay::BetaKopylov(), G4ParticleHPNBodyPhaseSpace::C(), G4QMDMeanField::Cal2BodyQuantities(), G4QMDCollision::CalFinalStateOfTheBinaryCollisionJQMD(), G4XAqmTotal::CrossSection(), G4NonEquilibriumEvaporator::deExcite(), G4DensityEffectCalculator::DEll(), G4DensityEffectCalculator::DeltaOnceSolved(), G4DensityEffectCalculator::DFRho(), G4DensityEffectCalculator::Ell(), G4NuclearAbrasionGeometry::F(), G4DensityEffectCalculator::FermiDeltaCalculation(), G4JTPolynomialSolver::FindRoots(), G4ScreeningMottCrossSection::FormFactor2UniformHelm(), G4DensityEffectCalculator::FRho(), G4Pow(), G4fissionEvent::G4SmpNEngCf252(), G4fissionEvent::G4SmpNuDistDataPu239(), G4fissionEvent::G4SmpNuDistDataPu239_241(), G4fissionEvent::G4SmpNuDistDataU232_234_236_238(), G4fissionEvent::G4SmpNuDistDataU233_235(), G4fissionEvent::G4SmpNuDistDataU235(), G4fissionEvent::G4SmpNuDistDataU238(), G4ComponentAntiNuclNuclearXS::GetAntiHadronNucleonElCrSc(), G4ComponentAntiNuclNuclearXS::GetAntiHadronNucleonTotCrSc(), G4InuclParamMomDst::GetMomentum(), G4HadronNucleonXsc::HadronNucleonXscNS(), G4Bessel::I0(), G4Bessel::I1(), G4HETCFragment::IntegrateEmissionProbability(), G4NuclearAbrasionGeometry::P(), G4PreCompoundIon::ProbabilityDistributionFunction(), G4PreCompoundNucleon::ProbabilityDistributionFunction(), G4NuclearRadii::RadiusECS(), G4InuclSpecialFunctions::randomInuclPowers(), G4QGSMSplitableHadron::SampleX(), G4FissionProductYieldDist::SetNubar(), and G4BigBanger::xProbability().

◆ powZ()

G4double G4Pow::powZ ( G4int  Z,
G4double  y 
) const
inline

◆ Z13()

G4double G4Pow::Z13 ( G4int  Z) const
inline

Definition at line 123 of file G4Pow.hh.

123{ return pz13[Z]; }

References pz13, and Z.

Referenced by G4GEMProbability::CalcAlphaParam(), G4NeutronEvaporationProbability::CalcAlphaParam(), G4StatMFMacroMultiplicity::CalcChemicalPotentialMu(), G4GEMCoulombBarrier::CalcCompoundRadius(), G4StatMFMicroCanonical::CalcEntropyOfCompoundNucleus(), G4StatMFMicroCanonical::CalcFreeInternalEnergy(), G4StatMFMicroPartition::CalcPartitionProbability(), G4GEMProbability::CalcProbability(), G4StatMFMacroCanonical::CalculateTemperature(), G4GEMProbabilityVI::ComputeTotalProbability(), G4Fancy3DNucleus::CoulombBarrier(), G4ionEffectiveCharge::EffectiveCharge(), G4CompetitiveFission::FissionKineticEnergy(), G4ScreeningMottCrossSection::FormFactor2UniformHelm(), G4StatMFMacroTemperature::FragsExcitEnergy(), G4InuclSpecialFunctions::G4cbrt(), G4GEMCoulombBarrier::G4GEMCoulombBarrier(), G4GEMProbabilityVI::G4GEMProbabilityVI(), G4IonisParamElm::G4IonisParamElm(), G4LevelManager::G4LevelManager(), G4NuclearFermiDensity::G4NuclearFermiDensity(), G4AntiProtonField::GetBarrier(), G4KaonMinusField::GetBarrier(), G4KaonPlusField::GetBarrier(), G4PionMinusField::GetBarrier(), G4PionPlusField::GetBarrier(), G4ProtonField::GetBarrier(), G4StatMFMicroPartition::GetCoulombEnergy(), G4IonsShenCrossSection::GetIsoCrossSection(), G4NuclearLevelData::GetLevelDensity(), G4StatMFMicroPartition::GetPartitionEnergy(), G4NistManager::GetZ13(), G4WentzelOKandVIxSection::InitialiseA(), G4VPreCompoundFragment::Initialize(), G4StatMFMacroCanonical::Initialize(), G4StatMFMicroCanonical::Initialize(), G4ICRU49NuclearStoppingModel::NuclearStoppingPower(), G4StatMFChannel::PlaceFragments(), G4NuclearRadii::Radius(), G4NuclearRadii::RadiusCB(), G4NuclearRadii::RadiusECS(), G4NuclearRadii::RadiusHNGG(), G4NuclearRadii::RadiusKNGG(), G4NuclearRadii::RadiusNNGG(), G4Generator2BS::SampleDirection(), G4AntiNuclElastic::SampleInvariantT(), G4HadronElastic::SampleInvariantT(), G4GEMChannel::SampleKineticEnergy(), G4IonCoulombCrossSection::SetScreenRSquare(), G4ScreeningMottCrossSection::SetupKinematic(), G4StatMFChannel::SolveEqOfMotion(), G4EvaporationProbability::TotalProbability(), and Z23().

◆ Z23()

G4double G4Pow::Z23 ( G4int  Z) const
inline

Definition at line 125 of file G4Pow.hh.

126{
127 G4double x = Z13(Z);
128 return x * x;
129}
G4double Z13(G4int Z) const
Definition: G4Pow.hh:123

References Z, and Z13().

Referenced by G4FissionBarrier::BarashenkovFissionBarrier(), G4GEMProbability::CalcBetaParam(), G4NeutronEvaporationProbability::CalcBetaParam(), G4StatMFMacroMultiplicity::CalcChemicalPotentialMu(), G4StatMFMacroChemicalPotential::CalcChemicalPotentialNu(), G4StatMFMacroBiNucleon::CalcEnergy(), G4StatMFMacroMultiNucleon::CalcEnergy(), G4StatMFMacroTetraNucleon::CalcEnergy(), G4StatMFMacroTriNucleon::CalcEnergy(), G4StatMFMacroMultiNucleon::CalcEntropy(), G4StatMFFragment::CalcExcitationEnergy(), G4StatMFMacroBiNucleon::CalcMeanMultiplicity(), G4StatMFMacroMultiNucleon::CalcMeanMultiplicity(), G4StatMFMacroTetraNucleon::CalcMeanMultiplicity(), G4StatMFMacroTriNucleon::CalcMeanMultiplicity(), G4StatMFMicroPartition::CalcPartitionProbability(), G4StatMFMacroMultiNucleon::CalcZARatio(), G4StatMFMacroCanonical::ChooseZ(), G4UrbanAdjointMscModel::ComputeCrossSectionPerAtom(), G4UrbanMscModel::ComputeCrossSectionPerAtom(), G4AtimaFluctuations::Dispersion(), G4InuclSpecialFunctions::FermiEnergy(), G4NuclearShellModelDensity::G4NuclearShellModelDensity(), G4AntiNuclElastic::GetcosTeta1(), G4StatMFMicroPartition::GetCoulombEnergy(), G4LindhardSorensenData::GetDeltaL(), G4StatMFFragment::GetEnergy(), G4StatMFMicroPartition::GetPartitionEnergy(), G4StatMFMacroCanonical::Initialize(), G4StatMFMicroCanonical::Initialize(), G4HETCFragment::IntegrateEmissionProbability(), G4HadronElastic::SampleInvariantT(), and G4WentzelOKandVIxSection::SetupTarget().

Field Documentation

◆ ener

G4DataVector G4Pow::ener
private

Definition at line 108 of file G4Pow.hh.

Referenced by G4Pow(), and logX().

◆ fact

G4DataVector G4Pow::fact
private

Definition at line 115 of file G4Pow.hh.

Referenced by factorial(), and G4Pow().

◆ fexp

G4DataVector G4Pow::fexp
private

Definition at line 114 of file G4Pow.hh.

Referenced by expA(), and G4Pow().

◆ fpInstance

G4Pow * G4Pow::fpInstance = nullptr
staticprivate

Definition at line 98 of file G4Pow.hh.

Referenced by GetInstance().

◆ logen

G4DataVector G4Pow::logen
private

Definition at line 109 of file G4Pow.hh.

Referenced by G4Pow(), and logX().

◆ logfact

G4DataVector G4Pow::logfact
private

Definition at line 116 of file G4Pow.hh.

Referenced by G4Pow(), and logfactorial().

◆ lowa13

G4DataVector G4Pow::lowa13
private

Definition at line 111 of file G4Pow.hh.

Referenced by A13Low(), and G4Pow().

◆ lz

G4DataVector G4Pow::lz
private

Definition at line 112 of file G4Pow.hh.

Referenced by G4Pow(), log10A(), log10Z(), logBase(), logZ(), and powZ().

◆ lz2

G4DataVector G4Pow::lz2
private

Definition at line 113 of file G4Pow.hh.

Referenced by G4Pow(), and logBase().

◆ max2

const G4int G4Pow::max2 = 5
private

Definition at line 101 of file G4Pow.hh.

Referenced by G4Pow(), and logBase().

◆ maxA

G4double G4Pow::maxA
private

Definition at line 103 of file G4Pow.hh.

Referenced by A13High(), G4Pow(), logBase(), and logX().

◆ maxA2

G4double G4Pow::maxA2
private

Definition at line 105 of file G4Pow.hh.

Referenced by G4Pow(), and logBase().

◆ maxAexp

G4double G4Pow::maxAexp
private

Definition at line 106 of file G4Pow.hh.

Referenced by expA(), and G4Pow().

◆ maxLowA

G4double G4Pow::maxLowA
private

Definition at line 104 of file G4Pow.hh.

Referenced by A13(), and G4Pow().

◆ onethird

const G4double G4Pow::onethird = 1.0 / 3.0
private

Definition at line 100 of file G4Pow.hh.

Referenced by A13High(), A13Low(), expA(), G4Pow(), and logBase().

◆ pz13

G4DataVector G4Pow::pz13
private

Definition at line 110 of file G4Pow.hh.

Referenced by A13High(), G4Pow(), and Z13().


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