API Reference

The package splits into a high-level interface — the four classes you instantiate directly — and a computational core of module-level functions that implement the forward model, photometry, isochrone fit and helpers. The core is documented here so the model in Theoretical background can be traced straight to the code that evaluates it.

High-level interface

Fitter

class astroARIADNE.fitter.Fitter[source]

Bases: object

The Fitter class handles the fitting routines and parameter estimation.

Examples

The fitter isn’t instantiaded with any arguments, rather you instantiate a Fitter object, then you set up the configurations and finally you initialize the object by running the initialize method. >>> f = Fitter() >>> f.star = s # s must be a valid Star object. >>> f.initialize()

star

The Star object to fit.

Type:

Star

setup

Sampler configuration. The first element selects the engine ('dynesty'); the remaining elements are [nlive, dlogz, bound, sample, threads, dynamic, walks].

Type:

list

out_folder

Folder where results and plots are written. Defaults to the star name.

Type:

str

verbose

Program verbosity. Default True.

Type:

bool

norm

If True, fit a single normalisation constant instead of radius and distance; the radius is recovered afterwards from the Gaia parallax.

Type:

bool

grid

Atmosphere grid for a single-grid fit (e.g. 'phoenix', 'btsettl', 'ck04', 'kurucz', 'coelho', 'bosz', 'sphinx', 'tlusty').

Type:

str

bma

If True, run Bayesian Model Averaging across several grids. This loads every model interpolator and fits all of them, so it is slower.

Type:

bool

models

The grids to average over when bma is True.

Type:

list

sequential

If True, fit the BMA grids one after another. See also n_grid_jobs to fit grids concurrently. Default True.

Type:

bool

n_grid_jobs

Number of BMA grids fit concurrently, one core each. Default 1 (sequential). Each job loads its own grid, so size this to the machine’s core count and memory.

Type:

int

n_samples

Number of samples drawn from the BMA-averaged posterior. None uses the combined chain length.

Type:

int or None

walks

Number of rwalk MCMC steps per dynesty proposal. Default 25; also settable as setup[7].

Type:

int

av_law

Extinction law. Default 'fitzpatrick'.

Type:

str

estimate_logg

If True, treat log g as a derived quantity via the isochrone step. Default False.

Type:

bool

estimate_age

If True, run the MIST isochrone fit for age, mass and EEP after the SED fit. Default True.

Type:

bool

isochrone_dlogz

Evidence tolerance for the MIST isochrone age/mass fit (the dominant cost of a BMA run). Lower is more precise but slower. Default 0.01.

Type:

float

prior_setup

Manual prior specification consumed by create_priors_from_setup.

Type:

dict or None

experimental

Enable experimental code paths. Default False.

Type:

bool

property av_law

Select extinction law.

static bayesian_model_average(outputs, grids, norm, nsamples, c='white')[source]

Perform Bayesian Model Averaging.

Parameters:
  • outputs (array_like) – An array or list of output pickle files from the individual model fits.

  • grids (array_like) – An array or list of the model grids used to model.

  • norm (bool) – Flag indicating if the normalization factor was fit for.

  • nsamples (int) – The number of samples to sample from the averaged distribution.

  • c (str, optional) – Termcolor color.

property bma

Bayesian Model Averaging (BMA).

Set to True if BMA is wanted. This loads every model grid interpolator and fits an SED to all of them, so the runtime will be slower!

colors = ['red', 'green', 'blue', 'yellow', 'grey', 'magenta', 'cyan', 'white']
create_priors_from_setup()[source]

Create priors from the manual setup.

static dynesty_results(results)[source]

Extract posterior samples, global evidence and its error.

fit()[source]

Run fitting routine.

fit_bma()[source]

Perform the fit with different models and the average the output.

Only works with dynesty.

fit_dynesty(out_file=None)[source]

Run dynesty.

fit_multinest(out_file=None)[source]

Run MultiNest.

get_ndim()[source]

Calculate number of dimensions.

property grid

Model grid selected.

initialize()[source]

Initialize the fitter.

To be run only after every input is added. This function calculates the number of dimensions, runs the prior creation, creates output directory, initializes coordinators and sets up global variables.

