Developer documentation

This page contains info for people interested in contributing code to Korg. If you have questions, do not hesitate to ask.

Julia development

If you've never developed a package in Julia, here are some tips.

  • Once you have a local copy of Korg, you can make it importable locally by deving it from the Julia environment in which you would like to run it: Pkg.dev("/path/to/Korg"). If you run Pkg.dev("Korg") instead, Julia will automatically clone the repo to ~/.julia/dev/Korg. (See the Pkg documentation for details.)
  • When working on your local copy of Korg, Revise is easiest way to make and test changes without constantly restarting your Julia session.
  • To run the test suite locally, start a Julia session in the Korg root directory and run ]activate . then test. This will automatically run test/runtests.jl in the test environment.
  • Documentation is generated with Documentor. Continuous integration on github will ensure that the documentation is generated without errors, but it won't catch all formatting problems. If you wish to generate documentation locally in order to check that everything is as expected, start a Julia session in Korg/docs, activate the test environment (]activate .), and run instantiate. This downloads and installs the docs dependences, and will only have to be run once. To generate documentation, run julia --project make.jl on the command line (the --project flag activates the local environment). The generated docs can be served from docs/build.

Code guidelines

  • Try to be explicit about units throughout the code, particularly when not using CGS.
  • Whenever possible, calculations should be precise up to a factor of $10^{-3}$. When it's easy and inexpensive, they should be precise to $10^{-5}$ or better. There are exceptions to this (e.g. the Voigt function), but ideally they will eventually be vanquished.
  • Ensure types are generic enough to support dual numbers and autodifferentiation.
  • Limit lines to 100 characters.
  • Use ASCII characters in the names of functions that are part of the public API.
  • Unless they will never be called elsewhere, provide docstrings describing the inputs, assumptions and outputs of any functions you write.

Continuum absorption

Steps for implementing new continuum sources of absorption:

  • Define a helper function that computes a single absorption coefficient (in units of cm⁻²). The function should accept ν (in Hz) and T (in K) as the first and second arguments, respectively. The convention is for it should share a name with the corresponding public function, but have an underscore pre-appended (e.g. we define _H_I_bf to help implement H_I_bf).
  • The public function is the function constructed and returned by Korg.ContinuumAbsorption.bounds_checked_absorption that wraps the above helper function. This wrapper function implements bounds-checking for ν and T and supports the keyword arguments described in Continuum Absorption Kwargs.
  • Add a docstring describing the new function. At the very least, please describe any non-standard arguments and include a reference in the docstring to the source where the function was taken from.
  • Add a line to doc/src/API.md under the Continuum absorption heading to render the docstring of your new function.
  • Add a line to total_continuum_absorption that calls the new public function for absorption.

The first two steps may not apply for sources that don't directly depend on ν and T (e.g. absorption from scattering).

Manifest.toml and Project.toml

Read more about these files here. These files are used by the julia package manager. Project.toml records dependencies, and you'll notice that the test and docs directories have their own Project.tomls for test and documentation-specific dependencies. Manifest.toml records exact package versions used when a package was run. It enables someone else to reproduce the exact environment and results later. There are a few directories containing scripts that generate Korg's data files which have there own Projects.tomls and Manifest.tomls, for example data/bf_cross-sections/.

Where to put data

If you are adding data to Korg, data/README provides an overview of the options and how to decide between them.

Complete API

Here are all the documented methods in Korg.

Korg.LineMethod
Line(wl::F, log_gf::F, species::Species, E_lower::F, 
     gamma_rad::Union{F, Missing}=missing, gamma_stark::Union{F, Missing}=missing, 
     vdw::Union{F, Tuple{F, F}, Missing}, missing) where F <: Real

Arguments:

  • wl: wavelength (Assumed to be in cm if < 1, otherwise in Å)
  • log_gf: (log base 10) oscillator strength (unitless)
  • species: the Species associated with the line
  • E_lower: The energy (excitiation potential) of the lower energy level (eV)

Optional Arguments (these override default recipes):

  • gamma_rad: Fundemental width

  • gamma_stark: per-perturber Stark broadening width at 10,000 K (s⁻¹).

  • vdW: If this is present, it may may be

    • log10(γ_vdW), assumed if negative
    • 0, corresponding to no vdW broadening
    • A fudge factor for the Unsoeld approximation, assumed if between 0 and 20
    • The ABO parameters as packed float (assumed if >= 20) or a Tuple, (σ, α).

    This behavior is intended to mirror that of Turbospectrum as closely as possible.

See approximate_gammas for more information on the default recipes for gamma_stark and vdW.

Note the the "gamma" values here are FWHM, not HWHM, of the Lorenztian component of the line profile, and are in units of s⁻¹.

source
Korg.MolecularCrossSectionMethod
MolecularCrossSection(linelist, wls; cutoff_alpha=1e-30, log_temp_vals=3:0.025:5, verbose=true)

Precompute the molecular absorption cross section for a given linelist and set of wavelengths. The MolecularCrossSection object can be passed to synthesize and potentially speed up the calculation significantly. At present, Korg only supports precomputed cross-sections created by this function.

Arguments

  • linelist: A vector of Line objects representing the molecular linelist. These must be of the same species.
  • wls: A vector of wavelength ranges (in Å) at which to precompute the cross section. These must match the wavelengths used for any subsequent synthesis exactly.

Keyword Arguments

  • cutoff_alpha (default: 1e-30): The value of the single-line absorption coefficient (in cm^-1) at which to truncate the profile.
  • log_temp_vals (default: 3:0.025:5): The log10 of the temperatures at which to precompute the cross-section.
  • verbose (default: true): Whether to print progress information.
Tip

