Cloudy
Spectral Synthesis Code for Astrophysics
Loading...
Searching...
No Matches
physconst_template.h File Reference
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Macros

#define pow2(a)
 
#define pow3(a)
 

Functions

 NEW_CONSTANT (EE, 2.718281828459045235360287)
 
 NEW_CONSTANT (EULER, 0.577215664901532860606512090082)
 
 NEW_CONSTANT (EXPEULER2, 0.74930600128844902360587)
 
 NEW_CONSTANT (PI, 3.141592653589793238462643)
 
 NEW_CONSTANT (PI2, 6.283185307179586476925287)
 
 NEW_CONSTANT (PI4, 12.56637061435917295385057)
 
 NEW_CONSTANT (PI8, 25.13274122871834590770115)
 
 NEW_CONSTANT (SQRT2, 1.414213562373095048801689)
 
 NEW_CONSTANT (SQRTPI, 1.772453850905516027298167)
 
 NEW_CONSTANT (SQRTPIBY2, 1.253314137315500251207883)
 
 NEW_CONSTANT (LN_TWO, 0.6931471805599453094172321)
 
 NEW_CONSTANT (LN_TEN, 2.302585092994045684017991)
 
 NEW_CONSTANT (LOG10_E, 0.4342944819032518276511289)
 
 NEW_CONSTANT (OPTDEP2EXTIN, 1.085736204758129569127822)
 
 NEW_CONSTANT (RADIAN, 57.29577951308232087679815)
 
 NEW_CONSTANT (SOLAR_MASS, 1.98841e33)
 
 NEW_CONSTANT (SOLAR_LUMINOSITY, 3.828e33)
 
 NEW_CONSTANT (MBOL_ZERO_LUMINOSITY, 3.0128e35)
 
 NEW_CONSTANT (AU, 1.49597870700e13)
 
 NEW_CONSTANT (ELECTRON_MASS, 9.1093837139e-28)
 
 NEW_CONSTANT (ELECTRON_MASS_U, 5.485799090441e-4)
 
 NEW_CONSTANT (PROTON_MASS, 1.67262192595e-24)
 
 NEW_CONSTANT (PROTON_MASS_U, 1.0072764665789)
 
 NEW_CONSTANT (ALPHA_MASS_U, 4.001506179129)
 
 NEW_CONSTANT (BOLTZMANN, 1.380649e-16)
 
 NEW_CONSTANT (SPEEDLIGHT, 2.99792458e10)
 
 NEW_CONSTANT (HPLANCK, 6.62607015e-27)
 
 NEW_CONSTANT (AVOGADRO, 6.02214076e23)
 
 NEW_CONSTANT (GRAV_CONST, 6.67430e-8)
 
 NEW_CONSTANT (ELEM_CHARGE, 1.602176634e-19)
 
 NEW_CONSTANT (RYD_INF, 1.0973731568157e5)
 
 NEW_CONSTANT (HIONPOT, 0.999466508345)
 
 NEW_CONSTANT (HE2IONPOT, 3.99963199547)
 
 NEW_CONSTANT (HMINUSIONPOT, 0.055432956)
 
 NEW_CONSTANT (ATOMIC_MASS_UNIT, 1.66053906892e-24)
 
 NEW_CONSTANT (MOL_MASS_CONST, AVOGADRO *ATOMIC_MASS_UNIT)
 
 NEW_CONSTANT (AS1RAD, RADIAN *3600.)
 
 NEW_CONSTANT (SQAS1SR, pow2(AS1RAD))
 
 NEW_CONSTANT (SQAS_SKY, PI4 *SQAS1SR)
 
 NEW_CONSTANT (PARSEC, AU *AS1RAD)
 
 NEW_CONSTANT (MEGAPARSEC, 1.e6 *PARSEC)
 
 NEW_CONSTANT (H_BAR, HPLANCK/(2.*PI))
 
 NEW_CONSTANT (ELEM_CHARGE_ESU, ELEM_CHARGE *SPEEDLIGHT/10.)
 
 NEW_CONSTANT (ELECTRIC_CONST, 1.e11/(PI4 *pow2(SPEEDLIGHT)))
 
 NEW_CONSTANT (HION_LTE_POP, pow2(HPLANCK)/(PI2 *BOLTZMANN *ELECTRON_MASS))
 
 NEW_CONSTANT (SAHA, 4.141330276355497e-16)
 
 NEW_CONSTANT (ERG1CM, HPLANCK *SPEEDLIGHT)
 
 NEW_CONSTANT (T1CM, HPLANCK *SPEEDLIGHT/BOLTZMANN)
 
 NEW_CONSTANT (KJMOL1CM, ERG1CM *AVOGADRO/1e10)
 
 NEW_CONSTANT (H_RYD_FACTOR, 1./(1.+ELECTRON_MASS_U/PROTON_MASS_U))
 
 NEW_CONSTANT (HE_RYD_FACTOR, 1./(1.+ELECTRON_MASS_U/ALPHA_MASS_U))
 
 NEW_CONSTANT (WAVNRYD, 1./RYD_INF)
 
 NEW_CONSTANT (RYDLAM, 1.e8/RYD_INF)
 
 NEW_CONSTANT (EN1RYD, HPLANCK *SPEEDLIGHT *RYD_INF)
 
 NEW_CONSTANT (TE1RYD, HPLANCK *SPEEDLIGHT *RYD_INF/BOLTZMANN)
 
 NEW_CONSTANT (EVDEGK, ELEM_CHARGE *1.e7/BOLTZMANN)
 
 NEW_CONSTANT (EVRYD, HPLANCK *SPEEDLIGHT *RYD_INF/ELEM_CHARGE *1.e-7)
 
 NEW_CONSTANT (EN1EV, EN1RYD/EVRYD)
 
 NEW_CONSTANT (FR1RYD, SPEEDLIGHT *RYD_INF)
 
 NEW_CONSTANT (HNU3C2, 2.*HPLANCK *SPEEDLIGHT *pow3(RYD_INF))
 
 NEW_CONSTANT (FR1RYDHYD, SPEEDLIGHT *RYD_INF *HIONPOT)
 
 NEW_CONSTANT (HBAReV, H_BAR/EN1EV)
 
 NEW_CONSTANT (RYDLAMHYD, RYDLAM/HIONPOT)
 
 NEW_CONSTANT (STEFAN_BOLTZ, pow2(PI *pow2(BOLTZMANN))/(60.*pow3(H_BAR) *pow2(SPEEDLIGHT)))
 
 NEW_CONSTANT (FREQ_1EV, SPEEDLIGHT *RYD_INF/EVRYD)
 
 NEW_CONSTANT (FINE_STRUCTURE, pow2(ELEM_CHARGE_ESU)/SPEEDLIGHT/H_BAR)
 
 NEW_CONSTANT (FINE_STRUCTURE2, pow2(FINE_STRUCTURE))
 
 NEW_CONSTANT (BOHR_RADIUS_CM, FINE_STRUCTURE/(PI4 *RYD_INF))
 
 NEW_CONSTANT (TWO_PHOT_CONST, 9.*pow3(FINE_STRUCTURE2) *FR1RYD/2048.)
 
 NEW_CONSTANT (COLL_CONST, SAHA *BOLTZMANN/HPLANCK)
 
 NEW_CONSTANT (MILNE_CONST, 4.123475589581428e+11)
 
 NEW_CONSTANT (TRANS_PROB_CONST, PI4 *HPLANCK *FINE_STRUCTURE/ELECTRON_MASS)
 
 NEW_CONSTANT (ABSOR_COEFF_CONST, SQRTPI *pow2(ELEM_CHARGE_ESU)/(ELECTRON_MASS *SPEEDLIGHT))
 
 NEW_CONSTANT (SIGMA_THOMSON, PI8/3.*pow2(FINE_STRUCTURE *H_BAR/(ELECTRON_MASS *SPEEDLIGHT)))
 
 NEW_CONSTANT (HC_ERG_ANG, HPLANCK *SPEEDLIGHT *1e8)
 
 NEW_CONSTANT (JEANS, PI *BOLTZMANN/(GRAV_CONST *ATOMIC_MASS_UNIT))
 
 NEW_CONSTANT (FREE_FREE_EMIS, 1.032526505912028e-11)
 
 NEW_CONSTANT (FREE_FREE_ABS, 1.036997355477227e-38)
 

