This page documents the Korg API for users of the package. Low-level functions that only people digging into the internals of Korg will be interested in can be found in the Developer documentation.

Top-level functions

If you are synthesize a spectrum Korg, these are the functions you will call. These functions are exported, so if you do using Korg, you can call them unqualified (i.e. synthesize instead of Korg.synthesize).

Korg.synthesizeFunction
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.read_linelistFunction
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_atmosphereFunction
read_model_atmosphere(filename)

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

source
Korg.interpolate_marcsFunction
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.format_A_XFunction
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

Fitting

Korg.Fit.ews_to_abundancesFunction
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

Secondary functions

You don't need use these to synthesize spectra, but they might be relevant depending on what you are doing.

Postprocessing

These are used to transform observational or synthetic spectra.

Korg.apply_LSFFunction
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.compute_LSF_matrixFunction
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.air_to_vacuumFunction
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.vacuum_to_airFunction
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

Line absorption

These functions can be used to directly compute line opacities.

Korg.line_absorption!Function
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_profileFunction
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.hydrogen_line_absorption!Function
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.voigt_hjertingFunction
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

Continuum absorption

These function can be used to directly compute continuum absorption coefficients.

Continuum Absorption Kwargs

All bound-free and free-free absorption functions (other than absorption_coef_bf_TOPBase) support the following keyword arguments:

  • error_oobounds::Bool: specifies the function's behavior when it encounters invalid ν or T values. When false (the default), the function assumes that the absorption coefficient for these invalid ν or T values is 0.0. Otherwise, a DomainError is thrown when invalid ν or T values are encountered.
  • out_α::Union{Nothing,AbstractVector}: When this is nothing (the default case), the function will simply allocate a new vector to store the output continuum absorption coefficients. Alternatively, this can be a vector (with the same length as the function's ν argument). In this case, the function will directly add the computed continuum absorption to the elements of this vector, in-place (the vector is also returned by the function).
Korg.ContinuumAbsorption.total_continuum_absorptionFunction
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.H_I_bfFunction
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.Hminus_bfFunction
Hminus_bf(ν, T, nH_I_div_partition, ne; kwargs...)

Compute the H⁻ bound-free linear absorption coefficient α, $\alpha_\nu = \sigma_{bf}(H^-) n(H⁻) (1 - \exp \left( \frac{-h\nu}{k T}\right))$

Arguments

  • ν::AbstractVector{<:Real}: sorted frequency vector in Hz
  • T: temperature in K
  • nH_I_div_partition: the total number density of H I divided by its partition function.
  • ne: the electron number density

This uses cross-sections from McLaughlin 2017. For a description of the kwargs, see Continuum Absorption Kwargs.

Notes

This function assumes that n(H⁻) ≪ n(H I) + n(H II). The number density of n(H⁻) isn't precomputed as part of Korg's molecular equlibrium, it's computed here instead.

Note

The McLaughlin tabulated values (downloaded as a text file) contained a stray line out of monotonic order. A version of the data file with the line deleted is saved at data/McLaughlin2017Hminusbf.dat for archival purposes. (Korg doesn't read this file, it reads data/McLaughlin2017Hminusbf.h5.)

source
Korg.ContinuumAbsorption.Hminus_ffFunction
Hminus_ff(ν, T, nH_I_div_partition, ne; kwargs...)

Compute the H⁻ free-free linear absorption coefficient α

The naming scheme for free-free absorption is counter-inutitive. This actually refers to the reaction: photon + e⁻ + H I -> e⁻ + H I.

Arguments

  • ν::AbstractVector{<:Real}: sorted frequency vector in Hz
  • T: temperature in K
  • nH_I_div_partition::Flt: the total number density of H I divided by its partition function.
  • ne: the number density of free electrons.

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

Notes

This is based on Table 1, in Bell & Berrington (1987), which tabulates values for "the H⁻ absorption coefficient", K (including the correction for stimulated emission). This quantity is in units of cm^4/dyn, and must be multiplied by the electron partial pressure and the ground-state neutral hydrogen number density to obtain a linear absorption coefficent, α.

The stipulation that the hydrogen should be ground-state only is based on the beginning of Section 2 in Bell and Berrington (1987), or alternately, Section 5.3 from Kurucz (1970). When Gray (2005) refers to this, it implicitly assumes that n(H I, n = 1) ≈ n(H I). Note that n(H I, n = 1) = n(H I)*gₙ₌₁/U(T)*exp(-Eₙ₌₁/(k*T)) = n(H I) * 2/U(T)*exp(0) = n(H I) * 2/U(T). Since U(T) ≈ 2 up to fairly large temperatures, this is not unreasonable.

source
Korg.ContinuumAbsorption.H2plus_bf_and_ffFunction
H2plus_bf_and_ff(ν, T, nH_I_div_partition, n_HII; kwargs...)

Compute the combined H₂⁺ bound-free and free-free linear absorption coefficient α using the tables from Stancil 1994. Note that this refers to interactions between free-protons and neutral hydrogen, not electrons and doubly ionized H2, i.e. the b-f interaction is photodissociation, not photoionization.

Arguments

  • ν::AbstractVector{<:Real}: sorted frequency vector in Hz
  • T: temperature in K
  • nH_I: the total number density of H I divided by its partition function.
  • nH_II: the number density of H II (not of H₂⁺).

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

Notes

This computes the H₂⁺ number density from those of H I and H II, since ionized molecules are not included in the molecular equlibrium calculation. Stancil provides approximate equilibrium constants, K, for the molecule, but the Barklem and Collet values used elsewhere by Korg may be more reliable. Once ionized molecules are fully supported, those values should be used instead. While cross sections are tabulated down to only 3150 K, the cross sections could be linearly interpolated 1000 K or so lower, if reliable K values are availble.

Because n(H₂⁺) is computed on the fly, the "cross-sections" used (internally) by this function have units of cm^-5, since they must be multiplied by n(H I) and n(H II) to obtain absorption coefficients.

Stancil, Bates (1952), and Gray (2005) all use n(H I), but Kurucz (1970) notes in section 5.2 that the photodisociation of H₂⁺ produces a ground state H atom and that we should therefore use n(H I, n=1) instead. This difference causes a descrepancy of a fraction of a percent (at most) for T ≤ 1.2e4 K. Here we use n(H I).

source
Korg.ContinuumAbsorption.Heminus_ffFunction
Heminus_ff(ν, T, nHe_I_div_partition, ne; kwargs...)

Compute the He⁻ free-free opacity κ.

The naming scheme for free-free absorption is counter-inutitive. This actually refers to the reaction: photon + e⁻ + He I -> e⁻ + He I.

Arguments

  • ν::AbstractVector{<:Real}: sorted frequency vector in Hz
  • T: temperature in K
  • nHe_I_div_partition: the total number density of H I divided by its partition function.
  • ne: the number density of free electrons.

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

Notes

This uses the tabulated values from John (1994). The quantity K is the same used by Bell and Berrington (1987). See Hminus_ff for an explanation.

According to John (1994), improved calculations are unlikely to alter the tabulated data for λ > 1e4Å, "by more than about 2%." The errors introduced by the approximations for 5.06e3 Å ≤ λ ≤ 1e4 Å "are expected to be well below 10%."

source
Korg.ContinuumAbsorption.electron_scatteringFunction
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.rayleighFunction
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.positive_ion_ff_absorption!Function
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

Statistical mechanics

These functions can be used to calculate the number densities of all species in a given atmospheric layer, or other statmech calculations.

Korg.saha_ion_weightsFunction
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.chemical_equilibriumFunction
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

Misc

Korg.blackbodyFunction
blackbody(T, λ)

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

source
Korg.prune_linelistFunction
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.merge_close_linesFunction
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