BASTA main functions

The main BASTA programme

The constants module

Collection of all constants used in BASTA

class constants.distanceranges[source]

Limits or ranges of different surveys

class constants.extinction[source]

Reddening law coefficients of the form Az = Rz*E(B-V). the coefficients are from Table 6 of Schlafly & Finkbeiner (2011) where available. The entries are in a polynomial format for Rz defined as: Rz = a0 + T4*(a1 + a2*T4) + a3*FeH with T4 = 1e-4*Teff. They are kept like this for backward compatibility reasons with Casagrande & VandenBerg (2014).

Coefficients were extracted from the following references: G19: Green et al. 2019 SF11: Schlafly & Finkbeiner 2011 SD18: Sanders & Das 2018 CV14: Casagrande & Vandenberg 2014 CV18: Casagrande & Vandenberg 2018 Y13: Yuan et al. 2013

We aim for homogeneity and prioritise those of SF11, and for systems not available in that compilation we use SD18 and CV14/18.

class constants.freqtypes[source]

Different possibilities of fitting frequencies, for global access

class constants.metallicityranges[source]

Limits in metallictity for colors

class constants.parameters[source]

All the different parameters in the form: (name, unit, pname, remark, color)

  • Note some parameters are only available for certain tracks.

  • Color is for the Kiel diagram

The list is available in table format in the documentation.

exclude_params()[source]

Takes a list of input parameters (or a single parameter) as strings and returns the entire params list, except for the params given as input.

get_keys()[source]

Takes a list of input parameters (or a single parameter) as strings and returns the correspding units, names shown on a plot and remarks for the params.

class constants.photsys[source]

Available photometric systems and mapping to internal names

class constants.statdata[source]

Constant values for statistics, to ensure consistensy across code

Contains

quantileslist

Median, lower and upper percentiles of Bayesian posterior distributions to draw

nsamplesint

Number of samples to draw when sampling

nsigmafloat

Fractional standard deviation used for smoothing

class constants.sydsun[source]

Default solar values from the SYD asteroseismic pipeline.

The distance module

The asteroseismic module

Fitting of frequencies and frequency-dependent properties.

freq_fit.BG14(joinkeys, join, scalnu, method='l1', onlyl0=False)[source]

Ball & Gizon frequency correction

Correcting stellar oscillation frequencies for near-surface effects. This routine follows the approach from Warrick H. Ball and Laurent Gizon. “A new correction of stellar oscillation frequencies for near-surface effects.” Astronomy & Astrophysics 568 (2014): A123.

The correction has the form:

\[d\nu = \frac{a \left(\frac{\nu}{\nu_{scal}}\right)^{-1} + b \left(\frac{\nu}{\nu_{scal}}\right)^3}{I}\]

where \(a\) and \(b\) are the found parameters, \(\nu_{scal}\) is a scaling frequency used to non-dimensionalize the frequencies \(\nu\), and \(I\) is the mode inertia for each mode.

Parameters:
  • joinkeys (int array) – Array containing the mapping of model mode to the observed mode. Model mode (l=joinkeys[i, 0], n=joinkeys[i, 1]) is mapped to observed mode (l=joinkeys[i, 0], n=joinkeys[i, 2]).

  • join (array) – Array containing the model frequency, model inertia, observed frequency, uncertainty in observed frequency for a pair of mapped modes. Model mode (l=joinkeys[i, 0], n=joinkeys[i, 1]) has frequency join[i, 0], inertia join[i, 1]. Observed mode (l=joinkeys[i, 0], n=joinkeys[i, 2]) has frequency join[i, 2] and uncertainty join[i, 3].

  • scalnu (float) – A scaling frequency used purely to non-dimensionalize the frequencies.

  • method (string, optional) – The name of the fitting method used. It can be ‘ransac’ (default), ‘l1’, or ‘l2’.

  • onlyl0 (bool) – Flag that if True only computes the correction based on the l=0 modes. This can be usefull if not many modes have been observed.