Detailed Description

* physical constants used by Cloudy, mostly taken from
* >>refer       phys    const   CODATA 2022, https://physics.nist.gov/cuu/Constants/
* <BR><BR>
* NB - these are all printed with the "print constants" command, 
* which is in parse_print.cpp. Using the NEW_CONSTANT macro guarantees
* that any constant added here is printed there as well 

Macro Definition Documentation

◆ pow2

#define pow2 ( a)
Value:
((a)*(a))

Referenced by anomal(), t_gaunt::brems_opac(), Bruggeman(), C2_PR78(), C6cs123(), Ca20cs123(), chi2_func(), cnewton(), ColliderList::ColliderList(), collision_strength_VF01(), CollisSuppres(), conorm(), CoolEvaluate(), CS_l_mixing_PS64(), CS_l_mixing_PS64_expI(), CS_PercivalRichards78(), CSresolver(), DebyeDeriv(), dense_parametric_wind(), Drude(), eeBremsCooling(), eeBremsSpectrum(), exp10(), exp10f(), exp10i(), FastVoigtH(), FastVoigtH(), Fe26cs123(), FillExtraLymanLine(), t_gaunt::gauntff(), gauss_legendre(), get_iso_statw(), GetProbDistr_HighLimit(), GrainMakeDiffuse(), h21_t_ge_10(), h21_t_lt_10(), Hydcs123(), hydro_Fujimoto_deexcit(), hydro_Lebedev_deexcit(), iso_cascade(), iso_collide(), iso_prt_pops(), iso_radiative_recomb_effective(), iso_solve(), iso_state_lifetime(), phymir_state< X, Y, NP, NSTR >::lgConvergedRestart(), linfit(), mie_auxiliary(), mie_cs(), mie_find_slope(), mie_integrate(), mie_read_opc(), mie_read_szd(), mole_effects(), Ne10cs123(), t_ran::normal(), OpacityAddTotal(), Atom_LevelN::operator()(), t_gaunt::p_gauntff_vec_sub(), phymir_state< X, Y, NP, NSTR >::p_phygrm(), t_wavl::p_RefIndex(), t_gaunt::p_setup_brems(), phymir_state< X, Y, NP, NSTR >::p_setup_next_hyperblock(), t_mesh::p_SetupMesh(), t_ran::p_ZigTailNormal(), pah1_fun(), pah2_fun(), pah3_fun(), qheat(), RebinSingleCell(), Recomb_Seaton59(), RT_continuum(), RT_diffuse(), RT_line_all(), S2Aul(), SanityCheckGaunt(), save_opacity(), sinpar(), size_distr(), Stognienko(), t_physconst::t_physconst(), uderiv(), y1psa(), and ZoneStart().

