Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Data Structures | Namespaces | Functions | Variables
G4Log.hh File Reference
#include <limits>
#include <stdint.h>
#include "G4Types.hh"

Go to the source code of this file.

Data Structures

union  G4LogConsts::ieee754
 

Namespaces

 G4LogConsts
 

Functions

G4double G4LogConsts::get_log_px (const G4double x)
 
G4double G4LogConsts::get_log_qx (const G4double x)
 
uint64_t G4LogConsts::dp2uint64 (G4double x)
 
G4double G4LogConsts::uint642dp (uint64_t ll)
 
G4float G4LogConsts::uint322sp (G4int x)
 
uint32_t G4LogConsts::sp2uint32 (G4float x)
 
G4double G4LogConsts::getMantExponent (const G4double x, G4double &fe)
 Like frexp but vectorising and the exponent is a double. More...
 
G4float G4LogConsts::getMantExponentf (const G4float x, G4float &fe)
 Like frexp but vectorising and the exponent is a float. More...
 
G4double G4Log (G4double x)
 
G4float G4LogConsts::get_log_poly (const G4float x)
 
G4float G4Logf (G4float x)
 
void logv (const uint32_t size, G4double const *__restrict__ iarray, G4double *__restrict__ oarray)
 
void G4Logv (const uint32_t size, G4double const *__restrict__ iarray, G4double *__restrict__ oarray)
 
void logfv (const uint32_t size, G4float const *__restrict__ iarray, G4float *__restrict__ oarray)
 
void G4Logfv (const uint32_t size, G4float const *__restrict__ iarray, G4float *__restrict__ oarray)
 

Variables

const G4double G4LogConsts::LOG_UPPER_LIMIT = 1e307
 
const G4double G4LogConsts::LOG_LOWER_LIMIT = 0
 
const G4double G4LogConsts::SQRTH = 0.70710678118654752440
 
const G4float G4LogConsts::MAXNUMF = 3.4028234663852885981170418348451692544e38f
 
const G4float G4LogConsts::LOGF_UPPER_LIMIT = MAXNUMF
 
const G4float G4LogConsts::LOGF_LOWER_LIMIT = 0
 
const G4float G4LogConsts::PX1logf = 7.0376836292E-2f
 
const G4float G4LogConsts::PX2logf = -1.1514610310E-1f
 
const G4float G4LogConsts::PX3logf = 1.1676998740E-1f
 
const G4float G4LogConsts::PX4logf = -1.2420140846E-1f
 
const G4float G4LogConsts::PX5logf = 1.4249322787E-1f
 
const G4float G4LogConsts::PX6logf = -1.6668057665E-1f
 
const G4float G4LogConsts::PX7logf = 2.0000714765E-1f
 
const G4float G4LogConsts::PX8logf = -2.4999993993E-1f
 
const G4float G4LogConsts::PX9logf = 3.3333331174E-1f
 
const G4float G4LogConsts::SQRTHF = 0.707106781186547524f
 

Function Documentation

G4double G4Log ( G4double  x)
inline

Definition at line 227 of file G4Log.hh.

References fe, G4LogConsts::get_log_px(), G4LogConsts::get_log_qx(), G4LogConsts::getMantExponent(), G4LogConsts::LOG_LOWER_LIMIT, G4LogConsts::LOG_UPPER_LIMIT, G4LogConsts::SQRTH, and test::x.

