StellarLimbDarkening

class exotic_ld.StellarLimbDarkening(M_H=None, Teff=None, logg=None, ld_model='mps1', ld_data_path='', interpolate_type='nearest', custom_wavelengths=None, custom_mus=None, custom_stellar_model=None, verbose=False)[source]

Bases: object

Stellar limb darkening class.

Compute the limb-darkening coefficients for a specified stellar model, instrument throughput, and limb-darkening law.

Parameters:
  • M_H (float) – Stellar metallicity [dex].

  • Teff (int) – Stellar effective temperature [kelvin].

  • logg (float) – Stellar log(g) [dex].

  • ld_model (string) – Choose between ‘kurucz’, ‘stagger’, ‘mps1’, ‘mps2’, or ‘custom’. kurucz are 1D stellar models, can be referenced as ‘1D’. stagger are 3D stellar models, can be referenced as ‘3D’. mps1 are the MPS-ATLAS set 1 models. mps2 are the MPS-ATLAS set 2 models. If custom, must also provide custom_wavelengths, custom_mus, and custom_stellar_model.

  • ld_data_path (string) – Path to ExoTiC-LD_data directory downloaded from Zenodo. These data include the stellar models and instrument throughputs.

  • interpolate_type (string) – Choose between ‘nearest’ and ‘trilinear’.

  • custom_wavelengths (numpy.ndarray, shape (n,)) – If ld_model=’custom’, pass the wavelengths of you stellar model in angstroms.

  • custom_mus (numpy.ndarray, shape (m,)) – If ld_model=’custom’, pass the mu values at which you stellar model is defined.

  • custom_stellar_model (numpy.ndarray, shape (n, m)) – If ld_model=’custom’, pass the specific intensity of your stellar for each wavelength and mu value. Note specific intensity must be in units of [n_photons / s / cm^2 / Angstrom / steradian].

  • verbose (boolean) – Print information during calculation. Default: False.

Examples

>>> from exotic_ld import StellarLimbDarkening
>>> sld = StellarLimbDarkening(
        M_H=0.1, Teff=6045, logg=4.2, ld_model='mps1',
        ld_data_path='path/to/ExoTiC-LD_data')
>>> c1, c2 = sld.compute_quadratic_ld_coeffs(
        wavelength_range=np.array([20000., 30000.]),
        mode='JWST_NIRSpec_prism')
compute_3_parameter_non_linear_ld_coeffs(wavelength_range, mode, custom_wavelengths=None, custom_throughput=None, mu_min=0.1, return_sigmas=False)[source]

Compute the three-parameter non-linear limb-darkening coefficients.

Parameters:
  • wavelength_range (array_like, (start, end)) – Wavelength range over which to compute the limb-darkening coefficients. Wavelengths must be given in angstroms and the values must fall within the supported range of the corresponding instrument mode.

  • mode (string) –

    Instrument mode that defines the throughput. Modes supported for Hubble:

    ’HST_STIS_G430L’, ‘HST_STIS_G750L’, ‘HST_WFC3_G280p1’, ‘HST_WFC3_G280n1’, ‘HST_WFC3_G102’, ‘HST_WFC3_G141’.

    Modes supported for JWST:

    ’JWST_NIRSpec_Prism’, ‘JWST_NIRSpec_G395H’, ‘JWST_NIRSpec_G395M’, ‘JWST_NIRSpec_G235H’, ‘JWST_NIRSpec_G235M’, ‘JWST_NIRSpec_G140H-f100’, ‘JWST_NIRSpec_G140M-f100’, ‘JWST_NIRSpec_G140H-f070’, ‘JWST_NIRSpec_G140M-f070’, ‘JWST_NIRISS_SOSSo1’, ‘JWST_NIRISS_SOSSo2’, ‘JWST_NIRCam_F322W2’, ‘JWST_NIRCam_F444’, ‘JWST_MIRI_LRS’.

    Modes for photometry:

    ’Spitzer_IRAC_Ch1’, ‘Spitzer_IRAC_Ch2’, ‘TESS’.

    Alternatively, use ‘custom’ mode. In this case the custom wavelength and custom throughput must also be specified.

  • custom_wavelengths (array_like, optional) – Wavelengths corresponding to custom_throughput [angstroms].

  • custom_throughput (array_like, optional) – Throughputs corresponding to custom_wavelengths.

  • mu_min (float) – Minimum value of mu to include in the fitting process.

  • return_sigmas (boolean) – Return the uncertainties, or standard deviations, of each fitted limb-darkening coefficient. Default: False.

Returns:

  • if return_sigmas == False

    (c1, c2, c3)tuple

    Limb-darkening coefficients for the three-parameter non-linear law.

  • else

    ((c1, c2, c3), (c1_sigma, c2_sigma, c3_sigma))tuple of tuples

    Limb-darkening coefficients for the three-parameter non-linear law and uncertainties on each coefficient.