◆ pow3

Function Documentation

◆ NEW_CONSTANT() [1/78]

NEW_CONSTANT ( ABSOR_COEFF_CONST ,
SQRTPI * pow2ELEM_CHARGE_ESU)/(ELECTRON_MASS *SPEEDLIGHT )

This is the atomic physics constant that enters the calculation of the absorption coefficient given the energy separation between two levels, and the oscillator strength. Numerical value: 0.01497361235435992 Units: esu^2 / (g cm/s) = cm^2 / s The absorption coefficient is defined as in Eq. (9-33) of Mihalas, Stellar Atmospheres, 2nd Edition

◆ NEW_CONSTANT() [2/78]

NEW_CONSTANT ( ALPHA_MASS_U ,
4. 001506179129 )

mass of the alpha particle, in u

◆ NEW_CONSTANT() [3/78]

NEW_CONSTANT ( AS1RAD ,
RADIAN * 3600. )

number of arcsec in 1 radian, 206264.806

◆ NEW_CONSTANT() [4/78]

NEW_CONSTANT ( ATOMIC_MASS_UNIT ,
1.66053906892e- 24 )

atomic mass unit, gram

◆ NEW_CONSTANT() [5/78]

NEW_CONSTANT ( AU ,
1. 49597870700e13 )

astronomical unit, cm, nearly the length of the semimajor axis of the Earth's elliptical orbit around the sun