static load_interpolator(model)[source]

Load a DFInterpolator.

property models

Models to be used in BMA.

static multinest_results(out_folder, ndim)[source]

Extract posterior samples, global evidence and its error.

property n_samples

Set number of samples for BMA.

property norm

Bool to decide if a normalization constant will be fitted.

Set as True to not fit for radius and distance and fit for a normalization constant. After the fit a radius is calculated using the Gaia parallax.

property out_folder

Output folder.

If none is provided the default will be the starname.

save(out_file, results=None)[source]

Save multinest/dynesty output and relevant information.

Saves a dictionary as a pickle file. The dictionary contains the following:

lnZ: The global evidence. lnZerr: The global evidence error. posterior_samples: A dictionary containing the samples of each

parameter (even if it’s fixed), the evidence, log likelihood, the prior, and the posterior for each set of sampled parameters.

fixed: An array with the fixed parameter values. coordinator : An array with the status of each parameter (1 for fixed

0 for free)

best_fit: The best fit is chosen to be the median of each sample.

It also includes the log likelihood of the best fit.

star: The Star object containing the information of the star (name,

magnitudes, fluxes, coordinates, etc)

engine: The fitting engine used (i.e. MultiNest or Dynesty)

Also creates a log file with the best fit parameters and 1 sigma error bars.

save_bma(avgd)[source]

Save BMA output and relevant information.

Saves a dictionary as a pickle file. The dictionary contains the following:

lnZ : The global evidences. posterior_samples : A dictionary containing the samples of each

parameter (even if it’s fixed), the evidence, log likelihood, the prior, and the posterior for each set of sampled parameters.

fixed : An array with the fixed parameter values. coordinator : An array with the status of each parameter (1 for fixed

0 for free)

best_fitThe best fit is chosen to be the median of each sample.

It also includes the log likelihood of the best fit.

starThe Star object containing the information of the star (name,

magnitudes, fluxes, coordinates, etc)

Also creates a log file with the best fit parameters, 1 sigma error bars and 3 sigma CIs.

property sequential

Set to True to make BMA sequentially instead of parallel.

property setup

Set up options.

show_priors(return_string=False)[source]

Display priors in statistical notation.

Parameters:

return_string (bool, optional) – If True, return the formatted string instead of printing. Default: False

Returns:

Formatted prior information if return_string=True, else None

Return type:

str or None

property star

Star to fit for.

to_dict()[source]

Return BMA results as a structured dict matching the ecosystem spec.

Returns a dictionary with keys posterior, observed_data, and sample_stats, each containing numpy arrays ready for conversion to arviz.InferenceData.

Raises:

RuntimeError – If the fitter has not been run yet (self.out is not populated).

to_netcdf(path)[source]

Export BMA results as arviz DataTree in netCDF4 format.

Produces an arviz.InferenceData-compatible .nc file following the ecosystem output spec.

Parameters:

path (str) – Output file path (e.g. 'ariadne_result.nc').

property verbose

Program verbosity. Default is True.

Star

class astroARIADNE.star.Star(starname, ra, dec, g_id=None, plx=None, plx_e=None, rad=None, rad_e=None, temp=None, temp_e=None, lum=None, lum_e=None, dist=None, dist_e=None, Av=None, Av_e=None, offline=False, mag_dict=None, verbose=True, ignore=None, dustmap='SFD')[source]

Bases: object

Object that holds stellar magnitudes and other relevant information.