Referenced by G4EnergyLossForExtrapolator::AverageScatteringAngle(), G4VEnergyLossProcess::BuildLambdaTable(), G4PairProductionRelModel::CalcLPMFunctions(), G4KleinNishinaModel::ComputeCrossSectionPerAtom(), G4KleinNishinaCompton::ComputeCrossSectionPerAtom(), G4BetheHeitlerModel::ComputeCrossSectionPerAtom(), G4UrbanMscModel::ComputeCrossSectionPerAtom(), G4eeToTwoGammaModel::ComputeCrossSectionPerElectron(), G4MollerBhabhaModel::ComputeCrossSectionPerElectron(), G4BraggIonModel::ComputeCrossSectionPerElectron(), G4ICRU73QOModel::ComputeCrossSectionPerElectron(), G4BraggModel::ComputeCrossSectionPerElectron(), G4MuBetheBlochModel::ComputeCrossSectionPerElectron(), G4BetheBlochModel::ComputeCrossSectionPerElectron(), G4MollerBhabhaModel::ComputeDEDXPerVolume(), G4BraggModel::ComputeDEDXPerVolume(), G4BraggIonModel::ComputeDEDXPerVolume(), G4ICRU73QOModel::ComputeDEDXPerVolume(), G4MuBetheBlochModel::ComputeDEDXPerVolume(), G4BetheBlochModel::ComputeDEDXPerVolume(), G4hBremsstrahlungModel::ComputeDMicroscopicCrossSection(), G4hPairProductionModel::ComputeDMicroscopicCrossSection(), G4MuBremsstrahlungModel::ComputeDMicroscopicCrossSection(), G4MuPairProductionModel::ComputeDMicroscopicCrossSection(), G4LivermoreBremsstrahlungModel::ComputeDXSectionPerAtom(), G4SeltzerBergerModel::ComputeDXSectionPerAtom(), G4GoudsmitSaundersonMscModel::ComputeGeomPathLength(), G4UrbanMscModel::ComputeGeomPathLength(), G4MuBremsstrahlungModel::ComputeMicroscopicCrossSection(), G4MuPairProductionModel::ComputeMicroscopicCrossSection(), G4UrbanMscModel::ComputeTheta0(), G4WentzelOKandVIxSection::ComputeTransportCrossSectionPerAtom(), G4WentzelVIRelXSection::ComputeTransportCrossSectionPerAtom(), G4UrbanMscModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTrueStepLength(), G4WentzelVIRelModel::ComputeTrueStepLength(), G4GoudsmitSaundersonMscModel::ComputeTrueStepLength(), G4UrbanMscModel::ComputeTrueStepLength(), G4MuPairProductionModel::ComputMuPairLoss(), G4EquilibriumEvaporator::deExcite(), G4ionEffectiveCharge::EffectiveCharge(), G4NucleiModel::fillZoneRadii(), G4InuclSpecialFunctions::G4cbrt(), G4F20GEMProbability::G4F20GEMProbability(), G4O17GEMProbability::G4O17GEMProbability(), G4PhysicsLnVector::G4PhysicsLnVector(), G4PhysicsLogVector::G4PhysicsLogVector(), G4NucleiModel::generateInteractionLength(), G4PreCompoundProton::GetOpt2(), G4mplIonisationWithDeltaModel::Initialise(), G4mplIonisationModel::Initialise(), G4MuPairProductionModel::Initialise(), G4VEmModel::InitialiseElementSelectors(), G4NucleiModel::initializeCascad(), G4SynchrotronRadiation::InvSynFracInt(), G4Pow::logX(), G4DipBustGenerator::PolarAngle(), G4VEmProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4InuclSpecialFunctions::randomGauss(), G4PhysicsLnVector::Retrieve(), G4PhysicsLogVector::Retrieve(), G4DipBustGenerator::SampleDirection(), G4DeltaAngle::SampleDirection(), G4RayleighAngularGenerator::SampleDirection(), G4ModifiedTsai::SampleDirection(), G4UniversalFluctuation::SampleFluctuations(), G4WentzelVIRelModel::SampleScattering(), G4WentzelVIModel::SampleScattering(), G4GoudsmitSaundersonMscModel::SampleScattering(), G4LivermoreGammaConversionModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreBremsstrahlungModel::SampleSecondaries(), G4SeltzerBergerModel::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4BetheHeitlerModel::SampleSecondaries(), G4PairProductionRelModel::SampleSecondaries(), G4eBremsstrahlungRelModel::SampleSecondaries(), G4eBremParametrizedModel::SampleSecondaries(), G4eeToTwoGammaModel::SampleSecondaries(), G4MuPairProductionModel::SampleSecondaries(), G4MuBremsstrahlungModel::SampleSecondaries(), G4MuBetheBlochModel::SampleSecondaries(), G4GoudsmitSaundersonTable::SampleTheta(), G4PhysicsLnVector::ScaleVector(), G4PhysicsLogVector::ScaleVector(), G4DiscreteGammaTransition::SelectGamma(), G4IonisParamMat::SetMeanExcitationEnergy(), G4EnergyLossForExtrapolator::TrueStepLength(), and G4NucleiModel::zoneIntegralWoodsSaxon().

