The AIMS program

A module which contains the main program for AIMS as well as various classes which intervene when calculating the priors and likelihood function:

  • Distribution: a class which represents a probability distribution
  • Prior_list: a class with a list of priors
  • Mode: a class used to represent observed modes
  • Combination: a class used to represent frequency combinations
  • Likelihood: a class used to represent the likelihood function
  • Probability: a class which groups the priors and likelihood function together

This module relies on the emcee package to apply an MCMC algorithm which will return a representative sample of models for a given set of seismic an classic constraints.

Warning

In various places in this module, for instance in the Prior_list and Likelihood classes, various methods return what is described as a \(\chi^2\) value. Technically, these are not \(\chi^2\) values, but rather \(-\chi^2/2\), i.e. the argument of the exponential function which intervenes in the Gaussian probability distribution.

class AIMS.Combination[source]

A class which contains indices and coefficients which intervene in:

  • linear combinations of frequencies
  • frequency ratios
add_den(j, coeff)[source]

Append the given index and coefficient to the list of denominator indices and coefficients.

Parameters:
  • j (int) – index of the mode
  • coeff (float) – coefficient used in the frequency combination
add_num(j, coeff)[source]

Append the given index and coefficient to the list of numerator indices and coefficients.

Parameters:
  • j (int) – index of the mode
  • coeff (float) – coefficient used in the frequency combination
den = None

Value of the denomenator in a frequency ratio.

den_coeff = None

Coefficients in the denominator of a frequency ratio, otherwise empty.

den_index = None

Indices in the denominator of a frequency ratio, otherwise empty.

num = None

Value of the frequency combination or numerator in a frequency ratio.

num_coeff = None

Coefficients in a linear combination or numerator of a frequency ratio.

num_index = None

Indices in a linear combination or numerator of a frequency ratio.

print_me()[source]

Print frequency combination.

value = None

Value of the frequency combination or ratio.

class AIMS.Distribution(_type, _values)[source]

A class which represents a probability distribution, and can yield its value for a given input parameter, or provide a random realisation.

Note

Derived from a class originally written by G. Davies.

Parameters:
  • _type (string) – type of probability function (current options include “Gaussian”, “Truncated_gaussian”, “Uniform”, “IMF1”, “IMF2”, “Uninformative”)
  • _values (list of floats) – list of parameters relevant to the probability function
error_bar

Returns an error bar based on the distribution. This does not necessarily correspond to the one-sigma value but rather to what is the most convenient value.

Returns:the error bar
Return type:float
mean

Returns the mean value of the probability distribution.

Returns:the mean value of the probability distribution
Return type:float
nparams

Return the number of relevant parameters for a given distribution.

Returns:the number of relevant parameters
Return type:int
print_me()[source]

Print type and parameters of probability distribution.

re_centre(value)[source]

Re-centre the probability distribution around the input value.

Parameters:value (float) – new value around which to centre the distribution
re_normalise(value)[source]

Re-normalise the probability distribution so that its characteristic width corresponds to the input value.

Parameters:value (float) – new value around for the chacteristic width
realisation(size=None)[source]

Return random values which statistically follow the probability distribution.

Parameters:size (int or tuple of ints) – shape of random variates
Returns:a set of random realisations
Return type:float
to_string()[source]

Produce nice string representation of the distribution.

Returns:nice string representation of the distribution
Return type:string
type = None

Type of probability function (“Uniform”, “Gaussian”, “Truncated_gaussian”, “IMF1”, “IMF2”, or “Uninformative”)

values = None

List of parameters relevant to probability function

class AIMS.Likelihood[source]

A class which described the likelihood function and allows users to evaluate it.

add_combinations(num_list, den_list=[], target_ell=None)[source]

This finds the indices of modes which intervene in a frequency combination or ratio, as specified by the mandatory and optional arguments. These indices, the relevant coefficients, the numerator, the denominator, and the resultant value of the combination are stored in the combinations variable.

Parameters:
  • num_list (list of (int,int,float)) – list of relative mode identifications and coefficients used to define a frequency combination or the numerator of a frequency ratio. This list contains tuples of the form (delta n, delta l, coeff).
  • den_list (list of (int,int,float)) – list of relative mode identifications and coefficients used to define the denominator of a frequency ratio. If absent, then, it is assumed that a linear combination of frequencies is represented. The form is the same as for num_list.
  • target_ell (int) – this is used to impose a specific l value on the first selected mode.