Parameters:
  • starname (str) – The name of the object. If ra and dec aren’t provided nor is a list of magnitudes with associated uncertainties provided, the search for stellar magnitudes will be done using the object’s name instead.

  • ra (float) – RA coordinate of the object in degrees.

  • dec (float) – DEC coordinate of the object in degrees.

  • g_id (int, optional) – The Gaia DR2 identifier.

  • plx (float, optional) – The parallax of the star in case no internet connection is available or if no parallax can be found on Gaia DR2.

  • plx_e (float, optional) – The error on the parallax.

  • rad (float, optional) – The radius of the star in case no internet connection is available or if no radius can be found on Gaia DR2.

  • rad_e (float, optional) – The error on the stellar radius.

  • temp (float, optional) – The effective temperature of the star in case no internet connection is available or if no effective temperature can be found on Gaia DR2.

  • temp_e (float, optional) – The error on the effective temperature.

  • lum (float, optional) – The stellar luminosity in case no internet connection is available or if no luminosity can be found on Gaia DR2.

  • lum_e (float, optional) – The error on the stellar luminosity.

  • dist (float, optional) – The distance in parsec.

  • dist_e (float, optional) – The error on the distance.

  • mag_dict (dictionary, optional) – A dictionary with the filter names as keys (names must correspond to those in the filter_names attribute) and with a tuple containing the magnitude and error for that filter as the value. Provide in case no internet connection is available.

  • offline (bool) – If False it overrides the coordinate search entirely.

  • verbose (bool, optional) – Set to False to suppress printed outputs.

  • ignore (list, optional) – A list with the catalogs to ignore for whatever reason.

full_grid

The full grid of fluxes.

Type:

ndarray

teff

The effective temperature axis of the flux grid.

Type:

ndarray

logg

The gravity axis of the flux grid

Type:

ndarray

z

If fixed_z is False, then z is the metallicity axis of the flux grid. Otherwise z has the same value as fixed_z

Type:

ndarray, float

starname

The name of the object.

Type:

str

ra

RA coordinate of the object.

Type:

float

dec

DEC coordinate of the object.

Type:

float

wave

An array containing the wavelengths associated to the different filters retrieved.

Type:

ndarray

flux

An array containing the fluxes of the different retrieved magnitudes.

Type:

ndarray

See class docstring.

add_mag(mag, err, filt)[source]

Add an individual photometry point to the SED.

calculate_distance()[source]

Calculate distance using parallax in solar radii.

clip_outlier_magnitudes(sigma=5.0, err_floor=0.05, warn_only=False)[source]

Iterative SED outlier rejection using synthetic photometry.

Uses per-filter zero-point fluxes (Vega or AB as appropriate) so that magnitudes across different photometric systems are compared on a consistent scale.

Seeds the trusted set from Gaia anchors (BP/G/RP), then iteratively promotes the best-fitting untrusted filter if its blackbody residual is within sigma. Filters never promoted are removed from the SED.

Falls back to iterative worst-outlier removal when no Gaia filters are present.

Parameters:
  • sigma (float) – Rejection threshold in units of max(photometric_error, err_floor).

  • err_floor (float) – Minimum assumed error (mag) to avoid over-rejecting precise points.

Returns:

Names of removed filters.

Return type:

list[str]

colors = ['red', 'green', 'blue', 'yellow', 'grey', 'magenta', 'cyan', 'white']
dustmaps = {'Bayestar': dustmaps.bayestar.BayestarQuery, 'Lenz': dustmaps.lenz2017.Lenz2017Query, 'Planck13': dustmaps.planck.PlanckQuery, 'Planck16': dustmaps.planck.PlanckGNILCQuery, 'SFD': dustmaps.sfd.SFDQuery}
estimate_logg(out='.')[source]

Estimate logg values from MIST isochrones.

filter_names = array(['GALEX_FUV', 'GALEX_NUV', 'STROMGREN_u', 'SkyMapper_u', 'SDSS_u',        'GROUND_JOHNSON_U', 'SkyMapper_v', 'STROMGREN_v', 'TYCHO_B_MvB',        'GROUND_JOHNSON_B', 'STROMGREN_b', 'SDSS_g', 'PS1_g',        'SkyMapper_g', 'GaiaDR2v2_BP', 'TYCHO_V_MvB', 'STROMGREN_y',        'GROUND_JOHNSON_V', 'SkyMapper_r', 'SDSS_r', 'PS1_r', 'PS1_w',        'KEPLER_Kp', 'GaiaDR2v2_G', 'GROUND_COUSINS_R', 'NGTS_I', 'SDSS_i',        'PS1_i', 'SkyMapper_i', 'GaiaDR2v2_RP', 'GROUND_COUSINS_I', 'TESS',        'PS1_z', 'SDSS_z', 'SkyMapper_z', 'PS1_y', '2MASS_J', '2MASS_H',        '2MASS_Ks', 'WISE_RSR_W1', 'SPITZER_IRAC_36', 'SPITZER_IRAC_45',        'WISE_RSR_W2', 'WISE_RSR_W3', 'WISE_RSR_W4', 'HERSCHEL_PACS_BLUE',        'HERSCHEL_PACS_GREEN', 'HERSCHEL_PACS_RED'], dtype='<U19')
load_grid(model)[source]

