API Reference

This is the API reference for the rvspecfit package.

class rvspecfit.vel_fit.ParamMapper(specParams, paramDict0, fixParam, vsiniMapper, fitVsini=True)

This class constructs the dictionary with human readable parameters out of the vector, as well as taking into account which params are fixed

forward(p0)

Convert a parameter vector into param dictionary Argument order in the vector velocity [vsini] stellar parameters vsini is optional stellar parameters are in the order given by specParams and the ones with fixParam are excluded

Parameters

p0: array

Vector with parameter values

Returns

ret: dictionary

The dict with all the parameter values

rvspecfit.vel_fit.chisq_func(p, args)

The function computes the chi-square of the fit This function is used for minimization

Parameters

p: array_like

Vector of parameters to be fitted (velocity, vsini, stellar params)

args: dict

Dictionary containing specdata, config, and other fitting arguments

Returns

chisq: float

Chi-square value for the given parameters

rvspecfit.vel_fit.firstguess(specdata, options=None, config=None, resolParams=None, vsinigrid=(None, 10, 100), paramsgrid=None)

Compute the starting point template parameters and radial velocity by brute force looping over small grid of templates and a grid of radial velocities. This can be useful before starting a maximum-likelihood fit using vel_fit.process() or for initializing MCMC samples.

Parameters

specdata: list of Specdata

Spectroscopic dataset

options: dict

Optional dictionary of options

config: dict

Optional dictionary config

resolParams: dict

dictionary indexed by the spectroscopic setup with resolution matrices

paramsgrid: dictionary

(optional) dictionary of template parameters to iterate over The default value is

paramsgrid = {
    'logg': [1, 2, 3, 4, 5],
    'teff': [3000, 5000, 8000, 10000],
    'feh': [-2, -1, 0],
    'alpha': [0]
}
vsinigrid: tuple

(optional) list of vsinis to consider

Returns

bestpar: dict

Dictionary of best parameters

rvspecfit.vel_fit.get_hess_inv(param_names)

The inverse hessian is approximately the errors^2 Here we set it up.

rvspecfit.vel_fit.hess_func(p, pdict, args)

The function computes the 0.5*chi-square and takes as input vector of parameters (not transformed) pdict is for fixedparameters specParams is the list of names of parameters that we are varying

rvspecfit.vel_fit.process(specdata, paramDict0, fixParam=None, options=None, config=None, resolParams=None, priors=None)

Process spectra by doing maximum likelihood fit to spectral data

Parameters

specdata: list of SpecData

List of spectroscopic datasets to fit

paramDict0: dict

Dictionary of parameters to start from

fixParam: tuple

Tuple of parameters that will be fixed

options: dict

Dictionary of options

config: dict

Configuration dictionary

priors: dict (optional)

Extra dictionary with Normal priors on paramaters I.e. {‘teff’:(5000,10)} for N(5000,10) prior on the effective temperature

resolParams: dictionary

dictionary of matrices for the resolution of current spectrum.

Returns

ret: dict

Dictionary of parameters Keys are ‘yfit’ – the list of best-fit models ‘vel’, ‘vel_kurtosis’, ‘vel_skewness’, ‘vel_err’ velocity and its constraints ‘param’ – parameter values ‘param_err’ – parameter uncertaintoes ‘chisq’ – -2*log(likelihood) value (does not equal to chisq, because of marginalization) ‘chisq_array’ – array of proper chi-squares of the mode

Examples

>>> ret = process(specdata, {'logg':10, 'teff':30, 'alpha':0, 'feh':-1,
    'vsini':0}, fixParam = ('feh','vsini'),
            config=config, resolParams = None)
class rvspecfit.spec_fit.LRUDict(N)

Simple implementation of LRU dictionary

class rvspecfit.spec_fit.SpecData(name, lam, spec, espec, badmask=None, resolution=None, dtype=<class 'numpy.float64'>)

Class describing a single spectrocopic dataset

rvspecfit.spec_fit.compute_vsini_kernel(R, eps=0.6)

Computes the discrete convolution kernel for a given broadening R (in pixels).

The signal is assumed to be piecewise linear (triangular basis). We calculate the weights W_k = Integral[ Lambda(k - R*x) * K(x) ] dx

Parameters

Rfloat

Broadening width in pixels: (vsini / c) / step_log_lambda

epsfloat

Limb darkening coefficient.

Returns

kernelnumpy array

The discrete kernel weights.