Returns:

  • cormod (array) – Array containing the corrected model frequencies and the unchanged mode inertias in the same format as the osc-arrays.

  • coeffs (array) – Array containing the coefficients a and b to the found correction.

freq_fit.HK08(joinkeys, join, nuref, bcor)[source]

Kjeldsen frequency correction

Correcting stellar oscillation frequencies for near-surface effects following the approach in Hans Kjeldsen, Timothy R. Bedding, and Jørgen Christensen-Dalsgaard. “Correcting stellar oscillation frequencies for near-surface effects.” The Astrophysical Journal Letters 683.2 (2008): L175.

The correction has the form:

\[r\nu{corr} - \nu{model} = \frac{a}{Q}\left(\frac{\nu{model}}{\nu_{max}}\right)^b .\]
Parameters:
  • joinkeys (int array) – Array containing the mapping of model mode to the observed mode. Model mode (l=joinkeys[i, 0], n=joinkeys[i, 1]) is mapped to observed mode (l=joinkeys[i, 0], n=joinkeys[i, 2]).

  • join (array) – Array containing the model frequency, model inertia, observed frequency, uncertainty in observed frequency for a pair of mapped modes. Model mode (l=joinkeys[i, 0], n=joinkeys[i, 1]) has frequency join[i, 0], inertia join[i, 1]. Observed mode (l=joinkeys[i, 0], n=joinkeys[i, 2]) has frequency join[i, 2] and uncertainty join[i, 3].

  • nuref (float) – The reference frequency used in the Kjeldsen correction. Normally the reference frequency used is numax of the observed star or numax of the Sun, but it shouldn’t make a difference.

  • bcor (float) – The exponent in the Kjeldsen correction.

Returns:

  • cormod (array) – Array containing the corrected model frequencies and the unchanged mode inertias in the same format as the osc-arrays.

  • coeffs (array) – Array containing the coefficients in the found correction.

freq_fit.apply_BG14(modkey, mod, coeffs, scalnu)[source]

Applies the BG14 frequency correction to frequency for a given set of coefficients

Parameters:
  • modkey (array) – Array containing the angular degrees and radial orders of mod.

  • mod (array) – Array containing the modes in the model.

  • coeffs (list) – List containing the two coefficients used in the BG14 correction.

  • scalnu (float) – A scaling frequency used purely to non-dimensionalize the frequencies.

Returns:

corosc (array) – Array containing the corrected model frequencies.

freq_fit.apply_HK08(modkey, mod, coeffs, scalnu)[source]

Applies the HK08 frequency correction to frequency for a given set of coefficients

Parameters:
  • modkey (array) – Array containing the angular degrees and radial orders of mod.

  • mod (array) – Array containing the modes in the model.

  • coeffs (list) – List containing [rpar, acor, bcor] used in the HK08 correction.

  • scalnu (float) – A scaling frequency used purely to non-dimensionalize the frequencies.

Returns:

corosc (array) – Array containing the corrected model frequencies.

freq_fit.apply_cubicBG14(modkey, mod, coeffs, scalnu)[source]

Applies the BG14 cubic frequency correction to frequency for a given set of coefficients

Parameters:
  • modkey (array) – Array containing the angular degrees and radial orders of mod.

  • mod (array) – Array containing the modes in the model.

  • coeffs (list) – List containing the two coefficients used in the BG14 correction.

  • scalnu (float) – A scaling frequency used purely to non-dimensionalize the frequencies.

Returns:

corosc (array) – Array containing the corrected model frequencies.

freq_fit.calc_join(mod, modkey, obs, obskey, obsintervals=None, dnu=None)[source]

This functions maps the observed modes to the model modes.