Load the model grid for interpolation.

print_mags(c=None)[source]

Pretty print of magnitudes and errors.

ra_dec_to_deg(ra, dec)[source]

Transform ra, dec from selected unit to degrees.

remove_mag(filt)[source]

Remove an individual photometry point.

save_mags(out)[source]

Save the used magnitudes in a file.

SEDPlotter

class astroARIADNE.plotter.SEDPlotter(input_files, out_folder, pdf=False, model=None, settings=None, method='averaged', save_model=False, ir_excess=False)[source]

Bases: object

Artist class for all things SED.

Parameters:
  • input_files (str) – Directory containing the code’s output files.

  • out_folder (type) – Directory where to put the output plots.

  • pdf (type) – Set to True to output plots in pdf.

  • png (type) – Set to True to output plots in png.

  • model (type) –

    Set to override the SED model that’s going to be plotted. Possible values are:

    • Phoenix

    • BTSettl

    • NextGen

    • CK04 (Castelli & Kurucz 2004)

    • Kurucz (Kurucz 1993)

Examples

Examples should be written in doctest format, and should illustrate how to use the function/class. >>>

chain_out

Output directory for chain plot.

Type:

str

like_out

Output directory for likelihood plot.

Type:

str

post_out

Output directory for posteriors plot.

Type:

str

moddir

Directory wheere the SED models are located.

Type:

type

out

SED fitting routine output.

Type:

dict

engine

Selected fitting engine.

Type:

str

star

The fitted Star object.

Type:

Star

coordinator

Array coordinating fixed parameters.

Type:

array_like

fixed

Array coordinating fixed parameters.

Type:

array_like

norm

norm is set to True if a normalization constant is fitted instead of radius + distance.

Type:

bool

grid

Selected model grid.

Type:

str

av_law

Exticntion law chosen for the fit.

Type:

function

order

Array coordinating parameter order.

Type:

array_like

interpolator

Interpolator function.

Type:

function

theta

Best fit parameter vector

Type:

array_like

See class docstring.

SED(ax)[source]

Plot the SED model.

clean()[source]

Close opened figures.

fetch_Phoenix()[source]

Fetch correct Phoenixv2 SED file.

The directory containing the Phoenix spectra must be called PHOENIXv2 Within PHOENIXv2 there should be the wavelength file called WAVE_PHOENIX-ACES-AGSS-COND-2011.fits and several folders called Z[-/+]X.X where X.X are the metallicities (e.g. Z-0.0, Z+1.0, etc)

fetch_btcond()[source]

Fetch correct BT-COND SED file.

The directory containing the BT-COND spectra must be called BTCOND. Within BTCOND there should be yet another directory called CIFIST2011, within BTCOND/CIFIST2011 there should be the SED fits files with the following naming convention:

lteTTT-G.G[-/+]Z.Za+0.0..BT-Cond.CIFIST2011.fits

where TTT are the first 3 digits of the effective temperature if it’s a number over 10000, else it’s the first 2 digit prepended by a 0. G.G is the log g and Z.Z the metallicity.

fetch_btnextgen()[source]

Fetch correct BT-NextGen SED file.

The directory containing the BT-NextGen spectra must be called BTNextGen. Within BTNextGen there should be yet another directory called AGSS2009, within BTNextGen/AGSS2009 there should be the SED fits files with the following naming convention:

lteTTT-G.G[-/+]Z.Za+0.0..BT-NextGen.AGSS2009.fits

where TTT are the first 3 digits of the effective temperature if it’s a number over 10000, else it’s the first 2 digit prepended by a 0. G.G is the log g and Z.Z the metallicity.

fetch_btsettl()[source]

Fetch correct BT-Settl SED file.