rvspecfit.spec_fit.construct_resol_mat(lam, resol=None, width=None)

Construct a sparse resolution matrix from a resolution number R

Parameters

lam: numpy

The wavelength vector

resol: real/numpy

The resolution value (R=lambda/delta lambda)

width: real

The Gaussian width of the kernel in angstrom (cannot be specified together with resol)

Returns

mat: scipy.sparse matrix

The matrix describing the resolution convolution operation

rvspecfit.spec_fit.convolve_resol(spec, resol_matrix)

Convolve the spectrum with the resolution matrix

Parameters

spec: numpy

The spectrum array

resol_matrix: ResolMatrix object

The resolution matrix object

Returns

spec: numpy

The spectrum array

rvspecfit.spec_fit.convolve_vsini(lam_templ, templ, vsini, eps=0.6)

Convolve the spectrum with the stellar rotation velocity kernel using analytic integration of the piecewise linear spectrum.

Robust for both low (sub-pixel) and high vsini.

Parameters

lam_templ: numpy

The wavelength vector (MUST be spaced logarithmically)

templ: numpy

The spectrum vector

vsini: float

The Vsini velocity in km/s

eps: float

Limb darkening coefficient (default 0.6)

Returns

spec: numpy

The convolved spectrum

rvspecfit.spec_fit.evalRV(interpol, vel, lams)

Evaluate the spectrum interpolator at a given velocity and given wavelengths

Parameters

interpol: scipy.interpolate object

Template interpolator

vel: real

Radial velocity

lams: numpy

Wavelength array

Returns

spec: numpy

Evaluated spectrum

rvspecfit.spec_fit.find_best(specdata, vel_grid, params_list, rot_params=None, resol_params=None, options=None, config=None, quadratic=True)

Find the best fit template and velocity from a grid

Parameters

specdata: list of SpecData

Spectroscopic dataset

vel_grid: ndarray

Array of velocities

params_list: list

Array of parameters to explore

rot_params: tuple

Parameters of rotation (or None)

resol_params: tuple

Parameters of resolution convolution (or None)

options: dict

Dictionary of options

config: dict

Dictionary of the configuration

quadratic: bool

if True we try to do velocity quadratic interpolation near the max

Returns

ret: dict

Dictionary with measured parameters. The keys are best_vel – velocity vel_err – velocity error best_param – best param kurtosis – kurtosis of velocity distribution skewness – skewness of velocity distribution probs – vector of ‘posterior’ probabilities over the grid

rvspecfit.spec_fit.getCurTempl(spec_setup, atm_param, rot_params, config)

Get the spectrum in the given setup with given atmospheric parameters and given config

Parameters

spec_setup: string

The name of the spectroscopic setup

atm_param: tuple

The atmospheric parameters

rot_params: tuple

The parameters of stellar rotation models (could be None)

config: dict

The configuration dictionary

Returns

outside: float

Flag indicating if parameters are outside the template grid

lam: numpy array

Wavelength array of the template

spec: numpy array

The template spectrum

templ_tag: int

Random tag for caching purposes

log_step: bool

Whether the wavelength grid is logarithmic

rvspecfit.spec_fit.getRVInterpol(lam_templ, templ, log_step=True)

Produce the spectrum interpolator to evaluate the spectrum at arbitrary wavelengths

Parameters

lam_templ: numpy

Wavelength array

templ: numpy

Spectral array

Returns

interpol: scipy.interpolate object

The object that can be used to evaluate template at any wavelength

rvspecfit.spec_fit.get_basis(specdata, npoly, rbf=True)

Get the precomputed polynomials for the continuum for a given specdata

Parameters

specdata: SpecData objects

The spectroscopic dataset objects

npoly: integer

The degree of polynomial to use

rbf: bool

Use the RBF basis instead of monomial basis

Returns

polys: numpy array(npolys, Nwave)

The array of continuum polynomials

rvspecfit.spec_fit.get_chisq(specdata, vel, atm_params, rot_params=None, resol_params=None, options=None, config=None, cache=None, full_output=False, fast_interp=False, espec_systematic=None, outside_penalty=True)
Find the chi-square of the dataset at a given velocity

atmospheric parameters, rotation parameters and resolution parameters

Parameters

specdata: spec_fit.SpecData

The object with the data to be fitted

vel: real

The radial velocity

atm_params: tuple

The tuple with parameters of the star

rot_parameters: tuple

The tuple with parameters of the rotation (can be None)

resol_parameters: dictionary