Parameters:
  • mod (array) – Array containing the modes in the model.

  • modkey (array) – Array containing the angular degrees and radial orders of mod

  • obs (array) – Array containing the modes in the observed data

  • obskey (array) – Array containing the angular degrees and radial orders of obs

  • obsintervals (array, optional) – Array containing the endpoints of the intervals used in the frequency fitting routine in freq_fit.calc_join. As it is the same in all iterations for the observed frequencies, when computing chi2, this is computed in su.prepare_obs once and given as an argument in order to save time. If obsintervals is None, it will be computed.

  • dnu (float, optional) – Large frequency separation in microhertz used for the computation of obsintervals. If left None, the large frequency separation will be computed as the median difference between consecutive l=0 modes.

Returns:

joins (list or None) – List containing joinkeys and join.

  • joinkeysint array

    Array containing the mapping of model mode to the observed mode. Model mode (l=joinkeys[i, 0], n=joinkeys[i, 1]) is mapped to observed mode (l=joinkeys[i, 0], n=joinkeys[i, 2]).

  • joinarray

    Array containing the model frequency, model inertia, observed frequency, uncertainty in observed frequency for a pair of mapped modes. Model mode (l=joinkeys[i, 0], n=joinkeys[i, 1]) has frequency join[i, 0], inertia join[i, 1]. Observed mode (l=joinkeys[i, 0], n=joinkeys[i, 2]) has frequency join[i, 2] and uncertainty join[i, 3].

freq_fit.compute_dnu_wfit(obskey, obs, numax)[source]

Compute large frequency separation weighted around numax, the same way as dnufit. Coefficients based on White et al. 2011.

Parameters:
  • obskey (array) – Array containing the angular degrees and radial orders of obs

  • obs (array) – Individual frequencies and uncertainties.

  • numax (scalar) – Frequency of maximum power

Returns:

  • dnu (scalar) – Large frequency separation obtained by fitting the radial mode observed frequencies.

  • dnu_err (scalar) – Uncertainty on dnudata.

freq_fit.compute_epsilondiff(osckey, osc, avgdnu, sequence='e012', nsorting=True, extrapolation=False, nrealisations=20000, debug=False)[source]

Compute epsilon differences and covariances.

From Roxburgh 2016: * Eq. 1: Epsilon(n,l) * Eq. 4: EpsilonDifference(l=0,l=(1,2))

Epsilon differences are independent of surface phase shift/outer layers when the epsilons are evaluated at the same frequency. It therefore relies on splining from epsilons at the observed frequencies of the given degree and order to the frequency of the compared/subtracted epsilon. See function “compute_epsilondiffseqs” for further clarification.

For MonteCarlo sampling of the covariances, it is replicated from the covariance determination of frequency ratios in BASTA, (sec 4.1.3 of Aguirre Børsen-Koch et al. 2022). A number of realisations of the epsilon differences are drawn from random Gaussian distributions of the individual frequencies within their uncertainty.

Parameters:
  • osckey (array) – Array containing the angular degrees and radial orders of the modes.

  • osc (array) – Array containing the modes (and inertias).

  • avgdnu (float) – Average value of the large frequency separation.

  • sequence (str, optional) – Similar to ratios, what sequence of epsilon differences to be computed. Can be e01, e02 or e012 for a combination of the two first.

  • nsorting (bool, optional) – If True (default), the sequences are sorted by n-value of the frequencies. If False, the entire 01 sequence is followed by the 02 sequence.

  • extrapolation (bool, optional) – If False (default), modes outside the range of the l=0 modes are discarded to avoid extrapolation.

  • nrealisations (int or bool, optional) – If int: number of realisations used for MC-sampling the covariances If bool: Whether to use MC (True) or analytic (False) deternubation of covariances. If True, nrealisations of 20,000 is used.

  • debug (bool, optional) – Print additional output and make plots for debugging (incl. a plot of the correlation matrix).

Returns:

  • epsdiff (array) – Array containing the modes in the observed data.

  • epsdiff_cov (array) – Covariances matrix.