add_constraint(constraint)[source]

Add a supplementary constraint to the list of constraints.

Parameters:constraint ((string, Distribution)) – supplementary constraint
add_dnu_constraint(l_targets=[0])[source]

Add the large frequency separation as a contraint. The coefficients are obtained via a least-squares approach. The approach taken here has two advantages:

  1. Correlations between the large frequency separation and other seismic constraints will be taken into account.
  2. The same modes will be used in the same way, both for the observations and the models.
Parameters:l_targets (list of int) – specifies for which l values the large frequency separation is to be calculated. If None is supplied, all modes will be used.

Note

This uses an analytical approach and is therefore the prefered method.

add_dnu_constraint_matrix(l_targets=[0])[source]

Add the large frequency separation as a contraint. The coefficients are obtained via a least-squares approach. The approach taken here has two advantages:

  1. Correlations between the large frequency separation and other seismic constraints will be taken into account.
  2. The same modes will be used in the same way, both for the observations and the models.
Parameters:l_targets (list of int) – specifies for which l values the large frequency separation is to be calculated. If None is supplied, all modes will be used.

Note

This uses a matrix approach and is therefore not the prefered method.

add_nu_min_constraint(target_ell=0, min_n=False)[source]

Add the minimun frequency/mode of a specific ell value as a seismic constraint. Typically, such constraints are used as an “anchor” when combined with constraints based on frequency ratios.

Parameters:
  • target_ell (int) – ell value of the minimum frequency/mode
  • min_n (boolean) – if False, look for minimum observational frequency. If True, look for minimum radial order.
add_seismic_constraint(string)[source]

Add seismic contraints based on the keyword given in string.

Parameters:string (string) –

keyword which specifies the type of constraint to be added. Current options include:

  • nu: individual frequencies
  • nu0: individual frequencies (radial modes only)
  • nu_min0: radial mode with minimum frequency
  • r02: \(r_{02}\) frequency ratios
  • r01: \(r_{01}\) frequency ratios
  • r10: \(r_{10}\) frequency ratios
  • dnu: individual large frequency separations (using all modes)
  • dnu0: individual large frequency separations (using radial modes only)
  • d2nu: second differences (using all modes)
  • avg_dnu: average large frequency separation (using all modes)
  • avg_dnu0: average large frequency separation (using radial modes only)
apply_constraints(my_model)[source]

Calculate a \(\chi^2\) value for the set of constraints (excluding seismic constraints based on mode frequencies).

Parameters:my_model (model.Model) – model for which the \(\chi^2\) value is being calculated
Returns:the \(\chi^2\) value deduced from classic constraints
Return type:float
assign_n(my_model)[source]

Assign the radial orders based on proximity to theoretical frequencies from an input model.

Parameters:my_model (model.Model) – input model
classic_weight = None

Absolute weight to be applied to classic constraints (incl. nu_max constraint).

clear_seismic_constraints()[source]

This clears the seismic constraints. Specifically, the list of seismic combinations, and associated covariance matrix and its inverse are reinitialised.

coeff = None

3D float array with the coefficients for each frequency combination. The indices are:

  1. The index of the term
  2. The type of term (0 = num, 1 = den)
  3. The index of the frequency combination
combinations = None

This contains indices and coefficients to frequency combinations and frequency ratios.

compare_frequency_combinations(my_model, mode_map, a=[])[source]

This finds a \(\chi^2\) value based on a comparison of frequencies combinations, as defined in the combinations variable.

Parameters:
  • my_model (model.Model) – model for which the \(\chi^2\) value is being calculated
  • mode_map (list of int) – a mapping which relates observed modes to theoretical ones
  • a (array-like) – parameters of surface correction terms
Returns:

the \(\chi^2\) value for the seismic constraints

Return type:

float

Note

I’m assuming none of the modes are missing (i.e. that mode_map doesn’t contain the value -1)

constraints = None

List of constraints which intervene in the likelihood function.

cov = None

Covariance matrix which intervenes when calculating frequency combinations.