compute_4_parameter_non_linear_ld_coeffs(wavelength_range, mode, custom_wavelengths=None, custom_throughput=None, mu_min=0.1, return_sigmas=False)[source]

Compute the four-parameter non-linear limb-darkening coefficients.

Parameters:
  • wavelength_range (array_like, (start, end)) – Wavelength range over which to compute the limb-darkening coefficients. Wavelengths must be given in angstroms and the values must fall within the supported range of the corresponding instrument mode.

  • mode (string) –

    Instrument mode that defines the throughput. Modes supported for Hubble:

    ’HST_STIS_G430L’, ‘HST_STIS_G750L’, ‘HST_WFC3_G280p1’, ‘HST_WFC3_G280n1’, ‘HST_WFC3_G102’, ‘HST_WFC3_G141’.

    Modes supported for JWST:

    ’JWST_NIRSpec_Prism’, ‘JWST_NIRSpec_G395H’, ‘JWST_NIRSpec_G395M’, ‘JWST_NIRSpec_G235H’, ‘JWST_NIRSpec_G235M’, ‘JWST_NIRSpec_G140H-f100’, ‘JWST_NIRSpec_G140M-f100’, ‘JWST_NIRSpec_G140H-f070’, ‘JWST_NIRSpec_G140M-f070’, ‘JWST_NIRISS_SOSSo1’, ‘JWST_NIRISS_SOSSo2’, ‘JWST_NIRCam_F322W2’, ‘JWST_NIRCam_F444’, ‘JWST_MIRI_LRS’.

    Modes for photometry:

    ’Spitzer_IRAC_Ch1’, ‘Spitzer_IRAC_Ch2’, ‘TESS’.

    Alternatively, use ‘custom’ mode. In this case the custom wavelength and custom throughput must also be specified.

  • custom_wavelengths (array_like, optional) – Wavelengths corresponding to custom_throughput [angstroms].

  • custom_throughput (array_like, optional) – Throughputs corresponding to custom_wavelengths.

  • mu_min (float) – Minimum value of mu to include in the fitting process.

  • return_sigmas (boolean) – Return the uncertainties, or standard deviations, of each fitted limb-darkening coefficient. Default: False.

Returns:

  • if return_sigmas == False

    (c1, c2, c3, c4)tuple

    Limb-darkening coefficients for the three-parameter non-linear law.

  • else

    ((c1, c2, c3, c4), (c1_sigma, c2_sigma, c3_sigma, c4_sigma)) :

    tuple of tuples Limb-darkening coefficients for the four-parameter non-linear law and uncertainties on each coefficient.

compute_linear_ld_coeffs(wavelength_range, mode, custom_wavelengths=None, custom_throughput=None, mu_min=0.1, return_sigmas=False)[source]

Compute the linear limb-darkening coefficients.

Parameters:
  • wavelength_range (array_like, (start, end)) – Wavelength range over which to compute the limb-darkening coefficients. Wavelengths must be given in angstroms and the values must fall within the supported range of the corresponding instrument mode.

  • mode (string) –

    Instrument mode that defines the throughput. Modes supported for Hubble:

    ’HST_STIS_G430L’, ‘HST_STIS_G750L’, ‘HST_WFC3_G280p1’, ‘HST_WFC3_G280n1’, ‘HST_WFC3_G102’, ‘HST_WFC3_G141’.

    Modes supported for JWST:

    ’JWST_NIRSpec_Prism’, ‘JWST_NIRSpec_G395H’, ‘JWST_NIRSpec_G395M’, ‘JWST_NIRSpec_G235H’, ‘JWST_NIRSpec_G235M’, ‘JWST_NIRSpec_G140H-f100’, ‘JWST_NIRSpec_G140M-f100’, ‘JWST_NIRSpec_G140H-f070’, ‘JWST_NIRSpec_G140M-f070’, ‘JWST_NIRISS_SOSSo1’, ‘JWST_NIRISS_SOSSo2’, ‘JWST_NIRCam_F322W2’, ‘JWST_NIRCam_F444’, ‘JWST_MIRI_LRS’.

    Modes for photometry:

    ’Spitzer_IRAC_Ch1’, ‘Spitzer_IRAC_Ch2’, ‘TESS’.

    Alternatively, use ‘custom’ mode. In this case the custom wavelength and custom throughput must also be specified.

  • custom_wavelengths (array_like, optional) – Wavelengths corresponding to custom_throughput [angstroms].

  • custom_throughput (array_like, optional) – Throughputs corresponding to custom_wavelengths.

  • mu_min (float) – Minimum value of mu to include in the fitting process.

  • return_sigmas (boolean) – Return the uncertainties, or standard deviations, of each fitted limb-darkening coefficient. Default: False.

Returns:

  • if return_sigmas == False

    (c1, )tuple

    Limb-darkening coefficients for the linear law.

  • else

    ((c1, ), (c1_sigma, ))tuple of tuples

    Limb-darkening coefficients for the linear law and uncertainties on each coefficient.