The dictionary with resollution matrices ResolMatrix (can be None) The keys are names of spectral configurations

options: dict

The dictionary with fitting options (such as npoly for the degree of the polynomial)

config: dict (optional)

The configuration objection

fast_interp: bool

If true, use the nearest neighbor interpolation

cache: dict (optional)

The cache object, to preserve info between calls

full_output: bool

If full_output is set more info is returned

espec_systematic: dict or float

This will be added in quadrature to the error vector when computing logl. If it is a dict it must be indexed by the spec setup otherwise

this constant will be used for all spectra.

outside_penalty: bool

if true the chi^2 will be penalize for being outside of the grid

Returns

ret: float or dictionary

If full_output is False, ret is float = -2*log(L) of the whole data If full_output is True ret is a dictionary with the following keys chisq – this is the -2*log(L) of the whole dataset logl – this is the log(L) of the whole dataset chisq_array – this is the array of chi-squares (proper ones) for each of the fited spectra redchisq_array – this is the array of reduced chi-squares models – array of best fit models raw_models – array of models not corrected by the polynomial

rvspecfit.spec_fit.get_chisq0(spec, templ, polys, get_coeffs=False, espec=None)

Get the chi-square values for the vector of velocities and grid of templates after marginalizing over continuum If espec is not provided we assume data and template are already normalized Importantly the returned chi-square is not the true chi-square, but instead the -2*log(L)

Parameters

spec: numpy

Spectrum array

templ: numpy

The template

polys: numpy

The continuum polynomials

get_coeffs: boolean (optional)

If true return the coefficients of polynomials

espec: numpy (optional)

If specified, this is the error vector. If not specified, then it is assumed that spectrum and template are already divided by the uncertainty

Returns

chisq: real

Chi-square of the data

coeffs: numpy

The polynomial coefficients (optional)

rvspecfit.spec_fit.get_chisq_continuum(specdata, options=None)

Fit the spectrum with continuum only

Parameters

specdata: list of Specdata

Input spectra

options: dict

Dictionary of options (npoly option is required)

Returns

ret: list

Array of chi-squares

rvspecfit.spec_fit.get_poly_basis(lam, npoly, rbf=True)

get polynomials for the grid of wavelength if rbf is equal true then the first 3 terms will be still polynomials, the rest will be Gaussian rbf

class rvspecfit.fitter_ccf.CCFCache

Singleton caching CCF information

rvspecfit.fitter_ccf.fit(specdata, config)

Process the data by doing cross-correlation with templates

Parameters

specdata: list of SpecData objects

The list of data that needs to be fitted from different spectral setups.

config: dict

The configuration dictionary

Returns

results: dict

The dictionary with results such as best template parameters, best velocity, best vsini.

rvspecfit.fitter_ccf.get_ccf_info(spec_setup, config)

Returns the CCF info from the pickled file for a given spectroscopic setup

Parameters

spec_setup: string

The spectroscopic setup needed

config: dict

The dictionary with the config

Returns

d: dict

The dictionary with the CCF Information as saved by the make_ccf code

rvspecfit.make_ccf.ccf_executor(spec_setup, ccfconf, prefix=None, oprefix=None, every=10, vsinis=None, revision='', nthreads=1, cmdline='')

Prepare the FFT transformations for the CCF

Parameters

spec_setup: string

The name of the spectroscopic spec_setup

ccfconf: CCFConfig

The CCF configuration object

prefix: string

The input directory where the templates are located

oprefix: string

The output directory

every: integer (optional)

Produce FFTs of every N-th spectrum

vsinis: list (optional)

Produce FFTS of the templates with Vsini from the list. Could be None (it means no rotation will be added)

revision: str (optional)

The revision of the files/run that will be tagged in the pickle file

Returns

Nothing

rvspecfit.make_ccf.get_ccf_config(logl0=None, logl1=None, npoints=None, splinestep=1000, maxcontpts=20)

Configure the cross-correlation

Parameters

logl0: float

The natural logarithm of the wavelength of the beginning of the CCF

logl1: float

The natural logarithm of the wavelength of the end of the CCF

npoints: integer

The number of points in the cross correlation functions

splinestep: float, optional

The stepsize in km/s that determines the smoothness of the continuum fit

maxcontpts: integer, optional

The maximum number of points used for the continuum fit

rvspecfit.make_ccf.get_continuum(lam0, spec0, espec0, ccfconf=None)

Determine the continuum of the spectrum by fitting a spline