The default values of vmic_vals, log_temp_vals, and cutoffalphawere chosen to ensure that lines in the APOGEE linelist ([getAPOGEEDR17linelist`](@ref)) could be accurately reproduced (better than 10^-3 everywhere). You should verify that they yield acceptable accuracy for other applications by comparing spectra synthesize with and without precomputing the molecular cross-section.

source
Korg.PlanarAtmosphereMethod
PlanarAtmosphere(atm::ShellAtmosphere)

Construct a planar atmosphere with the data from a shell atmosphere. Mostly useful for testing.

source
Korg.ShellAtmosphereMethod
ShellAtmosphere(atm::PlanarAtmosphere, R)

Construct a shell atmosphere with the data from a planar atmosphere and an outer radius. Mostly useful for testing.

source
Korg.SpeciesType

Represents an atom or molecule (a Formula) with a particular number of electrons (regardless of their configuration).

source
Korg.SpeciesMethod
Species(code)

Parse the "species code" in many of the forms in which it is often specified and return an object representing the species. code can be either a string or a float.

Examples

  • "H I" -> H I
  • "H 1" -> H I
  • "H 1" -> H I
  • "H_1" -> H I
  • "H.I" -> H I
  • "H 2" -> H II
  • "H2" -> H₂
  • "H" -> H I
  • "01.00" → H I
  • "02.01" → He II
  • "02.1000" → He II
  • "0608" → CO I
Note

To parse at compile time, use the species string macro, i.e. species"H I". This is important in hot inner loops.

Warning

MOOG codes which include isotopic information will not be parsed correctly by this function, though read_linelist handles them correctly.

source
Korg._get_multi_X_HMethod

Given a vector of abundances, A_X, get [I+J+K/H], where Zs = [I,J,K] is a vector of atomic numbers. This is used to calculate, for example, [α/H] and [metals/H].

source
Korg.air_to_vacuumMethod
air_to_vacuum(λ; cgs=λ<1)

Convert λ from an air to vacuum. λ is assumed to be in Å if it is ⩾ 1, in cm otherwise. Formula from Birch and Downs (1994) via the VALD website.

See also: vacuum_to_air.

source
Korg.apply_LSFMethod
apply_LSF(flux, wls, R; window_size=4)

Applies a gaussian line spread function the the spectrum with flux vector flux and wavelength vector wls with constant spectral resolution, $R = \lambda/\Delta\lambda$, where $\Delta\lambda$ is the LSF FWHM. The window_size argument specifies how far out to extend the convolution kernel in standard deviations.

For the best match to data, your wavelength range should extend a couple $\Delta\lambda$ outside the region you are going to compare.

If you are convolving many spectra defined on the same wavelenths to observational resolution, you will get much better performance using compute_LSF_matrix.

Keyword Arguments

  • window_size (default: 4): how far out to extend the convolution kernel in units of the LSF width (σ, not HWHM)
  • renormalize_edge (default: true): whether or not to renormalize the LSF at the edge of the wl range. This doen't matter as long as synth_wls extends to large and small enough wavelengths.
Warning
  • This is a naive, slow implementation. Do not use it when performance matters.

  • apply_LSF will have weird behavior if your wavelength grid is not locally linearly-spaced. It is intended to be run on a fine wavelength grid, then downsampled to the observational (or otherwise desired) grid.

source
Korg.apply_rotationMethod
apply_rotation(flux, wls, vsini, ε=0.6)

Given a spectrum flux sampled at wavelengths wls for a non-rotating star, compute the spectrum that would emerge given projected rotational velocity vsini and linear limb-darkening coefficient ε: $I(\mu) = I(1) (1 - \varepsilon + arepsilon \mu))$. See, for example, Gray equation 18.14.

source
Korg.approximate_gammasMethod
approximate_gammas(wl, species, E_lower; ionization_energies=Korg.ionization_energies)

A simplified form of the Unsoeld (1955) approximation for van der Waals broadening and the Cowley 1971 approximation for Stark broadening, evaluated at 10,000 K. Used for atomic lines with no vdW and stark broadening info in the linelist.

Returns (γ_stark, log10(γ_vdW)) in Hz, where these are the per-perturber quantities. For autoionizing lines (those for which Eupper > χ), Returns 0.0 for γvdW. Note the the "gamma" values here are FWHM, not HWHM, of the Lorenztian component of the line profile.

In the calculation of n*², uses the approximation that $\overbar{r^2} = 5/2 {n^*}^4 / Z^2$ which neglects the dependence on the angular momentum quantum number, l, in the the form given by Warner 1967 (the earliest english work reporting the Unsoeld result).

source
Korg.approximate_radiative_gammaMethod
approximate_radiative_gamma(wl, log_gf)

Approximate radiate broadening parameter. When using this, make sure that log_gf is the true value (not adjusted for isotopic abundance).

source
Korg.autodiffable_convMethod
autodiffable_conv(f, g)

Compute the convolution of two vectors, f and g. This is a wrapped version of DSP.conv that handles ForwardDiff.Dual types via the chain rule.

source
Korg.blackbodyMethod
blackbody(T, λ)

The value of the Planck blackbody function for temperature T at wavelength λ [cm].

source
Korg.bracket_line_interpolatorFunction
bracket_line_interpolator(m, λ₀, T, nₑ, ξ, λmin, λmax; kwargs...)

This routine numerically convolves the two components of the Brackett line stark profile (quasistatic/Holtsmark and impact) and the doppler profile, if necessary. It returns a pair containing the interpolator and the distance from the line center at which it is defined.

Arguments

  • m: the principle quantum number of the upper level
  • λ₀: the line center in Å
  • T: the temperature [K]
  • nₑ: the electron number density [cm^-3]
  • ξ: the microturbulence [cm/s]
  • λmin: the minimum wavelength at which the profile should be computed (used to avoid calculations outside the required region)
  • λmax: the mxinimum wavelength at which the profile should be computed

Keyword Arguments

  • n_wavelength_points (default=201): the number of wavelengths at which to sample the profiles quasistatic profiles to be convolved.
  • window_size (default=5): the size of the wavelength range over which the profiles should be

calculated, in units of the characteristic profile width

source
Korg.brackett_line_stark_profilesMethod
brackett_line_stark_profiles(m, λs, λ₀, T, nₑ)

Stark-broadened line profile (specialized to Brackett series). Translated and heavily adapted from HLINOP.f by Barklem, who adapted it from Peterson and Kurucz. Mostly follows Griem 1960, and Griem 1967. Ions and distant electrons have E fields which can be treated quasi-statically, leading to a Holtsmark broadening profile.

Returns a pair of vectors containing the impact and quasistatic profiles, respectively.

Arguents:

  • m: the upper level of the transition
  • λs: the wavelengths at which to calculate the profile [cm]
  • λ₀: the line center [cm]
  • T: temperature [K]
  • nₑ: electron number density [cm^-3]
source
Korg.brackett_oscillator_strengthMethod
brackett_oscillator_strength(n, m)

The oscillator strength of the transition from the 4th to the mth energy level of hydrogen. Adapted from HLINOP.f by Peterson and Kurucz.

Comparison to the values in Goldwire 1968 indicates that this is accurate to 10^-4 for the Brackett series.

source
Korg.chemical_equilibriumMethod
chemical_equilibrium(T, nₜ, nₑ, absolute_abundances, ionization_energies, 
                     partition_fns, log_equilibrium_constants; x0=nothing)

Iteratively solve for the number density of each species. Returns a pair containing the electron number density and Dict mapping species to number densities.

arguments:

  • the temperature, T, in K
  • the total number density nₜ
  • the electron number density nₑ
  • A Dict of absolute_abundances, NX/Ntotal
  • a Dict of ionization energies, ionization_energies. The keys of act as a list of all atoms.
  • a Dict of partition functions, partition_fns
  • a Dict of log molecular equilibrium constants, log_equilibrium_constants, in partial pressure form. The keys of equilibrium_constants act as a list of all molecules.

keyword arguments:

  • x0 (default: nothing) is an initial guess for the solution (in the format internal to chemical_equilibrium). If not supplied, a good guess is computed by neglecting molecules.
  • electron_number_density_warn_threshold (default: 0.25) is the fractional difference between the calculated electron number density and the model atmosphere electron number density at which a warning is issued.

The system of equations is specified with the number densities of the neutral atoms as free parameters. Each equation specifies the conservation of a particular species, e.g. (simplified)

n(O) = n(CO) + n(OH) + n(O I) + n(O II) + n(O III).

In this equation:

  • n(O), the number density of oxygen atoms in any form comes absolute_abundances and the total

number density (supplied later)

  • n(O I) is a free parameter. The numerical solver is varying this to satisfy the system of equations.
  • n(O II), and n(O III) come from the Saha (ionization) equation given n(O I)
  • n(CO) and n(OH) come from the molecular equilibrium constants K, which are precomputed over a range of temperatures.

Equilibrium constants are defined in terms of partial pressures, so e.g.

K(OH)  ==  (p(O) p(H)) / p(OH)  ==  (n(O) n(H)) / n(OH)) kT
source
Korg.compute_LSF_matrixMethod
compute_LSF_matrix(synth_wls, obs_wls, R; kwargs...)

Construct a sparse matrix, which when multiplied with a flux vector defined over wavelenths synth_wls, applies a gaussian line spead function (LSF) and resamples to the wavelenths obswls.

Arguments

  • synth_wls: the synthesis wavelengths. This should be a range, a vector of ranges, or a vector of wavelength values which is linearly spaced to within a tollerance of the step_tolerance kwarg.
  • obs_wls: the wavelengths of the observed spectrum
  • R: the resolving power, $R = \lambda/\Delta\lambda$

Keyword Arguments

  • window_size (default: 4): how far out to extend the convolution kernel in units of the LSF width (σ, not HWHM)
  • verbose (default: true): whether or not to emit warnings and information to stdout/stderr.
  • renormalize_edge (default: true): whether or not to renormalize the LSF at the edge of the wl range. This doen't matter as long as synth_wls extends to large and small enough wavelengths.
  • step_tolerance: the maximum difference between adjacent wavelengths in synth_wls for them to be considered linearly spaced. This is only used if synth_wls is a vector of wavelengths rather than a range or vector or ranges.

For the best match to data, your wavelength range should extend a couple $\Delta\lambda$ outside the region you are going to compare.

Korg.apply_LSF can apply an LSF to a single flux vector efficiently. This function is relatively slow, but one the LSF matrix is constructed, convolving spectra to observational resolution via matrix multiplication is fast.

source
Korg.containedMethod
contained(value, interval)

Returns whether value is contained by interval.

Examples

julia> contained(0.5, Interval(1.0,10.0))
false
julia> contained(5.0, Interval(1.0,10.0))
true
source
Korg.contained_sliceMethod
contained_slice(vals, interval)

Returns a range of indices denoting the elements of vals (which are assumed to be sorted in increasing order) that are contained by interval. When no entries are contained by interval, this returns (1,0) (which is a valid empty slice).

source
Korg.doppler_widthMethod
doppler_width(λ₀ T, m, ξ)

The standard deviation of of the doppler-broadening profile. In standard spectroscopy texts, the Doppler width often refers to σ√2, but this is σ

source
Korg.exponential_integral_1Method
exponential_integral_1(x)

Compute the first exponential integral, E1(x). This is a rough approximation lifted from Kurucz's VCSE1F. Used in brackett_line_profile.

source
Korg.format_A_XMethod
format_A_X(default_metals_H, default_alpha_H, abundances; kwargs... )

Returns a 92 element vector containing abundances in $A(X)$ ($\log_{10}(X/H) + 12$) format for elements from hydrogen to uranium.

Arguments

You can specify abundance with these positional arguments. All are optional, but if default_alpha_H is provided, default_metals_H must be as well.

  • default_metals_H (default: 0), i.e. [metals/H] is the $\log_{10}$ solar-relative abundance of elements heavier than He. It is overridden by default_alpha and abundances on a per-element basis.
  • default_alpha_H (default: same as default_metals_H), i.e. [alpha/H] is the $\log_{10}$ solar-relative abundance of the alpha elements (defined to be O, Ne, Mg, Si, S, Ar, Ca, and Ti). It is overridden by abundances on a per-element basis.
  • abundances is a Dict mapping atomic numbers or symbols to [$X$/H] abundances. (Set solar_relative=false to use $A(X)$ abundances instead.) These override default_metals_H. This is the only way to specify an abundance of He that is non-solar.

Keyword arguments

  • solar_relative (default: true): When true, interpret abundances as being in [$X$/H] ($\log_{10}$ solar-relative) format. When false, interpret them as $A(X)$ abundances, i.e. $A(x) = \log_{10}(n_X/n_\mathrm{H}) + 12$, where $n_X$ is the number density of $X$. Note that abundances not specified default to the solar value still depend on the solar value, as they are set according to default_metals_H and default_alpha_H.
  • solar_abundances (default: Korg.asplund_2020_solar_abundances) is the set of solar abundances to use, as a vector indexed by atomic number. Korg.asplund_2009_solar_abundances and Korg.grevesse_2007_solar_abundances are also provided for convenience.
source
Korg.get_APOGEE_DR17_linelistMethod
get_APOGEE_DR17_linelist(; include_water=true)

The APOGEE DR 17 linelist. It ranges from roughly 15,000 Å to 17,000 Å. It is nearly the same at the DR 16 linelist described in Smith+ 2021.

source
Korg.get_VALD_solar_linelistMethod
get_VALD_solar_linelist()

Get a VALD "extract stellar" linelist produced at solar parameters, with the "threshold" value set to 0.01. It was downloaded on 2021-05-20. It is intended to be used for quick tests only.

source
Korg.get_alpha_HMethod
get_alpha_H(A_X; solar_abundances=default_solar_abundances)

Calculate [α/H] given a vector, A_X of absolute abundances, $A(X) = \log_{10}(n_α/n_\mathrm{H})$. Here, the alpha elements are defined to be O, Ne, Mg, Si, S, Ar, Ca, Ti. See also get_metals_H.

Keyword Arguments

  • solar_abundances (default: Korg.asplund_2020_solar_abundances) is the set of solar abundances to use, as a vector indexed by atomic number. Korg.asplund_2009_solar_abundances, Korg.grevesse_2007_solar_abundances, and Korg.magg_2022_solar_abundances are also provided for convenience.
source
Korg.get_atomsMethod
get_atoms(x)

Returns an array view containing the atomic number of each atom that makes up the formula or species x. E.g. get_atoms(Korg.species"H2O") yields [1, 1, 8].

source
Korg.get_electron_number_densitiesMethod

getelectronnumberdensities(atm::ModelAtmosphere) = [l.electronnumber_density for l in atm.layers]

This is a convienince functions for making plots, etc. Note that it doesn't access quantities in a memory-efficient order.

source
Korg.get_gas_pressuresMethod
get_gas_pressures(atm::ModelAtmosphere)

This is a convienince functions for making plots, etc. Note that it doesn't access quantities in a memory-efficient order.

source
Korg.get_log_nKMethod
get_log_nK(mol, T, log_equilibrium_constants)

Given a molecule, mol, a temperature, T, and a dictionary of log equilibrium constants in partial pressure form, return the base-10 log equilibrium constant in number density form, i.e. log10(nK) where nK = n(A)n(B)/n(AB).

source
Korg.get_metals_HMethod
get_metals_H(A_X; solar_abundances=default_solar_abundances, ignore_alpha=true)

Calculate [metals/H] given a vector, A_X of absolute abundances, $A(X) = \log_{10}(n_M/n_\mathrm{H})$. See also get_alpha_H.

Keyword Arguments

  • solar_abundances (default: Korg.asplund_2020_solar_abundances) is the set of solar abundances to use, as a vector indexed by atomic number. Korg.asplund_2009_solar_abundances, Korg.grevesse_2007_solar_abundances, and Korg.magg_2022_solar_abundances are also provided for convenience.
  • ignore_alpha (default: true): Whether or not to ignore the alpha elements when calculating [metals/H]. If true, [metals/H] is calculated using all elements heavier than He. If false, the alpha elements (here defined as C, O, Ne, Mg, Si, S, Ar, Ca, Ti) are ignored.
source
Korg.get_number_densitiesMethod
get_number_densities(atm::ModelAtmosphere) = [l.number_density for l in atm.layers]

This is a convienince functions for making plots, etc. Note that it doesn't access quantities in a memory-efficient order.

source
Korg.get_tau_5000sMethod
get_tau_5000s(atm::ModelAtmosphere) = [l.tau_5000 for l in atm.layers]

This is a convienince functions for making plots, etc. Note that it doesn't access quantities in a memory-efficient order.

source
Korg.get_tempsMethod

get_temps(atm::ModelAtmosphere) = [l.temp for l in atm.layers]

This is a convienince functions for making plots, etc. Note that it doesn't access quantities in a memory-efficient order.

source
Korg.get_zsMethod

get_zs(atm::ModelAtmosphere) = [l.z for l in atm.layers]

This is a convienince functions for making plots, etc. Note that it doesn't access quantities in a memory-efficient order.

source
Korg.greim_1960_KnmMethod
greim_1960_Knm(n, m)

Knm constants as defined by Griem 1960 for the long range Holtsmark profile. This function includes only the values for Brackett lines.

$K_{nm} = C_{nm} 2 π c / λ^2$ where $C_{nm} F = Δω$ and $F$ is the ion field. See Griem 1960 EQs 7 and 12. This works out to $K_{nm} = λ/F$.

source
Korg.holtsmark_profileMethod
holtsmark_profile(β, P)

Calculates the Holtsmark profile for broadening of hydrogen lines by quasistatic charged particles. Adapted from SOFBET in HLINOP by Peterson and Kurucz. Draws heavily from Griem 1960.

source
Korg.hummer_mihalas_U_HMethod
hummer_mihalas_U_H(T, nH, nHe, ne)

!!!note This is experimental, and not used by Korg for spectral synthesis.

Calculate the partition function of neutral hydrogen using the occupation probability formalism from Hummer and Mihalas 1988. See hummer_mihalas_w for details.

source
Korg.hummer_mihalas_wMethod
hummer_mihalas_w(T, n_eff, nH, nHe, ne; use_hubeny_generalization=false)

Calculate the correction, w, to the occupation fraction of a hydrogen energy level using the occupation probability formalism from Hummer and Mihalas 1988, optionally with the generalization by Hubeny+ 1994. (Sometimes Daeppen+ 1987 is cited instead, but H&M seems to be where the theory originated. Presumably it was delayed in publication.)

The expression for w is in equation 4.71 of H&M. K, the QM correction used in defined in equation 4.24. Note that H&M's "N"s are numbers (not number densities), and their "V" is volume. These quantities apear only in the form N/V, so we use the number densities instead.

This is based partially on Paul Barklem and Kjell Eriksson's WCALC fortran routine (part of HBOP.f), which is used by (at least) Turbospectrum and SME. As in that routine, we do consider hydrogen and helium as the relevant neutral species, and assume them to be in the ground state. All ions are assumed to have charge 1. Unlike that routine, the generalization to the formalism from Hubeny+ 1994 is turned off by default because I haven't closely checked it. The difference effects the chargedterm only, and temperature is only used when `usehubeny_generalizationis set totrue`.

source
Korg.hydrogen_line_absorption!Method
hydrogen_line_absorption!(αs, λs, T, nₑ, nH_I, UH_I, ξ, window_size; kwargs...)

Calculate contribution to the the absorption coefficient, αs, from hydrogen lines in units of cm^-1, at wavelengths λs (a vector of ranges).

Uses profiles from Stehlé & Hutcheon (1999), which include Stark and Doppler broadening. For Halpha, Hbeta, and Hgamma, the p-d approximated profiles from Barklem, Piskunovet, and O'Mara 2000 are added to the absortion coefficient. This "convolution by summation" is inexact, but true convolution is expensive.

Arguments:

  • T: temperature [K]
  • nₑ: electron number density [cm^-3]
  • nH_I: neutral hydrogen number density [cm^-3]
  • UH_I: the value of the neutral hydrogen partition function
  • ξ: microturbulent velocity [cm/s]. This is only applied to Hα-Hγ. Other hydrogen lines profiles are dominated by stark broadening, and the stark broadened profiles are pre-convolved with a doppler profile.
  • window_size: the max distance from each line center [cm] at which to calculate the Stark and self-broadening profiles absorption for Hα-Hγ (those dominated by self-broadening).

Keyword arguments:

  • stark_profiles (default: Korg._hline_stark_profiles): tables from which to interpolate Stark profiles
  • use_MHD: whether or not to use the Mihalas-Daeppen-Hummer formalism to adjust the occupation probabilities of each hydrogen orbital for plasma effects. Default: true.
source
Korg.interpolate_marcsMethod
interpolate_marcs(Teff, logg, A_X; kwargs...)
interpolate_marcs(Teff, logg, m_H=0, alpha_m=0, C_m=0; kwargs...)

Returns a model atmosphere computed by interpolating models from MARCS ((Gustafsson+ 2008)[https://ui.adsabs.harvard.edu/abs/2008A&A...486..951G/abstract]). Along with Teff and logg, the atmosphere is specified by m_H, alpha_m, and C_m, which can be automatically determined from an A_X abundance vector (the recommended method, see format_A_X). Note that the MARCS atmosphere models were constructed with the Grevesse+ 2007 solar abundances (Korg.grevesse_2007_solar_abundances). This is handled automatically when A_X is provided.

interpolate_marcs uses three different interpolation schemes for different stellar parameter regimes. In the standard case the model atmosphere grid is the one generated for SDSS, transformed and linearly interpolated. For cool dwarfs (Teff ≤ 4000 K, logg ≥ 3.5), the grid is resampled onto unchanging tau_5000 values and interpolated with a cubic spline. For low-metallicity models (-5 ≤ m_H < -2.5), a grid of standard composition (i.e. fixed alpha and C) atmospheres is used. (The microturbulence is 1km/s for dwarfs and 2km/s for giants and the mass for spherical models is 1 solar mass.) The interpolation method is the same as in the standard case. See Wheeler+ 2024 for more details and a discussion of errors introduced by model atmosphere interpolation. (Note that the cubic scheme for cool dwarfs is referred to as not-yet-implemented in the paper but is now available.)

keyword arguments

  • spherical: whether or not to return a ShellAtmosphere (as opposed to a PlanarAtmosphere). By default true when logg < 3.5.
  • solar_abundances: (default: grevesse_2007_solar_abundances) The solar abundances to use when A_X is provided instead of M_H, alpha_M, and C_M. The default is chosen to match that of the atmosphere grid, and if you change it you are likely trying to do something else.
  • clamp_abundances: (default: false) allowed only when specifying A_X. Whether or not to clamp the abundance parameters to be within range to avoid throwing an out of bounds error.
  • perturb_at_grid_values (default: true): whether or not to add or a subtract a very small number to each parameter which is exactly at a grid value. This prevents null derivatives, which can cause problems for minimizers.
  • resampled_cubic_for_cool_dwarfs (default: true): whether or not to used specialized method for cool dwarfs.
  • archives: A tuple containing the atmosphere grids to use. For testing purposes.
source
Korg.inverse_gaussian_densityMethod
inverse_gaussian_density(ρ, σ)

Calculate the inverse of a (0-centered) Gaussian PDF with standard deviation σ, i.e. the value of x for which ρ = exp(-0.5 x^2/σ^2}) / √[2π], which is given by σ √[-2 log (√[2π]σρ)]. Returns 0 when ρ is larger than any value taken on by the PDF.

See also: inverse_lorentz_density.

source
Korg.inverse_lorentz_densityMethod
inverse_lorentz_density(ρ, γ)

Calculate the inverse of a (0-centered) Lorentz PDF with width γ, i.e. the value of x for which ρ = 1 / (π γ (1 + x^2/γ^2)), which is given by √[γ/(πρ) - γ^2]. Returns 0 when ρ is larger than any value taken on by the PDF.

See also: inverse_gaussian_density.

source
Korg.lazy_multilinear_interpolationMethod
lazy_multilinear_interpolation(params, nodes, grid; kwargs...)

This function is does multidimensional linear interpolation on a grid where the first two dimensions represent the values being interpolated. In other words, it interpolates matrices. It is written to minimize the number of reads from the grid, and access things in an efficient order. It is much faster than Interpolations.linear_interpolation when applied to memory-mapped grids. At present it is used only for model atmosphere and departure coefficient interpolation.

See also: interpolate_marcs, Korg.CubicSplines.CubicSpline

source
Korg.line_absorption!Method
line_absorption(linelist, λs, temp, nₑ, n_densities, partition_fns, ξ
               ; α_cntm=nothing, cutoff_threshold=1e-3, window_size=20.0*1e-8)

Calculate the opacity coefficient, α, in units of cm^-1 from all lines in linelist, at wavelengths λs [cm^-1].

other arguments:

  • temp the temerature in K (as a vector, for multiple layers, if you like)
  • n_densities, a Dict mapping species to absolute number density in cm^-3 (as a vector, if temp is a vector).
  • partition_fns, a Dict containing the partition function of each species
  • ξ is the microturbulent velocity in cm/s (n.b. NOT km/s)
  • α_cntm is as a callable returning the continuum opacity as a function of wavelength. The window within which a line is calculated will extend to the wavelength at which the Lorentz wings or Doppler core of the line are at cutoff_threshold * α_cntm[line.wl], whichever is greater.

Keyword Arguments

  • cuttoff_threshold (default: 3e-4): see α_cntm
  • verbose (default: false): if true, show a progress bar.
source
Korg.line_profileMethod
line_profile(λ₀, σ, γ, amplitude, λ)

A voigt profile centered on λ₀ with Doppler width σ (NOT √[2] σ, as the "Doppler width" is often defined) and Lorentz HWHM γ evaluated at λ (cm). Returns values in units of cm^-1.

source
Korg.load_atomic_partition_functionsFunction
load_atomic_partition_functions()

Loads saved tabulated values for atomic partition functions from disk. Returns a dictionary mapping species to interpolators over log(T).

source
Korg.load_exomol_partition_functionsMethod
load_exomol_partition_functions()

Loads the exomol partition functions for polyatomic molecules from the HDF5 archive. Returns a dictionary mapping species to interpolators over log(T).

source
Korg.merge_close_linesMethod
merge_close_lines(linelist; merge_distance=0.2)

Produce a list of the species and wavelengths of the lines in linelist, merging lines of the same species that are within merge_distance (default: 0.2 Å). This is useful for labeling lines in a plot after running prune_linelist.

Arguments

  • linelist: the linelist (a vector of Line objects)

Keyword Arguments

  • merge_distance=0.2: The maximum distance in Å between lines of the same species to be merged into a single entry

Returns

A vector of tuples (wl, wl_low, wl_high, species) where wl is gf-weighted wavelength of each set of merged lines (Å), wl_low and wl_high are their highest and lowest wavelength, and species is a string (not a Korg.Species) identifying the species of the line. These will be in wavelength order.

source
Korg.move_boundsMethod
move_bounds(λs, lb0, ub0, λ₀, window_size)

Using lb0 and ub0 as initial guesses, return the indices of λs, (lb, ub) corresponding to λ₀±window_size.λs` can either be a

  • Vector, in which case lb and ub are used as first guesses
  • Range, in which case lb and ub are ignored and the bounds are computed directly
  • Vector of Ranges, in which case lb and ub are also ignored. In this case, the bounds returned are indices into the concatenation of the ranges.

If the windows is outside the range of λs, the lb will be greater than ub.

Warning

Be careful with the values returned by this function. Indices into arrays and ranges work differently. For example, if lb = 1 and ub = 0 (as happens if the window is entirely below λs) then λs[lb:ub] will return an empty array if λs is an array and a non-empty array if λs is a range. Make sure to explicitly check that lb <= ub when appropriate.

source
Korg.n_atomsMethod
n_atoms(x)

The number of atoms in the Korg.Species or Korg.Formula x.

source
Korg.prune_linelistMethod
prune_linelist(atm, linelist, A_X, wls...; threshold=1.0, sort=true, synthesis_kwargs...)

Return the vector containing the strongest lines in linelist, (optionally) sorted by approximate equivalent width.

Arguments

  • atm: the atmosphere model
  • linelist: the linelist (a vector of Line objects)
  • A_X: the abundance of each element (see format_A_X)
  • wls...: the wavelength ranges to synthesize over. These are specified the same way as the wls for synthesize.

Keyword Arguments

  • threshold=0.1: The threshold ratio in the line center absorption to the continuum absorption computed at the photosphere for a line to be included in the returned list. 0.1 is a reasonable default for getting a sense of what might be measurable in a high-res, high-quality spectrum, but it should not be used to create a linelist for synthesis.
  • sort_by_EW=true: If true, the returned linelist will be sorted by approximate reduced equivalent width. If false, the linelist will be in wavelength order. Leaving the list in wavelength order is much faster, but sorting by strength is useful for visualizing the strongest lines.
  • verbose=true: If true, a progress bar will be displayed while measuring the EWs.

All other kwargs are passed to internal calls to synthesize.

Caution

While this function can be used to prune a linelist for synthesis, the default behavior too aggressive for this purpose. Set a much lower threshold (e.g. threshold=1e-4) and use sort=false if you are pruning the linelist to speedup synthesis. Note that Korg will dynamically choose which lines to include even if you use a large linelist (see the line_cutoff_threshold keyword argument to synthesize).

See also merge_close_lines if you are using this for plotting.

source
Korg.read_Barklem_Collet_tableMethod
function read_Barklem_Collet_table(fname; transform=identity)

Constructs a Dict holding tables containing partition function or equilibrium constant values across ln(temperature). Applies transform (which you can use to, e.g. change units) to each example.

source
Korg.read_linelistMethod
read_linelist(filename; format="vald", isotopic_abundances=Korg.isotopic_abundances)

Parse a linelist file, returning a vector of Lines.

The format keyword argument can be used to specify one of these linelist formats (default: "vald"):

  • "vald" for a VALD linelist. These can be either "short" or "long" format, "extract all" or "extract stellar". Air wavelengths will automatically be converted into vacuum wavelengths, and energy levels will be automatically converted from cm$^{-1}$ to eV.
  • "kurucz" for an atomic or molecular Kurucz linelist (format=kurucz_vac if it uses vacuum wavelengths; be warned that Korg will not assume that wavelengths are vacuum below 2000 Å),
  • "moog" for a MOOG linelist (doesn't support broadening parameters or dissociation energies).
  • "turbospectrum" for a Turbospectrum linelist in air wavelengths. Note that Korg doesn't make use of the (optional) orbital angular momentum quantum number, l, for the upper or lower levels, so it won't fall back on generic ABO recipes when the ABO parameters are not available. Korg's interpretation of the fdamp parameter is also slightly different from Turbospectrum's. See the documentation of the vdW parameter of Line for details. Korg will error if encounters an Unsoeld fudge factor, which it does not support.
  • "turbospectrum_vac" for a Turbospectrum linelist in vacuum wavelengths.

For VALD and Turbospectrum linelists with isotope information available, Korg will scale log gf values by isotopic abundance (unless VALD has already pre-scaled them), using isotopic abundances from NIST ([Korg.isotopicabundances]). To use custom isotopic abundances, just pass `isotopicabundances` with the same structure: a dict mapping atomic number to a dict mapping from atomic weight to abundance.

Be warned that for linelists which are pre-scaled for isotopic abundance, the estimation of radiative broadening from log(gf) is not accurate.

source
Korg.read_model_atmosphereMethod
read_model_atmosphere(filename)

Parse the provided model atmosphere file in MARCS ".mod" format. Returns either a PlanarAtmosphere or a ShellAtmosphere.

source
Korg.rectifyMethod
rectify(flux, wls; bandwidth=50, q=0.95, wl_step=1.0)

Rectify the spectrum with flux vector flux and wavelengths wls by dividing out a moving q-quantile with window size bandwidth. wl_step controls the size of the grid that the moving quantile is calculated on and interpolated from. Setting wl_step to 0 results in the exact calculation with no interpolation, but note that this is very slow when the wls is sampled for synthesis (~0.01 Å).

Experiments on real spectra show an agreement between the interpolated rectified spectrum and the "exact" one (with default values) at the 3 × 10^-4 level.

Warning

This function should not be applied to data with observational error, as taking a quantile will bias the rectification relative to the noiseless case. It is intended as a fast way to compute nice-looking rectified theoretical spectra.

source
Korg.saha_ion_weightsMethod
saha_ion_weights(T, nₑ, atom, ionization_energies, partition_functions)

Returns (wII, wIII), where wII is the ratio of singly ionized to neutral atoms of a given element, and wIII is the ration of doubly ionized to neutral atoms.

arguments:

  • temperature T [K]
  • electron number density nₑ [cm^-3]
  • atom, the atomic number of the element
  • ionization_energies is a collection indexed by integers (e.g. a Vector) mapping elements' atomic numbers to their first three ionization energies
  • partition_funcs is a Dict mapping species to their partition functions
source
Korg.scaled_vdWFunction
scaled_vdW(vdW, m, T)

The vdW broadening gamma scaled acording to its temperature dependence, using either simple scaling or ABO. See Anstee & O'Mara (1995) or https://www.astro.uu.se/~barklem/howto.html for the definition of the ABO γ.

vdW should be either γ_vdW evaluated at 10,000 K, or tuple containing the ABO params (σ, α). The species mass, m, is ignored in the former case.

source
Korg.setup_ionization_energiesFunction
setup_ionization_energies([filename])

Parses the table of ionization energies and returns it as a dictionary mapping elements to their ionization energies, [χ₁, χ₂, χ₃] in eV.

source
Korg.setup_partition_funcs_and_equilibrium_constantsMethod
setup_partition_funcs_and_equilibrium_constants()

Returns two dictionaries. One holding the default partition functions, and one holding the default log10 equilibrium constants.

Default partition functions

The partition functions are custom (calculated from NIST levels) for atoms, from Barklem & Collet 2016 for diatomic molecules, and from exomol for polyatomic molecules. For each molecule, we include only the most abundant isotopologue.

Note than none of these partition functions include plasma effects, e.g. via the Mihalas Hummer Daeppen occupation probability formalism. They are for isolated species. This can lead to a couple percent error for neutral alkalis and to greater errors for hydrogen in some atmospheres, particularly those of hot stars.

Default equilibrium constants

Molecules have equilibrium constants in addition to partition functions. For the diatomics, these are provided by Barklem and Collet, which extensively discusses the dissociation energies. For polyatomics, we calculate these ourselves, using atomization energies calculated from the enthalpies of formation at 0K from NIST's CCCDB.

Korg's equilibrium constants are in terms of partial pressures, since that's what Barklem and Collet provide.

source
Korg.sigma_lineMethod
sigma_line(wl)

The cross-section (divided by gf) at wavelength wl in Ångstroms of a transition for which the product of the degeneracy and oscillator strength is 10^log_gf.

source
Korg.synthesizeMethod
synthesize(atm, linelist, A_X, λ_start, λ_stop, [λ_step=0.01]; kwargs... )
synthesize(atm, linelist, A_X, wavelength_ranges; kwargs... )

Compute a synthetic spectrum.

Arguments

  • atm: the model atmosphere (see read_model_atmosphere)
  • linelist: A vector of Lines (see read_linelist, get_APOGEE_DR17_linelist, and get_VALD_solar_linelist).
  • A_X: a vector containing the A(X) abundances (log(X/H) + 12) for elements from hydrogen to uranium. (see format_A_X)
  • λ_start: the lower bound (in Å) of the region you wish to synthesize.
  • λ_stop: the upper bound (in Å) of the region you wish to synthesize.
  • λ_step (default: 0.01): the (approximate) step size to take (in Å).

If you provide a vector of wavelength ranges (or a single range) in place of λ_start and λ_stop, the spectrum will be synthesized over each range with minimal overhead. The ranges can be any Julia AbstractRange, for example: [5000:0.01:5010, 6000:0.03:6005]. they must be sorted and non-overlapping.

Returns

A named tuple with keys:

  • flux: the output spectrum
  • cntm: the continuum at each wavelength
  • alpha: the linear absorption coefficient at each wavelength and atmospheric layer a Matrix of size (layers x wavelengths)
  • number_densities: A dictionary mapping Species to vectors of number densities at each atmospheric layer
  • electron_number_density: the electron number density at each atmospheric layer
  • wavelengths: The vacuum wavelengths (in Å) over which the synthesis was performed. If air_wavelengths=true this will not be the same as the input wavelengths.
  • subspectra: A vector of ranges which can be used to index into flux to extract the spectrum for each range provided in wavelength_ranges. If you use the standard λ_start, λ_stop, λ_step arguments, this will be a vector containing only one range.

Example

to synthesize a spectrum between 5000 Å and 5100 Å, with all metal abundances set to 0.5 dex less than the solar value except carbon, except carbon, which we set to [C/H]=-0.25:

atm = read_model_atmosphere("path/to/atmosphere.mod")
linelist = read_linelist("path/to/linelist.vald")
A_X = format_A_X(-0.5, Dict("C" => -0.25))
solution = synthesize(atm, linelist, A_X, 5000, 5100)

Optional arguments:

  • vmic (default: 0) is the microturbulent velocity, $\xi$, in km/s.
  • air_wavelengths (default: false): Whether or not the input wavelengths are air wavelengths to be converted to vacuum wavelengths by Korg. The conversion will not be exact, so that the wavelength range can internally be represented by an evenly-spaced range. If the approximation error is greater than wavelength_conversion_warn_threshold, an error will be thrown. (To do wavelength conversions yourself, see air_to_vacuum and vacuum_to_air.)
  • wavelength_conversion_warn_threshold (default: 1e-4): see air_wavelengths. (In Å.)
  • line_buffer (default: 10): the farthest (in Å) any line can be from the provided wavelength range before it is discarded. If the edge of your window is near a strong line, you may have to turn this up.
  • cntm_step (default 1): the distance (in Å) between point at which the continuum opacity is calculated.
  • hydrogen_lines (default: true): whether or not to include H lines in the synthesis.
  • use_MHD_for_hydrogen_lines (default: true): whether or not to use the MHD occupation probability formalism for hydrogen lines. (MHD is always used for hydrogen bound-free absorption.)
  • hydrogen_line_window_size (default: 150): the maximum distance (in Å) from each hydrogen line center at which to calculate its contribution to the total absorption coefficient.
  • n_mu_points (default: 20): the number of μ values at which to calculate the surface flux when doing transfer in spherical geometry (when atm is a ShellAtmosphere). 20 points is sufficient for accuracy at the 10^-3 level.
  • line_cutoff_threshold (default: 3e-4): the fraction of the continuum absorption coefficient at which line profiles are truncated. This has major performance impacts, since line absorption calculations dominate more syntheses. Turn it down for more precision at the expense of runtime. The default value should effect final spectra below the 10^-3 level.
  • electron_number_density_warn_threshold (default: 1.0): if the relative difference between the calculated electron number density and the input electron number density is greater than this value, a warning is printed. Set to Inf to suppress this warning.
  • return_cntm (default: true): whether or not to return the continuum at each wavelength. If this is false, solution.cntm will be nothing.
  • ionization_energies, a Dict mapping Species to their first three ionization energies, defaults to Korg.ionization_energies.
  • partition_funcs, a Dict mapping Species to partition functions (in terms of ln(T)). Defaults to data from Barklem & Collet 2016, Korg.default_partition_funcs.
  • equilibrium_constants, a Dict mapping Species representing diatomic molecules to the base-10 log of their molecular equilibrium constants in partial pressure form. Defaults to data from Barklem and Collet 2016, Korg.default_log_equilibrium_constants.
  • bezier_radiative_transfer (default: false): Use the radiative transfer scheme. This is for testing purposes only.
  • molecular_cross_sections (default: []): A vector of precomputed molecular cross-sections. See MolecularCrossSection for how to generate these.
  • use_chemical_equilibrium_from (default: nothing): Takes another solution returned by synthesize. When provided, the chemical equilibrium solution will be taken from this object, rather than being recomputed. This is physically self-consistent only when the abundances, A_X, and model atmosphere, atm, are unchanged.
  • verbose (default: false): Whether or not to print information about progress, etc.
source
Korg.translational_UMethod
translational_U(m, T)

The (possibly inverse) contribution to the partition function from the free movement of a particle. Used in the Saha equation.

arguments

  • m is the particle mass
  • T is the temperature in K
source
Korg.vacuum_to_airMethod
vacuum_to_air(λ; cgs=λ<1)

convert λ from a vacuum to air. λ is assumed to be in Å if it is ⩾ 1, in cm otherwise. Formula from Birch and Downs (1994) via the VALD website.

See also: air_to_vacuum.

source
Korg.voigt_hjertingMethod
voigt_hjerting(α, v)

The Hjerting function, $H$, somtimes called the Voigt-Hjerting function. $H$ is defined as H(α, v) = ∫^∞_∞ exp(-y^2) / ((u-y)^2 + α^2) dy (see e.g. the unnumbered equation after Gray equation 11.47). It is equal to the ratio of the absorption coefficient to the value of the absorption coefficient obtained at the line center with only Doppler broadening.

If x = λ-λ₀, Δλ_D = σ√2 is the Doppler width, and Δλ_L = 4πγ is the Lorentz width,

voigt(x|Δλ_D, Δλ_L) = H(Δλ_L/(4πΔλ_D), x/Δλ_D) / (Δλ_D√π)
                    = H(γ/(σ√2), x/(σ√2)) / (σ√(2π))

Approximation from Hunger 1965.

source
Korg.λ_to_ν_boundMethod
λ_to_ν_bound(λ_bound)

Converts a λ Inverval (in cm) to an equivalent ν Interval (in Hz), correctly accounting for tricky floating point details at the bounds.

source
Korg.FitModule

Functions for fitting to data.

Warning

This submodule is in beta. It's API may change.

source
Korg.Fit.calculate_multilocal_masks_and_rangesMethod
calculate_multilocal_masks_and_ranges(obs_bounds_inds, obs_wls, synthesis_wls)

Given a vector of target synthesis ranges in the observed spectrum, return the masks, etc required.

Arguments: - windows: a vector of pairs of wavelength lower and upper bounds. - obs_wls: the wavelengths of the observed spectrum - synthesis_wls: the wavelengths of the synthesis spectrum - wl_buffer: the number of Å to add to each side of the synthesis range

Returns: - obs_wl_mask: a bitmask for obs_wls which selects the observed wavelengths - synthesis_wl_mask: a bitmask for synthesis_wls which selects the synthesis wavelengths needed to generated the masked observed spectrum. - multi_synth_wls: The vector of ranges to pass to Korg.synthesize.

source
Korg.Fit.ews_to_abundancesMethod
ews_to_abundances(atm, linelist, A_X, measured_EWs; kwargs... )

Compute per-line abundances on the linear part of the curve of growth given a model atmosphere and a list of lines with equivalent widths.

Arguments:

  • atm: the model atmosphere (see Korg.read_model_atmosphere and Korg.interpolate_marcs).
  • linelist: A vector of Korg.Lines (see Korg.read_linelist). The lines must be sorted by wavelength.
  • A_X: a vector containing the A(X) abundances (log(nX/nH) + 12) for elements from hydrogen to uranium (see Korg.format_A_X). All syntheses are done with these abundances, so if the resulting abundances deviate significantly from these, you may wish to iterate.
  • measured_EWs: a vector of equivalent widths (in mÅ)

Returns

A vector of abundances (A(X) = log10(n_X/n_H) + 12 format) for each line in linelist.

Optional arguments:

  • wl_step (default: 0.01) is the resolution in Å at which to synthesize the spectrum around each line.
  • ew_window_size (default: 2): the farthest (in Å) to consider equivalent width contributions for each line. It's very important that this is large enough to include each line entirely.
  • blend_warn_threshold (default: 0.01) is the minimum absorption between two lines allowed before triggering a warning that they may be blended.

All other keyword arguments are passed to Korg.synthesize when synthesizing each line.

source
Korg.Fit.ews_to_stellar_parametersFunction
ews_to_stellar_parameters(linelist, measured_EWs, [measured_EW_err]; kwargs...)

Find stellar parameters from equivalent widths the "old fashioned" way. This function finds the values of $T_\mathrm{eff}$, $\log g$, $v_{mic}$, and [m/H] which satisfy the following conditions (using a Newton-Raphson solver):

  • The slope of the abundances of neutral lines with respect to lower excitation potential is zero.
  • The difference between the mean abundances of neutral and ionized lines is zero.
  • The slope of the abundances of neutral lines with respect to reduced equivalent width is zero.
  • The difference between the mean abundances of all lines and the model-atmosphere input [m/H] is zero.

Here the "slope" refers to the slope of a linear fit to the abundances of the lines in question.

Arguments:

  • linelist: A vector of Korg.Line objects (see Korg.read_linelist). The lines must be sorted by wavelength.
  • measured_EWs: a vector of equivalent widths (in mÅ).
  • measured_EW_err (optional): the uncertainty in measured_EWs. If not specified, all lines are assumed to have the same uncertainty. These uncertainties are used when evaluating the equations above, and are propagated to provide uncertainties in the resulting parameters.

Returns:

A tuple containing:

  • the best-fit parameters: [Teff, logg, vmic, [m/H]] as a vector
  • the statistical uncertainties in the parameters, propagated from the uncertainties in the equivalent widths. This is zero when EW uncertainties are not specified.
  • the systematic uncertainties in the parameters, estimated from the scatter the abundances computed from each line not accounted for by the EW uncertainties.

Keyword arguments:

  • Teff0 (default: 5000.0) is the starting guess for Teff
  • logg0 (default: 3.5) is the starting guess for logg
  • vmic0 (default: 1.0) is the starting guess for vmic. Note that this must be nonzero in order to avoid null derivatives. Very small values are fine.
  • m_H0 (default: 0.0) is the starting guess for [m/H]
  • tolerances (default: [1e-3, 1e-3, 1e-3, 1e-3]) is the tolerance for the residuals each equation listed above. The solver stops when all residuals are less than the corresponding tolerance.
  • max_step_sizes (default: [1000.0, 1.0, 0.3, 0.5]) is the maximum step size to take in each parameter direction. This is used to prevent the solver from taking too large of a step and missing the solution. Be particularly cautious with the vmic (third) parameter, as the unadjusted step size is often too large.
  • parameter_ranges (default: [(2800.0, 8000.0), (-0.5, 5.5), (1e-3, 10.0), (-2.5, 1.0)]) is the allowed range for each parameter. This is used to prevent the solver from wandering into unphysical parameter space, or outside the range of the MARCS grid supported by Korg.interpolate_marcs. The default ranges $T_\mathrm{eff}$, $\log g$, and [m/H] are the widest supported by the MARCS grid. Note that vmic must be nonzero in order to avoid null derivatives.
  • callback: is a function which is called at each step of the optimization. It is passed three arguments:
    • the current values of the parameters
    • the residuals of each equation being solved
    • the abundances of each line computed with the current parameters.
    You can pass a callback function, to e.g. make a plot of the residuals at each step.
  • max_iterations (default: 30) is the maximum number of iterations to allow before stopping the optimization.
source
Korg.Fit.fit_spectrumFunction
fit_spectrum(obs_wls, obs_flux, obs_err, linelist, initial_guesses, fixed_params; kwargs...)

Find the parameters and abundances that best match a rectified observed spectrum.

Arguments:

  • obs_wls: the wavelengths of the observed spectrum in Å
  • obs_flux: the rectified flux of the observed spectrum
  • obs_err: uncertainty in flux
  • linelist: a linelist to use for the synthesis
  • initial_guesses: a NamedTuple specifying initial guesses for the parameters to be fit. See "Specifying parameters" below.
  • fixed_params: a NamedTuple specifying parameters to be held fixed. See "Specifying parameters" below.

initial_guesses and fixed_params can also be specified as Dicts instead of NamedTuples, which is more convenient when calling Korg from python.

Specifying parameters

Parameters are specified as named tuples or dictionaries. Named tuples look like this: (Teff=5000, logg=4.5, m_H=0.0). Single-element named tuples require a semicolon: (; Teff=5000).

Required parameters

Teff and logg must be specified in either initial_guesses or fixed_params.

Optional Parameters

These can be specified in either initial_guesses or fixed_params, but if they are not default values are used.

  • m_H: the metallicity of the star, in dex. Default: 0.0
  • alpha_H: the alpha enhancement of the star, in dex. Default: m_H. Note that, because of the parameter range supported by Korg.interpolate_marcs, only values within ±1 of m_H are supported.
  • vmic: the microturbulence velocity, in km/s. Default: 1.0
  • vsini: the projected rotational velocity of the star, in km/s. Default: 0.0. See Korg.apply_rotation for details.
  • epsilon: the linear limb-darkening coefficient. Default: 0.6. Used for applying rotational broadening only. See Korg.apply_rotation for details.
  • Individual elements, e.g. Na, specify the solar-relative ([X/H]) abundance of that element.
  • cntm_offset: a constant offset to the continuum. Default: 0.0 (no correction).
  • cntm_slope: the slope of a linear correction applied to the continuum. Default: 0.0 (no correction). The linear correction will be zero at the central wavelength. While it's possible to fit for a continuum slope with a fixed offset, it is highly disrecommended.
Tip

If you are doing more than a few fits, you will save a lot of time by precomputing the LSF matrix and synthesis wavelengths. See the keyword arguments below for how to do that.

Keyword arguments

  • windows (optional) is a vector of wavelength pairs, each of which specifies a wavelength "window" to synthesize and contribute to the total χ². If not specified, the entire spectrum is used. Overlapping windows are automatically merged.
  • LSF_matrix (optional) is a matrix which maps the synthesized spectrum to the observed spectrum. If not specified, it is calculated using Korg.compute_LSF_matrix. Computing the LSF matrix can be expensive, so you may want to precompute it if you are fitting many spectra with the same LSF.
  • synthesis_wls: a superset of the wavelengths to synthesize, as a range. If not specified, wavelengths spanning the first and last windows are used. If you pass in a precomputed LSF matrix, you must make sure that the synthesis wavelengths match it.
  • wl_buffer is the number of Å to add to each side of the synthesis range for each window.
  • precision specifies the tolerance for the solver to accept a solution. The solver operates on transformed parameters, so precision doesn't translate straightforwardly to Teff, logg, etc, but the default value, 1e-4, provides a theoretical worst-case tolerance of about 0.15 K in Teff, 0.0002 in logg, 0.0001 in m_H, and 0.0004 in detailed abundances. In practice the precision achieved by the optimizer is about 10x bigger than this.

Any additional keyword arguments will be passed to Korg.synthesize when synthesizing the spectra for the fit.

Returns

A NamedTuple with the following fields:

  • best_fit_params: the best-fit parameters
  • best_fit_flux: the best-fit flux, with LSF applied, resampled, and rectified.
  • obs_wl_mask: a bitmask for obs_wls which selects the wavelengths used in the fit (i.e. those in the windows)
  • solver_result: the result object from Optim.jl
  • trace: a vector of NamedTuples, each of which contains the parameters at each step of the optimization. This is empty for single parameter fits, because the underlying solver doesn't supply it.
  • covariance: a pair (params, Σ) where params is vector of parameter name (providing an order), and Σ is an estimate of the covariance matrix of the parameters. It is the approximate inverse hessian of the log likelihood at the best-fit parameter calculated by the BGFS algorithm, and should be interpreted with caution. For single-parameter fits (which are done with Nelder-Mead), this is not provided.
Tip

This function takes a long time to compile the first time it is called. Compilation performance is significantly better on Julia 1.10 than previous versions, so if you are using an older version of Julia, you may want to upgrade.

source
Korg.Fit.merge_boundsMethod

Sort a vector of lower-bound, upper-bound pairs and merge overlapping ranges. Used by fitspectrum and ewstostellarparameters.

Returns a pair containing:

  • a vector of merged bounds
  • a vector of vectors of indices of the original bounds which were merged into each merged bound
source
Korg.Fit.synthetic_spectrumMethod

Synthesize a spectrum, returning the flux, with LSF applied, resampled, and rectified. This is an internal function oned by fitting routines. See Korg.synthesize to synthesize spectra as a Korg user.

source
Korg.Fit.unscaleMethod

Unscale each parameter so that it lives on the appropriate range instead of (-∞, ∞).

source
Korg.Fit.validate_paramsMethod

Validate fitting parameters, and insert default values when needed. Used by fit_spectrum.

these can be specified in either initialguesses or fixedparams, but if they are not, these values are inserted into fixed_params

source
Korg.CubicSplines.CubicSplineMethod
CubicSpline(xs, ys; extrapolate=false)

Construct a interpolant using xs and ys as the knot coordinates. Assumes xs is sorted. Apply this object as a function to interpolate at any x value in the domain. If extrapolate is false, x values outside [xs[1], xs[end]] throw errors, if extrapolate is true, the interpolant uses flat extrapolation, i.e. it returns the extreme value.

source
Korg.CubicSplines.cumulative_integral!Method
cumulative_integral!(out, A, t1, t2)

Given a curve described by the spine, A, Calculates the integral from t1 to t for all t = t1, t2, and all spline knots in between. So if t1 is A.t[1] and t2 is A.t[end], out should have the same length as A.t.

source
Korg.ContinuumAbsorption.H_I_bfMethod
H_I_bf(νs, T, nH, nHe, ne, invU_H; n_max_MHD=6, use_hubeny_generalization=false, 
       taper=false, use_MHD_for_Lyman=false)

The bound-free linear absorption coefficient contributed by all energy states of a neutral Hydrogen atom. Even though the Mihalas-Hummer-Daeppen (MHD) occupation probability formalism is not used in Korg when computing the hydrogen partition function, it is used here. That means that series limit jumps (e.g. the Balmer jump), are "rounded off", as photons with less than the "classical" ionization energy can ionize if the upper level is dissolved into the continuum.

Required Arguments

  • νs: sorted frequency vector in Hz
  • T: temperature in K
  • nH_I: the total number density of neutral Hydrogen (in cm⁻³)
  • nHe_I: the total number density of neutral Helium (in cm⁻³)
  • ne: the number density of electrons (in cm⁻³)
  • invU_H: The inverse of the neutral hydrogen partition function (neglecting contributions from the MHD formalism)

For n=1 through n=n_max_MHD (default: 6), the cross-sections are computed using Nahar 2021. These are modified using the MHD formalism to account for level dissolution. For larger n, the cross-sections are calculated with a simple analytic formula (see simple_hydrogen_bf_cross_section).

Because MHD level dissolution applied to the the Lyman series limit leads to inflated cross-sections in the visible, we don't use MHD for bf absorption from n=1. This can be overridden by setting use_MHD_for_Lyman=true, in which case you will also want to set taper=true, which the same tapering of the cross-section as HBOP to fix the problem.

The use_hubeny_generalization keyword argument enables the generalization of the MHD from Hubeny 1994. It is experimental and switched off by default.

source
Korg.ContinuumAbsorption._ndens_HminusFunction
_ndens_Hminus(nH_I_div_partition, ne, T, ion_energy = _H⁻_ion_energy)

Compute the number density of H⁻ (implements eqn 5.10 of Kurucz 1970). This is an application of the saha equation where the "ground state" is H⁻ and the "first ionization state" is H I. The partition function of H⁻ is 1 at all temperatures.

source
Korg.ContinuumAbsorption.bounds_checked_absorptionMethod
bounds_checked_absorption(func; ν_bound, temp_bound)

Constructs a wrapped function that implements bounds checking and extrapolation

Parameters

  • func: a function that has a signature f(ν::Real, T::Real, args...)::Real, where ν is frequency (in Hz) and T is temperature (in K)
  • ν_bound::Interval: Interval of frequencies (in Hz) over which func is valid.
  • temp_bound::Interval: Interval of temperatures (in K) over which func is valid.

The resulting function will have this signature:

wrapped_func(ν::AbstractVector{<:Real}, T::Real, args...; kwargs...)

Wrapped Function Parameters

  • ν::AbstractVector{<:Real}: sorted vector of frequencies (in Hz)
  • T::Real: temperature (in K)
  • args...: function-specific arguments

For a description for the kwargs..., see Continuum Absorption Kwargs.

source
Korg.ContinuumAbsorption.electron_scatteringMethod
electron_scattering(nₑ)

Compute the linear absorption coefficient, α, from scattering off of free electrons. This has no wavelength dependence. It assumes isotropic scattering. (See, e.g. Gray p 160.)

Arguments

  • nₑ::F: number density of free electrons (in cgs)
source
Korg.ContinuumAbsorption.gaunt_ff_vanHoofMethod
gaunt_ff_vanHoof(log_u, log_γ2)

computes the thermally averaged, non-relativistic free-free gaunt factor by interpolating the table provided by van Hoof et al. (2014).

Arguments

  • log_u: Equal to log₁₀(u) = log₁₀(hν/(kTₑ))
  • log_γ2: Equal to log₁₀(γ²) = log₁₀(RydbergZ²/(kTₑ))

Rydberg is the "infinite mass unit of energy" and Tₑ is the temperature of free electrons (for our purposes, we assume that free electrons are in thermal equilibrium with ions and neutral species).

Notes

van Hoof et al. (2014) computed the associated data table with a non-relativistic approach, which is invalid at very high temperatures. They conclude (from comparisons with a different paper) that their "results should be accurate up to electron temperatures of roughly 100 MK". This is more than adequate for stellar atmospheres. In van Hoof et al. (2015), they find that relativistic effects introduce a ∼0.75% at 100MK, for Z = 1 (when Z > 1, the change is smaller).

This function currently uses linear interpolation. However, van Hoof et al. (2014) provides an implementation of a third-order Lagrange scheme, which "reaches a relative precision better than 1.5e-4 everywhere." The C and Fortran implementations of this scheme can be found here, and are copyrighted by a BSD-style license.

Earlier variants of this function used less-accurate data from section 5.1 of Kurucz (1970) that extended over a smaller interval of data. That table was originally derived from a figure in Karsas and Latter (1961) and it's now used for testing purposes.

source
Korg.ContinuumAbsorption.hydrogenic_ff_absorptionMethod
hydrogenic_ff_absorption(ν, T, Z, ni, ne)

computes the free-free linear absorption coefficient for a hydrogenic species

The naming convention for free-free absorption is counter-intuitive. A free-free interaction is named as though the species interacting with the free electron had one more bound electron (in other words it's named as though the free-electron and ion were bound together). In practice, this means that ni should refer to:

  • the number density of H II if computing the H I free-free absorption
  • the number density of He III if computing the He II free-free absorption
  • the number density of Li IV if computing the Li III free-free absorption

Arguments

  • Z::Integer: the charge of the ion. For example, this is 1 for ionized H.
  • ni: the number density of the ion species in cm⁻³.
  • ne: the number density of free electrons.
  • ν: frequency in Hz
  • T: temperature in K

Note

This approach was adopted from equation 5.8 from section 5.1 of Kurucz (1970). Comparison against equation 5.18b of Rybicki & Lightman (2004), reveals that the equation in Kurucz (1970) omits the dependence on ρ. According to Rybicki & Lightman (2004) the free-free absorption coefficient (corrected for stimulated emission) is:

    α = coef * Z² * ne * ni * (1 - exp(-hplanck*ν/(kboltz*T))) * g_ff / (sqrt(T) * ν³)

Note that the g_ff is the free-free gaunt factor and coef is ∼3.7e8 (a more exact coefficient can be computed from eqn 5.18a).

With this in mind, equation 5.8 of Kurucz (1970) should actually read

    κ = ne * n(H II) * F_ν(T) * (1 - exp(-hplanck*ν/(kboltz*T))) / ρ

where Fν(T) = coef * Z² * gff / (sqrt(T) * ν³).

See gaunt_ff_vanHoof for details about where our gaunt factor data comes from. For simplicity, we enforce temperature and ν bounds constraints that don't include all of the available gaunt factor data. In practice, this should never be a concern for stellar spectroscopy.

source
Korg.ContinuumAbsorption.metal_bf_absorption!Method
metal_bf_absorption!(α, νs, T, number_densities)

Adds to α the contributions of bf metal opacities. Uses precomputed tables from TOPBase for Li, Be, B, C, N, O, F, Ne, Na, Mg, Al, Si, S, Ar, and Ca. Uses tables from NORAD for Fe. For these elements, tables have been precomputed for the neutral and singly ionized species assuming and LTE distribution of energy levels. See Korg/data/metal_bf_cross-sections/ for the scripts which generate the tables.

Cross sections were computed for 100 K < T < 100,000 K and frequencies corresponding to 500 Å < λ < 30,000 Å. Outside of either of those ranges, flat extrapolation is used (i.e. the extreme value is returned).

source
Korg.ContinuumAbsorption.positive_ion_ff_absorption!Method
positive_ion_ff_absorption!(α_out::Vector{Real}, ν::Real, T::Real, number_densities::Dict, 
                            ne::Real, departure_coefficients=Peach1970.departure_coefficients)

Computes the linear absorption coefficient (in cm⁻¹) for all free-free interactions involving positively charged atomic species. Uses the provided departure coefficients when they are available, and the uncorrected hydrogenic approximation when they are not.

Arguments

  • ν: frequency in Hz
  • T: temperature in K
  • number_densities is a Dict mapping each Species to its number density.
  • ne: the number density of free electrons.
  • departure_coefficients (optional, defaults to ContinuumAbsorption.Peach1970.departure_coefficients: a dictionary mapping species to the departure coefficients for the ff process it participates in (e.g. species"C II" maps to the C III ff departure coefficients–see note below). Departure coefficients should be callables taking temperature and (photon energy / RydbergH / Zeff^2).
Note

A free-free interaction is named as though the species interacting with the free electron had one more bound electron (in other words it's named as though the free-electron and ion were bound together). To make matters more confusing, some sources present the free-free linear absorption coefficient's formula in terms of the number density of the species in the interaction's name (by including the Saha equation in the formula).

To compute the absorption for a given free-free interaction, this function accesses elements from the dictionary arguments that are associated with the species that participates in the interaction. For example:

  • Si I ff absorption: uses the number density of Si II (number_densities[species"Si II"]) and checks departure_coefficients[species"Si II"] for departure coefficients.
  • Si II ff absorption: uses the number density of Si III (number_densities[species"Si III"]) and checks departure_coefficients[species"Si III"] for departure coefficients.
source
Korg.ContinuumAbsorption.rayleighMethod
rayleigh(λs, nH_I, nHe_II, nH2)

Absorption coefficient from Rayleigh scattering by neutral H, He, and H2. Formulations for H and He are via Colgan+ 2016. Formulation for H2 from Dalgarno and Williams 1962.

The Dalgarno and Williams H2 is applicable redward of 1300 Å. Since Rayleigh scattering breaks down when the particle size to wavelength ratio gets large, we that all frequencies passed to this function be equivalent to 1300 Å or greater.

The formulations for H is adapted from Lee 2005, which states that it is applicable redward of Lyman alpha. See Colgan 2016 for details on He.

source
Korg.ContinuumAbsorption.simple_hydrogen_bf_cross_sectionMethod
simple_hydrogen_bf_cross_section(n::Integer, ν::Real)

Calculate the H I bf cross section in megabarns using a very simple approximation. See, for example, Kurucz 1970 equation 5.5 (though see note below). This implementation is used to extrapolate the cross-section past the ionization energy of an unperturbed hydrogen atom, as is required to take level dissolution into account with the MHD formalism.

Equation 5.5 of Kurucz had a typo in it. In the numerator of the fraction that is multiplied by the entire polynomial, Z² should be Z⁴. This was discovered during comparisons with data from the Opacity Project, and can be confirmed by looking at eqn 5.6 of Kurucz (it uses Z⁴ instead of Z²) or by comparison against equation 10.54 of Rybicki & Lightman.

source
Korg.ContinuumAbsorption.total_continuum_absorptionMethod
total_continuum_absorption(νs, T, nₑ, number_densities, partition_funcs; error_oobounds)

The total continuum linear absoprtion coefficient, α, at many frequencies, ν.

Arguments

  • νs are frequencies in Hz
  • T is temperature in K
  • nₑ is the electron number density in cm^-3
  • number_densities is a Dict mapping each Species to its number density
  • partition_funcs is a Dict mapping each Species to its partition function (e.g. Korg.partition_funcs)
  • error_oobounds::Bool specifies the behavior of most continnum absorption sources when passed frequencies or temperature values that are out of bounds for their implementation. When false (the default), those absorption sources are ignored at those values. Otherwise, an error is thrown.
Note

For efficiency reasons, νs must be sorted. While this function technically supports any sorted AbstractVector, it is most effient when passed an AbstractRange.

source
Korg.ContinuumAbsorption.Peach1970.departure_coefficientsMethod
Peach1970.departure_coefficients()

This module contains interpolators of the tabulated ff departure coeffients from Peach+ 1970, which we use to correct the hydrogenic ff absorption coefficient for H I ff, C I ff, Si I ff, and Mg I ff. It contains a dictionary (returned by departure_coefficients()), which maps Species to interpolator objects. Crucially, the dictionary is indexed by the species which actually participates in the interaction, not the one after which the interaction is named.

Outside the regime in which Peach 1970 provides data, the interpolators return 0, falling back to the hydreogenic approximation.

The species for which we use corrections are the same species which get corrected in MARCS/Turbospectrum (see Table 1 of Gustafsson+ 2008). The choices seem are largely motivated by which species have departure terms at normal stellar atmosphere conditions and which species are most abundant in the sun. For C II ff, we include only the contribution from the ¹S parent term, even though (in contrast to other speices) information is available for the ³Pᵒ term as well.

The free-free absorption coefficient (including stimulated emission) is given by:

$\alpha_{ m ff} = \alpha_{ m hydrogenic, ff}(\nu, T, n_i, n_e; Z) (1 + D(T, \sigma))$,

where

  • $\alpha_{ m hydrogenic, ff}(\nu, T, n_i, n_e; Z)$ should include the correction for stimulated emission.
  • $n_i$ is the number density fo the ion species that participates in the interation, not the species the interaction is named after.
  • $n_e$ is the number density of free electrons.
  • $D(T, \sigma)$ is specified as the departure arg, and is expected to interpolate over

the tabulated values specified in Table III of Peach (1970).

  • σ denotes the energy of the photon in units of RydbergH*Zeff²

It might not be immediately obvious how the above equation relates to the equations presented in Peach (1970). Peach describes the calculation for $k_ u^F$, the free-free absorption coefficient (uncorrected for stimulated emission) per particle of the species that the interaction is named after. In other words, he computes:

$k_ u^F = \alpha_{ m ff}/n_{i-1} \left(1 - e^\frac{-h\nu}{k T}\right)^{-1}$,

where $n_{i-1}$ is the number density of the species that the interaction is named after. $k_ u^F$ can directly be computed, under LTE, from just $u$, $T$, and $n_{i-1}$ (the Saha Equation relates $\alpha_{ m ff}$'s dependence on $n_e$ and $n_i$ to $n_{i-1}$ and $T$. Gray (2005) follows a similar convention when describing free-free absorption.

Warning

The tabulated data in this module was taken from Peach 1970 using OCR software, and may contain mis-read values, although they produce reasonable behavior and there are no obvious problems.

source
Korg.RadiativeTransfer.MoogStyleTransfer.radiative_transferFunction
radiative_transfer(atm::ModelAtmosphere, α, S, α_ref, mu_grid)

Returns the astrophysical flux at each wavelength.

inputs:

  • atm: the model atmosphere.
  • α: a matrix (atmospheric layers × wavelengths) containing the absoprtion coefficient
  • S: the source fuction as a matrix of the same shape.
  • α_ref: the continuum absoprtion coefficient at the reference wavelength (5000 Å). Used to rescale the total absorption to match the model atmosphere. This value should be calculated by Korg.
  • mu_grid: (required if atm is a ShellAtmosphere) the values of μ at which to calculate the surface intensity, which is integrated to obtain the astrophysical flux.

Keyword arguments:

  • n_nu_points: the number of rays to use when integrating over I_surface(μ) to obtain the astrophysical flux. This doesn't do anything if atm is a PlanarAtmosphere.
source
Korg.RadiativeTransfer.MoogStyleTransfer.ray_transfer_integralMethod
ray_transfer_integral(τ, S)

Compute exactly the solution to the transfer integral obtained be linearly interpolating the source function, S across optical depths τ, without approximating the factor of exp(-τ).

This breaks the integral into the sum of integrals of the form $\int (m\tau + b) \exp(-\tau)$ d\tau$ , which is equal to $ -\exp(-\tau) (m*\tau + b + m)$.

source
Korg.RadiativeTransfer.MoogStyleTransfer.spherical_transferMethod
spherical_transfer(α, S, τ_ref, α_ref, radii, μ_surface_grid)

Perform radiative transfer along rays emerging at the μ values in μ_surface_grid in a spherically symmetric atmosphere resolved at radii radii [cm]. See radiative_transfer for an explantion of the arguments. Note that radii should be in decreasing order.

Returns (flux, intensity), where flux is the astrophysical flux, and intensity, a matrix of shape (mu values × wavelengths), is the surface intensity as a function of μ.

source
Korg.RadiativeTransfer.BezierTransferModule

Functions to solve the radiative transfer formal solution using the de la Cruz Rodríguez and Piskunov method.

We would like to migrate to this transfer implementation as the default eventually, but it doens't seem to perform much better than RadiativeTransfer.MoogStyleTransfer` in practice. The treatment of both the path-length (ds) and what happens at the midplane of the star are more careful in this implementation, but it's not yet as fast and it doesn't produce results any closer to the MARCS flux files (even when adjusting the absorption coeffs to match the MARCS numbers).

source
Korg.RadiativeTransfer.BezierTransfer.planar_transferMethod
planar_transfer(α, S, z, n_μ_points)

Perform radiative transfer in a planar atmosphere. See radiative_transfer for an explantion of the arguments.

Returns (flux, intensity), where flux is the astrophysical flux, and intensity, a matrix of shape (mu values × wavelengths), is the surface intensity as a function of μ.

source
Korg.RadiativeTransfer.BezierTransfer.radiative_transferMethod
radiative_transfer(atm::ModelAtmosphere, α, S, n_μ_points)

Returns the astrophysical flux at each wavelength.

inputs:

  • atm: the model atmosphere.
  • α: a matrix (atmospheric layers × wavelengths) containing the absoprtion coefficient
  • S: the source fuction as a matrix of the same shape. rescale the total absorption to match the model atmosphere. This value should be calculated by Korg.
  • n_μ_points: the number of quadrature points to use when integrating over I_surface(μ) to obtain the astrophysical flux.
source
Korg.RadiativeTransfer.BezierTransfer.spherical_transferMethod
spherical_transfer(α, S, τ_ref, α_ref, radii, μ_surface_grid)

Perform radiative transfer along rays emerging in a spherically symmetric atmosphere resolved at radii radii [cm]. See radiative_transfer for an explantion of the arguments. Note that radii should be in decreasing order.

Returns (flux, intensity), where flux is the astrophysical flux, and intensity, a matrix of shape (mu values × wavelengths), is the surface intensity as a function of μ.

source