compute_quadratic_ld_coeffs(wavelength_range, mode, custom_wavelengths=None, custom_throughput=None, mu_min=0.1, return_sigmas=False)[source]

Compute the quadratic limb-darkening coefficients.

Parameters:
  • wavelength_range (array_like, (start, end)) – Wavelength range over which to compute the limb-darkening coefficients. Wavelengths must be given in angstroms and the values must fall within the supported range of the corresponding instrument mode.

  • mode (string) –

    Instrument mode that defines the throughput. Modes supported for Hubble:

    ’HST_STIS_G430L’, ‘HST_STIS_G750L’, ‘HST_WFC3_G280p1’, ‘HST_WFC3_G280n1’, ‘HST_WFC3_G102’, ‘HST_WFC3_G141’.

    Modes supported for JWST:

    ’JWST_NIRSpec_Prism’, ‘JWST_NIRSpec_G395H’, ‘JWST_NIRSpec_G395M’, ‘JWST_NIRSpec_G235H’, ‘JWST_NIRSpec_G235M’, ‘JWST_NIRSpec_G140H-f100’, ‘JWST_NIRSpec_G140M-f100’, ‘JWST_NIRSpec_G140H-f070’, ‘JWST_NIRSpec_G140M-f070’, ‘JWST_NIRISS_SOSSo1’, ‘JWST_NIRISS_SOSSo2’, ‘JWST_NIRCam_F322W2’, ‘JWST_NIRCam_F444’, ‘JWST_MIRI_LRS’.

    Modes for photometry:

    ’Spitzer_IRAC_Ch1’, ‘Spitzer_IRAC_Ch2’, ‘TESS’.

    Alternatively, use ‘custom’ mode. In this case the custom wavelength and custom throughput must also be specified.

  • custom_wavelengths (array_like, optional) – Wavelengths corresponding to custom_throughput [angstroms].

  • custom_throughput (array_like, optional) – Throughputs corresponding to custom_wavelengths.

  • mu_min (float) – Minimum value of mu to include in the fitting process.

  • return_sigmas (boolean) – Return the uncertainties, or standard deviations, of each fitted limb-darkening coefficient. Default: False.

Returns:

  • if return_sigmas == False

    (c1, c2)tuple

    Limb-darkening coefficients for the quadratic law.

  • else

    ((c1, c2), (c1_sigma, c2_sigma))tuple of tuples

    Limb-darkening coefficients for the quadratic law and uncertainties on each coefficient.

compute_squareroot_ld_coeffs(wavelength_range, mode, custom_wavelengths=None, custom_throughput=None, mu_min=0.1, return_sigmas=False)[source]

Compute the square root limb-darkening coefficients.

Parameters:
  • wavelength_range (array_like, (start, end)) – Wavelength range over which to compute the limb-darkening coefficients. Wavelengths must be given in angstroms and the values must fall within the supported range of the corresponding instrument mode.

  • mode (string) –

    Instrument mode that defines the throughput. Modes supported for Hubble:

    ’HST_STIS_G430L’, ‘HST_STIS_G750L’, ‘HST_WFC3_G280p1’, ‘HST_WFC3_G280n1’, ‘HST_WFC3_G102’, ‘HST_WFC3_G141’.

    Modes supported for JWST:

    ’JWST_NIRSpec_Prism’, ‘JWST_NIRSpec_G395H’, ‘JWST_NIRSpec_G395M’, ‘JWST_NIRSpec_G235H’, ‘JWST_NIRSpec_G235M’, ‘JWST_NIRSpec_G140H-f100’, ‘JWST_NIRSpec_G140M-f100’, ‘JWST_NIRSpec_G140H-f070’, ‘JWST_NIRSpec_G140M-f070’, ‘JWST_NIRISS_SOSSo1’, ‘JWST_NIRISS_SOSSo2’, ‘JWST_NIRCam_F322W2’, ‘JWST_NIRCam_F444’, ‘JWST_MIRI_LRS’.

    Modes for photometry:

    ’Spitzer_IRAC_Ch1’, ‘Spitzer_IRAC_Ch2’, ‘TESS’.

    Alternatively, use ‘custom’ mode. In this case the custom wavelength and custom throughput must also be specified.

  • custom_wavelengths (array_like, optional) – Wavelengths corresponding to custom_throughput [angstroms].

  • custom_throughput (array_like, optional) – Throughputs corresponding to custom_wavelengths.

  • mu_min (float) – Minimum value of mu to include in the fitting process.

  • return_sigmas (boolean) – Return the uncertainties, or standard deviations, of each fitted limb-darkening coefficient. Default: False.

Returns:

  • if return_sigmas == False

    (c1, c2)tuple

    Limb-darkening coefficients for the square root law.

  • else

    ((c1, c2), (c1_sigma, c2_sigma))tuple of tuples

    Limb-darkening coefficients for the square root law and uncertainties on each coefficient.