Parameters

lam0: numpy array

The wavelength vector

spec0: numpy array

The spectral vector

espec0: numpy array

The vector of spectral uncertainties

ccfconf: CCFConfig object

The CCF configuration object

Returns

cont: numpy array

The continuum vector

rvspecfit.make_ccf.interleave_bits(X)

Take the 2d array with shape nsamp,ndim with values between 0 and 1 and convert into integer z-curve number (Morton curve)

rvspecfit.make_ccf.interp_masker(lam, spec, badmask)

Fill the gaps spectrum by interpolating across a badmask. The gaps are filled by linear interpolation. The edges are just using the value of the closest valid pixel.

Parameters

lam: numpy array

The array of wavelengths of pixels

spec: numpy array

The spectrum array

badmask: boolean array

The array identifying bad pixels

Returns

spec: numpy array

The array with bad pixels interpolated away

rvspecfit.make_ccf.preprocess_data(lam, spec0, espec, ccfconf=None, badmask=None, maxerr=10)

Preprocess data in the same manner as the template spectra, normalize by the continuum, apodize and pad

Parameters

lam: numpy array

The wavelength vector

spec0: numpy array

The input spectrum vector

espec0: numpy array

The error-vector of the spectrum

ccfconf: CCFConfig object

The CCF configuration

badmask: Numpy array(boolean), optional

The optional mask for the CCF

maxerr: integer

The maximum value of error to be masked in units of median(error)

Returns

cap_spec: numpy array

The processed apodized/normalized/padded spectrum

rvspecfit.make_ccf.preprocess_model(logl, lammodel, model0, vsini=None, ccfconf=None)

Take the input template model and return prepared for FFT vectors. That includes padding, apodizing and normalizing by continuum

Parameters

logl: numpy array

The array of log wavelengths on which we want to outputted spectra

lammodel: numpy array

The wavelength array of the model

model0: numpy array

The initial model spectrum

vsini: float, optional

The stellar rotation, vsini

ccfconf: CCFConfig object, required

The CCF configuration object

Returns

xlogl: Numpy array

The log-wavelength of the resulting processed model0

cpa_model: Numpy array

The continuum normalized/subtracted, apodized and padded spectrum

rvspecfit.make_ccf.preprocess_model_list(lammodels, models, params, ccfconf, vsinis=None, nthreads=1)

Apply preprocessing to the array of models

Parameters

lammodels: numpy array

The array of wavelengths of the models (assumed to be the same for all models)

models: numpy array

The 2D array of modell with the shape [number_of_models, len_of_model]

params: numpy array

The 2D array of template parameters (i.e. stellar atmospheric parameters) with the shape [number_of_models,length_of_parameter_vector]

ccfconf: CCFConfig object

CCF configuration

vsinis: list of floats

The list of possible Vsini values to convolve model spectra with Could be None

Returns

ret: tuple

FILL/CHECK ME 1) log wavelenghts 2) processed spectra, 3) spectral params 4) list of vsini

rvspecfit.make_interpol.extract_spectrum(param, dbfile, prefix, wavefile, normalize=True, log_spec=True)

Extract a spectrum of a given parameters then apply the resolution smearing and divide by the continuum

Parameters

param: dict

The dictionary of key value pairs of parameters

dbfile: string

Path to the sqlite database

prefix: string

Prefix to the data files

wavefile: string

Path to the file with wavelengths

normalize: boolean

Normalize the spectrum by a linear continuum

log_spec: boolean

If True, take the logarithm of the spectrum

Returns

spec: numpy array

The processed spectrum (logarithmic if log_spec=True)

lognorm: float

The logarithm of the normalization factor

rvspecfit.make_interpol.get_line_continuum(lam, spec)

Determine the extremely simple linear in log continuum to remove away continuum trends in templates

Parameters

lam: numpy array

Wavelength vector

spec: numpy array

spectrum

Returns

cont: numpy array

Continuum model

rvspecfit.make_interpol.process_all(setupInfo, parnames=('teff', 'logg', 'feh', 'alpha'), dbfile='/tmp/files.db', oprefix='psavs/', prefix=None, wavefile=None, air=False, resolution0=None, normalize=True, revision='', cmdline='', nthreads=8, log_parameters=None)

Process the whole library of spectra and prepare the pickle file with arrays of convolved spectra, wavelength arrays, transformed parameters

Parameters

setupInfo: string

The name of spectral configuration

parnames: list of strings

