Package ivs :: Package sed :: Module fit
[hide private]
[frames] | no frames]

Module fit

source code

Fit SED models to observed data using various approaches.

Functions [hide private]
    PCA functions
N-D numpy array,(1Darray,1Darray,1Darray)
get_PCA_grid(colors, res=3, teffrange=(-np.inf,np.inf), loggrange=(-np.inf,np.inf), ebvrange=(-np.inf,np.inf), **kwargs)
Transform a flux grid to a colour grid.
source code
N-D array, N-D array, (1D array, 1D array)
get_PCA(A)
Find the principal components of a color grid.
source code
n-tuple
calibrate_PCA(T, pars, function='linear')
Define the interpretations of the principal components.
source code
list of lists
get_PCA_parameters(obsT, calib, P, means, stds, e_obsT=None, mc=None)
Derive fundamental parameters of a sample of observations given a PCA.
source code
    Grid search
float,float,float
stat_chi2(meas, e_meas, colors, syn, full_output=False, **kwargs)
Calculate Chi2 and compute angular diameter.
source code
 
generate_grid_single_pix(photbands, points=None, clear_memory=True, **kwargs)
Generate a grid of parameters.
source code
 
generate_grid_pix(photbands, points=None, clear_memory=False, **kwargs)
Generate a grid of parameters for 1 or more stars.
source code
record array
generate_grid_single(photbands, teffrange=(-inf,inf), loggrange=(-inf,inf), ebvrange=(-inf,inf), zrange=(-inf,inf), points=None, res=None, clear_memory=True, **kwargs)
Generate grid points at which to fit an interpolated grid of SEDs.
source code
record array
generate_grid(photbands, teffrange=((-inf,inf),(-inf,inf)), loggrange=((-inf,inf),(-inf,inf)), ebvrange=(-inf,inf), zrange=((-inf,inf),(-inf,inf)), radiusrange=((1,1),(0.1,10.)), grids=None, points=None, res=None, clear_memory=False, type='single', primary_hottest=False, **kwargs)
Generate grid points at which to fit an interpolated grid of SEDs.
source code
    Fitting: grid search
array
igrid_search_pix(meas, e_meas, photbands, constraints={}, **kwargs)
Run over gridpoints and evaluate model model_func via stat_func.
source code
4/5X1d array
igrid_search(meas, e_meas, photbands, *args, **kwargs)
Run over gridpoints and evaluate model model_func via stat_func.
source code
    Fitting: minimizer
 
create_parameter_dict(**pars) source code
 
_get_info_from_minimizer(minimizers, startpars, photbands, meas, e_meas, **fitkws) source code
 
_iminimize_model(varlist, x, *args, **kws) source code
 
_iminimize_residuals(synth, meas, weights=None, **kwargs) source code
 
iminimize(meas, e_meas, photbands, points=None, return_minimizer=False, **kwargs)
minimizer based on the sigproc.fit lmfit minimizer.
source code
 
calculate_iminimize_CI(meas, e_meas, photbands, CI_limit=0.66, **kwargs)
Calculate the confidence intervalls for every parameter.
source code
 
calculate_iminimize_CI2D(meas, e_meas, photbands, xpar, ypar, df=None, **kwargs)
Calculate a 2D confidence grid for the given parameters.
source code
 
iminimize2(meas, e_meas, photbands, *args, **kwargs) source code
Variables [hide private]
  logger = logging.getLogger("SED.FIT")
Function Details [hide private]

get_PCA_grid(colors, res=3, teffrange=(-np.inf,np.inf), loggrange=(-np.inf,np.inf), ebvrange=(-np.inf,np.inf), **kwargs)

source code 

Transform a flux grid to a colour grid.

The resulting grid is actually log10(flux ratio) instead of flux ratio! This works better in the PCA.

Extra keyword arguments can be used to set the atmosphere model.

Parameters:
  • colors (list of strings) - list of desired colours (['GENEVA.U-B','STROMGREN.M1'...])
  • res (integer) - resolution in E(B-V)
  • teffrange (tuple) - range of Teffs to use
  • loggrange (tuple) - range of loggs to use
  • ebvrange (tuple) - range of E(B-V)s to use