create_combination_arrays()[source]

Create array form of frequency combinations to be used with a fortran based routine for calculating the seismic chi^2 value.

create_mode_arrays()[source]

Create arrays with mode parameters (n, l, freq), which can be interfaced with fortran methods more easily.

dfvalues = None

Array with the error bars on the observed frequencies

evaluate(my_model)[source]

Calculate ln of likelihood function (i.e. a \(\chi^2\) value) for a given model.

Parameters:my_model (model.Model) – model for which the \(\chi^2\) value is being calculated
Returns:the \(\chi^2\) value, and optionally the optimal surface amplitudes (depending on the value of AIMS_configure.surface_option)
Return type:float, np.array (optional)

Note

This avoids model interpolation and can be used to gain time.

find_covariance()[source]

This prepares the covariance matrix and its inverse based on the frequency combinations in combinations.

Warning

This method should be called after all of the methods which add to the list of frequency combinations.

find_l_list(l_targets)[source]

Find a list of l values with the following properties:

  • each l value only occurs once
  • each l value given in the parameter l_targets is in the result l list, except if there is 1 or less modes with this l value
  • if the parameter l_targets is None, look for all l values with 2 or more modes associated with them
Parameters:l_targets (list of int) – input list of l values
Returns:new list of l values with the above properties
Return type:list of int
find_map(my_model, use_n)[source]

This finds a map which indicates the correspondance between observed modes and theoretical modes from my_model.

Parameters:
  • my_model – model for which the \(\chi^2\) value is being calculated
  • use_n (boolean) – specify whether to use the radial order when finding the map from observed modes to theoretical modes. If False, the map is based on frequency proximity.
Returns:

the correspondance between observed and theoretical modes from the above model, and the number of observed modes which weren’t mapped onto theoretical modes

Return type:

list of int, int

Note

  • a value of -1 is used to indicate that no theoretical mode corresponds to a particular observed mode.
  • only zero or one observed mode is allowed to correspond to a theoretical mode
find_vec(a_combination)[source]

This finds a set of coefficients which intervene when constructing the coviance matrix for frequency combinations.

Parameters:a_combination (Combination) – variable which specifies the frequency combination.
Returns:the above set of coefficients
Return type:np.array
find_weights()[source]

Find absolute weights for seismic and classic constraints based on options in AIMS_configure.py.

find_weights_new()[source]

Find absolute weights for seismic and classic constraints based on options in AIMS_configure.py.

fvalues = None

Array with the observed frequencies

get_optimal_surface_amplitudes(my_model, mode_map)[source]

Find optimal surface correction amplitude, for the surface correction specified by surface_option.

Parameters:
  • my_model (model.Model) – the model for which we’re finding the surface correction amplitude
  • mode_map (list of int) – a mapping which relates observed modes to theoretical ones
Returns:

optimal surface correction amplitudes

Return type:

np.array

guess_dnu(with_n=False)[source]

Guess the large frequency separation based on the radial modes.

Parameters:with_n (boolean) – specifies whether to use the n values already stored with each mode, when calculating the large frequency separation.
Returns:the large frequency separation
Return type:float
guess_n()[source]

Guess the radial order of the observed pulsations modes.

This method uses the large frequency separation, as calculated with guess_dnu(), to estimate the radial orders. These orders are subsequently adjusted to avoid multiple modes with the same identification. The resultant radial orders could be off by a constant offset, but this is not too problematic when computing frequency combinations or ratios.

indices = None

3D int array with the mode indices for each frequency combination. The indices are:

  1. The index of the term
  2. The type of term (0 = num, 1 = den)
  3. The index of the frequency combination
invcov = None

Inverse of covariance matrix, Likelihood.cov.

is_outside(params)[source]

Test to see if the given set of parameters lies outside the grid of models. This is done by evaluate the probability and seeing if the result indicates this.

Parameters:params (array-like) – input set of parameters
Returns:True if the set of parameters corresponds to a point outside the grid.
Return type:boolean
lvalues = None

Array with the l values of the observed modes

modes = None

List of pulsation modes (of type Mode).

ncoeff = None

2D int array with the number of terms for each frequency combination. The indices are:

  1. The type of term (0 = num, 1 = den)
  2. The index of the frequency combination