228 {
229  const G4double original_x = x;
230 
231  /* separate mantissa from exponent */
232  G4double fe;
234 
235  // blending
236  x > G4LogConsts::SQRTH? fe+=1. : x+=x ;
237  x -= 1.0;
238 
239  /* rational form */
241 
242  //for the final formula
243  const G4double x2 = x*x;
244  px *= x;
245  px *= x2;
246 
247  const G4double qx = G4LogConsts::get_log_qx(x);
248 
249  G4double res = px / qx ;
250 
251  res -= fe * 2.121944400546905827679e-4;
252  res -= 0.5 * x2 ;
253 
254  res = x + res;
255  res += fe * 0.693359375;
256 
257  if (original_x > G4LogConsts::LOG_UPPER_LIMIT)
258  res = std::numeric_limits<G4double>::infinity();
259  if (original_x < G4LogConsts::LOG_LOWER_LIMIT) // THIS IS NAN!
260  res = - std::numeric_limits<G4double>::quiet_NaN();
261 
262  return res;
263 }
G4double get_log_qx(const G4double x)
Definition: G4Log.hh:123
const G4double SQRTH
Definition: G4Log.hh:79
G4double getMantExponent(const G4double x, G4double &fe)
Like frexp but vectorising and the exponent is a double.
Definition: G4Log.hh:186
G4double get_log_px(const G4double x)
Definition: G4Log.hh:100
const G4double LOG_LOWER_LIMIT
Definition: G4Log.hh:77
const G4double LOG_UPPER_LIMIT
Definition: G4Log.hh:76
double G4double
Definition: G4Types.hh:76
G4fissionEvent * fe
G4float G4Logf ( G4float  x)
inline

Definition at line 308 of file G4Log.hh.

References fe, G4LogConsts::get_log_poly(), G4LogConsts::getMantExponentf(), G4LogConsts::LOGF_LOWER_LIMIT, G4LogConsts::LOGF_UPPER_LIMIT, G4LogConsts::SQRTHF, and test::x.

309 {
310  const G4float original_x = x;
311 
312  G4float fe;
314 
315  x > G4LogConsts::SQRTHF? fe+=1.f : x+=x ;
316  x -= 1.0f;
317 
318  const G4float x2 = x*x;
319 
321  res *= x2*x;
322 
323  res += -2.12194440e-4f * fe;
324  res += -0.5f * x2;
325 
326  res= x + res;
327 
328  res += 0.693359375f * fe;
329 
330  if (original_x > G4LogConsts::LOGF_UPPER_LIMIT)
331  res = std::numeric_limits<G4float>::infinity();
332  if (original_x < G4LogConsts::LOGF_LOWER_LIMIT)
333  res = -std::numeric_limits<G4float>::quiet_NaN();
334 
335  return res;
336 }
G4float getMantExponentf(const G4float x, G4float &fe)
Like frexp but vectorising and the exponent is a float.
Definition: G4Log.hh:210
const G4float LOGF_UPPER_LIMIT
Definition: G4Log.hh:269
float G4float
Definition: G4Types.hh:77
const G4float LOGF_LOWER_LIMIT
Definition: G4Log.hh:270
const G4float SQRTHF
Definition: G4Log.hh:303
G4float get_log_poly(const G4float x)
Definition: G4Log.hh:282
G4fissionEvent * fe
void G4Logfv ( const uint32_t  size,
G4float const *__restrict__  iarray,
G4float *__restrict__  oarray 
)
void G4Logv ( const uint32_t  size,
G4double const *__restrict__  iarray,
G4double *__restrict__  oarray 
)
void logfv ( const uint32_t  size,
G4float const *__restrict__  iarray,
G4float *__restrict__  oarray 
)
void logv ( const uint32_t  size,
G4double const *__restrict__  iarray,
G4double *__restrict__  oarray 
)