The directory containing the BT-Settl spectra must be called BTSettl Within BTSettl there should be yet another directory called AGSS2009, within BTSettl/AGSS2009 there should be the SED fits files with the following naming convention:

lteTTT-G.G[-/+]Z.Za+0.0.BT-Settl.AGSS2009.fits

where TTT are the first 3 digits of the effective temperature if it’s a number over 10000, else it’s the first 2 digit prepended by a 0. G.G is the log g and Z.Z the metallicity.

fetch_ck04()[source]

Fetch correct Castelli-Kurucz 2004 SED file.

The directory containing the Castelli-Kurucz spectra must be called Castelli_Kurucz. Within Castelli_Kurucz there should be a group of directories called ck[pm]ZZ where ZZ is the metalicity without the dot. Within each directory there are fits files named:

ck[pm]ZZ_TTTT.fits

where ZZ is metalicity as previous and TTTT is the effective temperature.

fetch_coelho()[source]

Fetch correct Coelho 2014 SED file.

The directory containing the Coelho spectra must be called Coelho14. Within Coelho14 there should be a group of files called t[0X]XXXX_g[+-]Y.Y_[mp]ZZp0[14]_sed.fits where X is the temperature, Y is the logg and Z the metallicity.

fetch_kurucz()[source]

Fetch correct Kurucz 1993 SED file.

The directory containing the Kurucz spectra must be called Kurucz. Within Kurucz there should be a group of directories called k[pm]ZZ where ZZ is the metalicity without the dot. Within each directory there are fits files named:

k[pm]ZZ_TTTT.fits

where ZZ is metalicity as previous and TTTT is the effective temperature.

plot_SED()[source]

Create the plot of the SED.

plot_SED_no_model(s=None)[source]

Plot raw photometry.

plot_bma_HR(nsamp)[source]

Plot HR diagram for the star.

plot_bma_hist()[source]

Plot histograms.

plot_corner()[source]

Make corner plot.

plot_hist()[source]
plot_like()[source]

Plot Likelihoods.

plot_post()[source]

Plot posteriors.

plot_trace()[source]

Plot SED chains.

Librarian

The photometry/astrometry retrieval layer. Librarian resolves a target and gathers catalogue photometry; adapt_librarian bridges its result to the Star/Fitter parameter conventions.

class astroARIADNE.librarian.Librarian(ra: float, dec: float, *, gaia_id: int | None = None, radius: float = 3.0, ignore: Sequence[str] = (), verbose: bool = True, feh_source: str = 'hypatia', hypatia_statistic: str = 'median', hypatia_uncertainty: str = 'std', hypatia_solarnorm: str = 'asplund09', feh_floor: float = 0.1)[source]

Bases: object

Automated photometry and astrometry retrieval for isochrone fitting.

Uses Gaia DR3 best_neighbour tables for catalog crossmatching, with VizieR XMatch fallback when Gaia TAP is unavailable.

Parameters:
  • ra (float) – Right ascension and declination in degrees.

  • dec (float) – Right ascension and declination in degrees.

  • gaia_id (int, optional) – Gaia DR3 source_id (skips cone search if provided).

  • radius (float) – Search radius in arcseconds (default 3).

  • ignore (sequence of str) – Catalog names to skip (e.g. [“SDSS”, “APASS”]).

  • verbose (bool) – Print retrieval summary.

property age: tuple[float, float] | None

(age_Gyr, error_Gyr) from FLAME.

property distance: tuple[float, float] | None

(distance_pc, error_pc) from Bailer-Jones EDR3.

property gaia_id: int | None
property luminosity: tuple[float, float] | None

(L_Lsun, error_Lsun) from FLAME.

property magnitudes: dict[str, tuple[float, float]]

(mag, error)}.

Type:

{‘pyphot_filter_name’

property mass: tuple[float, float] | None

(M_Msun, error_Msun) from FLAME.

property parallax: tuple[float, float] | None

(parallax_mas, error_mas) with Lindegren+21 correction.

property radius: tuple[float, float] | None

(R_Rsun, error_Rsun) from FLAME (5x inflated error).

property rave_params: dict | None

RAVE DR6 spectroscopic params, or None. Backward-compatible.

property spectroscopic_params: dict | None

Spectroscopic params dict with source tag, or None.

Keys: teff, teff_err, logg, logg_err, feh, feh_err, source.

property teff: tuple[float, float] | None

(Teff_K, error_K) from Gaia GSP-Phot.

property tic_id: int | None
astroARIADNE.librarian._adapter.adapt_librarian(lib) AdaptedStar[source]

Convert a new-style Librarian into ARIADNE’s Star contract.

Parameters:

lib (object) – Anything exposing the new Librarian interface: magnitudes dict and properties parallax, distance, radius, teff, luminosity (each (value, error) or None), plus gaia_id, tic_id, _kic_id, rave_params, spectroscopic_params.

Returns:

Filled arrays + scalars + pass-through identifiers/params.

Return type:

AdaptedStar

Computational core

sed_library — model, likelihood and priors

The forward model and Gaussian likelihood evaluated on every nested-sampling proposal, plus the unit-cube prior transforms. See Theoretical background for the equations these implement.

sed_library contain the model, prior and likelihood to be used.

astroARIADNE.sed_library.build_params(theta, flux, flux_e, filts, coordinator, fixed, use_norm)[source]

Build the parameter vector that goes into the model.

astroARIADNE.sed_library.fast_loglik(res, ers)[source]
astroARIADNE.sed_library.get_interpolated_flux(temp, logg, z, filts, interpolators)[source]

Interpolate the grid of fluxes in a given teff, logg and z.

Parameters:
  • temp (float) – The effective temperature.

  • logg (float) – The superficial gravity.

  • z (float) – The metallicity.

  • filts (str) – The desired filter.

Returns:

flux – The interpolated flux at temp, logg, z for filter filt.

Return type:

float

astroARIADNE.sed_library.get_residuals(theta, flux, flux_er, wave, filts, interpolators, use_norm, av_law)[source]

Calculate residuals of the model.

astroARIADNE.sed_library.log_likelihood(theta, flux, flux_er, wave, filts, interpolators, use_norm, av_law)[source]

Calculate log likelihood of the model.

astroARIADNE.sed_library.model_grid(theta, filts, wave, interpolators, use_norm, av_law)[source]

Return the model grid in the selected filters.

Parameters:

thetaarray_like

The parameters of the fit: teff, logg, z, radius, distance

starStar

The Star object containing all relevant information regarding the star.

interpolatorsdict

A dictionary with the interpolated grid.

use_normbool

False for a full fit (including radius and distance). True to fit for a normalization constant instead.

returns:

model – A dictionary whose keys are the filters and the values are the interpolated fluxes

rtype:

dict

astroARIADNE.sed_library.prior_transform_dynesty(u, flux, flux_er, filts, prior_dict, coordinator, use_norm)[source]

Transform the prior from the unit cube to the parameter space.

astroARIADNE.sed_library.prior_transform_multinest(u, flux, flux_er, filts, prior_dict, coordinator, use_norm)[source]

Transform the prior from the unit cube to the parameter space.

phot_utils — photometric conversions

Magnitude/flux conversions, zero-points and bandpass lookups used to turn catalogue magnitudes into the flux densities the likelihood compares against.

phot_utils module for SED fitting.

This module contains useful functions in order to obtain fluxes from different broadband filters. It also has functions to convert to different units of flux

It uses the module pyphot to get the fluxes and bandpasses of different broadband filters.

astroARIADNE.phot_utils.convert_f_lambda_to_f_nu(f, l)[source]

Convert flux from erg s-1 cm-2 lambda-1 to erg s-1 cm-2 Hz-1.

astroARIADNE.phot_utils.convert_f_nu_to_f_lambda(f, l)[source]

Convert flux from erf s-1 cm-2 Hz-1 to erg s-1 cm-2 lambda-1.

astroARIADNE.phot_utils.convert_jansky_to_ergs(j)[source]

Convert flux from jansky to erg s-1 cm-2 Hz-1.

astroARIADNE.phot_utils.convert_jansky_to_ergs_lambda(j, l)[source]

Convert flux from jansky to erg s-2 cm-2 lambda-1 in the units of l.

astroARIADNE.phot_utils.extract_info(magnitudes, errors, filters)[source]

Extract the flux information for a Star.

astroARIADNE.phot_utils.flux_to_mag(flux, flux_err, band)[source]

Convert from flux to magnitude.

The flux is expected to be in the units of erg s-1 cm-2 um-1

astroARIADNE.phot_utils.get_bandpass(band)[source]

Get the bandpass of a specific filter in um.

astroARIADNE.phot_utils.get_effective_wavelength(band)[source]

Get central wavelength of a specific filter in um.

astroARIADNE.phot_utils.get_zero_flux(band)[source]

Look for the filter information in the pyphot library of filters.

astroARIADNE.phot_utils.mag_to_flux(mag, mag_err, band)[source]

Convert from magnitude to flux.

mag_to_flux performs the conversion from magnitude to flux in erg s-1 cm-2 um-1.

The band parameter is a string representing the filter used and it must match exactly the name in pyphots filter database

If the filter is from PanSTARRS or SDSS, then the magnitude is in the AB system. Else it’s in the Vega system.

astroARIADNE.phot_utils.mag_to_flux_AB(mag, mag_err)[source]

Calculate flux in erg s-1 cm-2 Hz-1.

isochrone — MIST age/mass estimation

The post-fit step that maps the sampled (Teff, logg, [Fe/H], radius) posterior onto MIST tracks for age, mass and EEP.

Estimate logg using MIST isochrones.

astroARIADNE.isochrone.estimate(bands, params, logg=True, out_folder='.', threads=1, dlogz=0.01)[source]

Estimate logg (or age/mass/eep) using MIST isochrones.

Parameters:
  • threads (int) – Number of processes for the nested-sampling pool. >1 parallelises the isochrone fit (the dominant cost of a BMA run).

  • dlogz (float) – Evidence tolerance for the isochrone nested sampling. Lower is more precise but slower.

astroARIADNE.isochrone.get_isochrone(logage, feh)[source]

Retrieve isochrone for given age and feh.

astroARIADNE.isochrone.loglike(theta)[source]
astroARIADNE.isochrone.prior_transform(u)[source]

utils — helpers

Credibility intervals, KDE summaries, and the formatted output routines.

Various utilities used throughout the code.

Here go various utilities that don’t belong directly in any class, photometry utils module nor or SED model module.

astroARIADNE.utils.create_dir(path)[source]

Create a directory.

astroARIADNE.utils.credibility_interval(post, alpha=1.0)[source]

Calculate bayesian credibility interval.

Parameters:

postarray_like

The posterior sample over which to calculate the bayesian credibility interval.

alphafloat, optional

Confidence level.

Returns:

medfloat

Median of the posterior.

lowfloat

Lower part of the credibility interval.

upfloat

Upper part of the credibility interval.

astroARIADNE.utils.display_routine(engine, live_points, dlogz, ndim, bound=None, sample=None, nthreads=None, dynamic=None)[source]

Display program information.

What is displayed is: Program name Program author Star selected Algorithm used (i.e. Multinest or Dynesty) Setup used (i.e. Live points, dlogz tolerance)

astroARIADNE.utils.display_star_fin(star, c)[source]

Display stellar information.

astroARIADNE.utils.display_star_init(star, c)[source]

Display stellar information.

astroARIADNE.utils.end(coordinator, elapsed_time, out_folder, engine, use_norm)[source]

Display end of run information.

astroARIADNE.utils.estimate_cdf(distribution, hdr=False)[source]

Estimate the CDF of a distribution.

astroARIADNE.utils.estimate_pdf(distribution)[source]

Estimates the PDF of a distribution using a gaussian KDE.

Parameters:

distribution (array_like) – The distribution.

Returns:

  • xx (array_like) – The x values of the PDF.

  • pdf (array_like) – The estimated PDF.

astroARIADNE.utils.execution_time(start)[source]

Calculate run execution time.

astroARIADNE.utils.get_noise_name(filt)[source]

Retrieve parameter name for white noise.

astroARIADNE.utils.norm_fit(x, mu, sigma, A)[source]

Gaussian function.

astroARIADNE.utils.out_filler(samp, logdat, param, name, out, fmt='f', fixed=False, method='averaged')[source]

Fill up the output file.