nvalues = None

Array with the n values of the observed modes

read_constraints(filename, factor=1.0)[source]

Read a file with pulsation data and constraints.

Parameters:
  • filename (string) – name of file with pulsation data.
  • factor (float) – multiplicative factor for pulsation frequencies. Can be used for conversions.
seismic_weight = None

Absolute weight to be applied to seismic constraints

sort_modes()[source]

Sort the modes. The ordering will depend on the value of use_n from the AIMS_configure.py file.

values = None

1D float array with the value for each frequency combination

class AIMS.Mode(_n, _l, _freq, _dfreq)[source]

A class which describes an observed pulsation mode.

Parameters:
  • _n (int) – radial order of observed mode
  • _l (int) – harmonic degree of observed mode.
  • _freq (float) – pulsation frequency (in \(\mathrm{\mu Hz}\)).
  • _dfreq (float) – error bar on pulsation frequency (in \(\mathrm{\mu Hz}\)).

Warning

Negative values are not accepted for _l, _freq, or _dfreq.

dfreq = None

Error bar on pulsation frequency (in \(\mathrm{\mu Hz}\)).

freq = None

Pulsation frequency (in \(\mathrm{\mu Hz}\)).

l = None

Harmonic degree of observed mode.

match(a_mode)[source]

Check to see if input mode has the same (n,l) values as the current mode.

Parameters:a_mode (Mode) – input mode which is being compared with current mode.
Returns:True if the input mode has the same (n,l) values as the current mode.
Return type:boolean
n = None

Radial order of observed mode.

class AIMS.Prior_list[source]

A class which contains a list of priors as well as convenient methods for adding priors and for evaluating them.

add_prior(aPrior)[source]

Add a prior to the list.

Parameters:aPrior (Distribution) – prior which is to be added to the list.
priors = None

A list of probability distributions which correspond to priors.

realisation(size=None)[source]

Return an array with realisations for each prior. The last dimension will correspond to the different priors.

Parameters:size (int or tuple of ints) – shape of random variates (for each prior)
Returns:a set of realisations
Return type:numpy float array
class AIMS.Probability(_priors, _likelihood)[source]

A class which combines the priors and likelihood function, and allows the the user to evalute ln of the product of these.

Parameters:
  • _priors (Prior_list) – input set of priors
  • _likelihood (Likelihood) – input likelihood function
evaluate(my_model)[source]

Evalulate the ln of the product of the priors and likelihood function, i.e. the probability, for a given model, to within an additive constant.

Parameters:my_model (model.Model) – input model
Returns:the ln of the probability
Return type:float

Note

This avoids model interpolation and can be used to gain time.

likelihood = None

The likelihood function.

priors = None

The set of priors.

AIMS.accepted_parameters = []

list of parameters associated with accepted models

AIMS.append_osm_parameter(config_osm, name, value, step, rate, bounds)[source]

Add a parameter in xlm format in the file with the classic constraints for OSM.

Parameters:
  • config_osm (lxml.etree._Element) – XLM etree element to which to add the parameter
  • name (string) – name of the parameter
  • value (float) – value of the parameter
  • step (float) – parameter step (this intervenes when numerically calculating derivatives with respect to this parameter)
  • rate (float) – parameter rate (this corresponds to a tolerance on this parameter)
  • bounds (float tuple) – bounds on the parameter
AIMS.append_osm_surface_effects(modes_osm, name, numax, values)[source]

Add a method with which to calculate surface effects to the OSM contraint file.

Parameters:
  • modes_osm (lxml.etree._Element) – XML element to which to add the surface effects method
  • name (string) – name of the method
  • numax (float) – value of numax
  • values (float tuple) – values which intervene in the method
AIMS.best_MCMC_model = None

best model from the MCMC run

AIMS.best_MCMC_params = None

parameters for the model best_MCMC_model

AIMS.best_MCMC_result = -1e+300

ln(probability) result for the model best_MCMC_model

AIMS.best_age_range = 0.0

Age range on track with best_model_model

AIMS.best_grid_model = None

best model from a scan of the entire grid

AIMS.best_grid_params = None

parameters for the model best_grid_model