◆ NEW_CONSTANT() [6/78]

NEW_CONSTANT ( AVOGADRO ,
6. 02214076e23 )

Avogadro constant, exact

◆ NEW_CONSTANT() [7/78]

NEW_CONSTANT ( BOHR_RADIUS_CM ,
FINE_STRUCTURE/ PI4 *RYD_INF )

Bohr radius in cm, 5.29177249e-9

◆ NEW_CONSTANT() [8/78]

NEW_CONSTANT ( BOLTZMANN ,
1.380649e- 16 )

this is the Boltzmann factor, erg/K, exact

◆ NEW_CONSTANT() [9/78]

NEW_CONSTANT ( COLL_CONST ,
SAHA *BOLTZMANN/ HPLANCK )

this is the square of the value roughly equal to 8.629e-6 that appears in converting collision strengths to rates. The constant is h^2/((2PI*me)^3/2 * k^1/2).

◆ NEW_CONSTANT() [10/78]

NEW_CONSTANT ( EE ,
2. 718281828459045235360287 )

the number e

◆ NEW_CONSTANT() [11/78]

NEW_CONSTANT ( ELECTRIC_CONST ,
1.e11/ PI4 *pow2(SPEEDLIGHT) )

electric constant, in F/m, 8.854e-12

◆ NEW_CONSTANT() [12/78]

NEW_CONSTANT ( ELECTRON_MASS ,
9.1093837139e- 28 )

electron mass, gram

◆ NEW_CONSTANT() [13/78]

NEW_CONSTANT ( ELECTRON_MASS_U ,
5.485799090441e- 4 )

electron mass, in u

◆ NEW_CONSTANT() [14/78]

NEW_CONSTANT ( ELEM_CHARGE ,
1.602176634e- 19 )

elementary charge, in C (SI units), exact, to use this must convert to cgs

◆ NEW_CONSTANT() [15/78]

NEW_CONSTANT ( ELEM_CHARGE_ESU ,
ELEM_CHARGE *SPEEDLIGHT/ 10. )

elementary charge, in ESU, 4.8032e-10

◆ NEW_CONSTANT() [16/78]

NEW_CONSTANT ( EN1EV ,
EN1RYD/ EVRYD )

ergs per eV, 1.602176e-012

◆ NEW_CONSTANT() [17/78]

NEW_CONSTANT ( EN1RYD ,
HPLANCK *SPEEDLIGHT * RYD_INF )

ergs per inf mass Ryd, 2.180e-11

◆ NEW_CONSTANT() [18/78]

NEW_CONSTANT ( ERG1CM ,
HPLANCK * SPEEDLIGHT )

number of ergs per wavenumber, 1.9864e-16

◆ NEW_CONSTANT() [19/78]

NEW_CONSTANT ( EULER ,
0. 577215664901532860606512090082 )

the Euler constant (aka Euler-Mascheroni constant or gamma)

◆ NEW_CONSTANT() [20/78]

NEW_CONSTANT ( EVDEGK ,
ELEM_CHARGE *1.e7/ BOLTZMANN )

Kelvins per eV, 1.1604e4

◆ NEW_CONSTANT() [21/78]

NEW_CONSTANT ( EVRYD ,
HPLANCK *SPEEDLIGHT *RYD_INF/ELEM_CHARGE *1.e- 7 )

eV per inf mass Ryd, 13.606

◆ NEW_CONSTANT() [22/78]

NEW_CONSTANT ( EXPEULER2 ,
0. 74930600128844902360587 )

the Euler constant parameter exp(-gamma/2)

◆ NEW_CONSTANT() [23/78]

NEW_CONSTANT ( FINE_STRUCTURE ,
pow2(ELEM_CHARGE_ESU)/SPEEDLIGHT/ H_BAR )

