PhysData.jl
Luna.PhysData.N_A
— ConstantAvogadro constant
Luna.PhysData.amg
— ConstantAmagat (Loschmidt constant)
Luna.PhysData.atm
— ConstantPressure in Pascal at standard conditions (atmospheric pressure)
Luna.PhysData.au_Efield
— ConstantAtomic unit of electric field
Luna.PhysData.au_energy
— ConstantAtomic unit of energy
Luna.PhysData.au_polarisability
— ConstantAtomic unit of electric polarisability
Luna.PhysData.au_time
— ConstantAtomic unit of time
Luna.PhysData.bar
— ConstantPressure in Pascal of 1 bar
Luna.PhysData.c
— ConstantSpeed of light
Luna.PhysData.e_ratio
— ConstantRatio of electron charge squared to electron mass (for plasma)
Luna.PhysData.electron
— ConstantElectron charge
Luna.PhysData.k_B
— ConstantBoltzmann constant
Luna.PhysData.m_e
— ConstantElectron mass
Luna.PhysData.m_u
— ConstantAtomic mass unit
Luna.PhysData.roomtemp
— ConstantRoom temperature in Kelvin (ca 20 deg C)
Luna.PhysData.ħ
— ConstantReduced Planck's constant
Luna.PhysData.ε_0
— ConstantPermittivity of vacuum
Luna.PhysData.μ_0
— ConstantPermeability of vacuum
Luna.PhysData.Cnl_ADK
— MethodCnl_ADK(material)
Return the value of Cₙₗ from the ADK paper for the material
. For material
S other than noble gases, this returns missing
.
Reference: Ammosov, M. V., Delone, N. B. & Krainov, V. P. Tunnel Ionization Of Complex Atoms And Atomic Ions In Electromagnetic Field. Soviet Physics JETP 64, 1191–1194 (1986).
Luna.PhysData.data_metal
— MethodLookup tables for complex refractive indices of metals.
Luna.PhysData.density
— Functiondensity(material::Symbol, P=1.0, T=roomtemp)
For a gas material
, return the number density [m^-3] at pressure P
[bar] and temperature T
[K]. For a glass, this simply returns 1.0.
Luna.PhysData.densityspline
— Methoddensityspline(gas; Pmax, Pmin=0, N=2^10, T=roomtemp)
Create a CSpline
interpolant for the density of the gas
between pressures Pmin
and Pmax
at temperature T
. The spline is created using N
samples.
Luna.PhysData.dispersion
— Functiondispersion(order, material, λ, P=1.0, T=roomtemp; lookup=nothing)
Calculate the dispersion of order order
of a given material
at a wavelength λ
.
For gases the pressure P
(default:atmosphere) and the temperature T
(default: room temp) can also be specified. lookup::Bool
determines whether a lookup table or a Sellmeier expansion is used for the refractive index (default is material dependent).
Examples
julia> dispersion(2, :BK7, 400e-9) * 1e30 * 1e-3 # convert to fs^2/mm
122.03632107303108
Luna.PhysData.dispersion_func
— Functiondispersion_func(order, material, P=1, T=roomtemp; lookup=nothing)
Get a function to calculate dispersion. Arguments are the same as for dispersion
.
Luna.PhysData.dispersion_func
— Methoddispersion_func(order, n)
Get a function that calculates dispersion of order order
for a refractive index given by n(λ)
.
Luna.PhysData.fresnel
— Methodfresnel(n2, θi; n1=1.0)
Calcualte reflection coefficients from Fresnel's equations.
Luna.PhysData.ionisation_potential
— Methodionisation_potential(material; unit=:SI)
Return the first ionisation potential of the material
in a specific unit (default: SI). Possible units are :SI
, :atomic
and :eV
.
Luna.PhysData.lookup_glass
— Methodlookup_glass(material::Symbol)
Create a CSpline
interpolant for look-up-table values of the refractive index.
Luna.PhysData.lookup_metal
— Methodlookup_metal(material::Symbol)
Create a CSpline
interpolant for look-up-table values of the refractive index.
Luna.PhysData.lookup_mirror
— Methodlookup_mirror(type)
Create a CSpline
interpolant for the complex-valued reflectivity of a mirror of type
.
Luna.PhysData.polarisability_difference
— Methodpolarisability_difference(material; unit=:SI)
Return the difference in polarisability between the ground state and the ion for the material
. unit
can be :SI
or :atomic
Reference: Wang, K. et al. Static dipole polarizabilities of atoms and ions from Z=1 to 20 calculated within a single theoretical scheme. Eur. Phys. J. D 75, 46 (2021).
Luna.PhysData.pressure
— Functionpressure(gas, density, T=roomtemp)
Calculate the pressure in bar of the gas
at number density density
and temperature T
.
Luna.PhysData.process_mirror_data
— Methodprocess_mirror_data(λR, R, λGDD, GDD, λ0, λmin, λmax; fitorder=5, windowwidth=20e-9)
Process reflectivity and group-delay dispersion data for a mirror and create a transfer function for a frequency-domain electric field representing the mirror.
Arguments:
λR
: wavelength samples for reflectivity in SI units (m)R
: mirror reflectivity (between 0 and 1)λGDD
: wavelength samples for GDD in SI units (m)GDD
: GDD in SI units (s²)λ0
: central wavelength (used to remove any overall group delay)λmin
,λmax
: bounds of the wavelength region to apply the transfer function over
Keyword arguments
fitorder
: order of polynomial fit to use in removing overall group delay (default: 5)windowwidth
: wavelength width of the smoothing region outside(λmin, λmax)
for the window in SI units (default: 20e-9, i.e. 20 nm)
Luna.PhysData.quantum_numbers
— Methodquantum_numbers(material)
Return the quantum numbers of the material
for use in the PPT ionisation rate.
Luna.PhysData.raman_parameters
— Methodraman_parameters(material)
Get the Raman parameters for material
.
Fields
Fields in the returned named tuple must include:
kind::Symbol
: one of:molecular
or:intermediate
or:normedsdo
If kind == :molecular
then the following must also be specified:
rotation::Symbol
: only:nonrigid
or:none
supported at present.vibration::Symbol
: only:sdo
or:none
supported at present.
If rotation == :nonrigid
then the following must also be specified:
B::Real
: the rotational constant [1/m]Δα::Real
: molecular polarizability anisotropy [m^3]qJodd::Integer
: nuclear spin parameter for oddJ
qJeven::Integer
: nuclear spin parameter for evenJ
D::Real=0.0
: centrifugal constant [1/m]
Along with one of:
τ2r::Real
: coherence time [s]Bρr::Real
: density dependent broadening coefficient [Hz/amagat]
If both τ2r
and Bρr
are specified, then Bρr
takes precedence. If Bρr
is specified then we also need:
Aρr::Real
: self diffusion coefficient [Hz amagat]
If vibration == :sdo
then the following must also be specified:
Ωv::Real
: vibrational frequency [rad/s]dαdQ::Real
: isotropic averaged polarizability derivative [m^2]μ::Real
: reduced molecular mass [kg]
Along with one of:
τ2v::Real
: coherence time [s]Bρv::Real
: density dependent broadening coefficient [Hz/amagat]
If both τ2v
and Bρv
are specified, then Bρv
takes precedence. If Bρv
is specified then we also need:
Aρv::Real
: self diffusion coefficient [Hz amagat]
And can also add (if necessary) a constant offset:
Cv::Real
: constant linewidth offset [Hz]
If kind == :intermediate
then the following must be specified
ωi::Vector{Real}
[rad/s], central angular freqenciesAi::Vector{Real}
, amplitudesΓi::Vector{Real}
[rad/s], Gaussian widthsγi::Vector{Real}
[rad/s], Lorentzian widths
References
[1] Phys. Rev. A, 94, 023816 (2016) [2] Phys. Rev. A, 85, 043820 (2012) [3] Phys. Rev. A, 92, 063828 (2015) [4] Journal of Raman Spectroscopy 2, 133 (1974) [5] J. Phys. Chem., 91, 41 (1987) [6] Applied Spectroscopy 23, 211 (1969) [7] Phys. Rev. A, 34, 3, 1944 (1986) [8] Can. J. Phys., 44, 4, 797 (1966) [9] G. V. MIKHAtLOV, SOVIET PHYSICS JETP, vol. 36, no. 9, (1959). [10] Phys. Rev. A, 33, 5, 3113 (1986) [11] IEEE Journal of Quantum Electronics 1986, 22 (2), 332–336. https://doi.org/10.1109/JQE.1986.1072945. [12] Phys. Rev. Lett. 1998, 81 (6), 1215–1218. https://doi.org/10.1103/PhysRevLett.81.1215. [13] Optics Communications 1987, 64 (4), 393–397. https://doi.org/10.1016/0030-4018(87)90258-6. [14] Science Advances 2020, 6 (34), eabb5375. https://doi.org/10.1126/sciadv.abb5375. [15] Long, The Raman Effect; John Wiley & Sons, Ltd, 2002; [16] IEEE Journal of Quantum Electronics, vol. 24, no. 10, pp. 2076–2080, Oct. 1988, doi: 10.1109/3.8545. [17] Journal of Raman Spectroscopy, vol. 22, no. 11, pp. 607–611, 1991, doi: 10.1002/jrs.1250221103. [18] Hollenbeck and Cantrell, JOSA B 19, 2886-2892 (2002). https://doi.org/10.1364/JOSAB.19.002886
Luna.PhysData.ref_index
— Functionref_index(material, λ, P=1.0, T=roomtemp; lookup=nothing)
Get refractive index for any material at wavelength given in SI units.
Luna.PhysData.ref_index_fun
— Functionref_index_fun(material::Symbol, P=1.0, T=roomtemp; lookup=nothing)
Get function which returns refractive index.
Luna.PhysData.ref_index_fun
— Methodref_index_fun(gases, T=roomtemp)
Get function which returns ref index for mixture as function of wavelength and densities.
Luna.PhysData.ref_index_fun
— Methodref_index_fun(gases, P, T=roomtemp; lookup=nothing)
Get function which returns refractive index for gas mixture. gases
is a Tuple
of gas identifiers (Symbol
s) and P
is a Tuple
of equal length containing pressures.
Luna.PhysData.sellmeier_crystal
— Methodsellmeier_crystal(material, axis)
Sellmeier for crystals. Returns function of wavelength in μm which in turn returns the refractive index directly. Possible values for axis
depend on the type of crystal.
Luna.PhysData.sellmeier_gas
— Methodsellmeier_gas(material::Symbol)
Return function for linear polarisability γ, i.e. susceptibility of a single particle, calculated from Sellmeier expansions.
Luna.PhysData.sellmeier_glass
— Methodsellmeier_glass(material::Symbol)
Sellmeier for glasses. Returns function of wavelength in μm which in turn returns the refractive index directly
Luna.PhysData.wlfreq
— Methodwlfreq(ωλ)
Change from ω (angular frequency) to λ (wavelength) and vice versa
Luna.PhysData.ΔλΔω
— MethodΔλΔω(Δλ, λ)
Convert Δλ (wavelength bandwidth) at λ (central wavelength) to Δω (angular frequency bandwidth)
Luna.PhysData.γ3_gas
— Methodγ3_gas(material::Symbol; source=nothing)
Calculate single-molecule third-order hyperpolarisability of a gas at given wavelength(s) and at room temperature. If source
== :Bishop
: Uses reference values to calculate γ If source
== :Lehmeier
(default): Uses scaling factors to calculate χ3 at 1 bar and scales by density to get to a single molecule i.e. the hyperpolarisability
References: [1] Journal of Chemical Physics, AIP, 91, 3549-3551 (1989) [2] Chemical Reviews, 94, 3-29 (1994) [3] Optics Communications, 56(1), 67–72 (1985) [4] Phys. Rev. A, vol. 42, 2578 (1990) [5] Optics Letters Vol. 40, No. 24 (2015)) [6] Phys. Rev. A 2012, 85 (4), 043820. https://doi.org/10.1103/PhysRevA.85.043820. [7] Phys. Rev. A, 32, no. 6, 3454, (1985), doi: 10.1103/PhysRevA.32.3454.
Luna.PhysData.γ_Börzsönyi
— Methodγ_Börzsönyi(B1, C1, B2, C2)
Sellmeier expansion for linear susceptibility from Applied Optics 47, 27, 4856 (2008) at room temperature and atmospheric pressure
Luna.PhysData.γ_JCT
— Methodγ_JCT(B1, C1, B2, C2, B3, C3)
Adapted Sellmeier expansion for helium made to fit high frequency data Phys. Rev. A 92, 033821 (2015)
Luna.PhysData.γ_Peck
— Methodγ_Peck(B1, C1, B2, C2, dens)
Sellmeier expansion for linear susceptibility from J. Opt. Soc. Am. 67, 1550 (1977)
Luna.PhysData.γ_QuanfuHe
— Methodγ_QuanfuHe(A, B, C, dens)
Sellmeier expansion for CH4, SF6 and N2O from Atmospheric Chemistry and Physics 2021, 21 (19), 14927–14940. https://doi.org/10.5194/acp-21-14927-2021.
Luna.PhysData.γ_Zhang
— Methodγ_Zhang(A, B, C, dens)
Sellmeier expansion for Oxygen from Applied Optics 50, 35, 6484 (2011)
Luna.PhysData.χ1
— Functionχ1(gas::Symbol, λ, P=1.0, T=roomtemp)
Calculate χ1 at wavelength λ in SI units, pressure P in bar and temperature T in Kelvin. Gases only.
Luna.PhysData.χ1_fun
— Methodχ1_fun(gas::Symbol)
Get function to return χ1 (linear susceptibility) for gases as a function of wavelength in SI units, pressure in bar, and temperature in Kelvin.