AIMS.best_grid_result = -1e+300

ln(probability) result for the model best_grid_model

AIMS.check_configuration()[source]

Test the values of the variables in check_configuration to make sure they’re acceptable. If an unacceptable value is found, then this will stop AIMS and explain what variable has an erroneous value.

AIMS.echelle_diagram(my_model, my_params, model_name)[source]

Write text file with caracteristics of input model.

Parameters:
  • my_model (model.Model) – model for which we’re writing a text file
  • my_params (array-like) – parameters of the model
  • model_name (string) – name used to describe this model. This is also used when naming the text file.
AIMS.find_a_blob(params)[source]

Find a blob (i.e. supplementary output parameters) for a given set of parameters (for one model). The blob also includes the log(P) value as a first entry.

Parameters:params (array-like) – input set of parameters
Returns:list of supplementary output parameters
Return type:list of floats
AIMS.find_best_model()[source]

Scan through grid of models to find “best” model for a given probability function (i.e. the product of priors and a likelihood function).

AIMS.find_best_model_in_track(ntrack)[source]

Scan through an evolutionary track to find “best” model for prob, the probability function (i.e. the product of priors and a likelihood function).

Parameters:ntrack (int) – number of the evolutionary track
Returns:the ln(probability) value, and the “best” model
Return type:(float, model.Model)
AIMS.find_blobs(samples)[source]

Find blobs (i.e. supplementary output parameters) from a set of samples (i.e. for multiple models).

Parameters:samples (list/array of array-like) – input set of samples
Returns:set of supplementary output parameters
Return type:np.array
AIMS.grid = None

grid of models

AIMS.grid_params_MCMC = ()

parameters used in the MCMC run (excluding surface correction parameters)

AIMS.grid_params_MCMC_with_surf = ()

parameters used in the MCMC run (including surface correction parameters)

AIMS.init_walkers()[source]

Initialise the walkers used in emcee.

Returns:array of starting parameters
Return type:np.array
AIMS.interpolation_tests(filename)[source]

Carry out various interpolation tests and write results in binary format to file.

Parameters:filename (string) – name of file in which to write test results.

Note

The contents of this file may be plotted using methods from plot_interpolation_test.py.

AIMS.load_binary_data(filename)[source]

Read a binary file with a grid of models.

Parameters:filename (string) – name of file with grid in binary format
Returns:the grid of models
Return type:model.Model_grid
AIMS.log0 = -1e+300

a large negative value used to represent ln(0)

AIMS.my_map = None

pointer to the map function (either the parallel or sequential versions)

AIMS.ndims = 0

number of dimensions for MCMC parameters (includes nsurf)

AIMS.nsurf = 0

number of surface term parameters

AIMS.output_folder = None

folder in which to write the results

AIMS.plot_distrib_iter(samples, labels, folder)[source]

Plot individual distribution of walkers as a function of iterations.

Parameters:
  • samples (np.array) – samples from the emcee run
  • labels (list of strings) – labels for the different dimensions in parameters space
  • folder (string) – specify name of file in which to save plots of walkers.

Warning

This method must be applied before the samples are reshaped, and information on individual walkers lost.

AIMS.plot_histograms(samples, names, fancy_names, truths=None)[source]

Plot a histogram based on a set of samples.

Parameters:
  • samples (np.array) – samples form the emcee run
  • names (list of strings) – names of the quantities represented by the samples. This will be used when naming the file with the histogram
  • fancy_names (list of strings) – name of the quantities represented by the samples. This will be used as the x-axis label in the histogram.
  • truths (list of floats) – reference values (typically the true values or some other important values) to be added to the histograms as a vertical line
AIMS.plot_walkers(samples, labels, filename, nw=3)[source]

Plot individual walkers.

Parameters:
  • samples (np.array) – samples from the emcee run
  • labels (list of strings) – labels for the different dimensions in parameters space
  • filename (string) – specify name of file in which to save plots of walkers.
  • nw (int) – number of walkers to be plotted

Warning

This method must be applied before the samples are reshaped, and information on individual walkers lost.

AIMS.pool = None

pool from which to carry out parallel computations

AIMS.prob = None

Probability type object that represents the probability function which includes the likelihood and priors