Returns: N-D numpy array,(1Darray,1Darray,1Darray)
log10(colorgrid), (teffs,loggs,ebvs)

get_PCA(A)

source code 

Find the principal components of a color grid.

Parameters:
  • A - color grid obtained via get_PCA_grid. @type A. numpy N-D array
Returns: N-D array, N-D array, (1D array, 1D array)
PCA loadings, PCA scores, (column means, column standard devs)

calibrate_PCA(T, pars, function='linear')

source code 

Define the interpretations of the principal components.

T are the PCA scores, pars a list of axes (e.g. teff, logg and ebv).

Parameters:
  • T (N-D array) - PCA scores
  • pars (list of 1D arrays) - real-life axes
Returns: n-tuple
Rbf interpolating function

get_PCA_parameters(obsT, calib, P, means, stds, e_obsT=None, mc=None)

source code 

Derive fundamental parameters of a sample of observations given a PCA.

Monte Carlo simulations are only available when you give 1 target.

Parameters:
  • obsT (list of lists) - observed colours of the sample
  • calib (list of Rbf functions) - PCA calibration obtained via calibrate_PCA
  • P (N-D array) - PCA loadings
  • means (numpy array) - means of PCA grid
  • stds (numpy array) - standard deviations of PCA grid
  • mc (integer) - number of MonteCarlo simulations
Returns: list of lists
list of fundamental parameters

stat_chi2(meas, e_meas, colors, syn, full_output=False, **kwargs)

source code 

Calculate Chi2 and compute angular diameter.

Colors and absolute fluxes are used to compute the Chi2, only absolute fluxes are used to compute angular diameter. If no absolute fluxes are given, the angular diameter is set to 0.

Parameters:
  • meas (1D array) - array of measurements
  • e_meas (1D array) - array containing measurements errors
  • colors (1D boolean array) - boolean array separating colors (True) from absolute fluxes (False)
  • syn (1D array) - synthetic fluxes and colors
  • full_output (boolean) - set to True if you want individual chisq
Returns: float,float,float
chi-square, scale, e_scale

generate_grid_pix(photbands, points=None, clear_memory=False, **kwargs)

source code 

Generate a grid of parameters for 1 or more stars. Based on the generate_grid_single_pix method. The radius of the components is based on the masses if given, otherwise on the radrange arguments. If masses are given the radrange arguments are ignored, meaning that the returned radius can be outside those limits. If only 1 component is provided no radius will be returned.

This function can handle multiple components in the same way as a single component. parameter ranges are provided as <parname><component>range=(...,...). fx: teffrange = (5000, 10000), loggrange = (3.5, 4.5), teff2range = (10000, 20000), logg2range = (5.5, 6.0)

For the first (or only) component both no component number (teffrange) or 1 as component number (teff1range) can be used.

Returns a dictionary with for each parameter that has a range provided, an array of values. In the case of multiple components, the radius will be returned even is no radius ranges are provided.

generate_grid_single(photbands, teffrange=(-inf,inf), loggrange=(-inf,inf), ebvrange=(-inf,inf), zrange=(-inf,inf), points=None, res=None, clear_memory=True, **kwargs)

source code 

Generate grid points at which to fit an interpolated grid of SEDs.

If points=None, the points are chosen on the predefined grid points. Otherwise, points grid points will be generated, uniformly distributed between the ranges defined by teffrange, loggrange and ebvrange. If you set the resolution to 2, one out of every two points will be selected.

Extra keyword arguments can be used to give more details on the atmosphere models to use.

Colors are automatically detected.

You can fix one parameter e.g. via setting teffrange=(10000,10000).

>>> photbands = ['GENEVA.G','GENEVA.B-V']

Start a figure:

>>> p = pl.figure()
>>> rows,cols = 2,4

On the grid points, but only one in every 100 points (otherwise we have over a million points):