the fine-structure constant a= 2pi e^2/hc 7.297 352 533 x 10-3

◆ NEW_CONSTANT() [24/78]

NEW_CONSTANT ( FINE_STRUCTURE2 ,
pow2(FINE_STRUCTURE)  )

the square of the fine-structure constant

◆ NEW_CONSTANT() [25/78]

NEW_CONSTANT ( FR1RYD ,
SPEEDLIGHT * RYD_INF )

frequency of one Ryd for infinite mass nuclei, 3.289842e15

◆ NEW_CONSTANT() [26/78]

NEW_CONSTANT ( FR1RYDHYD ,
SPEEDLIGHT *RYD_INF * HIONPOT )

frequency of ionization potential of H (not inf mass), 3.288087e15 - never used

◆ NEW_CONSTANT() [27/78]

NEW_CONSTANT ( FREE_FREE_ABS ,
1.036997355477227e- 38 )

◆ NEW_CONSTANT() [28/78]

NEW_CONSTANT ( FREE_FREE_EMIS ,
1.032526505912028e- 11 )

◆ NEW_CONSTANT() [29/78]

NEW_CONSTANT ( FREQ_1EV ,
SPEEDLIGHT *RYD_INF/ EVRYD )

the frequency of one eV, 2.418e14

◆ NEW_CONSTANT() [30/78]

NEW_CONSTANT ( GRAV_CONST ,
6.67430e- 8 )

Gravitational constant, cm^3/g/s^2

◆ NEW_CONSTANT() [31/78]

NEW_CONSTANT ( H_BAR ,
HPLANCK/ 2.*PI )

h/2pi, 1.05457e-27

◆ NEW_CONSTANT() [32/78]

NEW_CONSTANT ( H_RYD_FACTOR ,
1./ 1.+ELECTRON_MASS_U/PROTON_MASS_U )

ratio of Rydberg constants R_H/R_inf, 0.999456

◆ NEW_CONSTANT() [33/78]

NEW_CONSTANT ( HBAReV ,
H_BAR/ EN1EV )

H_BAR in eV sec, 6.582e-16

◆ NEW_CONSTANT() [34/78]

NEW_CONSTANT ( HC_ERG_ANG ,
HPLANCK *SPEEDLIGHT * 1e8 )

Express hc in erg Angstroms, approximately equal to 1.9864e-8 Used to get h*nu from wavelength in Angstroms

◆ NEW_CONSTANT() [35/78]

NEW_CONSTANT ( HE2IONPOT ,
3. 99963199547 )

ionization potential of real He^+ ion, in inf mass ryd, based on CODATA 2006, uncertainty 13e-11, calculated by Peter van Hoof

◆ NEW_CONSTANT() [36/78]

NEW_CONSTANT ( HE_RYD_FACTOR ,
1./ 1.+ELECTRON_MASS_U/ALPHA_MASS_U )

ratio of Rydberg constants R_He/R_inf, 0.999863

◆ NEW_CONSTANT() [37/78]

NEW_CONSTANT ( HION_LTE_POP ,
pow2(HPLANCK)/(PI2 *BOLTZMANN *ELECTRON_MASS)  )

this is the factor that appears in front of Boltzmann factor to get LTE level populations for hydrogenic ions. It is given in the first parts of section 5 of part 2 of hazy, and is actually ( planck^2 / (2 pi m_e k ) )^3/2, but cannot evaluate powers here, must raise this to 3/2 when used, HION_LTE_POP, 5.556e-11 cm^2 K

◆ NEW_CONSTANT() [38/78]

NEW_CONSTANT ( HIONPOT ,
0. 999466508345 )

ionization potential of real hydrogen atom, in inf mass ryd, based on CODATA 2006, uncertainty 10e-12, calculated by Peter van Hoof

◆ NEW_CONSTANT() [39/78]

NEW_CONSTANT ( HMINUSIONPOT ,
0. 055432956 )

ionization potential of H- in inf mass ryd,

>refer phys const T. Anderson, 2004, Phys. Reports, 394, 157