AIMS.rejected_parameters = []

list of parameters associated with rejected models

AIMS.run_emcee()[source]

Run the emcee program.

Returns:the emcee sampler for the MCMC run
AIMS.statistical_model = None

model corresponding to statistical parameters

AIMS.statistical_params = None

parameters for the model statistical_model

AIMS.statistical_result = -1e+300

ln(probability) result for the model statistical_model

AIMS.string_to_title(string)[source]

Create fancy title from string.

Parameters:string (string) – string from which the title is created.
Returns:the fancy string title
Return type:string
AIMS.threshold = -1e+290

threshold for “accepted” models. Needs to be greater than log0

AIMS.tight_ball_distributions = None

Prior_list type object with the distributions for the initial tight ball

AIMS.write_LEGACY_summary(filename, KIC, labels, samples)[source]

Write a one line summary of the statistical properties based on a sequence of realisations to a file. The format matches that of the LEGACY project.

The results include:

  • average values for each variable (statistical mean)
  • error bars for each variable (standard mean deviation)
Parameters:
  • filename (string) – name of file in which to write the statistical properties
  • KIC (string) – KIC number of the star
  • labels (list of strings) – names of relevant variables
  • samples (np.array) – samples for which statistical properties are calculated
AIMS.write_binary_data(infile, outfile)[source]

Read an ascii file with a grid of models, and write corresponding binary file.

Parameters:
  • infile (string) – input ascii file name
  • outfile (string) – output binary file name
AIMS.write_combinations(filename, samples)[source]

Produce a list of linear combinations of grid models (based on interpolation) corresponding to the provided model parameters.

Parameters:
  • filename (string) – name of the file to which to write the model combinations
  • samples (np.array) – set of model parameters for which we would like to obtain the grid models and interpolation coefficients
AIMS.write_list_file(filename)[source]

Write list file from which to generate binary grid. Various filters can be included to reduce the number of models.

Note

This code is intended for developpers not first time users.

AIMS.write_model(my_model, my_params, my_result, model_name, extended=False)[source]

Write text file with caracteristics of input model.

Parameters:
  • my_model (model.Model) – model for which we’re writing a text file
  • my_params (array-like) – parameters of the model
  • my_result (float) – ln(P) value obtained for the model
  • model_name (string) – name used to describe this model. This is also used when naming the text file.
  • extended – if set to True, all of the theoretical modes are saved in the text file, including those not matched to observations
AIMS.write_osm_don(filename, my_model)[source]

Write file with choice of physical ingredients to be used by CESAM or CESTAM and OSM.

Parameters:
  • filename (string) – name of file which will contain the physical ingredients
  • my_model (model.Model) – model from which which is derived various physical constraints/settings

Note

Written by B. Herbert.

AIMS.write_osm_frequencies(filename, my_model)[source]

Write file with frequencies for Optimal Stellar Model (OSM), written by R. Samadi.

Parameters:
  • filename (string) – name of file which will contain the frequencies
  • my_model (model.Model) – model from which are derived the radial orders

Note

Written by B. Herbert.

AIMS.write_osm_xml(filename, my_params, my_model)[source]

Write file with classic constraints for OSM

Parameters:
  • filename (string) – name of file with classic constraints
  • my_model (model.Model) – model used in deriving some of the constraints

Note

Originally written by B. Herbert. Includes some modifications.

AIMS.write_readme(filename, elapsed_time)[source]

Write parameters relevant to this MCMC run.

Parameters:filename (string) – name of file in which to write the statistical properties
AIMS.write_samples(filename, labels, samples)[source]

Write raw samples to a file.

Parameters:
  • filename (string) – name of file in which to write the samples
  • labels (list of strings) – names of relevant variables (used to write a header)
  • samples (array-like) – samples for which statistical properties are calculated
AIMS.write_statistics(filename, labels, samples)[source]

Write statistical properties based on a sequence of realisations to a file. The results include:

  • average values for each variable (statistical mean)
  • error bars for each variable (standard mean deviation)
  • correlation matrix between the different variables
Parameters:
  • filename (string) – name of file in which to write the statistical properties
  • labels (list of strings) – names of relevant variables
  • samples (np.array) – samples for which statistical properties are calculated