The parameter names of spectra

log_parameters: integer positions of parameters that needs

to be log10() for interpolation. I.e. if the first parameter si teff and we want to perform interpolation in log(teff) space this needs to be [0]

air: boolean

Transform from vacuum to air

rvspecfit.make_nd.execute(spec_setup, prefix=None, regular=False, perturb=True, revision='', cmdline='')

Prepare the triangulation objects for the set of spectral data for a given spec_setup.

Parameters

spec_setup: string

The spectroscopic configuration

prefix: string

The prefix where the data are located and where the triangulation will be stored

perturb: boolean

Boolean flag whether to perturb a little bit the points before doing a triangulation. This prevents issues with degenerate vertices and stability of triangulation. Without perturbation find_simplex for example may revert to brute force search.

Returns

None

rvspecfit.make_nd.getedgevertices(vec)

Given a set of n-dimentional points, return vertices of an n-dimensional cube that fully encompass/surrounds the data, This is sort of the envelope around the data

Parameters

vec: numpy (Ndim, Npts)

The array of input points

Returns

vec: numpy (Ndim, Nretpts)

The returned array of surrounding points

class rvspecfit.read_grid.LogParamMapper(log_ids)

Class used to map stellar atmospheric parameters into more suitable space used for interpolation

forward(vec)

Map atmospheric parameters into parameters used in the grid for interpolation. That includes logarithming the teff

Parameters

vec: numpy array

The vector of atmospheric parameters i.e. Teff, logg, feh, alpha

Returns

ret: numpy array

The vector of transformed parameters used in interpolation

inverse(vec)

Map transformed parameters used in the grid for interpolation back into the atmospheric parameters. That includes exponentiating the log10(teff)

Parameters

vec: numpy array

The vector of transformed atmospheric parameters log(Teff), logg, feh, alpha

Returns

ret: numpy array

The vector of original atmospheric parameters.

rvspecfit.read_grid.gau_integrator(A, B, x1, x2, l1, l2, s)