◆ NEW_CONSTANT() [40/78]

NEW_CONSTANT ( HNU3C2 ,
2.*HPLANCK *SPEEDLIGHT * pow3RYD_INF )

2 h FR1RYD^3 / c^2 for infinite mass nucleus, 0.5250

◆ NEW_CONSTANT() [41/78]

NEW_CONSTANT ( HPLANCK ,
6.62607015e- 27 )

Planck's constant, exact

◆ NEW_CONSTANT() [42/78]

NEW_CONSTANT ( JEANS ,
PI *BOLTZMANN/ GRAV_CONST *ATOMIC_MASS_UNIT )

Jeans constant as in Hazy 3 – other constant factors are quoted

◆ NEW_CONSTANT() [43/78]

NEW_CONSTANT ( KJMOL1CM ,
ERG1CM *AVOGADRO/ 1e10 )

kJ/mol per unit wavenumber, 0.01196266

◆ NEW_CONSTANT() [44/78]

NEW_CONSTANT ( LN_TEN ,
2. 302585092994045684017991 )

ln(10)

◆ NEW_CONSTANT() [45/78]

NEW_CONSTANT ( LN_TWO ,
0. 6931471805599453094172321 )

ln(2)

◆ NEW_CONSTANT() [46/78]

NEW_CONSTANT ( LOG10_E ,
0. 4342944819032518276511289 )

log(e)

◆ NEW_CONSTANT() [47/78]

NEW_CONSTANT ( MBOL_ZERO_LUMINOSITY ,
3. 0128e35 )

luminosity of a star with absolute bolometric magnitude 0, in erg s-1

>refer phys const https://pdg.lbl.gov/2021/reviews/rpp2021-rev-astrophysical-constants.pdf

◆ NEW_CONSTANT() [48/78]

NEW_CONSTANT ( MEGAPARSEC ,
1.e6 * PARSEC )

megaparsec in cm, 3.085678e24

◆ NEW_CONSTANT() [49/78]

NEW_CONSTANT ( MILNE_CONST ,
4.123475589581428e+ 11 )

this is the square of the value roughly equal to 4.123e11 that appears in the integration of photoionization cross-sections to obtain recombination coefficients.

◆ NEW_CONSTANT() [50/78]

NEW_CONSTANT ( MOL_MASS_CONST ,
AVOGADRO * ATOMIC_MASS_UNIT )

molar mass constant, g/mol

◆ NEW_CONSTANT() [51/78]

NEW_CONSTANT ( OPTDEP2EXTIN ,
1. 085736204758129569127822 )

factor that converts optical depth into extinction in mags, 2.5 log e

◆ NEW_CONSTANT() [52/78]

NEW_CONSTANT ( PARSEC ,
AU * AS1RAD )

parsec in cm, 3.085678e18

◆ NEW_CONSTANT() [53/78]

NEW_CONSTANT ( PI ,
3. 141592653589793238462643 )

pi

◆ NEW_CONSTANT() [54/78]

NEW_CONSTANT ( PI2 ,
6. 283185307179586476925287 )

2*pi

◆ NEW_CONSTANT() [55/78]

NEW_CONSTANT ( PI4 ,
12. 56637061435917295385057 )

4*pi

◆ NEW_CONSTANT() [56/78]

NEW_CONSTANT ( PI8 ,
25. 13274122871834590770115 )

8*pi

◆ NEW_CONSTANT() [57/78]

NEW_CONSTANT ( PROTON_MASS ,
1.67262192595e- 24 )

proton mass, gram

◆ NEW_CONSTANT() [58/78]

NEW_CONSTANT ( PROTON_MASS_U ,
1. 0072764665789 )

proton mass, in u

◆ NEW_CONSTANT() [59/78]

NEW_CONSTANT ( RADIAN ,
57. 29577951308232087679815 )

180/pi

◆ NEW_CONSTANT() [60/78]

NEW_CONSTANT ( RYD_INF ,
1. 0973731568157e5 )

infinite mass rydberg constant, in cm^-1

◆ NEW_CONSTANT() [61/78]