freq_fit.compute_epsilondiffseqs(osckey, osc, avgdnu, sequence, nsorting=True)[source]

Computed epsilon differences, based on Roxburgh 2016 (eq. 1 and 4)

Epsilons E of frequency v with order n and degree l is determined as: E(n,l) = E(v(n,l)) = v(n,l)/dnu - n - l/2

From this, an epsilon is determined for each original frequncy. These are not independent on the surface layers, but their differences between different degrees are, if evaluated at the same frequency. Therefore, the epsilon differences dE of e.g. E(n,l=0) and E(n,l=2), dE(02) is determined from interpolating/splining the l=0 epsilon sequence SE0 and evaluating it at v(n,l=2), and subtracting the corresponding E(n,l=2). Therefore, the epsilon difference can be summarised as dE(0l) = SE0(v(n,l)) - E(n,l)

Parameters:
  • osckey (array) – Array containing the angular degrees and radial orders of the modes

  • osc (array) – Array containing the modes (and inertias)

  • avgdnu (float) – Average large frequency separation

  • sequence (str) – Similar to ratios, what sequence of epsilon differences to be computed. Can be 01, 02 or 012 for a combination of the two first.

  • nsorting (bool) – If True (default), the sequences are sorted by n-value of the frequencies. If False, the entire 01 sequence is followed by the 02 sequence.

Returns:

deps (array) – Array containing epsilon differences. First index correpsonds to: 0 - Epsilon differences 1 - Indentifying frequencies 2 - Identifying degree l 3 - Radial degree n of identifying l={1,2} mode

freq_fit.compute_ratios(obskey, obs, ratiotype, nrealisations=10000, threepoint=False)[source]

Routine to compute the ratios r02, r01 and r10 from oscillation frequencies, and return the desired ratio sequence, both individual and combined sequences, along with their covariance matrix and its inverse.

Developers note: We currently do not store identifying information of the sequence origin in the combined sequences. We rely on them being constructed/sorted identically when computed.

Parameters:
  • obskey (array) – Harmonic degrees, radial orders and radial orders of frequencies.

  • obs (array) – Frequencies and their error, following the structure of obs.

  • ratiotype (str) – Which ratio sequence to determine, see constants.freqtypes.rtypes for possible sequences.

  • nrealizations (int) – Number of realizations of the sampling for the computation of the covariance matrix.

  • threepoint (bool) – If True, use three point definition of r01 and r10 ratios instead of default five point definition.

Returns:

  • ratio (array) – Ratio requested from ratiotype.

  • ratio_cov (array) – Covariance matrix of the requested ratio.

freq_fit.compute_ratioseqs(obskey, obs, sequence, threepoint=False)[source]

Routine to compute the ratios r02, r01 and r10 from oscillation frequencies, and return the desired ratio sequence, both individual and combined sequences.

Developers note: We currently do not store identifying information of the sequence origin in the combined sequences. We rely on them being constructed/sorted identically when computed.

Parameters:
  • obskey (array) – Harmonic degrees, radial orders and radial orders of frequencies

  • obs (array) – Frequencies and their error, following the structure of obs

  • sequence (str) – Which ratio sequence to determine, see constants.freqtypes.rtypes for possible sequences.

  • threepoint (bool) – If True, use three point definition of r01 and r10 ratios instead of default five point definition.

Returns:

ratio (array) – Ratio requested from sequence. First index correspond to: 0 - Frequency ratios 1 - Defining/corresponding frequency 2 - Identifying integer (r01: 1, r02: 2, r10: 10) 3 - Identifying radial order n

freq_fit.cubicBG14(joinkeys, join, scalnu, method='l1', onlyl0=False)[source]

Ball & Gizon frequency correction

Correcting stellar oscillation frequencies for near-surface effects. This routine follows the approach from Warrick H. Ball and Laurent Gizon. “A new correction of stellar oscillation frequencies for near-surface effects.” Astronomy & Astrophysics 568 (2014): A123.