>>> teffs,loggs,ebvs,zs = generate_grid(photbands,res=100)
>>> p = pl.subplot(rows,cols,1)
>>> p = pl.scatter(teffs,loggs,c=ebvs,s=(zs+5)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlim(pl.xlim()[::-1]);p = pl.ylim(pl.ylim()[::-1])
>>> p = pl.xlabel('Teff');p = pl.ylabel('Logg')
>>> p = pl.subplot(rows,cols,1+cols)
>>> p = pl.scatter(ebvs,zs,c=teffs,s=(loggs+2)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlabel('E(B-V)');p = pl.ylabel('Z')

Randomly distributed over the grid's ranges:

>>> teffs,loggs,ebvs,zs = generate_grid(photbands,points=10000)
>>> p = pl.subplot(rows,cols,2)
>>> p = pl.scatter(teffs,loggs,c=ebvs,s=(zs+5)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlim(pl.xlim()[::-1]);p = pl.ylim(pl.ylim()[::-1])
>>> p = pl.xlabel('Teff');p = pl.ylabel('Logg')
>>> p = pl.subplot(rows,cols,2+cols)
>>> p = pl.scatter(ebvs,zs,c=teffs,s=(loggs+2)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlabel('E(B-V)');p = pl.ylabel('Z')

Confined to a small area in the grid's range:

>>> teffs,loggs,ebvs,zs = generate_grid(photbands,teffrange=(8000,10000),loggrange=(4.1,4.2),zrange=(0,inf),ebvrange=(1.,2),points=10000)
>>> p = pl.subplot(rows,cols,3)
>>> p = pl.scatter(teffs,loggs,c=ebvs,s=(zs+5)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlim(pl.xlim()[::-1]);p = pl.ylim(pl.ylim()[::-1])
>>> p = pl.xlabel('Teff');p = pl.ylabel('Logg')
>>> p = pl.subplot(rows,cols,3+cols)
>>> p = pl.scatter(ebvs,zs,c=teffs,s=(loggs+2)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlabel('E(B-V)');p = pl.ylabel('Z')

Confined to a small area in the grid's range with some parameters fixed:

>>> teffs,loggs,ebvs,zs = generate_grid(photbands,teffrange=(8765,8765),loggrange=(4.1,4.2),zrange=(0,0),ebvrange=(1,2),points=10000)
>>> p = pl.subplot(rows,cols,4)
>>> p = pl.scatter(teffs,loggs,c=ebvs,s=(zs+5)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlim(pl.xlim()[::-1]);p = pl.ylim(pl.ylim()[::-1])
>>> p = pl.xlabel('Teff');p = pl.ylabel('Logg')
>>> p = pl.subplot(rows,cols,4+cols)
>>> p = pl.scatter(ebvs,zs,c=teffs,s=(loggs+2)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlabel('E(B-V)');p = pl.ylabel('Z')

]include figure]]ivs_sed_fit_grids.png]

Parameters:
  • photbands (list of strings) - a list of photometric passbands, corresponding each measurement
  • teffrange (2-tuple) - range of temperatures to use
  • loggrange (2-tuple) - range of surface gravities to use
  • ebvrange (2-tuple) - range of reddenings to use
  • points (int) - points to sample (when None, predefined grid points are used)
  • res (int) - resolution of the original grid (the higher, the coarser)
  • clear_memory (boolean) - flag to clear memory from previously loaded SED tables. If you set it to False, you can easily get an overloaded memory!
Returns: record array
record array containing the searched grid, chi-squares and scale factors

generate_grid(photbands, teffrange=((-inf,inf),(-inf,inf)), loggrange=((-inf,inf),(-inf,inf)), ebvrange=(-inf,inf), zrange=((-inf,inf),(-inf,inf)), radiusrange=((1,1),(0.1,10.)), grids=None, points=None, res=None, clear_memory=False, type='single', primary_hottest=False, **kwargs)

source code 

Generate grid points at which to fit an interpolated grid of SEDs.

If points=None, the points are chosen on the predefined grid points. Otherwise, points grid points will be generated, uniformly distributed between the ranges defined by teffrange, loggrange and ebvrange. If you set the resolution to 2, one out of every two points will be selected.

Setting primary_hottest makes sure that, when doing a binary fit, that the the effective temperature of the first star is hotter than the second.

Extra keyword arguments can be used to give more details on the atmosphere models to use.

Colors are automatically detected.

You can fix one parameter e.g. via setting teffrange=(10000,10000).

>>> photbands = ['GENEVA.G','GENEVA.B-V']

Start a figure:

>>> p = pl.figure()
>>> rows,cols = 2,4

On the grid points, but only one in every 100 points (otherwise we have over a million points):

>>> teffs,loggs,ebvs,zs = generate_grid(photbands,res=100)
>>> p = pl.subplot(rows,cols,1)
>>> p = pl.scatter(teffs,loggs,c=ebvs,s=(zs+5)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlim(pl.xlim()[::-1]);p = pl.ylim(pl.ylim()[::-1])
>>> p = pl.xlabel('Teff');p = pl.ylabel('Logg')
>>> p = pl.subplot(rows,cols,1+cols)
>>> p = pl.scatter(ebvs,zs,c=teffs,s=(loggs+2)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlabel('E(B-V)');p = pl.ylabel('Z')

Randomly distributed over the grid's ranges:

>>> teffs,loggs,ebvs,zs = generate_grid(photbands,points=10000)
>>> p = pl.subplot(rows,cols,2)
>>> p = pl.scatter(teffs,loggs,c=ebvs,s=(zs+5)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlim(pl.xlim()[::-1]);p = pl.ylim(pl.ylim()[::-1])
>>> p = pl.xlabel('Teff');p = pl.ylabel('Logg')
>>> p = pl.subplot(rows,cols,2+cols)
>>> p = pl.scatter(ebvs,zs,c=teffs,s=(loggs+2)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlabel('E(B-V)');p = pl.ylabel('Z')

Confined to a small area in the grid's range:

>>> teffs,loggs,ebvs,zs = generate_grid(photbands,teffrange=(8000,10000),loggrange=(4.1,4.2),zrange=(0,inf),ebvrange=(1.,2),points=10000)
>>> p = pl.subplot(rows,cols,3)
>>> p = pl.scatter(teffs,loggs,c=ebvs,s=(zs+5)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlim(pl.xlim()[::-1]);p = pl.ylim(pl.ylim()[::-1])
>>> p = pl.xlabel('Teff');p = pl.ylabel('Logg')
>>> p = pl.subplot(rows,cols,3+cols)
>>> p = pl.scatter(ebvs,zs,c=teffs,s=(loggs+2)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlabel('E(B-V)');p = pl.ylabel('Z')

Confined to a small area in the grid's range with some parameters fixed:

>>> teffs,loggs,ebvs,zs = generate_grid(photbands,teffrange=(8765,8765),loggrange=(4.1,4.2),zrange=(0,0),ebvrange=(1,2),points=10000)
>>> p = pl.subplot(rows,cols,4)
>>> p = pl.scatter(teffs,loggs,c=ebvs,s=(zs+5)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlim(pl.xlim()[::-1]);p = pl.ylim(pl.ylim()[::-1])
>>> p = pl.xlabel('Teff');p = pl.ylabel('Logg')
>>> p = pl.subplot(rows,cols,4+cols)
>>> p = pl.scatter(ebvs,zs,c=teffs,s=(loggs+2)*10,edgecolors='none',cmap=pl.cm.spectral)
>>> p = pl.xlabel('E(B-V)');p = pl.ylabel('Z')

]include figure]]ivs_sed_fit_grids.png]

Parameters:
  • photbands (list of strings) - a list of photometric passbands, corresponding each measurement
  • teffrange (2-tuple) - range of temperatures to use
  • loggrange (2-tuple) - range of surface gravities to use
  • ebvrange (2-tuple) - range of reddenings to use
  • points (int) - points to sample (when None, predefined grid points are used)
  • res (int) - resolution of the original grid (the higher, the coarser)
  • clear_memory (boolean) - flag to clear memory from previously loaded SED tables. If you set it to False, you can easily get an overloaded memory!
Returns: record array
record array containing the searched grid, chi-squares and scale factors

igrid_search_pix(meas, e_meas, photbands, constraints={}, **kwargs)

source code 

Run over gridpoints and evaluate model model_func via stat_func.

The measurements are defined via meas, e_meas, photbands, colors and should be 1d arrays of equal length. colors should be a boolean array, photbands should be a string array.

The grid points are defined via args. args should be a tuple of 1x? dimensional arrays of equal length. For single stars, this is typically effective temperatures, loggs, reddenings and metallicities. For multiple systems, (at least some of the) previously mentioned parameters are typically doubled, and radius ratios are added. Remember to specify the model_func to match single or multiple systems.

At each grid point, the pre-calculated photometry will be retrieved via the keyword model_func and compared to the measurements via the function definded via stat_func. This function should be of the same form as stat_chi2.

Extra arguments are passed to parallel_gridsearch for parallelization and to {model_func} for further specification of grids etc.

The index array is returned to trace the results after parallelization.

Parameters:
  • meas (1D numpy array of floats) - the measurements that have to be compared with the models
  • e_meas (1D numpy array of floats) - errors on the measurements
  • photbands (1D numpy array of strings) - names of the photometric passbands
  • model_func (function) - function to translate parameters to synthetic (model) data
  • stat_func (function) - function to evaluate the fit
Returns: array
(chi squares, scale factors, error on scale factors, absolute luminosities (R=1Rsol)

igrid_search(meas, e_meas, photbands, *args, **kwargs)

source code 

Run over gridpoints and evaluate model model_func via stat_func.

The measurements are defined via meas, e_meas, photbands, colors and should be 1d arrays of equal length. colors should be a boolean array, photbands should be a string array.

The grid points are defined via args. args should be a tuple of 1x? dimensional arrays of equal length. For single stars, this is typically effective temperatures, loggs, reddenings and metallicities. For multiple systems, (at least some of the) previously mentioned parameters are typically doubled, and radius ratios are added. Remember to specify the model_func to match single or multiple systems.

At each grid point, the pre-calculated photometry will be retrieved via the keyword model_func and compared to the measurements via the function definded via stat_func. This function should be of the same form as stat_chi2.

Extra arguments are passed to parallel_gridsearch for parallelization and to {model_func} for further specification of grids etc.

The index array is returned to trace the results after parallelization.

Parameters:
  • meas (1D numpy array of floats) - the measurements that have to be compared with the models
  • e_meas (1D numpy array of floats) - errors on the measurements
  • photbands (1D numpy array of strings) - names of the photometric passbands
  • model_func (function) - function to translate parameters to synthetic (model) data
  • stat_func (function) - function to evaluate the fit
Returns: 4/5X1d array
(chi squares, scale factors, error on scale factors, absolute luminosities (R=1Rsol), index
Decorators:
  • @parallel_gridsearch
  • @make_parallel

iminimize(meas, e_meas, photbands, points=None, return_minimizer=False, **kwargs)

source code 

minimizer based on the sigproc.fit lmfit minimizer. provide the observed data, the fitting model, residual function and parameter information about the variables, and this function will return the best fit parameters together with extra information about the fit.

if the fitkws keyword is supplied, this dict will be made available to the model_func (fit model) during the fitting process. The order of the parameters will also be made available as the 'pnames' keyword.

calculate_iminimize_CI(meas, e_meas, photbands, CI_limit=0.66, **kwargs)

source code 

Calculate the confidence intervalls for every parameter. When ci fails, the boundaries are returned as ci.

calculate_iminimize_CI2D(meas, e_meas, photbands, xpar, ypar, df=None, **kwargs)

source code 

Calculate a 2D confidence grid for the given parameters. You can provide limits on the parameters yourself, if none are provided, the function will take care of it, but it might crash if the limits are outside the model range. returns: x vals, y vals, ci values