NEW_CONSTANT ( RYDLAM ,
1.e8/ RYD_INF )

Angstrom per infinite mass Ryd, 911.2671

◆ NEW_CONSTANT() [62/78]

NEW_CONSTANT ( RYDLAMHYD ,
RYDLAM/ HIONPOT )

wavelength (A) of ionization potential of Hydrogen, 911.7535 - never used

◆ NEW_CONSTANT() [63/78]

NEW_CONSTANT ( SAHA ,
4.141330276355497e- 16 )

SAHA is ( h^2/2/pi/m/k )^3/2, is correct constant for free electron SAHA, 4.14132e-16 cm^3 K^(3/2)

◆ NEW_CONSTANT() [64/78]

NEW_CONSTANT ( SIGMA_THOMSON ,
PI8/3.* pow2FINE_STRUCTURE *H_BAR/(ELECTRON_MASS *SPEEDLIGHT) )

Thomson cross-section, cm^2

◆ NEW_CONSTANT() [65/78]

NEW_CONSTANT ( SOLAR_LUMINOSITY ,
3. 828e33 )

◆ NEW_CONSTANT() [66/78]

NEW_CONSTANT ( SOLAR_MASS ,
1. 98841e33 )

◆ NEW_CONSTANT() [67/78]

NEW_CONSTANT ( SPEEDLIGHT ,
2. 99792458e10 )

speed of light, cm/s, exact

◆ NEW_CONSTANT() [68/78]

NEW_CONSTANT ( SQAS1SR ,
pow2(AS1RAD)  )

number of square arcsec in 1 steradian, 4.254517e10

◆ NEW_CONSTANT() [69/78]

NEW_CONSTANT ( SQAS_SKY ,
PI4 * SQAS1SR )

number of square arcsec in the whole sky, 5.3463838e11

◆ NEW_CONSTANT() [70/78]

NEW_CONSTANT ( SQRT2 ,
1. 414213562373095048801689 )

sqrt(2)

◆ NEW_CONSTANT() [71/78]

NEW_CONSTANT ( SQRTPI ,
1. 772453850905516027298167 )

sqrt(pi)

◆ NEW_CONSTANT() [72/78]

NEW_CONSTANT ( SQRTPIBY2 ,
1. 253314137315500251207883 )

sqrt(pi/2)

◆ NEW_CONSTANT() [73/78]

NEW_CONSTANT ( STEFAN_BOLTZ ,
pow2(PI *pow2(BOLTZMANN))/(60.*pow3(H_BAR) *pow2(SPEEDLIGHT))  )

Stefan-Boltzmann constant, 5.6704e-5

◆ NEW_CONSTANT() [74/78]

NEW_CONSTANT ( T1CM ,
HPLANCK *SPEEDLIGHT/ BOLTZMANN )

degrees kelvin per unit wavenumber, 1.4388

◆ NEW_CONSTANT() [75/78]

NEW_CONSTANT ( TE1RYD ,
HPLANCK *SPEEDLIGHT *RYD_INF/ BOLTZMANN )

the temperature of 1 Rydberg te1ryd is h/k is temp of 1 Rydberg, 1.579e5

◆ NEW_CONSTANT() [76/78]

NEW_CONSTANT ( TRANS_PROB_CONST ,
PI4 *HPLANCK *FINE_STRUCTURE/ ELECTRON_MASS )

This is the constant used in converting oscillator strengths to As. The formula is Aul, TRANS_PROB_CONST * f(u,l) * wavenumber^2. TRANS_PROB_CONST is 0.667025

◆ NEW_CONSTANT() [77/78]

NEW_CONSTANT ( TWO_PHOT_CONST ,
9.*pow3(FINE_STRUCTURE2) *FR1RYD/ 2048. )

the two photon constant as defined by Breit & Teller, as in equation 4 of Spitzer & Greenstein 51, 2.18313

◆ NEW_CONSTANT() [78/78]

NEW_CONSTANT ( WAVNRYD ,
1./ RYD_INF )

number of Ryd per wavenumber, 9.11267e-6