The correction has the form:

\[d\nu = \frac{b \left(\frac{\nu}{\nu_{scal}}\right)^3}{I}\]

where \(b\) are the found parameters, \(\nu_{scal}\) is a scaling frequency used to non-dimensionalize the frequencies \(\nu\), and \(I\) is the mode inertia for each mode.

Parameters:
  • joinkeys (int array) – Array containing the mapping of model mode to the observed mode. Model mode (l=joinkeys[i, 0], n=joinkeys[i, 1]) is mapped to observed mode (l=joinkeys[i, 0], n=joinkeys[i, 2]).

  • join (array) – Array containing the model frequency, model inertia, observed frequency, uncertainty in observed frequency for a pair of mapped modes. Model mode (l=joinkeys[i, 0], n=joinkeys[i, 1]) has frequency join[i, 0], inertia join[i, 1]. Observed mode (l=joinkeys[i, 0], n=joinkeys[i, 2]) has frequency join[i, 2] and uncertainty join[i, 3].

  • scalnu (float) – A scaling frequency used purely to non-dimensionalize the frequencies.

  • method (string, optional) – The name of the fitting method used. It can be ‘ransac’ (default), ‘l1’, or ‘l2’.

  • onlyl0 (bool) – Flag that if True only computes the correction based on the l=0 modes. This can be usefull if not many modes have been observed.

Returns:

  • cormod (array) – Array containing the corrected model frequencies and the unchanged mode inertias in the same format as the osc-arrays.

  • coeffs (array) – Array containing the coefficients a and b to the found correction.

freq_fit.make_intervals(osc, osckey, dnu=None)[source]

This function computes the interval bins used in the frequency mapping in freqfit.calc_join().

Parameters:
  • osc (array) – Array containing either model or observed modes.

  • osckey (array) – Array containing the angular degrees and radial orders of osc

  • dnu (float, optional) – Large frequency separation in microhertz. If left None, the large frequency separation will be computed as the median difference between consecutive l=0 modes.

Returns:

intervals (array) – Array containing the endpoints of the intervals used in the frequency fitting routine in freq_fit.calc_join`().

Auxiliary functions for glitch fitting

glitch_fit.compute_glitchseqs(osckey: array, osc: array, sequence: str, dnu: float, fitfreqs: dict, ac_depths: bool = False, debug: bool = False) array[source]

Routine to compute glitch parameters of given frequencies, based on the given method options.

Parameters:
  • osckey (numpy array) – Spherical degrees and radial orders of the frequencies to be used.

  • osc (numpy array) – Frequencies and corresponding uncertainties.

  • sequence (str) – Glitch sequence to be computed, see constants.freqtypes.glitches.

  • dnu (float) – Value of large frequency separation to be used in the computation.

  • fitfreqs (dict) – Dictionary containing frequency fitting options/controls.

  • ac_depts (bool or dict) – Acoustic depths used to search for glitch signatures. If not provided as a dict, they will be calculated as a simple estimate.

Returns:

glitchseq (numpy array) – Determined glitch parameters (and ratios) from the provided frequencies. If computation failed, the glitch parameters will be NaNs.

glitch_fit.compute_observed_glitches(osckey: array, osc: array, sequence: str, dnu: float, fitfreqs: dict, debug=False) tuple[array, array][source]

Routine to compute glitch parameters (and ratios) with full covariance matrix using MC sampling.

Parameters:
  • osckey (numpy array) – Spherical degrees and radial orders of the frequencies to be used.

  • osc (numpy array) – Frequencies and corresponding uncertainties.

  • sequence (str) – Glitch sequence to be computed, see constants.freqtypes.glitches.

  • dnu (float) – Value of large frequency separation to be used in the computation.

  • fitfreqs (dict) – Dictionary containing frequency fitting options/controls.

Returns:

  • glitchseq (numpy array) – Computed glitch parameters (and ratios) as median values from MC sampling

  • glitchseq_cov (numpy array) – Determined covariance matrix of glitch parameters (and ratios)

The stats module

Definition of priors

Define any prior functions here that you want to use in BASTA!!

The prior function must be of the form: PRIOR = PRIORFUN(LIBITEM, INDEX)

Any prior defined here can be used from an .xml input file.

priors.baldryglazebrook2003(libitem, index)[source]

Initial mass function from Baldry & Glazebrook (2003) https://ui.adsabs.harvard.edu/abs/2003ApJ…593..258B The global normalisation is not needed as we normalise later.

priors.chabrier2003(libitem, index)[source]

Initial mass function from Chabrier (2003) https://ui.adsabs.harvard.edu/abs/2003PASP..115..763C/abstract Note that this is in linear mass space, hence the (1/m)

priors.kennicutt1994(libitem, index)[source]

Initial mass function from Kennicutt et al. 1994 https://ui.adsabs.harvard.edu/abs/1994ApJ…435…22K The global normalisation is not needed as we normalise later.

priors.kroupa2001(libitem, index)[source]

Initial mass function from Kroupa (2001) https://ui.adsabs.harvard.edu/abs/2001MNRAS.322..231K https://ui.adsabs.harvard.edu/abs/2002Sci…295…82K The global normalisation is not needed as we normalise later.

priors.millerscalo1979(libitem, index)[source]

Initial mass function from Miller & Scalo (1979) https://ui.adsabs.harvard.edu/abs/1979ApJS…41..513M The global normalisation is not needed as we normalise later.

priors.salpeter1955(libitem, index)[source]

Initial mass function from Salpeter (1955) https://ui.adsabs.harvard.edu/abs/1955ApJ…121..161S

priors.scalo1998(libitem, index)[source]

Initial mass function from Scalo (1998) https://ui.adsabs.harvard.edu/abs/1998ASPC..142..201S The global normalisation is not needed as we normalise later.

Key statistics functions

class stats.Trackdnusurf(dnusurf)
dnusurf

Alias for field number 0

class stats.Trackglitchpar(AHe, dHe, tauHe)
AHe

Alias for field number 0

dHe

Alias for field number 1

tauHe

Alias for field number 2

class stats.Trackstats(index, logPDF, chi2)
chi2

Alias for field number 2

index

Alias for field number 0

logPDF

Alias for field number 1

stats.calc_key_stats(x, centroid, uncert, weights=None)[source]

Calculate and report the wanted format of centroid value and uncertainties.

Parameters:
  • x (list) – Parameter values of the models with non-zero likelihood.

  • centroid (str) – Options for centroid value definition, ‘median’ or ‘mean’.

  • uncert (str) – Options for reported uncertainty format, ‘quantiles’ or ‘std’ for standard deviation.

  • weights (list) – Weights needed to be applied to x before extracting statistics.

Returns:

  • xcen (float) – Centroid value

  • xm (float) – 16’th percentile if quantile selected, 1 sigma for standard deviation.

  • xp (float, None) – 84’th percentile if quantile selected, None for standard deviation.

stats.chi2_astero(obskey, obs, obsfreqmeta, obsfreqdata, obsintervals, libitem, ind, fitfreqs, warnings=True, shapewarn=0, debug=False, verbose=False)[source]

chi2_astero calculates the chi2 contributions from the asteroseismic fitting e.g., the frequency, ratio, or glitch fitting.

Parameters:
  • obskey (array) – Array containing the angular degrees and radial orders of obs

  • obs (array) – Array containing the modes in the observed data

  • obsfreqmeta (dict) – The requested information about which frequency products to fit or plot, unpacked for easier access later.

  • obsfreqdata (dict) – Combinations of frequencies read from the frequency input files. Contains ratios, glitches, epsilon differences and covariance matrices.

  • obsintervals (array) – Array containing the endpoints of the intervals used in the frequency fitting routine in :func:’freq_fit.calc_join’. As it is the same in all iterations for the observed frequencies, this is computed in su.prepare_obs once and given as an argument in order to save time and memory.

  • libitem (hdf5 group) – Contains the entire track of models being processed.

  • ind (int) – The index in libitem corresponding to the current model being processed.

  • fitfreqs (dict) – Contains all user inputted frequency fitting options.

  • warnings (bool) – If True, print something when it fails.

  • shapewarn (int) – If a mismatch in array dimensions of the fitted parameters, or range of frequencies is encountered, this is set to corresponding integer in order to warn the user at the end of the run.

  • debug (str) – Flag for print control.

  • verbose (str) – Flag for print control.

Returns:

  • chi2rut (float) – The calculated chi2 value for the given model to the observed data. If the calculated value is less than zero, chi2rut will be np.inf.

  • warnings (bool) – See ‘warnings’ above.

stats.chi_for_plot(selectedmodels)[source]

FOR VALIDATION PLOTTING!

Could this be combined with the functions above in a wrapper?

Return chi2 value of HLM and LCM.

stats.get_highest_likelihood(Grid, selectedmodels, inputparams)[source]

Find highest likelihood model and print info.

Parameters:
  • Grid (hdf5 object) – The already loaded grid, containing the tracks/isochrones.

  • selectedmodels (dict) – Contains information on all models with a non-zero likelihood.

  • inputparams (dict) – The standard bundle of all fitting information

Returns:

  • maxPDF_path (str) – Path to model with maximum likelihood.

  • maxPDF_ind (int) – Index of model with maximum likelihood.

stats.get_lowest_chi2(Grid, selectedmodels, inputparams)[source]

Find model with lowest chi2 value and print info.

Parameters:
  • Grid (hdf5 object) – The already loaded grid, containing the tracks/isochrones.

  • selectedmodels (dict) – Contains information on all models with a non-zero likelihood.

  • inputparams (dict) – The standard bundle of all fitting information

Returns:

  • minchi2_path (str) – Path to model with lowest chi2

  • minchi2_ind (int) – Index of model with lowest chi2

stats.lowest_chi2(selectedmodels)[source]

Find the index of the model with the lowest chi2 (only gaussian fitting parameters)

Parameters:

selectedmodels (dict) – Contains information on all models with a non-zero likelihood.

Returns:

  • minchi2_path (str) – Path to the model with the lowest chi2.

  • minchi2_ind (int) – Index of the model with the lowest chi2.

stats.most_likely(selectedmodels)[source]

Find the index of the model with the highest probability

Parameters:

selectedmodels (dict) – Contains information on all models with a non-zero likelihood.

Returns:

  • maxPDF_path (str) – The path to the model with the highest probability.

  • maxPDF_ind (int) – Index of the model with the highest probability.

stats.posterior(x, nonzeroprop, sampled_indices, nsigma=0.25)[source]

Compute posterior of x as a weighted histogram of x with weights y. The bin width is determined by the Freedman Diaconis Estimator The histogram is smoothed using a Gaussian kernel with a width of nsigma*sigma

Parameters:
  • x (array) – samples

  • y (array) – likelihood of samples

  • nonzeroprop (array) – indices of x with non-zero likelihood

  • sampled_indices – indices of weighted draws from y, used for bin estimation

  • nsigma (float) – fractional standard deviation used for smoothing

Returns:

like (function) – interpolated posterior

stats.priorlogPDF

alias of Trackstats

stats.quantile_1D(data, weights, quantile)[source]

Compute the weighted quantile of a 1D numpy array. The function is borrowed from the Python package wquantiles

Parameters:
  • data (ndarray) – Input array (one dimension).

  • weights (ndarray) – Array with the weights of the same size of data.

  • quantile (float) – Quantile to compute. It must have a value between 0 and 1.

Returns:

result (float) – The output value.