This computes the integral of (Ax+B)/sqrt(2pi)/s*exp(-1/2*(x-y)^2/s^2 for x=x1..x2 y=l1..l2

Here is the Mathematica code FortranForm[ Simplify[Integrate[(A*x + B)/Sqrt[2*Pi]/s* Exp[-1/2*(x - y)^2/s^2], {x, x1, x2}, {y, l1, l2}]]]

Parameters A: ndarray B: ndarray x1: ndarray x2: ndarray l1: ndarray l2: ndarray s: ndarray

Returns

ret: float

Integral of the function

rvspecfit.read_grid.get_spec(params, dbfile=None, prefix=None, wavefile=None)

Returns individual spectra for a given spectral parameters

Parameters

params: dict

The dictionary of parameters

dbfile: string

The pathname to the database sqlite file of templates

prefix: string

The prefix path to templates

wavefile: string

The filename of fits file with the wavelength vector

Returns

lam: ndarray

1-D array of wavelength spec: ndarray 1-D array of spectrum

Examples

> lam,spec=read_grid.get_spec(dict(logg=1,teff=5250,feh=-1,alpha=0.4))

rvspecfit.read_grid.make_rebinner(lam00, lam, resolution_function, resolution0=None, toair=True)

Make the sparse matrix that convolves a given spectrum to a given resolution and new wavelength grid

Parameters

lam00: array

The input wavelength grid of the templates

lam: array

The desired wavelength grid of the output

resolution_function: function

The function that when called with the wavelength as an argument will return the resolution of the desired spectra (R=l/dl) I.e. it could be just lambda x: 5000

toair: bool

Convert the input spectra into the air frame

resolution0: float

The resolution of input templates

Returns

The sparse matrix to perform interpolation

rvspecfit.read_grid.makedb(prefix='/PHOENIX-ACES-AGSS-COND-2011/', dbfile='files.db', keywords=None, mask=None, extra_params=None, name_metallicity=None, name_alpha=None)

Create an sqlite database of the templates

Parameters

prefix: str

The path to PHOENIX

dbfile: str

The output file with the sqlite DB

keywords: dict

The dictionary with the map of teff,logg,feh,alpha to keyword names in the headers

mask: string

The string how to identify spectral files, i.e. */*fits

extra_params: dict or None

The dictionary of variable name vs FITS name to read from spectral files

rvspecfit.read_grid.pix_integrator(x1, x2, l1, l2, s)

Integrate the flux within the pixel given the LSF We assume that the flux before LSF convolution is given by linear interpolation from x1, x2. The l1,l2 are scalar edges of the pixel within which we want to compute the flux. s is the sigma of the LSF The function returns two values of weights for y1,y2 where y1,y2 are the values of the non-convolved spectra at x1,x2

Parameters

x1: ndarray x2: ndarray l1: ndarray l2: ndarray s: ndarray

Returns

ret: tuple of ndarray

Two weights for the integral (c1,c2) i.e. if the fluxes at x1,x2 are f1, f2 the integral is c1*f1+c2*f2

rvspecfit.read_grid.rebin(lam0, spec0, newlam, resolution)

Rebin a given spectrum lam0, spec0 to the new wavelength and new resolution

Parameters

lam0: ndarray

1d numpy array with wavelengths of the template pixels

spec0: ndarray

1d numpy array with the template spectrum

newlam: ndarray

1d array with the wavelengths of the output spectrum

resolution: float

Resolution lam/dlam

Returns

spec: ndarray

Rebinned spectrum

Examples

>>> lam,spec=read_grid.get_spec(dict(logg=1,teff=5250,feh=-1,alpha=0.4))
>>> newlam = np.linspace(4000,9000,4000)
>>> newspec=read_grid.rebin(lam, spec, newlam, 1800)
rvspecfit.regularize_grid.converter(path, opath, smooth=0, min_feh=None, max_feh=None, step_feh=None, min_alpha=None, max_alpha=None, step_alpha=None, cmdline='')

Read the input spectrum pickle file and convert it into the file with gaps filled and smaller step sizes

rvspecfit.regularize_grid.findbestoverlaps(x, intervals)

find the interval where the value is closer to the center I.e. given intervals [0,10],[1,11],[2,12],[3,13],[4,14],[5,15],[6,16] and value of 8 it return [3,13]

rvspecfit.serializer.load_dict_from_hdf5(filename)

Loads the dictionary from an HDF5 file.

rvspecfit.serializer.recursively_load_dict_contents_from_group(h5file, path)

Recursively loads dictionary contents from HDF5 groups and datasets.

rvspecfit.serializer.recursively_save_dict_contents_to_group(h5file, path, dic, allow_pickle=False)

Recursively saves dictionary contents to HDF5 groups and datasets.

rvspecfit.serializer.save_dict_to_hdf5(filename, dic, allow_pickle=False)

Saves the provided dictionary to an HDF5 file.

rvspecfit.serializer.verify_data(original, loaded, path='/')

Recursively verify that two dictionaries are identical in both value and type.

class rvspecfit.spec_inter.SpecInterpolator(name, interper, extraper, lam, mapper, parnames, revision='', filename='', creation_soft_version='', log_step=None)

Spectrum interpolator object

eval(param0)

Evaluate the spectrum at a given parameter

outsideFlag(param0)

Check if the point is outside the interpolation grid

Parameters

param0: tuple

parameter vector

Returns

retfloat

> 0 if point outside the grid

rvspecfit.spec_inter.getInterpolator(HR, config, warmup_cache=False, cache=None)

Return the spectrum interpolation object for a given instrument setup HR and config. This function also checks the cache

Parameters

HR: string

Spectral configuration

config: dict

Configuration

cache: dict or None

Dictionary like object with the cache. If None, internal cache is used instead.

warmup_cache: bool

If True we read the whole file to warm up the OS cache

Returns

ret: SpecInterpolator

The spectral interpolator

rvspecfit.spec_inter.getSpecParams(setup, config)

Return the ordered list of spectral parameters of a given spectroscopic setup

Parameters

setup: str

Spectral configuration

config: dict

Configuration dictionary

Returns

ret: list

List of parameters names

class rvspecfit.spec_inter.interp_cache

Singleton object caching the interpolators

rvspecfit.utils.freezeDict(d)

Take the input object and if it is a dictionary, freeze it (i.e. return frozendict) If not, do nothing

Parameters

d: dict

Input dictionary

Returns

d: frozendict

Frozen input dictionary

rvspecfit.utils.get_default_config()

Create a default parameter config dictionary

Returns

ret: dict

Dictionary with config params

rvspecfit.utils.read_config(fname=None, override_options=None)

Read the configuration file and return the frozendict with it

Parameters

fname: string, optional

The path to the configuration file. If not given config.yaml in the current directory is used

override_options: dictionary, optional

Update the options

Returns

config: frozendict

The dictionary with the configuration from a file