full documentation to sedcreator¶
- class sedcreator.__init__.FitterContainer(models, master_dir='./', extinction_law='kmh', lambda_array=None, flux_array=None, err_flux_array=None, upper_limit_array=None, dist=None)¶
A class to store the results from the SedFluxer class
- get_model_info(keys=['mcore', 'sigma', 'mstar', 'theta_view'], tablename=None)¶
Gets the model information from the results of the SED fit. Automatically sorts the results by chisq. Columns units are given as astropy.units
- Parameters:
keys (str or array of str) – keys to be used to select the unique final table. It should be use either [‘mcore’,’sigma’,’mstar’] for the 432 physical models considering the best (by chisq) viewing angle or [‘mcore’,’sigma’,’mstar’,’theta_view’] for the 8640 models including the viewing angle Default is [‘mcore’,’sigma’,’mstar’,’theta_view’]
tablename (str, optional) – A path, or a Python file-like object. Note that fname is used verbatim, and there is no attempt to make the extension. Default is None. Note that one can choose the format of the table by changing the extension.
- Returns:
table
- Return type:
astropy.table with all the model information based on the given keys
- class sedcreator.__init__.FluxerContainer(data=None, flux_bkgsub=None, flux=None, fluc_error=None, bkg=None, central_coords=None, aper_rad=None, inner_annu=None, outer_annu=None, x_source=None, y_source=None, aper_rad_pixel=None, wcs_header=None, aperture=None, annulus_aperture=None, mask=None, flux_method=None)¶
A class to store the results from the SedFluxer class
- get_info()¶
Prints the results and useful information of the image based on the header and aperture photometry
- get_value()¶
Outputs the background-subtracted flux, the flux including the background, the flux error, and the background estimate
- Returns:
flux_bkgsub,flux,fluc_error,background
- Return type:
numpy array
- plot(figsize=(6, 4), cmap='gray', percent=100.0, stretch='log', colorbar=True, aperture_color='black', annulus_color='red', plot_mask=False, title=None, figname=None)¶
Plots the used image together with the aperture used for the flux and the annulus for the background subtraction.
- Parameters:
cmap (str) – Color map from matplolib to be used in the image. Default is ‘gray’
figsize (tuple) – specify the size of the figure. Default is (6,4)
stretch ({'linear', 'sqrt', 'power', 'log', 'asinh'}, optional) – The stretch function to apply to the image. The default is ‘log’.
percent (float, optional) – The percentage of the image values used to determine the pixel values of the minimum and maximum cut levels. The lower cut level will set at the (100 - percent) / 2 percentile, while the upper cut level will be set at the (100 + percent) / 2 percentile. The default is 100.0.
aperture_color (str) – Circular aperture color to show in the image. Default is ‘black’
annulus_color (str) – Annulus aperture color to show in the image. Default is ‘red’
plot_mask (bool, optional) – Plots the mask used in when getting the flux. Default is False
title (str, optional) – Set title for the image. Default is None
figname (str, optional) – A path, or a Python file-like object. Note that fname is used verbatim, and there is no attempt to make the extension. Default is None. Note that one can choose the format of the figure by changing the extension, e.g., figure.pdf would generate the figure in PDF format.
- class sedcreator.__init__.ModelPlotter(fittercontainer)¶
A class used to plot the SED model results
- plot2d(models, figsize=(10, 5), marker='s', markersize=50, cmap='rainbow_r', title=None, figname=None)¶
2D Plots of the SEDs results.
- Parameters:
models (astropy.table) – Table with the sed results.
figsize (tuple) – specify the size of the figure. Default is (10,5)
marker (str) – defines the marker style like matplotlib. Default is ‘s’
markersize (int) – defines the marker size. Default is 50
cmap (str) – defines the colormap like matplotlib of the 2D plot. Default is ‘rainbow_r’
title (str, optional) – Defines the title of the plot. Default is None.
figname (str, optional) – A path, or a Python file-like object. Note that fname is used verbatim, and there is no attempt to make the extension. Default is None. Note that one can choose the format of the figure by changing the extension, e.g., figure.pdf would generate the figure in PDF format.
- Return type:
plot the 2D results
- plot3d(models, figsize=(8, 8), markersize=200, cmap='rainbow_r', title=None, figname=None)¶
3D plot for the SED model results
- Parameters:
models (astropy.table) – Table with the sed results.
figsize (tuple) – specify the size of the figure. Default is (8,8)
markersize (int) – defines the marker size. Default is 200
cmap (str) – defines the colormap like matplotlib of the 2D plot. Default is ‘rainbow_r’
title (str, optional) – Defines the title of the plot. Default is None.
figname (str, optional) – A path, or a Python file-like object. Note that fname is used verbatim, and there is no attempt to make the extension. Default is None. Note that one can choose the format of the figure by changing the extension, e.g., figure.pdf would generate the figure in PDF format.
- Return type:
plot the 3D results
- plot_best_sed(models=None, figsize=(6, 4), xlim=[1.0, 1000.0], ylim=[1e-12, 1e-05], marker='k*', markersize=6, title=None, figname=None)¶
Plots the best SED model
- Parameters:
models (astropy.table) – Table with the sed results.
figsize (tuple) – specify the size of the figure. Default is (6,4)
xlim (array) – array to specify the limits in the x axis. Default is [1.0,1000.0]
ylim (array) – array to specify the limits in the y axis. Default is [1.0e-12,1.0e-5]
marker (str) – defines the marker style and color like matplotlib. Default is ‘k*’
markersize (int) – defines the marker size. Default is 6
title (str, optional) – Defines the title of the plot. Default is None.
figname (str, optional) – A path, or a Python file-like object. Note that fname is used verbatim, and there is no attempt to make the extension. Default is None. Note that one can choose the format of the figure by changing the extension, e.g., figure.pdf would generate the figure in PDF format.
- Return type:
plot of the best SED
- plot_multiple_seds(models=None, figsize=(6, 4), xlim=[1.0, 1000.0], ylim=[1e-12, 1e-05], marker='k*', markersize=6, cmap='rainbow_r', colorbar=True, title=None, figname=None)¶
Plots multiple SEDs.
- Parameters:
models (astropy.table) – Table with the sed results.
figsize (tuple) – specify the size of the figure. Default is (6,4)
xlim (array) – array to specify the limits in the x axis. Default is [1.0,1000.0]
ylim (array) – array to specify the limits in the y axis. Default is [1.0e-12,1.0e-5]
marker (str) – defines the marker style and color like matplotlib. Default is ‘k*’
markersize (int) – defines the marker size. Default is 6
cmap (str) – defines the colormap like matplotlib of the SEDs. Default is ‘rainbow_r’
colorbar (bool) – Turn on or off the colorbar of the plot. Default is True
title (str, optional) – Defines the title of the plot. Default is None.
figname (str, optional) – A path, or a Python file-like object. Note that fname is used verbatim, and there is no attempt to make the extension. Default is None. Note that one can choose the format of the figure by changing the extension, e.g., figure.pdf would generate the figure in PDF format.
- Return type:
plot of the multiple SEDs
- class sedcreator.__init__.PentagonPlot(fig, titles, labels, rect=None)¶
A class used to plot a pentagon plot
- plot¶
Creates the pentagon plot
- Type:
method
- class sedcreator.__init__.SedFitter(extc_law='kmh', lambda_array=None, flux_array=None, err_flux_array=None, upper_limit_array=None, filter_array=None)¶
A class used to fit the SED model grid
- add_filter(filter_name, instrument, lambda_array, response_array)¶
Adds an user defined filter to the database. It adds a txt file to the database with columns lambda_array and response_array. It also adds the fits file with the convolved fluxes for the use of the idl method.
- Parameters:
filter_name (str) – string to define the name of the filter. Following the convention LLNN, where L is the letter of the Telescope/Instrument and N is the number refering to the wavelength
instrument (str) – string to specify the name of the telescope_instrument
lambda_array (array) – Wavelength array in microns where the response of the filter is defined
response_array (array) – Response array that corresponds to instrumental response of the filter. The response array does not need to be normalised.
- Return type:
None
- add_square_filter(filter_name, instrument, filter_lambda, filter_width)¶
Adds square filter to the database given the central filter_lambda and the filter_width. It adds a txt file to the database with columns lambda_array and response_array. It also adds the fits file with the convolved fluxes for the use of the idl method.
- Parameters:
filter_name (str) – string to define the name of the filter. Following the convention LLNN, where L is the letter of the Telescope/Instrument and N is the number refering to the wavelength
instrument (str) – string to specify the name of the telescope_instrument
filter_lambda (float) – Wavelength array in microns where the response of the filter is defined
filter_width (float) – Response array that corresponds to instrumental response of the filter. The response array does not need to be normalised.
- Return type:
None
- chisq(flux_model, av)¶
chisq function as defined in Eq 4 of Z&T 2018. This equation is used to compute the chisq stored in the tables.
- chisq_to_minimize(av, flux_model)¶
chisq function to be used with scipy.optmize.minimize function. Note the slight change in sysntax from chisq() function above. In short it is exactly the same, but the parameters are called in a different order.
- get_average_model(models, number_of_models=5, chisq_cut=None, core_radius_cut=None, method=None, tablename=None)¶
Get the average model by means of calculating the geometric mean for all parameters except for av, theta_view, and theta_w_esc from which arithmetic mean is calculated
models must be the astropy table obtained from the SED results that has column names defined.
There are two methods for the average. The first just take the number of models specified (default 5) and the second takes into account a chisq cut and/or radius_cut. When set the chisq_cut, it will consider all models below that chisq_cut value. When set the core_radius_cut, it will consider all models with the model output core radius is smaller than the core_radius_cutvalue.
- Parameters:
models (table, astropy.table) – models to be averaged over. It should be an astropy table got from the get_model_info() function.
number_of_models (int) – Number of models to be used for the average for method 1.
chisq_cut (float) – Chisq cut value to consider all models below
core_radius_cut (float) – Core radius cut value to consider all models below
method (str, {'liu','moser'} optional) – Defines the method for backwards compatibilty with previous works. liu takes the best 5 models OR FEWER with the condition chi2<chis2min+5 and Rc<2Rap moser takes the best 10 models OR FEWER with the condition chis2<2*chi2min. The chis2 and core radius conditions have to be given as usual. The method only takes the best 5 (or 10) OR FEWER models satisfying the condition.
tablename (str, optional) – A path, or a Python file-like object. Note that fname is used verbatim, and there is no attempt to make the extension. Default is None. Note that one can choose the format of the figure by changing the extension.
- Returns:
table
- Return type:
astropy.table with all the model information based on the given keys
- get_extinction_laws_available()¶
print all the extinction laws available
- get_norm_extinction_law(extc_law='kmh')¶
Computes the normalised extinction law based on the extc_law file
- get_sed(mcore=10, sigma=0.1, mstar=0.5, theta_view=22.33, av=0, dist=1000, filter_array=None)¶
Retrive the SED for give physical parameters. Note that not all combinations of the parameters are present in the models database, an error will be raised, please try with other combination.
- Parameters:
mc (float) – mass of the core (Msun) for the SED to be retrieved. Allowed values are: [10.0,20.0,30.0,40.0,50.0,60.0,80.0,100.0,120.0,160.0,200.0,240.0,320.0,400.0,480.0]
sigma (float) – mass surface density of the clump (g/cm2) for the SED to be retrieved. Allowed values are: [0.1,0.316,1.0,3.16]
mstar (float) – mass of the star (Msun) for the SED to be retrieved. Allowed values are: [0.5,1.0,2.0,4.0,8.0,12.0,16.0,24.0,32.0,48.0,64.0,96.0,128.0,160.0]
theta_view (float) – viewing angle (deg) for the SED to be retrieved. Allowed values are: [12.84, 22.33, 28.96, 34.41, 39.19, 43.53, 47.55, 51.32, 54.90, 58.33, 61.64, 64.85, 67.98, 71.03, 74.04, 77.00, 79.92, 82.82, 85.70, 88.57]
Av (float) – visual extintion (mag) for the SED to be extincted. Allowed values are: >0
filter_array (str (or array of str)) – Name(s) of the filter(s) to convolve the SED. It convolves by the filter responde the model SED value. Default None.
- Returns:
lambda_model, flux_model_extincted ((array,array))
wavelength and flux for the selected physical parameters
If filter_array is not None, lambda_model, flux_model_extincted have the same length as filter_array.
- plot_filter(filter_name, figsize=(6, 4), legend=False, title=None, figname=None)¶
Plots a filter in the database. To see the available filter use SedFitter().print_default_filters
- Parameters:
filter_name (array of str) – array of strings with the name of the filter. E.g. if one filter is given [‘filter1’], if two filters (or more) are given [‘filter1’,’filter2’]
legend (bool, optional) – set the legend with the name of the filter. Default is False.
figname (str, optional) – A path, or a Python file-like object. Note that fname is used verbatim, and there is no attempt to make the extension. Default is None. Note that one can choose the format of the figure by changing the extension, e.g., figure.pdf would generate the figure in PDF format.
- Return type:
A figure with the normalised response for the given filter_name
- remove_filter(filter_name=None)¶
Removes a given filter from the database
- Parameters:
filter_name (str) – string with the name of the filter.
- sed_convolution(filter_wave_resp, lambda_array_model, sed_lambda_model, flux_model_extincted)¶
Performs the convolution of the extincted SED with the filter response
- sed_extinction(flux_model_log, av)¶
Extincts the flux_model given av
- sed_fit(dist, AV_min=0.0, AV_max=1000, method='minimize', progress=True, avopt=0)¶
Fits the SED observations to the Z&T18 set of models
- Parameters:
dist (float) – distance given in pc to the object.
AV_min (float) – Minimum visual extunction value to consider in the fit. The fit is perform in the range [AV_min,AV_max]. Default is 0.0
AV_max (float) – Maximum visual extunction value to consider in the fit. The fit is perform in the range [AV_min,AV_max]. Default is 1000.0
method ({'minimize', 'grid_search', 'idl'}) – Method to perform the fit. ‘minimize’ method uses the scipy.optimize minimize function over chisq_to_minimize() to find the best AV for each 8640 model. ‘grid_search’ performs a for loop over the AV_array = np.arange(AV_min,AV_max+1.0,1.0) and calculates the chi square as define in Z&T18 (See also De Buizer et al. 2017). ‘idl’ is a translation of the IDL version that also performs a grid search, it keeps the compatibility with the previous version using the same constants and fits files.
progress (bool) – progress bar either True or False.
avopt (int) – This is only relevant when choosing method = ‘idl’. It sets the visual extunction option. 0 is equally distributed AV point in the range of (0,AV_max) and 1 is define taking into account the mass surface density value. Default is 0 and should be the one used if wants to compare results with minimize or grid search.
- Returns:
FitterContainer
- Return type:
class with needed information to use the functions.
- table2latex(table, keys=['chisq', 'mcore', 'sigma', 'rcore', 'mstar', 'theta_view', 'av', 'massenv', 'theta_w_esc', 'mdotd', 'lbol_iso', 'lbol'], tablename=None)¶
Outputs the astropy table into a latex table given the keys properly formatted and with the correct units
- Parameters:
table (astropy.table) – Table to be transformed into latex
keys (str array) – Array with the keys that wants to be used. Default is keys=[‘chisq’,’mcore’,’sigma’,’rcore’,’mstar’,’theta_view’,’av’,’massenv’,’theta_w_esc’,’mdotd’,’lbol_iso’,’lbol’]
tablename (str) – A path, or a Python file-like object. Note that fname is used verbatim, and there is no attempt to make the extension. Default is None. Note that one can choose the format of the figure by changing the extension.
- Return type:
Latex formatted table in the path with the name tablename
- class sedcreator.__init__.SedFluxer(image)¶
A class used to measure flux
- get_data¶
Get the data to measure flux
- Type:
function
- get_flux(central_coords, aper_rad, inner_annu, outer_annu, mask=None)¶
Performs circular aperture photometry for a given image, specifying the central coordinates and the aperture radius. It returns the background subtracted flux and the flux including the background in units of Jy. Considering the header of the image, the function makes the units transformations.
Current units supported are (BUNIT): MJy/sr; MJY/SR; Jy/beam; mJy/beam; Jy/pix
If you are not sure about the units in your header, use get_raw_flux() function and make separately the transformation.
- Parameters:
central_coords (astropy.coordinates.SkyCoord) – Central coordinates where the apertures will be centred on the data image. This must be a SkyCoord statement.
aper_rad (float) – Aperture radius given in arcsec. The script will take care to convert it into pixels.
inner_annu (float) – Inner annulus radius given in arcsec. The script will take care to convert it into pixels. Usually this is set to 1*aper_rad.
outer_annu (float) – Outer annulus radius given in arcsec. The script will take care to convert it into pixels. Usually this is set to 2*aper_rad.
mask (2D array bool or int) – mask with the same shape of the image where the flux wants to be calculated. True values (or 1) will be ignore in the aperture_photometry.
- Returns:
array – with flux with background subtracted and flux including background in units of Jy
- Return type:
float,float
- get_optimal_aperture(central_coords, ap_inner=5.0, ap_outer=60.0, step_size=0.25, aper_increase=1.3, threshold=1.1, profile_plot=False)¶
Finds the optimal aperture based on the following method: EXPLAIN METHOD. In short: if one increases by 30% the aperture size, the flux does not increase by more than threshold %
- Parameters:
image (fits file) – Image for which we would like to calculate the “optimal” aperture radius
central_coords (astropy.coordinates.SkyCoord) – Central coordinates where the apertures will be centred on the data image. This must be a SkyCoord statement.
ap_inner (float) – Inner aperture radius given in arcsec. The script will take care to convert it into pixels. It the starting point to find the optimal aperture. Default is 5.0 arcsec.
ap_outer (float) – Outer aperture radius given in arcsec. The script will take care to convert it into pixels. It the ending point to find the optimal aperture. Default is 60.0 arcsec.
step_size (float) – Step size for sampling aperture radii. Default is 0.25arcsec
aper_increase (float) – It is the increase in aperture to meet the condition aper_increase*optimal radius <= flux at opt.rad. * threshold. Default is 1.30, i.e. 30% increase in aperture.
threshold (float) – It is the condition to be met such as, at the optimal aperture radius, 1.3*optimal radius <= flux at opt.rad. * threshold. Default is 1.08, i.e. 10% increase in flux.
profile_plot (bool) – Plots the flux profile versus the aperture radius size. Default is False
- Returns:
opt_rad – Optimal aperture radius defined by this method given in arcsec
- Return type:
float
- get_raw_flux(central_coords, aper_rad, inner_annu, outer_annu, mask=None)¶
Performs circular aperture photometry for a given image, specifying the central coordinates and the aperture radius. It returns the background subtracted flux and the flux including the background without any transformations of units. This functions is intended in case the user wants to apply a custom tranformation.
- Parameters:
central_coords (astropy.coordinates.SkyCoord) – Central coordinates where the apertures will be centred on the data image. This must be a SkyCoord statement.
aper_rad (float) – Aperture radius given in arcsec. The script will take care to convert it into pixels.
inner_annu (float) – Inner annulus radius given in arcsec. The script will take care to convert it into pixels. Usually this is set to 1*aper_rad.
outer_annu (float) – Outer annulus radius given in arcsec. The script will take care to convert it into pixels. Usually this is set to 2*aper_rad.
mask (2D array bool or int) – mask with the same shape of the image where the flux wants to be calculated. True values (or 1) will be ignore in the aperture_photometry.
- Returns:
array – with flux with background subtracted and flux including background without any units transformations
- Return type:
float,float